25 #ifndef _GEMA_INT2D_SHAPE_H_ 26 #define _GEMA_INT2D_SHAPE_H_ 51 virtual void naturalCoordLimits(
int coord,
double* min,
double* max)
const { Q_UNUSED(coord); *min = -1.0; *max = 1.0; }
60 elementCoord = srcEdgeCoord;
67 Q_UNUSED(face); Q_UNUSED(srcFaceCoord); Q_UNUSED(elementCoord);
return false;
75 assert(J.n_elem == 2);
78 double dy = (J.n_cols == 2) ? J(0, 1) : J(1, 0);
79 return sqrt(dx*dx + dy*dy);
84 GmMatrix& P,
bool transposed =
false)
const;
88 const GmMatrix& X,
bool transposed =
false)
const 90 Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
96 const GmMatrix& X,
bool transposed =
false)
const 98 Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
104 const GmMatrix& X,
bool transposed =
false)
const 106 Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
116 virtual bool gradientBasedCartesianToNatural(
const GmVector& coord,
const GmMatrix& X,
double tol,
117 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const 119 Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
153 GmVector& ncoord,
bool* inside)
const;
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmInt2DShape.h:140
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmInt2DShape.h:45
virtual bool translateFacePoint(int face, const GmVector &srcFaceCoord, GmVector &elementCoord) const
Given a "quad like" or "tri like" face coordinate (-1 to 1 pair for quad faces and 0 to 1 barycentric...
Definition: gmInt2DShape.h:65
virtual double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmInt2DShape.h:73
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmInt2DShape.h:54
virtual void naturalToCartesian(const GmVector &ncoord, const GmMatrix &X, GmVector &coord) const
Given a set of natural coordinates 'ncoord' and a matrix with node coordinates 'X',...
Definition: gmShape.cpp:106
virtual void jacobianAndPartials(const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed=false) const =0
This function does the same calculations as the jacobian() call, but also filling the extra parameter...
virtual void shapePartials(const GmVector &ncoord, GmMatrix &dN, bool transposed=false) const =0
Function used to calculate shape function partial derivatives with respect to its natural coordinates...
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:137
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:170
A 2D linear interface element with 4 nodes.
Definition: gmCellType.h:45
virtual bool translateEdgePoint(int edge, const GmVector &srcEdgeCoord, GmVector &elementCoord) const
Given a "bar like" edge coordinate (from -1 to 1), at the given edge, fills elementCoord with the equ...
Definition: gmInt2DShape.h:57
GmShape specialization for a 2D quadratic interface element with 8 nodes object 2 nodes (x,...
Definition: gmInt2DShape.h:208
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:186
Shape function handling base classe.
Definition: gmShape.h:37
virtual void shapeValues(const GmVector &ncoord, GmVector &N) const =0
Function used to evaluate the set of shape functions over a point defined by its natural coordinates.
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:212
GmShape specialization for a 2D linear interface element with 4 nodes object.
Definition: gmInt2DShape.h:133
GmShape specialization for Cohesive elements, containing common functions used to model fracture proc...
Definition: gmInt2DShape.h:41
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmInt2DShape.h:189
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmInt2DShape.h:48
#define GMC_API_EXPORT
Macro for controling if the class is being exported (GEMA_CORE_LIB defined) or imported (GEMA_CORE_LI...
Definition: gmCoreConfig.h:35
virtual bool geometryBasedCartesianToNatural(const GmVector &coord, const GmMatrix &X, double natTol, GmVector &ncoord, bool *inside) const
Virtual function that should be implemented if a shape can provide a geometric algorithm for converti...
Definition: gmShape.h:410
GmCellType
Mesh Cell types. Don't change type orders or add types without reading comments below.
Definition: gmCellType.h:30
A 2D linear interface element with 6 nodes (2 nodes of pore pressure dof at the middle)
Definition: gmCellType.h:46
virtual double shapeCartesianPartialsFromPJ(const GmMatrix &P, const GmMatrix &J, GmMatrix &dN, bool transposed=false) const =0
Alternative version of shapeCartesianPartialsXxxx() to calculate shape function partial derivatives w...
virtual double borderScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
Returns the scaling factor needed when calculating integrals over borders (edges or faces)....
Definition: gmInt2DShape.h:87
virtual double edgeScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
Returns the scaling factor needed when calculating integrals over element edges. This has the same ef...
Definition: gmInt2DShape.h:95
A 2D quadratic interface element with 6 nodes.
Definition: gmCellType.h:47
virtual void naturalCoordLimits(int coord, double *min, double *max) const
Fills min and max with the domain limits for the given natural coordinate (between 0 and numNaturalCo...
Definition: gmInt2DShape.h:51
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmInt2DShape.h:150
virtual void nodeNaturalCoord(int node, GmVector &coord) const =0
Fills the coord vector with the set of natural coordinates for the reference shape function node....
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
GmShape specialization for a 2D quadratic interface element with 6 nodes object.
Definition: gmInt2DShape.h:182
Declaration of the GmShape base class.
virtual double faceScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
Returns the scaling factor needed when calculating integrals over element faces. This has the same ef...
Definition: gmInt2DShape.h:103
GmShape specialization for a 2D coupled interface element with 6 nodes object with 4 nodes (x,...
Definition: gmInt2DShape.h:166
A 2D quadratic interface element with 8 nodes (2 nodes of pore pressure dof at the middle)
Definition: gmCellType.h:48