NeoPZ
tpzprism.h
Go to the documentation of this file.
1 
6 #ifndef PZTOPOLOGYTPZPRISM_H
7 #define PZTOPOLOGYTPZPRISM_H
8 
9 #include "pzfmatrix.h"
10 #include "pzstack.h"
11 #include "pztrnsform.h"
12 #include "pzeltype.h"
13 #include "pznumeric.h"
14 #include "pzaxestools.h"
15 #include "TPZTopologyUtils.h"
16 
17 class TPZIntPoints;
18 class TPZIntPrism3D;
20 
21 class TPZCompEl;
22 class TPZGeoEl;
23 class TPZCompMesh;
24 
26 namespace pztopology {
34  class TPZPrism : public TPZSavable {
35  public:
36  friend void pztopology::GetPermutation<TPZPrism>(const int permute, TPZVec<int> &permutation);
38  enum {NSides = 21, NCornerNodes = 6, Dimension = 3, NFaces = 5, NPermutations = 12};
39 
40  virtual int ClassId() const override;
41  void Read(TPZStream& buf, void* context) override;
42  void Write(TPZStream& buf, int withclassid) const override;
43 
44 
45 
48  }
49 
51  virtual ~TPZPrism() {
52  }
53 
58  static int SideDimension(int side);
59 
61  static void LowerDimensionSides(int side,TPZStack<int> &smallsides);
63  static void LowerDimensionSides(int side,TPZStack<int> &smallsides, int DimTarget);
64 
70  static void HigherDimensionSides(int side, TPZStack<int> &high);
71 
73  static int NSideNodes(int side);
75  static int SideNodeLocId(int side, int node);
76 
78  static int NumSides();
79 
81  static int NumSides(int dimension);
82 
84  static int NContainedSides(int side);
86  static int ContainedSideLocId(int side, int c);
87 
88 
90  static void Shape(TPZVec<REAL> &loc,TPZFMatrix<REAL> &phi,TPZFMatrix<REAL> &dphi){
91  TShape(loc, phi, dphi);
92  }
94  template<class T>
95  static void TShape(const TPZVec<T> &loc,TPZFMatrix<T> &phi,TPZFMatrix<T> &dphi);
104  template<class T>
105  static void BlendFactorForSide(const int &side, const TPZVec<T> &xi, T &blendFactor,
106  TPZVec<T> &corrFactorDxi);
113  static void CenterPoint(int side, TPZVec<REAL> &center);
114 
116  static bool IsInParametricDomain(const TPZVec<REAL> &pt, REAL tol = pztopology::gTolerance);
117 
118  #ifdef _AUTODIFF
119 
120  static bool IsInParametricDomain(const TPZVec<Fad<REAL>> &pt, REAL tol = pztopology::gTolerance){
121  TPZVec<REAL> xi(pt.size());
122  for(int i = 0; i < pt.size(); i++) xi[i]= pt[i].val();
123  return IsInParametricDomain(xi,tol);
124  }
125  #endif
126 
127 
129  static void RandomPoint(TPZVec<REAL> &pt);
130 
138  template<class T>
139  static bool CheckProjectionForSingularity(const int &side, const TPZVec<T> &xiInterior);
140 
141  template<class T>
142  static void MapToSide(int side, TPZVec<T> &InternalPar, TPZVec<T> &SidePar, TPZFMatrix<T> &JacToSide);
143 
144  static void ParametricDomainNodeCoord(int node, TPZVec<REAL> &nodeCoord);
145 
152  static MElementType Type();
153 
155  static MElementType Type(int side);
156 
168  static TPZTransform<> SideToSideTransform(int sidefrom, int sideto);
169 
175  static TPZTransform<> TransformSideToElement(int side);
181  static TPZTransform<> TransformElementToSide(int side);
182 
188  static int GetTransformId(TPZVec<int64_t> &id);
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 1.0L; }
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 
253  template <class TVar>
254  static void ComputeHCurlDirections(TPZFMatrix<TVar> &gradx, TPZFMatrix<TVar> &directions, const TPZVec<int> &transformationIds);
255 
256 
260  static int NBilinearSides();
261 
262  protected:
267  static int FaceNodes[5][4];
268 
270  static int SideNodes[9][2];
271 
273  static int ShapeFaceId[5][4];
274 
276  static int fPermutations[12][21];
279  };
280 
281 }
282 
283 #endif
static void CenterPoint(int side, TPZVec< REAL > &center)
Returns the barycentric coordinates in the master element space of the original element.
Definition: tpzprism.cpp:588
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
Definition: tpzprism.cpp:1923
static int bilinearounao[81]
Definition: tpzcube.cpp:405
static bool IsInParametricDomain(const TPZVec< REAL > &pt, REAL tol=pztopology::gTolerance)
Verifies if the parametric point pt is in the element parametric domain.
Definition: tpzprism.cpp:991
static void GetSideHDivPermutation(int transformationid, TPZVec< int > &permgather)
Identifies the permutation of the nodes needed to make neighbouring elements compatible in terms of o...
Definition: tpzprism.cpp:1420
static void ParametricDomainNodeCoord(int node, TPZVec< REAL > &nodeCoord)
Definition: tpzprism.cpp:1276
static MElementType Type()
Returns the type of the element as specified in file pzeltype.h.
Definition: tpzprism.cpp:896
virtual int ClassId() const override
Define the class id associated with the class.
Definition: tpzprism.cpp:1915
static void ComputeHDivDirections(TPZFMatrix< TVar > &gradx, TPZFMatrix< TVar > &directions)
Compute the directions of the HDiv vectors.
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
static void RandomPoint(TPZVec< REAL > &pt)
Generates a random point in the master domain.
Definition: tpzprism.cpp:1008
TPZIntPrism3D IntruleType
Typedef to numerical integration rule.
Definition: tpzprism.h:211
static TPZIntPoints * CreateSideIntegrationRule(int side, int order)
Create an integration rule over side.
Definition: tpzprism.cpp:878
static int GetTransformId(TPZVec< int64_t > &id)
Method which identifies the transformation based on the IDs of the corner nodes.
Definition: tpzprism.cpp:1340
static void MapToSide(int side, TPZVec< T > &InternalPar, TPZVec< T > &SidePar, TPZFMatrix< T > &JacToSide)
Definition: tpzprism.cpp:1085
static void GetSideHDivDirections(TPZVec< int > &sides, TPZVec< int > &dir, TPZVec< int > &bilinearounao)
Definition: tpzprism.cpp:1752
Definition: fad.h:54
static void HigherDimensionSides(int side, TPZStack< int > &high)
Returns all sides whose closure contains side.
Definition: tpzprism.cpp:536
Defines enum MElementType and contains the implementation of MElementType_NNodes(...) functions.
static int NContainedSides(int side)
Returns the number of nodes (not connectivities) associated with a side.
Definition: tpzprism.cpp:941
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
static int SideDimension(int side)
Returns the dimension of the side.
Definition: tpzprism.cpp:598
static bool CheckProjectionForSingularity(const int &side, const TPZVec< T > &xiInterior)
Definition: tpzprism.cpp:1019
static int NumSides()
Returns the number of connects of the element (21)
Definition: tpzprism.cpp:936
static REAL gTolerance
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
static int ShapeFaceId[5][4]
Ids of the shape face.
Definition: tpzprism.h:273
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
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: tpzprism.h:90
static void ComputeHCurlDirections(TPZFMatrix< TVar > &gradx, TPZFMatrix< TVar > &directions, const TPZVec< int > &transformationIds)
Definition: tpzprism.cpp:1788
static int SideNodeLocId(int side, int node)
Returns the local node number of the node "node" along side "side".
Definition: tpzprism.cpp:569
static const double tol
Definition: pzgeoprism.cpp:23
Groups all classes defining the structure of the master element.
Definition: PrismExtend.cpp:15
static int SideNodes[9][2]
Nodes over lines sides (1d)
Definition: tpzprism.h:270
static TPZTransform TransformSideToElement(int side)
Returns the transformation which transform a point from the side to the interior of the element...
Definition: tpzprism.cpp:749
Implements the graphical element for a prism using a degenerated cube element. Post processing...
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
Handles the numerical integration for three-dimensional problems using prism elements. Numerical Integration This cubature rule uses a cubature rule for one dimension (zeta) and a cubature rule for triangle (base).
Definition: pzquad.h:417
Contains TPZMatrixclass which implements full matrix (using column major representation).
static int ContainedSideLocId(int side, int c)
Returns the local connect number of the connect "c" along side "side".
Definition: tpzprism.cpp:951
void Read(TPZStream &buf, void *context) override
read objects from the stream
Definition: tpzprism.cpp:1919
static int NBilinearSides()
Definition: tpzprism.cpp:512
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: tpzprism.cpp:348
Defines the topology of a Prism. Topology Sides 0 to 7 are vertices, sides 7 to 14 are lines...
Definition: tpzprism.h:34
static void LowerDimensionSides(int side, TPZStack< int > &smallsides)
Get all sides with lower dimension on side.
Definition: tpzprism.cpp:518
static int fPermutations[12][21]
Valid permutations between nodes.
Definition: tpzprism.h:276
A simple stack.
virtual ~TPZPrism()
Default destructor.
Definition: tpzprism.h:51
static void ComputeDirections(int side, TPZFMatrix< REAL > &gradx, TPZFMatrix< REAL > &directions, TPZVec< int > &sidevectors)
Definition: tpzprism.cpp:1531
TPZGraphElPrismMapped GraphElType
Typedef to graphical element type.
Definition: tpzprism.h:213
static TPZTransform SideToSideTransform(int sidefrom, int sideto)
Returns the transformation which takes a point from the side sidefrom to the side sideto...
Definition: tpzprism.cpp:606
MElementType
Define the element types.
Definition: pzeltype.h:52
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
static int NSideNodes(int side)
Returns the number of nodes (not connectivities) associated with a side.
Definition: tpzprism.cpp:552
static TPZTransform TransformElementToSide(int side)
Returns the transformation which projects a point from the interior of the element to the side...
Definition: tpzprism.cpp:645
static void BlendFactorForSide(const int &side, const TPZVec< T > &xi, T &blendFactor, TPZVec< T > &corrFactorDxi)
Definition: tpzprism.cpp:380
Contains declaration of the TPZAxesTools class which implements verifications over axes...
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
Implements an affine transformation between points in parameter space. Topology Utility.
Definition: pzmganalysis.h:14
static int FaceNodes[5][4]
Nodes over quadrilateral sides (2d - faces).
Definition: tpzprism.h:267
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
Contains declaration of the TPZNumeric class which implements several methods to calculation.
static constexpr REAL RefElVolume()
Volume of the master element.
Definition: tpzprism.h:227
TPZPrism()
Default constructor.
Definition: tpzprism.h:47