NeoPZ
tpzline.h
Go to the documentation of this file.
1 
6 #ifndef PZTOPOLOGYTPZLINE_H
7 #define PZTOPOLOGYTPZLINE_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 
16 #ifdef _AUTODIFF
17 #include "fadType.h"
18 #endif
19 
20 #include "TPZTopologyUtils.h"
21 
22 class TPZIntPoints;
23 class TPZInt1d;
24 class TPZGraphEl1dd;
25 
26 class TPZCompEl;
27 class TPZGeoEl;
28 class TPZCompMesh;
29 
31 namespace pztopology {
38  class TPZLine : public TPZSavable {
39  public:
40  friend void pztopology::GetPermutation<TPZLine>(const int permute, TPZVec<int> &permutation);
42  enum {NCornerNodes = 2, NSides = 3, Dimension = 1, NFaces = 2, NPermutations = 2};
43 
44  public:
45  int ClassId() const override;
46  void Read(TPZStream &buf, void *context) override;
47  void Write(TPZStream &buf, int withclassid) const override;
48 
51  }
52 
54  virtual ~TPZLine() {
55  }
56 
58  static void Shape(TPZVec<REAL> &loc,TPZFMatrix<REAL> &phi,TPZFMatrix<REAL> &dphi){
59  TShape(loc, phi, dphi);
60  }
62  template<class T>
63  static void TShape(const TPZVec<T> &loc,TPZFMatrix<T> &phi,TPZFMatrix<T> &dphi);
64 
76  template<class T>
77  static void BlendFactorForSide(const int &side, const TPZVec<T> &xi, T &blendFactor,
78  TPZVec<T> &corrFactorDxi);
80  static int SideDimension(int side);
81 
83  static void LowerDimensionSides(int side,TPZStack<int> &smallsides);
85  static void LowerDimensionSides(int side,TPZStack<int> &smallsides, int DimTarget);
86 
92  static void HigherDimensionSides(int side, TPZStack<int> &high);
93 
95  static int NSideNodes(int side);
96 
98  static int SideNodeLocId(int side, int node);
99 
101  static int NumSides();
103  static int NumSides(int dimension);
104 
106  static int NContainedSides(int side);
107 
109  static int ContainedSideLocId(int side, int c);
110 
117  static void CenterPoint(int side, TPZVec<REAL> &center);
118 
120  static bool IsInParametricDomain(const TPZVec<REAL> &pt, REAL tol = pztopology::gTolerance);
121 
122  #ifdef _AUTODIFF
123 
124  static bool IsInParametricDomain(const TPZVec<Fad<REAL>> &pt, REAL tol = pztopology::gTolerance){
125  TPZVec<REAL> xi(pt.size());
126  for(int i = 0; i < pt.size(); i++) xi[i]= pt[i].val();
127  return IsInParametricDomain(xi,tol);
128  }
129  #endif
130 
131  static void RandomPoint(TPZVec<REAL> &pt);
132 
140  template<class T>
141  static bool CheckProjectionForSingularity(const int &side, const TPZVec<T> &xiInterior);
142 
143  template<class T>
144  static void MapToSide(int side, TPZVec<T> &InternalPar, TPZVec<T> &SidePar, TPZFMatrix<T> &JacToSide);
145 
146  static void ParametricDomainNodeCoord(int node, TPZVec<REAL> &nodeCoord);
147 
154  static MElementType Type();
156  static MElementType Type(int side);
157 
168  static TPZTransform<> SideToSideTransform(int sidefrom, int sideto);
169 
175  static TPZTransform<> TransformSideToElement(int side);
176 
182  static TPZTransform<> TransformElementToSide(int side);
183 
189  static int GetTransformId(TPZVec<int64_t> &id);
190 
191 
192 
200  static int GetTransformId(int side, TPZVec<int64_t> &id);
201 
212  static TPZIntPoints *CreateSideIntegrationRule(int side, int order);
213 
218 
228  static void GetSideHDivPermutation(int transformationid, TPZVec<int> &permgather);
229 
231  static constexpr REAL RefElVolume() { return 2.0; }
232 
233  /* Given side and gradx the method returns directions needed for Hdiv space */
234  static void ComputeDirections(int side, TPZFMatrix<REAL> &gradx, TPZFMatrix<REAL> &directions, TPZVec<int> &sidevectors);
236  static void GetSideHDivDirections(TPZVec<int> &sides, TPZVec<int> &dir, TPZVec<int> &bilinearounao, TPZVec<int> &sidevectors);
237 
239  template <class TVar>
240  static void ComputeHDivDirections(TPZFMatrix<TVar> &gradx, TPZFMatrix<TVar> &directions)
241  {
242  TVar detjac = TPZAxesTools<TVar>::ComputeDetjac(gradx);
243 
244  for (int i=0; i<3; i++) {
245  directions(i,0) = -gradx(i,0)/detjac;
246  directions(i,1) = gradx(i,0)/detjac;
247  directions(i,2) = gradx(i,0)/detjac;
248  }
249  }
250 
251 
252 
256  static int NBilinearSides();
257  protected:
259  static int fPermutations [2][3];
260 
261  };
262 
263 }
264 
265 #endif
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: tpzline.h:58
static int bilinearounao[81]
Definition: tpzcube.cpp:405
Contains the TPZInt1d, TPZIntTriang, TPZIntQuad, TPZIntCube3D, TPZIntTetra3D, TPZIntPyram3D and TPZIn...
static TPZTransform SideToSideTransform(int sidefrom, int sideto)
Returns the transformation which takes a point from the side sidefrom to the side sideto...
Definition: tpzline.cpp:253
static void LowerDimensionSides(int side, TPZStack< int > &smallsides)
Get all sides with lower dimension on side.
Definition: tpzline.cpp:136
To export a graphical one dimensional discontinuous element. Post processing.
Definition: pzgraphel1dd.h:22
static int NumSides()
Returns the number of connects for a set dimension.
Definition: tpzline.cpp:376
static bool IsInParametricDomain(const TPZVec< REAL > &pt, REAL tol=pztopology::gTolerance)
Verifies if the parametric point pt is in the element parametric domain.
Definition: tpzline.cpp:154
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
static int NBilinearSides()
Definition: tpzline.cpp:60
static int SideDimension(int side)
Returns the dimension of the side.
Definition: tpzline.cpp:245
Definition: fad.h:54
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 int GetTransformId(TPZVec< int64_t > &id)
Method which identifies the transformation based on the IDs of the corner nodes.
Definition: tpzline.cpp:412
void Read(TPZStream &buf, void *context) override
read objects from the stream
Definition: tpzline.cpp:544
static TVar ComputeDetjac(TPZFMatrix< TVar > &gradx)
Definition: pzaxestools.h:111
TPZInt1d IntruleType
Typedef to numerical integration rule.
Definition: tpzline.h:215
Handles the numerical integration for one-dimensional problems. Numerical Integration.
Definition: pzquad.h:36
TPZGraphEl1dd GraphElType
Typedef to graphical element type.
Definition: tpzline.h:217
static REAL gTolerance
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
static void HigherDimensionSides(int side, TPZStack< int > &high)
Returns all sides whose closure contains side.
Definition: tpzline.cpp:206
static void ComputeHDivDirections(TPZFMatrix< TVar > &gradx, TPZFMatrix< TVar > &directions)
Compute the directions of the HDiv vectors.
Definition: tpzline.h:240
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
static void MapToSide(int side, TPZVec< T > &InternalPar, TPZVec< T > &SidePar, TPZFMatrix< T > &JacToSide)
Definition: tpzline.cpp:170
static const double tol
Definition: pzgeoprism.cpp:23
Groups all classes defining the structure of the master element.
Definition: PrismExtend.cpp:15
static TPZTransform TransformElementToSide(int side)
Returns the transformation which transform a point from the interior of the element to the side...
Definition: tpzline.cpp:293
int ClassId() const override
Define the class id associated with the class.
Definition: tpzline.cpp:540
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).
static void RandomPoint(TPZVec< REAL > &pt)
Generates a random point in the master domain.
Definition: tpzline.cpp:236
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: tpzline.cpp:447
static int NSideNodes(int side)
Returns the number of nodes (not connectivities) associated with a side.
Definition: tpzline.cpp:120
static TPZIntPoints * CreateSideIntegrationRule(int side, int order)
Create an integration rule over side.
Definition: tpzline.cpp:344
static int fPermutations[2][3]
Valid permutations between nodes.
Definition: tpzline.h:259
static void BlendFactorForSide(const int &side, const TPZVec< T > &xi, T &blendFactor, TPZVec< T > &corrFactorDxi)
Definition: tpzline.cpp:78
static int ContainedSideLocId(int side, int c)
Returns the local connect number of the connect "c" along side "side".
Definition: tpzline.cpp:389
Defines the topology of a line element. Topology Sides 0 and 1 are vertices, side 2 is the line...
Definition: tpzline.h:38
static TPZTransform TransformSideToElement(int side)
Returns the transformation which transform a point from the side to the interior of the element...
Definition: tpzline.cpp:316
A simple stack.
TPZLine()
Default constructor.
Definition: tpzline.h:50
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: tpzline.cpp:69
static void GetSideHDivDirections(TPZVec< int > &sides, TPZVec< int > &dir, TPZVec< int > &bilinearounao)
Definition: tpzline.cpp:506
MElementType
Define the element types.
Definition: pzeltype.h:52
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
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
virtual ~TPZLine()
Default destructor.
Definition: tpzline.h:54
static bool CheckProjectionForSingularity(const int &side, const TPZVec< T > &xiInterior)
Definition: tpzline.cpp:165
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 int SideNodeLocId(int side, int node)
Returns the local node number of the node "node" along side "side".
Definition: tpzline.cpp:216
static void ParametricDomainNodeCoord(int node, TPZVec< REAL > &nodeCoord)
Definition: tpzline.cpp:180
static constexpr REAL RefElVolume()
Volume of the master element (measure of the element)
Definition: tpzline.h:231
static int NContainedSides(int side)
Returns the number of nodes (not connectivities) associated with a side.
Definition: tpzline.cpp:381
static void CenterPoint(int side, TPZVec< REAL > &center)
Returns the barycentric coordinates in the master element space of the original element.
Definition: tpzline.cpp:224
static void ComputeDirections(int side, TPZFMatrix< REAL > &gradx, TPZFMatrix< REAL > &directions, TPZVec< int > &sidevectors)
Definition: tpzline.cpp:481
static MElementType Type()
Returns the type of the element as specified in file pzeltype.h.
Definition: tpzline.cpp:357
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
Definition: tpzline.cpp:548