NeoPZ
tpztriangle.h
Go to the documentation of this file.
1 
6 #ifndef PZTOPOLOGYTPZTRIANGLE_H
7 #define PZTOPOLOGYTPZTRIANGLE_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 
16 #include "TPZTopologyUtils.h"
17 
18 class TPZIntPoints;
19 class TPZGraphElTd;
20 class TPZIntTriang;
22 
23 class TPZCompEl;
24 class TPZGeoEl;
25 class TPZCompMesh;
26 
28 namespace pztopology {
35  class TPZTriangle : public TPZSavable{
36  public:
37  friend void pztopology::GetPermutation<TPZTriangle>(const int permute, TPZVec<int> &permutation);
39  enum {NSides = 7, NCornerNodes= 3, Dimension = 2, NFaces = 3, NPermutations = 6};
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 ~TPZTriangle(){
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);
86  static int ContainedSideLocId(int side, int c);
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);
105 
112  static void CenterPoint(int side, TPZVec<REAL> &center);
113 
115  static bool IsInParametricDomain(const TPZVec<REAL> &pt, REAL tol = pztopology::gTolerance);
116 
117  #ifdef _AUTODIFF
118 
119  static bool IsInParametricDomain(const TPZVec<Fad<REAL>> &pt, REAL tol = pztopology::gTolerance){
120  TPZVec<REAL> xi(pt.size());
121  for(int i = 0; i < pt.size(); i++) xi[i]= pt[i].val();
122  return IsInParametricDomain(xi,tol);
123  }
124  #endif
125  #ifdef _AUTODIFF
126  template<typename T,
127  typename std::enable_if<std::is_same<T,Fad<REAL>>::value>::type* = nullptr>
128  static bool IsInParametricDomain(const TPZVec<T> &pt, REAL tol){
129  TPZVec<REAL> qsiReal(pt.size(),-1);
130  for(int i = 0; i < qsiReal.size(); i++) qsiReal[i] = pt[i].val();
131  return IsInParametricDomain(qsiReal,tol);
132  }
133  #endif
134 
135  static void RandomPoint(TPZVec<REAL> &pt);
136 
144  template<class T>
145  static bool CheckProjectionForSingularity(const int &side, const TPZVec<T> &xiInterior);
146 
147  template<class T>
148  static void MapToSide(int side, TPZVec<T> &InternalPar, TPZVec<T> &SidePar, TPZFMatrix<T> &JacToSide);
149 
150  static void ParametricDomainNodeCoord(int node, TPZVec<REAL> &nodeCoord);
151 
158  static MElementType Type();// { return ETriangle;}
159 
161  static MElementType Type(int side);
162 
174  static TPZTransform<> SideToSideTransform(int sidefrom, int sideto);
175 
181  static TPZTransform<> TransformSideToElement(int side);
187  static TPZTransform<> TransformElementToSide(int side);
188 
189 
195  static int GetTransformId(TPZVec<int64_t> &id);
196 
203  static int GetTransformId(int side, TPZVec<int64_t> &id);
204 
208  static void GetHDivGatherPermute(int transformid, TPZVec<int> &permute);
209 
220  static TPZIntPoints * CreateSideIntegrationRule(int side, int order);
221 
226 
236  static void GetSideHDivPermutation(int transformationid, TPZVec<int> &permgather);
237 
239  static constexpr REAL RefElVolume() { return 0.5L; }
240 
241  /* Given side and gradx the method returns directions needed for Hdiv space */
242 // static void ComputeDirections(int side, TPZFMatrix<REAL> &gradx, TPZFMatrix<REAL> &directions, TPZVec<int> &sidevectors);
243 
245  static void GetSideHDivDirections(TPZVec<int> &sides, TPZVec<int> &dir, TPZVec<int> &bilinearounao, TPZVec<int> &sidevectors);
246 
248  template <class TVar>
249  static void ComputeHDivDirections(TPZFMatrix<TVar> &gradx, TPZFMatrix<TVar> &directions);
250 
266  template <class TVar>
267  static void ComputeHCurlDirections(TPZFMatrix<TVar> &gradx, TPZFMatrix<TVar> &directions, const TPZVec<int> &transformationIds);
268 
276  template <class TVar>
277  static void ComputeHCurlFaceDirections(TPZVec<TVar> &v1, TPZVec<TVar> &v2, int transformationId);
278 
279 
283  static int NBilinearSides();
284 
285  static int SideNodes[3][2];
286  static int FaceNodes[1][3];
287  static REAL gTrans2dT[6][2][2] ;
288 
289 
290  protected:
292  static int fPermutations [6][7];
293  static REAL fTangentVectors [12][2];
294  };
295 
296 }
297 
298 #endif
static void GetSideHDivDirections(TPZVec< int > &sides, TPZVec< int > &dir, TPZVec< int > &bilinearounao)
static int NContainedSides(int side)
Returns the number of nodes (not connectivities) associated with a side.
static int bilinearounao[81]
Definition: tpzcube.cpp:405
static REAL gTrans2dT[6][2][2]
Definition: tpztriangle.h:287
static void ComputeHDivDirections(TPZFMatrix< TVar > &gradx, TPZFMatrix< TVar > &directions)
Compute the directions of the HDiv vectors.
static void ParametricDomainNodeCoord(int node, TPZVec< REAL > &nodeCoord)
static void ComputeHCurlDirections(TPZFMatrix< TVar > &gradx, TPZFMatrix< TVar > &directions, const TPZVec< int > &transformationIds)
Handles the numerical integration for two-dimensional problems using triangular elements. Numerical Integration.
Definition: pzquad.h:102
static bool IsInParametricDomain(const TPZVec< REAL > &pt, REAL tol=pztopology::gTolerance)
Verifies if the parametric point pt is in the element parametric domain.
static void GetHDivGatherPermute(int transformid, TPZVec< int > &permute)
return the vector which permutes the connects according to the transformation id
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: tpztriangle.cpp:27
static TPZTransform TransformSideToElement(int side)
Returns the transformation which transform a point from the side to the interior of the element...
static void BlendFactorForSide(const int &side, const TPZVec< T > &xi, T &blendFactor, TPZVec< T > &corrFactorDxi)
Definition: tpztriangle.cpp:53
static void LowerDimensionSides(int side, TPZStack< int > &smallsides)
Get all sides with lower dimension on side.
TPZGraphElT2dMapped GraphElType
Typedef to graphical element type.
Definition: tpztriangle.h:225
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
Defines the topology of a triangle element. Topology Sides 0 to 2 are vertices, sides 3 to 5 are line...
Definition: tpztriangle.h:35
static constexpr REAL RefElVolume()
Volume of the master element (measure)
Definition: tpztriangle.h:239
static TPZTransform TransformElementToSide(int side)
Returns the transformation which transform a point from the interior of the element to the side...
static void ComputeHCurlFaceDirections(TPZVec< TVar > &v1, TPZVec< TVar > &v2, int transformationId)
static int ContainedSideLocId(int side)
Returns the local connect along side "side" especial for hdivspace.
static TPZIntPoints * CreateSideIntegrationRule(int side, int order)
Create an integration rule over side.
TPZTriangle()
Default constructor.
Definition: tpztriangle.h:46
static void GetSideHDivPermutation(int transformationid, TPZVec< int > &permgather)
Identifies the permutation of the nodes needed to make neighbouring elements compatible in terms of o...
static REAL fTangentVectors[12][2]
Definition: tpztriangle.h:293
Definition: fad.h:54
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: tpztriangle.h:88
Defines enum MElementType and contains the implementation of MElementType_NNodes(...) functions.
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
static bool CheckProjectionForSingularity(const int &side, const TPZVec< T > &xiInterior)
static MElementType Type()
Returns the type of the element as specified in file pzeltype.h.
static int SideNodeLocId(int side, int node)
Returns the local node number of the node "node" along side "side".
static REAL gTolerance
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
Implements a graphical element for a triangle mapped into de quadrilateral element. Post processing.
static void HigherDimensionSides(int side, TPZStack< int > &high)
Returns all sides whose closure contains side.
void Read(TPZStream &buf, void *context) override
read objects from the stream
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
static int FaceNodes[1][3]
Definition: tpztriangle.h:286
static const double tol
Definition: pzgeoprism.cpp:23
Groups all classes defining the structure of the master element.
Definition: PrismExtend.cpp:15
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
static int SideDimension(int side)
Returns the dimension of the side.
Contains TPZMatrixclass which implements full matrix (using column major representation).
int ClassId() const override
Define the class id associated with the class.
To export a graphical discontinuous triangular element. Post processing.
Definition: pztrigraphd.h:16
static int GetTransformId(TPZVec< int64_t > &id)
Method which identifies the transformation based on the IDs of the corner nodes.
static void MapToSide(int side, TPZVec< T > &InternalPar, TPZVec< T > &SidePar, TPZFMatrix< T > &JacToSide)
static void RandomPoint(TPZVec< REAL > &pt)
Generates a random point in the master domain.
A simple stack.
static void CenterPoint(int side, TPZVec< REAL > &center)
Returns the barycentric coordinates in the master element space of the original element.
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
MElementType
Define the element types.
Definition: pzeltype.h:52
static int NSideNodes(int side)
Returns the number of nodes (not connectivities) associated with a side.
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
static int SideNodes[3][2]
Definition: tpztriangle.h:285
static int NumSides()
Returns the number of connects of the element (7)
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
static int fPermutations[6][7]
Valid permutations between nodes.
Definition: tpztriangle.h:292
Implements an affine transformation between points in parameter space. Topology Utility.
Definition: pzmganalysis.h:14
virtual ~TPZTriangle()
Default destructor.
Definition: tpztriangle.h:50
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.
TPZIntTriang IntruleType
Typedef to numerical integration rule.
Definition: tpztriangle.h:223
static int NBilinearSides()
static TPZTransform SideToSideTransform(int sidefrom, int sideto)
Returns the transformation which takes a point from the side sidefrom of the side sideto...