24 #ifndef _GEMA_PYRA_SHAPE_H_ 25 #define _GEMA_PYRA_SHAPE_H_ 50 virtual void naturalCoordLimits(
int coord,
double* min,
double* max)
const { Q_UNUSED(coord); *min = -1.0; *max = 1.0; }
71 const GmMatrix& X,
bool transposed =
false)
const 77 const GmMatrix& X,
bool transposed =
false)
const;
80 const GmMatrix& X,
bool transposed =
false)
const;
83 const int* fixedEdgeCoord(
int border)
const;
84 const double* fixedEdgeValue(
int border)
const;
85 int fixedFaceCoord(
int border)
const;
86 double fixedFaceValue(
int border)
const;
89 virtual bool gradientBasedCartesianToNatural(
const GmVector& coord,
const GmMatrix& X,
double tol,
90 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const virtual bool translateFacePoint(int face, const GmVector &srcFaceCoord, GmVector &elementCoord) const =0
Given a "quad like" or "tri like" face coordinate (-1 to 1 pair for quad faces and 0 to 1 barycentric...
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: gmPyraShape.h:50
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: gmPyraShape.h:70
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmPyraShape.h:111
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 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...
GmShape specialization for a Tetrahedron with 4 nodes object.
Definition: gmPyraShape.h:129
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...
virtual double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmPyraShape.h:61
Shape function handling base classe.
Definition: gmShape.h:37
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmPyraShape.h:53
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmPyraShape.h:47
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 jacIndependentNatCoord(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 inde...
Definition: gmShape.cpp:800
GmShape specialization for Pyramid, containing common functions used by concrete Pyra specializations...
Definition: gmPyraShape.h:40
virtual double faceScalingFactor(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 faces. This has the same ef...
A pyramid with 13 nodes - quadratic.
Definition: gmCellType.h:65
#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
GmCellType
Mesh Cell types. Don't change type orders or add types without reading comments below.
Definition: gmCellType.h:30
GmShape specialization for a Pyramid with 5 nodes object.
Definition: gmPyraShape.h:104
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 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: gmPyraShape.h:64
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmPyraShape.h:108
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmPyraShape.h:44
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
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmPyraShape.h:136
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.
A pyramid with 5 nodes - linear.
Definition: gmCellType.h:64
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmPyraShape.h:133
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...