NeoPZ
tpzpyramid.h
Go to the documentation of this file.
1 
6 #ifndef PZTOPOLOGYTPZPYRAMID_H
7 #define PZTOPOLOGYTPZPYRAMID_H
8 
9 #include "pzfmatrix.h"
10 #include "pzstack.h"
11 #include "pztrnsform.h"
12 #include "pzquad.h"
13 #include "pzeltype.h"
14 #include "pzaxestools.h"
15 #include "TPZTopologyUtils.h"
16 
17 class TPZIntPoints;
18 class TPZIntPyram3D;
20 
21 class TPZCompEl;
22 class TPZGeoEl;
23 class TPZCompMesh;
24 
26 namespace pztopology {
27 
35  class TPZPyramid : public TPZSavable{
36  public:
37  friend void pztopology::GetPermutation<TPZPyramid>(const int permute, TPZVec<int> &permutation);
39  enum {NSides = 19, NCornerNodes = 5, Dimension = 3, NFaces = 5, NPermutations = 8};
40 
41  int ClassId() const override;
42  void Read(TPZStream &buf, void *context) override;
43  void Write(TPZStream &buf, int withclassid) const override;
44 
47  }
48 
50  virtual ~TPZPyramid() {
51  }
52 
57  static int SideDimension(int side);
58 
60  static void LowerDimensionSides(int side,TPZStack<int> &smallsides);
62  static void LowerDimensionSides(int side,TPZStack<int> &smallsides, int DimTarget);
63 
69  static void HigherDimensionSides(int side, TPZStack<int> &high);
70 
72  static int NSideNodes(int side);
74  static int SideNodeLocId(int side, int node);
75 
77  static int NumSides();
79  static int NumSides(int dimension);
80 
82  static int NContainedSides(int side);
84  static int ContainedSideLocId(int side, int c);
85 
86 
88  static void Shape(TPZVec<REAL> &loc,TPZFMatrix<REAL> &phi,TPZFMatrix<REAL> &dphi){
89  TShape(loc, phi, dphi);
90  }
92  template<class T>
93  static void TShape(const TPZVec<T> &loc,TPZFMatrix<T> &phi,TPZFMatrix<T> &dphi);
102  template<class T>
103  static void BlendFactorForSide(const int &side, const TPZVec<T> &xi, T &blendFactor,
104  TPZVec<T> &corrFactorDxi);
111  static void CenterPoint(int side, TPZVec<REAL> &center);
112 
114  static bool IsInParametricDomain(const TPZVec<REAL> &pt, REAL tol = pztopology::gTolerance);
115 
116  #ifdef _AUTODIFF
117 
118  static bool IsInParametricDomain(const TPZVec<Fad<REAL>> &pt, REAL tol = pztopology::gTolerance){
119  TPZVec<REAL> xi(pt.size());
120  for(int i = 0; i < pt.size(); i++) xi[i]= pt[i].val();
121  return IsInParametricDomain(xi,tol);
122  }
123  #endif
124 
125  static void RandomPoint(TPZVec<REAL> &pt);
126 
134  template<class T>
135  static bool CheckProjectionForSingularity(const int &side, const TPZVec<T> &xiInterior);
136 
137  template<class T>
138  static void MapToSide(int side, TPZVec<T> &InternalPar, TPZVec<T> &SidePar, TPZFMatrix<T> &JacToSide);
139 
140  static void ParametricDomainNodeCoord(int node, TPZVec<REAL> &nodeCoord);
141 
148  static MElementType Type();
149 
151  static MElementType Type(int side);
152 
164  static TPZTransform<> SideToSideTransform(int sidefrom, int sideto);
165 
171  static TPZTransform<> TransformSideToElement(int side);
177  static TPZTransform<> TransformElementToSide(int side);
178 
185  {
186  DebugStop();
187  return -1;
188  }
189 
196  static int GetTransformId(int side, TPZVec<int64_t> &id);
197 
208  static TPZIntPoints * CreateSideIntegrationRule(int side, int order);
209 
214 
224  static void GetSideHDivPermutation(int transformationid, TPZVec<int> &permgather);
225 
227  static constexpr REAL RefElVolume() {return (4.L/3.L); }
228 
229  /* Given side and gradx the method returns directions needed for Hdiv space */
230  static void ComputeDirections(int side, TPZFMatrix<REAL> &gradx, TPZFMatrix<REAL> &directions, TPZVec<int> &sidevectors);
232  static void GetSideHDivDirections(TPZVec<int> &sides, TPZVec<int> &dir, TPZVec<int> &bilinearounao, TPZVec<int> &sidevectors);
233 
235  template <class TVar>
236  static void ComputeHDivDirections(TPZFMatrix<TVar> &gradx, TPZFMatrix<TVar> &directions);
237 
239  template <class TVar>
240  static void AdjustTopDirections(int ConstrainedFace,TPZFMatrix<TVar> &gradx, TPZFMatrix<TVar> &directions);
241 
242 
246  static int NBilinearSides();
247 
254  static void CornerShape(const TPZVec<REAL> &pt, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &dphi);
255 
256 
257 
259  static int FaceNodes[5][4];
260 
262  static int SideNodes[8][2];
263 
264  protected:
271  static int ShapeFaceId[5][4];
272 
274  static int fPermutations[8][19];
275 
278  };
279 
280 }
281 
282 #endif
static int FaceNodes[5][4]
Nodes over quadrilateral sides (2d - faces).
Definition: tpzpyramid.h:259
Defines the topology of a Pyramid element. Topology Sides 0 to 4 are vertices, sides 5 to 12 are line...
Definition: tpzpyramid.h:35
TPZPyramid()
Default constructor.
Definition: tpzpyramid.h:46
static void GetSideHDivDirections(TPZVec< int > &sides, TPZVec< int > &dir, TPZVec< int > &bilinearounao)
static int bilinearounao[81]
Definition: tpzcube.cpp:405
static int fPermutations[8][19]
Valid permutations between nodes.
Definition: tpzpyramid.h:274
static TPZTransform SideToSideTransform(int sidefrom, int sideto)
Returns the transformation which takes a point from the side sidefrom to the side sideto...
Definition: tpzpyramid.cpp:528
Contains the TPZInt1d, TPZIntTriang, TPZIntQuad, TPZIntCube3D, TPZIntTetra3D, TPZIntPyram3D and TPZIn...
static int ShapeFaceId[5][4]
Ids of the shape face.
Definition: tpzpyramid.h:271
static int ContainedSideLocId(int side, int c)
Returns the local connect number of the connect "c" along side "side".
Definition: tpzpyramid.cpp:923
static bool IsInParametricDomain(const TPZVec< REAL > &pt, REAL tol=pztopology::gTolerance)
Verifies if the parametric point pt is in the element parametric domain.
Definition: tpzpyramid.cpp:952
static constexpr REAL RefElVolume()
Volume of the master element.
Definition: tpzpyramid.h:227
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
static void BlendFactorForSide(const int &side, const TPZVec< T > &xi, T &blendFactor, TPZVec< T > &corrFactorDxi)
Definition: tpzpyramid.cpp:365
static int NSideNodes(int side)
Returns the number of nodes (not connectivities) associated with a side.
Definition: tpzpyramid.cpp:82
static void HigherDimensionSides(int side, TPZStack< int > &high)
Returns all sides whose closure contains side.
Definition: tpzpyramid.cpp:482
static void LowerDimensionSides(int side, TPZStack< int > &smallsides)
Get all sides with lower dimension on side.
Definition: tpzpyramid.cpp:464
static void MapToSide(int side, TPZVec< T > &InternalPar, TPZVec< T > &SidePar, TPZFMatrix< T > &JacToSide)
void Read(TPZStream &buf, void *context) override
read objects from the stream
Definition: fad.h:54
Implements the graphical element for a pyramid using a map to the cube element. Post processing...
Defines enum MElementType and contains the implementation of MElementType_NNodes(...) functions.
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
Handles the numerical integration for three-dimensional problems using pyramid elements. Numerical Integration.
Definition: pzquad.h:369
static void RandomPoint(TPZVec< REAL > &pt)
Generates a random point in the master domain.
Definition: tpzpyramid.cpp:973
static REAL gTolerance
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
static void ParametricDomainNodeCoord(int node, TPZVec< REAL > &nodeCoord)
static const double tol
Definition: pzgeoprism.cpp:23
Groups all classes defining the structure of the master element.
Definition: PrismExtend.cpp:15
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
static void AdjustTopDirections(int ConstrainedFace, TPZFMatrix< TVar > &gradx, TPZFMatrix< TVar > &directions)
Adjust the directions associated with the tip of the pyramid, considering that one of the faces is co...
static bool CheckProjectionForSingularity(const int &side, const TPZVec< T > &xiInterior)
Definition: tpzpyramid.cpp:985
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
Contains TPZMatrixclass which implements full matrix (using column major representation).
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
static void CornerShape(const TPZVec< REAL > &pt, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi)
Computes the corner shape functions of the element.
virtual ~TPZPyramid()
Default destructor.
Definition: tpzpyramid.h:50
static int SideDimension(int side)
Returns the dimension of the side.
Definition: tpzpyramid.cpp:520
int ClassId() const override
Define the class id associated with the class.
static int SideNodes[8][2]
Nodes over lines sides (1d)
Definition: tpzpyramid.h:262
static int NContainedSides(int side)
Returns the number of nodes (not connectivities) associated with a side.
Definition: tpzpyramid.cpp:913
static int GetTransformId(TPZVec< int64_t > &id)
Method which identifies the transformation based on the IDs of the corner nodes.
Definition: tpzpyramid.h:184
static void CenterPoint(int side, TPZVec< REAL > &center)
Returns the barycentric coordinates in the master element space of the original element.
Definition: tpzpyramid.cpp:510
static int NBilinearSides()
Definition: tpzpyramid.cpp:458
A simple stack.
static MElementType Type()
Returns the type of the element as specified in file pzeltype.h.
Definition: tpzpyramid.cpp:859
static void Shape(TPZVec< REAL > &loc, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi)
Compute the shape being used to construct the x mapping from local parametric coordinates.
Definition: tpzpyramid.h:88
static int SideNodeLocId(int side, int node)
Returns the local node number of the node "node" along side "side".
Definition: tpzpyramid.cpp:492
static void ComputeHDivDirections(TPZFMatrix< TVar > &gradx, TPZFMatrix< TVar > &directions)
Compute the directions of the HDiv vectors.
MElementType
Define the element types.
Definition: pzeltype.h:52
static TPZIntPoints * CreateSideIntegrationRule(int side, int order)
Create an integration rule over side.
Definition: tpzpyramid.cpp:842
static TPZTransform TransformElementToSide(int side)
Returns the transformation which projects a point from the interior of the element to the side...
Definition: tpzpyramid.cpp:563
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
static void TShape(const TPZVec< T > &loc, TPZFMatrix< T > &phi, TPZFMatrix< T > &dphi)
Compute the shape being used to construct the x mapping from local parametric coordinates.
Definition: tpzpyramid.cpp:302
Contains declaration of the TPZAxesTools class which implements verifications over axes...
TPZIntPyram3D IntruleType
Typedef to numerical integration rule.
Definition: tpzpyramid.h:211
Contains the TPZTransform<> class which implements an affine transformation between points in paramet...
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
TPZGraphElPyramidMapped GraphElType
Typedef to graphical element type.
Definition: tpzpyramid.h:213
static int NumSides()
Number of connects of the element (21)
Definition: tpzpyramid.cpp:897
static TPZTransform TransformSideToElement(int side)
Returns the transformation which transform a point from the side to the interior of the element...
Definition: tpzpyramid.cpp:711
Implements an affine transformation between points in parameter space. Topology Utility.
Definition: pzmganalysis.h:14
This class defines the interface to save and restore objects from TPZStream objects. Persistency.
Definition: TPZSavable.h:67
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
static void ComputeDirections(int side, TPZFMatrix< REAL > &gradx, TPZFMatrix< REAL > &directions, TPZVec< int > &sidevectors)
static void GetSideHDivPermutation(int transformationid, TPZVec< int > &permgather)
Identifies the permutation of the nodes needed to make neighbouring elements compatible in terms of o...