Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ElemFieldPattern.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef __Panzer_ElemFieldPattern_hpp__
12#define __Panzer_ElemFieldPattern_hpp__
13
14#include "PanzerDofMgr_config.hpp"
15
16#include <ostream>
17#include <vector>
18
20
21#include "Shards_CellTopology.hpp"
22
23namespace panzer {
24
29public:
30
32
33 ElemFieldPattern(const shards::CellTopology & ct);
34
36 virtual ~ElemFieldPattern() {}
37
38 Teuchos::RCP<panzer::FieldPattern> clone() const
39 {return Teuchos::rcp(new ElemFieldPattern(*this));}
40
42 void setCellTopology(const shards::CellTopology & ct);
43
53 virtual int getSubcellCount(int dim) const;
54
66 virtual const std::vector<int> & getSubcellIndices(int dim,int cellIndex) const;
67
78 virtual void getSubcellClosureIndices(int dim,int cellIndex,std::vector<int> & indices) const;
79
84 virtual int getDimension() const;
85
88 virtual shards::CellTopology getCellTopology() const
89 { return cellTopo_; }
90
91public:
92 shards::CellTopology cellTopo_;
93 std::vector<std::vector<int> > ElemIndices_;
94 std::vector<int> empty_;
95};
96
97}
98
99#endif
virtual ~ElemFieldPattern()
Do nothing destructor.
std::vector< std::vector< int > > ElemIndices_
Teuchos::RCP< panzer::FieldPattern > clone() const
virtual const std::vector< int > & getSubcellIndices(int dim, int cellIndex) const
virtual void getSubcellClosureIndices(int dim, int cellIndex, std::vector< int > &indices) const
virtual shards::CellTopology getCellTopology() const
void setCellTopology(const shards::CellTopology &ct)
Set the cell topology for this field pattern.
virtual int getSubcellCount(int dim) const