25 #ifndef _GEMA_TRI_SHAPE_H_ 26 #define _GEMA_TRI_SHAPE_H_ 51 virtual void naturalCoordLimits(
int coord,
double* min,
double* max)
const { Q_UNUSED(coord); *min = 0.0; *max = 1.0; }
61 Q_UNUSED(face); Q_UNUSED(srcFaceCoord); Q_UNUSED(elementCoord);
return false;
77 const GmMatrix& X,
bool transposed =
false)
const 84 const GmMatrix& X,
bool transposed =
false)
const 86 Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
91 static int fixedEdgeCoord(
int border);
96 virtual bool gradientBasedCartesianToNatural(
const GmVector& coord,
const GmMatrix& X,
double tol,
97 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const 126 const GmMatrix& X,
bool transposed =
false)
const;
133 GmVector& ncoord,
bool* inside)
const;
160 const GmMatrix& X,
bool transposed =
false)
const;
185 virtual bool gradientBasedCartesianToNatural(
const GmVector& coord,
const GmMatrix& X,
double tol,
186 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const 188 Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
213 virtual bool gradientBasedCartesianToNatural(
const GmVector& coord,
const GmMatrix& X,
double tol,
214 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const 216 Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:206
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: gmTriShape.h:83
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmTriShape.h:182
bool gradientBasedCartesianToNaturalTriTet(const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
An implementation of the gradient based cartesian to natural algorithm specific to Tri and Tet elemen...
Definition: gmShape.cpp:236
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmTriShape.h:118
virtual void jacobianAndPartials(const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed=false) const
This function does the same calculations as the jacobian() call, but also filling the extra parameter...
Definition: gmTriShape.h:70
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmTriShape.h:176
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: gmTriShape.h:76
GmShape specialization for a Triangle with 6 nodes object.
Definition: gmTriShape.h:145
A 3D quadratic triangle with 6 nodes.
Definition: gmCellType.h:44
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 double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmTriShape.h:67
A 3D linear triangle with only 3 nodes.
Definition: gmCellType.h:43
virtual bool translateEdgePoint(int edge, const GmVector &srcEdgeCoord, GmVector &elementCoord) const =0
Given a "bar like" edge coordinate (from -1 to 1), at the given edge, fills elementCoord with the equ...
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 void jacDependentNatCoord(const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed) const
A generic implementation of jacobianAndPartials() for the case where all natural coordinates are NOT ...
Definition: gmShape.cpp:822
A quadratic triangle with 6 nodes.
Definition: gmCellType.h:42
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmTriShape.h:48
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmTriShape.h:45
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmTriShape.h:152
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmTriShape.h:130
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:175
A linear triangle with only 3 nodes.
Definition: gmCellType.h:41
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmTriShape.h:207
#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
GmShape specialization for a 3D Triangle with 6 nodes object.
Definition: gmTriShape.h:202
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 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: gmTriShape.h:59
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:115
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:149
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
GmShape specialization for a Triangle with 3 nodes object.
Definition: gmTriShape.h:111
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
Declaration of the GmShape base class.
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmTriShape.h:54
GmShape specialization for Triangles, containing common functions used by concrete Tri specialization...
Definition: gmTriShape.h:41
GmShape specialization for a 3D Triangle with 3 nodes object.
Definition: gmTriShape.h:171
virtual double edgeScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const =0
Returns the scaling factor needed when calculating integrals over element edges. This has the same ef...
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: gmTriShape.h:51