25 #ifndef _GEMA_HEX_SHAPE_H_ 26 #define _GEMA_HEX_SHAPE_H_ 51 virtual void naturalCoordLimits(
int coord,
double* min,
double* max)
const { Q_UNUSED(coord); *min = -1.0; *max = 1.0; }
72 const GmMatrix& X,
bool transposed =
false)
const 78 const GmMatrix& X,
bool transposed =
false)
const;
81 const GmMatrix& X,
bool transposed =
false)
const;
84 const int* fixedEdgeCoord(
int border)
const;
85 const double* fixedEdgeValue(
int border)
const;
86 int fixedFaceCoord(
int border)
const;
87 double fixedFaceValue(
int border)
const;
90 virtual bool gradientBasedCartesianToNatural(
const GmVector& coord,
const GmMatrix& X,
double tol,
91 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...
An hexahedron(brick) with 20 nodes.
Definition: gmCellType.h:50
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: gmHexShape.h:65
GmShape specialization for a Hexahedron with 20 nodes object.
Definition: gmHexShape.h:128
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmHexShape.h:45
An hexahedron(brick) with 27 nodes.
Definition: gmCellType.h:51
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...
bool gradientBasedCartesianToNaturalQuadHex(const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
Given a set of cartesian coordinates 'coord' and a matrix with node coordinates 'X',...
Definition: gmShape.cpp:174
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 int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmHexShape.h:112
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: gmHexShape.h:51
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 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
A linear hexahedron (brick) with 8 nodes.
Definition: gmCellType.h:49
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: gmHexShape.h:71
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmHexShape.h:54
GmShape specialization for a Hexahedron with 27 nodes object.
Definition: gmHexShape.h:151
GmShape specialization for a Hexahedron with 8 nodes object.
Definition: gmHexShape.h:105
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...
GmShape specialization for Hexahedrons, containing common functions used by concrete Hex specializati...
Definition: gmHexShape.h:41
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmHexShape.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
GmCellType
Mesh Cell types. Don't change type orders or add types without reading comments below.
Definition: gmCellType.h:30
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 GmCellType elemType() const
Returns the type of this element.
Definition: gmHexShape.h:132
virtual double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmHexShape.h:62
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmHexShape.h:158
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmHexShape.h:109
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 GmCellType elemType() const
Returns the type of this element.
Definition: gmHexShape.h:155
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 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 int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmHexShape.h:135