Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_FieldPattern.cpp
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 #include "Panzer_FieldPattern.hpp"
12 
13 namespace panzer {
14 
16 
18 {
19  int count = 0;
20  int dim = getDimension();
21 
22  // compute number of IDs
23  for(int i=0;i<dim+1;i++) {
24  for(int sc=0;sc<getSubcellCount(i);sc++)
25  count += getSubcellIndices(i,sc).size();
26  }
27 
28  return count;
29 }
30 
31 void FieldPattern::print(std::ostream & os) const
32 {
33  int dim = getDimension()+1;
34  os << "FieldPattern: " << dim << " Subcell types" << std::endl;
35  for(int i=0;i<dim;i++) {
36  int subcells = getSubcellCount(i);
37  os << "FieldPattern: " << subcells << " subcells of type " << i << std::endl;
38 
39  for(int j=0;j<subcells;j++) {
40  const std::vector<int> & indices = getSubcellIndices(i,j);
41  os << "FieldPattern: subcell " << j << " = [ ";
42  for(std::size_t k=0;k<indices.size();k++)
43  os << indices[k] << " ";
44  os << "]" << std::endl;
45  }
46  }
47 }
48 
50 {
51  bool equal = true;
52 
53  // test same dimension
54  std::size_t dim = getDimension();
55  equal &= (dim==(std::size_t) fp.getDimension());
56 
57  // check sub cells
58  for(std::size_t d=0;d<dim;d++)
59  equal &= getSubcellCount(d)==fp.getSubcellCount(d);
60 
61  return equal;
62 }
63 
65 {
66  bool consistent = true;
67 
68  std::size_t dim = getDimension();
69  for(std::size_t d=0;d<dim+1;d++) {
70  int numSC = getSubcellCount(d);
71  std::size_t sz = getSubcellIndices(d,0).size();
72  for(int i=1;i<numSC;i++) {
73  consistent &= (sz==getSubcellIndices(d,i).size());
74  }
75  }
76 
77  return consistent;
78 }
79 
80 bool FieldPattern::equals(const FieldPattern & fp) const
81 {
82  // same geometry is required
83  if(not this->sameGeometry(fp))
84  return false;
85 
86  // check to make sure subcell indices are equal
87  int dimension = this->getDimension();
88  for(int d=0;d<dimension+1;d++) {
89  for(int sc=0;sc<this->getSubcellCount(d);sc++) {
90  const std::vector<int> & myVector = this->getSubcellIndices(d,sc);
91  const std::vector<int> & argVector = fp.getSubcellIndices(d,sc);
92 
93  // check size of vectors
94  if(myVector.size()!=argVector.size())
95  return false;
96 
97  // check content of vectors
98  bool eq = std::equal(myVector.begin(),myVector.end(),argVector.begin());
99  if(not eq)
100  return false;
101  }
102  }
103 
104  return true;
105 }
106 
107 std::ostream & operator<<(std::ostream & os,const FieldPattern & fp)
108 {
109  fp.print(os);
110  return os;
111 }
112 
113 }
virtual int getDimension() const =0
virtual bool sameGeometry(const FieldPattern &fp) const
virtual bool consistentSubcells() const
virtual int numberIds() const
virtual ~FieldPattern()=0
Do nothing destructor.
virtual void print(std::ostream &os) const
virtual bool equals(const FieldPattern &fp) const
std::ostream & operator<<(std::ostream &os, const AssemblyEngineInArgs &in)
virtual const std::vector< int > & getSubcellIndices(int dim, int cellIndex) const =0
virtual int getSubcellCount(int dim) const =0