NeoPZ
pzgeoelrefless.h
Go to the documentation of this file.
1 
6 #ifndef PZGEOELREFLESS_H
7 #define PZGEOELREFLESS_H
8 
9 #include "pzgeoel.h"
10 #include "pzgeoelside.h"
11 #include "pzgeom_utility.h"
12 
13 class TPZGeoElSide;
14 class TPZCompMesh;
15 class TPZGeoMesh;
16 class TPZCompEl;
17 template<class T,int N>
18 class TPZStack;
19 
31 template <class TGeo>
32 class TPZGeoElRefLess : public TPZGeoEl {
33  // int fSubElement;
34 protected:
35  TGeo fGeo;
36  // int fNodeIndexes[TGeo::NNodes];
38 public:
39 
40 virtual int ClassId() const override;
41 
42  virtual ~TPZGeoElRefLess();
44 
46  TPZGeoElRefLess(const TPZGeoElRefLess &gel);
47 
48  virtual TPZGeoEl *Clone(TPZGeoMesh &dest) const override
49  {
50  return new TPZGeoElRefLess(dest,*this);
51  }
52 
54  TPZGeoElRefLess(TPZGeoMesh &DestMesh, const TPZGeoElRefLess &cp);
55 
56  //virtual void ParametricDomainNodeCoord(int node, TPZVec<REAL> &nodeCoord);
57 
69  TPZGeoElRefLess(TPZGeoMesh &DestMesh,
70  const TPZGeoElRefLess &cp,
71  std::map<int64_t,int64_t> & gl2lcNdMap,
72  std::map<int64_t,int64_t> & gl2lcElMap);
73 
74  TPZGeoEl *ClonePatchEl(TPZGeoMesh &destmesh,std::map<int64_t,int64_t> & gl2lcNdMap,
75  std::map<int64_t,int64_t> & gl2lcElMap) const override
76  {
77  return new TPZGeoElRefLess(destmesh,*this,gl2lcNdMap, gl2lcElMap);
78  }
79 
80  TPZGeoElRefLess(int64_t id,TPZVec<int64_t> &nodeindexes,int matind,TPZGeoMesh &mesh);
81 
82  TPZGeoElRefLess(TGeo &Geo, int matind,TPZGeoMesh &mesh);
83 
84  TPZGeoElRefLess(TPZVec<int64_t> &nodeindices,int matind,TPZGeoMesh &mesh);
85 
86  TPZGeoElRefLess(TPZVec<int64_t> &nodeindices,int matind,TPZGeoMesh &mesh,int64_t &index);
87 
88  virtual void Read(TPZStream &str, void *context) override;
89 
90  virtual void Write(TPZStream &str, int withclassid) const override;
91 
92  virtual void Initialize() override
93  {
94  fGeo.Initialize(this);
95  }
96 
97 
98  static int main_refless();
99 
101  virtual void Divide(TPZVec < TPZGeoEl * > & pv) override {
102  DebugStop();
103  }
104 
106  virtual int HasSubElement() const override {return 0;}//fSubEl[0]!=0;}
107 
109  virtual TPZGeoElSide Neighbour(int side) override {
110 #ifdef PZDEBUG
111  if (fNeighbours[side] < 0 || fNeighbours[side] >= this->Mesh()->NElements()) {
112  DebugStop();
113  }
114 #endif
115  return TPZGeoElSide(fNeighbours[side],this->Mesh());
116  }
117 
119  virtual int64_t NeighbourIndex(int side) const override{
120 #ifdef PZDEBUG
121  if (fNeighbours[side] < 0 || fNeighbours[side] >= this->Mesh()->NElements()) {
122  DebugStop();
123  }
124 #endif
125  return this->fNeighbours[side].ElementIndex();
126  }
127 
128  virtual int64_t NodeIndex(int node) const override;
129 
130  void CornerCoordinates(TPZFMatrix<REAL> &coord) const;
131  //HDiv
132 
133 // virtual void Directions(int side, TPZVec<REAL> &pt, TPZFMatrix<REAL> &directions, TPZVec<int> &vectorsides) override;
134 
135  virtual void HDivDirectionsMaster(TPZFMatrix<REAL> &directions) override;
136 
137  virtual void HDivDirections(TPZVec<REAL> &pt, TPZFMatrix<REAL> &directions, int ConstrainedFace = -1) override;
138 
139 
140 #ifdef _AUTODIFF
141  virtual void HDivDirections(TPZVec<REAL> &pt, TPZFMatrix<Fad<REAL> > &directions, int ConstrainedFace = -1) override;
142 #endif
143 
144  //virtual void VecHdiv(TPZFMatrix<REAL> &normalvec ,TPZVec<int> &sidevector) override;
145 
146 
148  virtual void HDivPermutation(int side, TPZVec<int> &permutegather) override;
149 
150 
152  virtual void SetNeighbour(int side,const TPZGeoElSide &neighbour) override {
153  fNeighbours[side]=neighbour;
154  }
155 
156  virtual void Print(std::ostream &out) override
157  {
158  TPZGeoEl::Print(out);
159  out << "fGeo:\n";
160  fGeo.Print(out);
161  }
163  virtual void PrintTopologicalInfo(std::ostream &out) override
164  {
165  out << "Geo Element - fId " << fId << "\t Type " << fGeo.TypeName() << "\t";
166  fGeo.Print(out);
167  int i;
168  out << std::endl << "\t";
169  for (i = 0;i < NNodes();i++) out << NodePtr(i)->Id() << " ";
170  }
171 
172  virtual int64_t SideNodeIndex(int side,int node) const override;
173 
174  virtual int SideNodeLocIndex(int side,int node) const override;
175 
177  virtual void SetSideDefined(int side) override { fNeighbours[side] = TPZGeoElSide(this,side); }
178 
179  virtual void SetSubElement(int id, TPZGeoEl *el) override;
180 
181 
187  void GetPermutation(const int& i, TPZVec<int> &permutation) const override;
192  virtual TPZIntPoints * CreateSideIntegrationRule(int side, int order) override;
193 
195  virtual MElementType Type() const override {
196  return TGeo::Type();
197  }
198 
200  virtual MElementType Type(int side) const override {
201  return TGeo::Type(side);
202  }
203 
205  virtual std::string TypeName() const override
206  {
207  return fGeo.TypeName();
208  }
209 
211  virtual int NNodes() const override;
212 
214  virtual int NCornerNodes() const override;
215 
217  virtual int NSides() const override;
218 
223  virtual int NPermutations() const override;
224 
226  virtual int SideNodeLocId(int side, int node) const;
227 
229  virtual REAL RefElVolume() override;
230 
232  virtual int NSideNodes(int side) const override;
233 
235  virtual void MidSideNodeIndex(int side,int64_t &index) const override;
236 
240  virtual int SideIsUndefined(int side) override;
241 
243  virtual int NSubElements() const override;
244 
246  virtual int NSideSubElements(int side) const override;
247 
252  virtual TPZGeoEl *CreateBCGeoEl(int side, int bc) override;
253 
258  virtual TPZGeoEl *CreateBCGeoBlendEl(int side, int bc);
259 
261  virtual TPZGeoEl *CreateGeoElement(MElementType type,
262  TPZVec<int64_t>& nodeindexes,
263  int matid,
264  int64_t& index) override;
265 
267  virtual void SetNodeIndex(int i,int64_t nodeindex) override;
268 
273  virtual TPZTransform<> SideToSideTransform(int sidefrom,int sideto) override;
274 
276  virtual TPZGeoEl *SubElement(int is) const override;
277 
279  virtual int SideDimension(int side) const override;
280 
282  virtual int Dimension() const override;
283 
284  virtual TPZGeoElSide HigherDimensionSides(int side,int targetdimension);
285 
286  virtual void AllHigherDimensionSides(int side,int targetdimension,TPZStack<TPZGeoElSide> &elsides) override;
287 
288  virtual void LowerDimensionSides(int side,TPZStack<int> &smallsides) const override;
289 
292  virtual void BuildTransform(int side, TPZGeoEl *father,TPZTransform<> &t);
293 
294  virtual TPZTransform<> BuildTransform2(int side, TPZGeoEl *father,TPZTransform<> &t) override;
295 
297  virtual void X(TPZVec<REAL> &coordinate,TPZVec<REAL> &result) const override;
298 
300  virtual void GradX(TPZVec<REAL> &coordinate, TPZFMatrix<REAL> &gradx) const override;
301 
302 #ifdef _AUTODIFF
303 
304  virtual void X(TPZVec<Fad<REAL> > &coordinate,TPZVec<Fad<REAL> > &result) const override;
305 
307  virtual void GradX(TPZVec<Fad<REAL> > &coordinate, TPZFMatrix<Fad<REAL> > &gradx) const override;
308 #endif
309 
310  virtual bool IsLinearMapping( int side) const override;
311  virtual bool IsGeoBlendEl() const override;
312 
321  bool ResetBlendConnectivity(const int64_t &side, const int64_t &index) override;
322  TGeo &Geom() { return fGeo; }
323 
324  virtual TPZTransform<> GetTransform(int side,int son) override;
325 
327  virtual void CenterPoint(int side, TPZVec<REAL> &masscent) const override;
328 
329  virtual TPZGeoElSide Father2(int side) const override;
330 
331  virtual int FatherSide(int side, int son) override {
332  return side;
333  }
334 
335  virtual void GetSubElements2(int side, TPZStack<TPZGeoElSide> &subel) const override;
336 
337  virtual void ResetSubElements() override {
338  DebugStop();
339  }
340 
341  virtual void SetNeighbourInfo(int side, TPZGeoElSide &neigh, TPZTransform<> &trans) override
342  {
343  Geom().SetNeighbourInfo(side,neigh,trans);
344  }
345 
347  virtual void RandomPoint(TPZVec<REAL> &pt) override
348  {
349  Geom().RandomPoint(pt);
350  }
351 
353  virtual bool IsInParametricDomain(const TPZVec<REAL> &pt, REAL tol = 1e-6) override;
354 
362  virtual int ProjectInParametricDomain(TPZVec<REAL> &qsi, TPZVec<REAL> &qsiInDomain) override;
363 
369  virtual int ProjectBissectionInParametricDomain(TPZVec<REAL> &qsi, TPZVec<REAL> &qsiInDomain) override;
370 };
371 
372 template<class TGeo>
373 inline
375  const bool result = fGeo.IsInParametricDomain(pt,tol);
376  return result;
377 }
378 
379 template<class TGeo>
380 inline
382  const int side = ::ProjectInParametricDomain<TGeo>(pt, ptInDomain);
383  return side;
384 }
385 
386 template<class TGeo>
387 inline
389  const int side = ::ProjectBissectionInParametricDomain<TGeo>(pt, ptInDomain);
390  return side;
391 }
392 
393 template<class TGeo>
395  return Hash("TPZGeoElRefLess") ^ TPZGeoEl::ClassId() << 1 ^ TGeo().ClassId() << 2;
396 }
397 
398 //template<class TGeo>
399 //inline
400 //void TPZGeoElRefLess<TGeo>::ParametricDomainNodeCoord(int node, TPZVec<REAL> &nodeCoord)
401 //{
402 // fGeo.ParametricDomainNodeCoord(node, nodeCoord);
403 //}
404 
405 #include "pzgeoelrefless.h.h"
406 
407 #endif
virtual ~TPZGeoElRefLess()
virtual void SetNeighbour(int side, const TPZGeoElSide &neighbour) override
Fill in the data structure for the neighbouring information.
virtual TPZGeoEl * CreateGeoElement(MElementType type, TPZVec< int64_t > &nodeindexes, int matid, int64_t &index) override
Creates a geometric element according to the type of the father element.
virtual void SetNeighbourInfo(int side, TPZGeoElSide &neigh, TPZTransform<> &trans) override
virtual int NSubElements() const override
Returns the number of subelements of the element independent of the fact hether the element has alrea...
virtual int64_t NodeIndex(int node) const override
Returns the index of the ith node the index is the location of the node in the nodevector of the mesh...
virtual void Read(TPZStream &str, void *context) override
read objects from the stream
bool ResetBlendConnectivity(const int64_t &side, const int64_t &index) override
TPZGeoNode * NodePtr(int i) const
Returns a pointer to the ith node of the element.
Definition: pzgeoel.cpp:2566
Contains declaration of TPZGeoElSide class which represents an element and its side, and TPZGeoElSideIndex class which represents an TPZGeoElSide index.
virtual TPZTransform SideToSideTransform(int sidefrom, int sideto) override
compute the transformation between the master element space of one side of an element to the master e...
virtual void Initialize() override
virtual void Print(std::ostream &out) override
Print all relevant data of the element to cout.
virtual int NSideNodes(int side) const override
Returns the number of nodes for a particular side.
clarg::argBool bc("-bc", "binary checkpoints", false)
virtual void RandomPoint(TPZVec< REAL > &pt) override
Generates a random point in the master domain.
virtual void PrintTopologicalInfo(std::ostream &out) override
Prints topological information of: TGeo (TPZGeoCube, TPZGeoPrism, TPZGeoQuad, ...) ...
virtual void AllHigherDimensionSides(int side, int targetdimension, TPZStack< TPZGeoElSide > &elsides) override
TPZGeoElSideIndex fNeighbours[TGeo::NSides]
virtual int HasSubElement() const override
Returns 1 if the element has subelements along side.
virtual int ProjectInParametricDomain(TPZVec< REAL > &qsi, TPZVec< REAL > &qsiInDomain) override
Ortogonal projection from given qsi to a qsiInDomain (all in the element parametric domain) ...
virtual TPZTransform GetTransform(int side, int son) override
Utility class which represents an element with its side. The Geometric approximation classes Geometry...
Definition: pzgeoelside.h:83
Definition: fad.h:54
Implements the mapping between the master element and deformed element. Geometry. ...
TPZGeoMesh * Mesh() const
Returns the mesh to which the element belongs.
Definition: pzgeoel.h:220
virtual void LowerDimensionSides(int side, TPZStack< int > &smallsides) const override
virtual REAL RefElVolume() override
Volume of the master element.
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
int64_t ElementIndex() const
Definition: pzgeoelside.h:391
Utility class which represents an element index with its side. Geometry.
Definition: pzgeoelside.h:33
virtual void HDivPermutation(int side, TPZVec< int > &permutegather) override
Compute the permutation for an HDiv side.
virtual int FatherSide(int side, int son) override
virtual TPZTransform BuildTransform2(int side, TPZGeoEl *father, TPZTransform<> &t) override
int64_t fId
Traditional element number or element id.
Definition: pzgeoel.h:51
virtual TPZGeoEl * Clone(TPZGeoMesh &dest) const override
int ClassId() const override
Define the class id associated with the class.
Definition: pzgeoel.cpp:2562
virtual void ResetSubElements() override
Reset all subelements to NULL.
virtual void CenterPoint(int side, TPZVec< REAL > &masscent) const override
It returns the coordinates of the center of the side of the element.
static const double tol
Definition: pzgeoprism.cpp:23
virtual int NSides() const override
Returns the number of connectivities of the element.
virtual void SetSideDefined(int side) override
Flags the side as defined, this means no neighbouring element was found.
virtual int NSideSubElements(int side) const override
Returns the number of subelements of the same dimension of the element at the side.
virtual void GetSubElements2(int side, TPZStack< TPZGeoElSide > &subel) const override
This method will return a partition of the side of the current element as the union of sub elements/...
virtual bool IsGeoBlendEl() const override
TPZGeoElRefLess()
Contains the implementation of the TPZGeoElRefLess methods.
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
virtual int Dimension() const override
Returns the dimension of the element.
virtual int SideNodeLocIndex(int side, int node) const override
Returns the local index of a node on a side.
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
virtual int64_t SideNodeIndex(int side, int node) const override
Returns the index of the nodenum node of side.
virtual void X(TPZVec< REAL > &coordinate, TPZVec< REAL > &result) const override
Returns the coordinate in real space of the point coordinate in the master element space...
virtual int SideDimension(int side) const override
Return the dimension of side.
virtual void SetSubElement(int id, TPZGeoEl *el) override
Sets the subelement of index i.
virtual TPZGeoElSide Neighbour(int side) override
Returns a pointer to the neighbour and the neighbourside along side of the current element...
void CornerCoordinates(TPZFMatrix< REAL > &coord) const
Gets the corner node coordinates in coord.
virtual int NNodes() const override
Returns the number of nodes of the element.
virtual MElementType Type() const override
Returns the type of the element acording to the definition in pzeltype.h.
virtual int SideNodeLocId(int side, int node) const
Returns the local node number of the node "node" along side "side".
virtual int NCornerNodes() const override
Returns the number of corner nodes of the element.
virtual TPZIntPoints * CreateSideIntegrationRule(int side, int order) override
Creates an integration rule for the topology of the corresponding side and able to integrate a polyno...
void GetPermutation(const int &i, TPZVec< int > &permutation) const override
virtual void GradX(TPZVec< REAL > &coordinate, TPZFMatrix< REAL > &gradx) const override
Return the gradient of the transformation at the point.
virtual int ClassId() const override
Define the class id associated with the class.
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
TPZGeoEl * ClonePatchEl(TPZGeoMesh &destmesh, std::map< int64_t, int64_t > &gl2lcNdMap, std::map< int64_t, int64_t > &gl2lcElMap) const override
Creates a clone of this element into a new patch mesh.
virtual bool IsInParametricDomain(const TPZVec< REAL > &pt, REAL tol=1e-6) override
Verifies if the parametric point pt is in the element parametric domain.
virtual bool IsLinearMapping() const
Definition: pzgeoel.h:268
virtual int64_t NeighbourIndex(int side) const override
Returns the neighbour index for a given side.
virtual void HDivDirections(TPZVec< REAL > &pt, TPZFMatrix< REAL > &directions, int ConstrainedFace=-1) override
virtual void HDivDirectionsMaster(TPZFMatrix< REAL > &directions) override
virtual int ProjectBissectionInParametricDomain(TPZVec< REAL > &qsi, TPZVec< REAL > &qsiInDomain) override
Projection from given qsi to a qsiInDomain (in the element boundary) using bissection method from giv...
static int main_refless()
virtual void BuildTransform(int side, TPZGeoEl *father, TPZTransform<> &t)
Accumulates the transformation of the jacobian which maps the current master element space into the s...
virtual MElementType Type(int side) const override
Returns the type of the element acording to the definition in pzeltype.h.
virtual TPZGeoEl * CreateBCGeoBlendEl(int side, int bc)
Method which creates a blend geometrical boundary condition element based on the current geometric el...
This class implements a geometric mesh for the pz environment. Geometry.
Definition: pzgmesh.h:48
MElementType
Define the element types.
Definition: pzeltype.h:52
This class implements a stack object. Utility.
Definition: pzcheckmesh.h:14
virtual TPZGeoElSide HigherDimensionSides(int side, int targetdimension)
virtual void Divide(TPZVec< TPZGeoEl * > &pv) override
Divides the element and puts the resulting elements in the vector.
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
virtual TPZGeoEl * SubElement(int is) const override
Returns a pointer to the subelement is.
int Id() const
Returns the identity of the current node.
Definition: pzgnode.h:68
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
Implements an affine transformation between points in parameter space. Topology Utility.
Definition: pzmganalysis.h:14
virtual void Print(std::ostream &out=std::cout)
Print all relevant data of the element to cout.
Definition: pzgeoel.cpp:238
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
virtual TPZGeoEl * CreateBCGeoEl(int side, int bc) override
Method which creates a computational boundary condition element based on the current geometric elemen...
virtual void MidSideNodeIndex(int side, int64_t &index) const override
Returns the midside node index along a side of the element.
virtual void Write(TPZStream &str, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
virtual int SideIsUndefined(int side) override
Returns 1 if the side has not been defined by buildconnectivity.
virtual std::string TypeName() const override
Returns the type of the element as a string.
virtual int NPermutations() const override
virtual void SetNodeIndex(int i, int64_t nodeindex) override
Initializes the node i of the element.
virtual TPZGeoElSide Father2(int side) const override
Returns the father/side of the father which contains the side of the sub element. ...