Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_FieldAggPattern.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_FieldAggPattern_hpp__
12 #define __Panzer_FieldAggPattern_hpp__
13 
14 #include "Panzer_FieldPattern.hpp"
15 #include "Panzer_FieldType.hpp"
16 
17 #include <map>
18 #include <vector>
19 #include <tuple>
20 
21 #include "Teuchos_RCP.hpp"
22 #include "Teuchos_Tuple.hpp"
23 #include "Phalanx_KokkosDeviceTypes.hpp"
24 
25 namespace panzer {
26 
34 class FieldAggPattern : public FieldPattern {
35 public:
36  // The FEI format for prescribing fields seems quite general,
37  // so it is my belief that an abstract can be written to it.
38 
40 
43  FieldAggPattern(std::vector<std::tuple<int,panzer::FieldType,Teuchos::RCP<const FieldPattern> > > & patterns,
44  const Teuchos::RCP<const FieldPattern> & geomAggPattern=Teuchos::null);
45 
46  virtual ~FieldAggPattern() {}
47 
51 
56  virtual void buildPattern(const std::vector<std::tuple<int,panzer::FieldType,Teuchos::RCP<const FieldPattern> > > & patterns,
57  const Teuchos::RCP<const FieldPattern> & geomAggPattern=Teuchos::null);
58 
59  // functions from the abstract FieldPattern
60  virtual int getDimension() const;
61  virtual int getSubcellCount(int dimension) const;
62  virtual const std::vector<int> & getSubcellIndices(int dimension, int subcell) const;
63  virtual void getSubcellClosureIndices(int, int, std::vector<int> &) const;
64  virtual shards::CellTopology getCellTopology() const;
65 
67  virtual void print(std::ostream & os) const;
68 
75  virtual Teuchos::RCP<const FieldPattern> getFieldPattern(int fieldId) const;
76 
83  virtual FieldType getFieldType(int fieldId) const;
84 
86 
87 
89  const std::vector<int> & numFieldsPerId() const { return numFields_; }
90 
94  const std::vector<int> & fieldIds() const { return fieldIds_; }
95 
97 
99 
100 
118  const std::vector<int> & localOffsets(int fieldId) const;
119 
137  const PHX::View<const int*> localOffsetsKokkos(int fieldId) const;
138 
161  const std::pair<std::vector<int>,std::vector<int> > &
162  localOffsets_closure(int fieldId,int subcellDim,int subcellId) const;
163 
165 
166 protected:
167  typedef Teuchos::RCP<const FieldPattern> FPPtr; // for convenience
168 
173 
179  void buildFieldIdsVector();
180 
191  void mergeFieldPatterns(const FieldType& fieldType);
192 
196  void addAllPatternSubcellIndices(int dim,int sc,std::vector<int> & indices);
197 
200  void buildFieldPatternData();
201 
204  void localOffsets_build(int fieldId,std::vector<int> & offsets) const;
205 
206  std::size_t dimension_;
208 
211  std::vector<std::vector<std::vector<int> > > patternData_;
212 
213  // FEI patterns
214  std::vector<int> numFields_; // number of fields at each geometric ID
215  std::vector<int> fieldIds_; // field IDs at each geometric ID
216 
217  // storage for patterns
218  std::vector<std::tuple<int,panzer::FieldType,Teuchos::RCP<const FieldPattern> > > patterns_;
219  std::map<int,int> fieldIdToPatternIdx_;
220 
222  mutable std::map<int, std::vector<int> > fieldOffsets_;
223 
225  mutable std::map<int, PHX::View<int*> > fieldOffsetsKokkos_;
226 
227  struct LessThan
228  { bool operator()(const Teuchos::Tuple<int,3> & a,const Teuchos::Tuple<int,3> & b) const; };
229  mutable std::map<Teuchos::Tuple<int,3>, std::pair<std::vector<int>,std::vector<int> >,LessThan>
231 };
232 
233 }
234 
235 #endif
const std::pair< std::vector< int >, std::vector< int > > & localOffsets_closure(int fieldId, int subcellDim, int subcellId) const
std::map< int, std::vector< int > > fieldOffsets_
Stores the Field offsets for the fieldId key. Note that the key is the fieldId, not the index into th...
std::map< Teuchos::Tuple< int, 3 >, std::pair< std::vector< int >, std::vector< int > >, LessThan > fieldSubcellOffsets_closure_
const std::vector< int > & numFieldsPerId() const
Lenght of vector is number of Ids, value is how many ids per field.
virtual Teuchos::RCP< const FieldPattern > getFieldPattern(int fieldId) const
virtual Teuchos::RCP< const FieldPattern > getGeometricAggFieldPattern() const
std::map< int, int > fieldIdToPatternIdx_
virtual int getSubcellCount(int dimension) const
std::vector< std::vector< std::vector< int > > > patternData_
FieldType
The type of discretization to use for a field pattern.
virtual FieldType getFieldType(int fieldId) const
virtual int getDimension() const
Teuchos::RCP< const FieldPattern > geomAggPattern_
virtual void getSubcellClosureIndices(int, int, std::vector< int > &) const
virtual void buildPattern(const std::vector< std::tuple< int, panzer::FieldType, Teuchos::RCP< const FieldPattern > > > &patterns, const Teuchos::RCP< const FieldPattern > &geomAggPattern=Teuchos::null)
PHX::View< const int * > offsets
void addAllPatternSubcellIndices(int dim, int sc, std::vector< int > &indices)
void localOffsets_build(int fieldId, std::vector< int > &offsets) const
std::map< int, PHX::View< int * > > fieldOffsetsKokkos_
Stores the Field offsets for the fieldId key. Note that the key is the fieldId, not the index into th...
const std::vector< int > & fieldIds() const
virtual void print(std::ostream &os) const
Print this pattern.
std::vector< std::tuple< int, panzer::FieldType, Teuchos::RCP< const FieldPattern > > > patterns_
bool operator()(const Teuchos::Tuple< int, 3 > &a, const Teuchos::Tuple< int, 3 > &b) const
const PHX::View< const int * > localOffsetsKokkos(int fieldId) const
virtual const std::vector< int > & getSubcellIndices(int dimension, int subcell) const
Teuchos::RCP< const FieldPattern > FPPtr
virtual shards::CellTopology getCellTopology() const
void mergeFieldPatterns(const FieldType &fieldType)
const std::vector< int > & localOffsets(int fieldId) const