NeoPZ
TPZAgglomerateEl.h
Go to the documentation of this file.
1 
6 #ifndef AGGLOMERATEELEMHPP
7 #define AGGLOMERATEELEMHPP
8 
9 #include "pzreal.h"
10 #include "pzstack.h"
11 #include "pzeltype.h"
12 #include "TPZCompElDisc.h"
13 #include <iostream>
14 #include <map>
15 #include <set>
16 
17 #include <stdlib.h>
18 
19 class TPZGeoEl;
20 class TPZCompEl;
21 class TPZCompMesh;
22 class TPZGeoElSide;
24 struct TPZElementMatrix;
25 class TPZAgglomerateMesh;
26 
35 
36 private:
37 
45 
53 
60 
62  int fNFaces;
63 
66 
67 public:
68 
70  TPZAgglomerateElement(int nummat,int64_t &index,TPZCompMesh &aggcmesh,TPZCompMesh *finemesh);
71 
73 
76  void InitializeElement();
77 
79  void SetInnerRadius(REAL InnerRadius) override { fInnerRadius = InnerRadius;}
80 
83  REAL InnerRadius2() {return fInnerRadius;}
84 
87  REAL InnerRadius() override {
88  int64_t nsubel = this->NIndexes();
89  REAL value = 0.;
90  TPZCompElDisc * disc;
91  for (int64_t i = 0; i < nsubel; i++){
92  disc = dynamic_cast<TPZCompElDisc *>(this->SubElement(i) );
93 #ifdef PZDEBUG
94  if (!disc) {
95  PZError << "TPZAgglomerateElement::InnerRadius FineElement must be a TPZCompElDisc" << std::endl;
96  exit (-1);
97  }
98 #endif
99  value += disc->InnerRadius() * disc->VolumeOfEl();
100  }
101  value = value / this->VolumeOfEl();
102  return value;
103  }
104 
106  void SetNInterfaces(int nfaces) override {fNFaces = nfaces; }
107 
109  int NInterfaces() override {return fNFaces;}
110 
112  static void AddSubElementIndex(TPZCompMesh *aggcmesh,int64_t subel,int64_t destind);
113 
116 
118  MElementType Type() override {return EAgglomerate;}
119 
122 
124  virtual void AccumulateIntegrationRule(int degree, TPZStack<REAL> &point, TPZStack<REAL> &weight) override;
125 
127  virtual void AccumulateVertices(TPZStack<TPZGeoNode *> &nodes) override;
128 
130  void CenterPoint();
131 
133  virtual void CenterPoint(TPZVec<REAL> &center) override;
134 
136  REAL VolumeOfEl() override;
137 
140 
141  void CalcResidual(TPZElementMatrix &ef) override
142  {
143  std::cout << __PRETTY_FUNCTION__ << " is not implemented\n";
144  exit(-1);
145  }
146 
148  void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef) override;
149 
151  int64_t NIndexes() const { return fIndexes.NElements(); }
152 
155 
156  // /**os geométricos agrupados apontam para o computacional*/
157  /* void SetReference(){
158  int nindex = NIndexes(),i;
159  for(i=0;i<nindex;i++){
160  TPZCompEl *cel = SubElement(i);
161  int type = cel->Type();
162  //caso comp é aglomerado: chamada recursiva
163  if(type == EAgglomerate){//aglomerado
164  SetReference();
165  } else if(type == EDiscontinuous){//descontínuo
166  //o geométrico agrupado apontará para o atual computacional
167  cel->Reference()->SetReference(this);
168  }
169  }
170  }
171 
172  void SetReference2(int sub){
173 
174  TPZCompEl *cel = SubElement(sub);
175  int type = cel->Type();
176  //caso comp é aglomerado: chamada recursiva
177  if(type == EAgglomerate){//aglomerado
178  dynamic_cast<TPZAgglomerateElement *>(cel)->SetReference();
179  } else if(type == EDiscontinuous){//descontínuo
180  //o geométrico agrupado apontará para o atual computacional
181  cel->Reference()->SetReference(this);
182  }
183  } */
184 
187  virtual REAL LesserEdgeOfEl();
188 
190  TPZCompEl *SubElement(int64_t sub) const;
191 
192  REAL NormalizeConst() override;
193 
195  virtual int64_t CreateMidSideConnect() override;
196 
198  int Dimension() const override;
199 
201  virtual void Print(std::ostream & out = std::cout) const override;
202 
204  void ListOfDiscEl(TPZStack<TPZCompEl *> &elvec);
207 
209  int NSides();
210 
212  void CreateGraphicalElement(TPZGraphMesh &grmesh, int dimension) override;
213 
214  static void ListOfGroupings(TPZCompMesh *finemesh,TPZVec<int64_t> &accumlist,int nivel,int64_t &numaggl,int dim);
215 
216  // void FineSolution(TPZVec<REAL> &x,TPZCompElDisc *disc,TPZVec<REAL> &uh);
217 
218  void Print(TPZStack<int64_t> &listindex);
219 
220  void ProjectSolution(TPZFMatrix<STATE> &projectsol);
221 
222 
223  static TPZAgglomerateMesh *CreateAgglomerateMesh(TPZCompMesh *finemesh,TPZVec<int64_t> &accumlist,int64_t numaggl);
224 
225  static void ComputeNeighbours(TPZCompMesh *mesh, std::map<TPZCompElDisc *,std::set<TPZCompElDisc *> > &neighbours);
226 
228  public:
229 int ClassId() const override;
230 
231  /*@brief Save the element data to a stream */
232  void Write(TPZStream &buf, int withclassid) const override;
233 
235  void Read(TPZStream &buf, void *context) override;
236 
237 };
238 
239 #endif
void IndexesDiscSubEls(TPZStack< int64_t > &elvec)
Returns a vector of all indexes of the discontinuous elements in cluster.
Represents a graphical mesh used for post processing purposes. Post processing.
Definition: pzgraphmesh.h:34
void Read(TPZStream &buf, void *context) override
Read the element data from a stream.
virtual void Print(std::ostream &out=std::cout) const override
Prints the features of the element.
static TPZAgglomerateMesh * CreateAgglomerateMesh(TPZCompMesh *finemesh, TPZVec< int64_t > &accumlist, int64_t numaggl)
void SetInnerRadius(REAL InnerRadius) override
Sets the inner radius value.
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
void CenterPoint()
Computes the center of the mass to clustered elements.
void ProjectSolution(TPZFMatrix< STATE > &projectsol)
void SetNInterfaces(int nfaces) override
Sets element&#39;s number of interfaces.
void degree(int root, int adj_num, int adj_row[], int adj[], int mask[], int deg[], int *iccsze, int ls[], int node_num)
Definition: rcm.cpp:875
REAL InnerRadius() override
Returns the inner radius value.
Utility class which represents an element with its side. The Geometric approximation classes Geometry...
Definition: pzgeoelside.h:83
Defines enum MElementType and contains the implementation of MElementType_NNodes(...) functions.
TPZCompMesh * fMotherMesh
Mesh for the clusters which elements are part and from that the current mesh was obtained.
int Dimension() const override
It returns dimension from the elements.
TPZStack< int64_t > fIndexes
Indexes into the fine mesh of computational subelements in the clusters by the current.
virtual void AccumulateIntegrationRule(int degree, TPZStack< REAL > &point, TPZStack< REAL > &weight) override
Accumulates integration rule to deformed element.
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
virtual REAL LesserEdgeOfEl()
Computes a measure of the element.
REAL fInnerRadius
Stores the element&#39;s inner radius.
int ClassId() const override
Returns the unique identifier for reading/writing objects to streams.
virtual int64_t CreateMidSideConnect() override
It creates new conect that it associates the degrees of freedom of the element and returns its index...
virtual void AccumulateVertices(TPZStack< TPZGeoNode *> &nodes) override
Accumulate the vertices of the agglomerated elements.
static void AddSubElementIndex(TPZCompMesh *aggcmesh, int64_t subel, int64_t destind)
Insert the subelement index.
Contains declaration of TPZCompelDisc class which implements a computational element for discontinuou...
~TPZAgglomerateElement()
Destructor.
A simple stack.
Computes the contribution over an interface between two discontinuous elements. Computational Element...
void CalcResidual(TPZElementMatrix &ef) override
Only computes the element residual.
This class associates an element matrix with the coeficients of its contribution in the global stiffn...
Definition: pzelmat.h:30
static void ComputeNeighbours(TPZCompMesh *mesh, std::map< TPZCompElDisc *, std::set< TPZCompElDisc *> > &neighbours)
void CalcResidual(TPZFMatrix< REAL > &Rhs, TPZCompElDisc *el)
Computes the residual of the solution to father element from clustered subelements.
void Write(TPZStream &buf, int withclassid) const override
void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef) override
Assembles the differential equation to model over the element defined by clustered subelements...
void ListOfDiscEl(TPZStack< TPZCompEl *> &elvec)
Returns a vector of all discontinuous elements in cluster.
int fMaterialId
Material id of the agglomerated element.
Implements a mesh that contains agglomerated elements. Computational Mesh.
void CreateGraphicalElement(TPZGraphMesh &grmesh, int dimension) override
TPZCompEl * SubElement(int64_t sub) const
Returns the "sub" subelement.
MElementType
Define the element types.
Definition: pzeltype.h:52
REAL VolumeOfEl() override
Returns the volume of the geometric element referenced.
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
MElementType Type() override
Type of the element.
void InitializeElement()
Initialize the characteristics data of the clustered elements.
int64_t NIndexes() const
Returns the number of clustered subelements.
int NInterfaces() override
Retunrs the number of interfaces.
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
Contains the declaration of TPZFlopCounter class and TPZCounter struct.
This class implements a discontinuous element (for use with discontinuous Galerkin). Computational Element.
Definition: TPZCompElDisc.h:35
virtual REAL VolumeOfEl() override
Returns the volume of the geometric element associated.
REAL InnerRadius2()
Returns the inner radius value.
virtual REAL InnerRadius()
Returns the inner radius value.
int fNFaces
Stores the number of interfaces of the element.
TPZCompMesh * MotherMesh()
Returns father mesh.
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
int NSides()
Returns the number of sides. If all the volumes agglomerated have the same number, it returns this number, else it returns -1.
TPZGeoEl * CalculateReference()
Returns the geometric element to which this element references.
static void ListOfGroupings(TPZCompMesh *finemesh, TPZVec< int64_t > &accumlist, int nivel, int64_t &numaggl, int dim)
REAL NormalizeConst() override
Calculates the normalizing constant of the bases of the element.
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
Implements an agglomerated discontinuous element. Computational Element.