24 #ifndef _GEMA_SHAPE_H_ 25 #define _GEMA_SHAPE_H_ 45 virtual int numFunctions()
const = 0;
48 virtual int numNaturalCoord()
const = 0;
51 virtual int numCartesianCoord()
const = 0;
54 virtual void naturalCoordLimits(
int coord,
double* min,
double* max)
const = 0;
64 virtual void nodeNaturalCoord(
int node,
GmVector& coord)
const = 0;
71 virtual void naturalCenter(
GmVector& coord)
const = 0;
73 void fillNaturalCoordinates(
GmMatrix& coord,
bool transposed =
false)
const;
87 virtual bool translateEdgePoint(
int edge,
const GmVector& srcEdgeCoord,
GmVector& elementCoord)
const = 0;
102 virtual bool translateFacePoint(
int face,
const GmVector& srcFaceCoord,
GmVector& elementCoord)
const = 0;
132 return cartesianToNatural(coord, X, _preferGeometricConfig, _gradTolConfig, _maxGradIterConfig,
133 _outsideNatTolConfig, ncoord, inside);
170 int maxGradIter,
double outsideNatTol,
GmVector& ncoord,
bool* inside)
const 172 if(preferGeometric && hasGeometryBasedCartesianToNatural())
173 return geometryBasedCartesianToNatural(coord, X, outsideNatTol, ncoord, inside);
175 return gradientBasedCartesianToNatural(coord, X, gradTol, maxGradIter, outsideNatTol, ncoord, inside);
223 virtual void shapePartials(
const GmVector& ncoord,
GmMatrix& dN,
bool transposed =
false)
const = 0;
225 virtual double shapeCartesianPartialsFromCoord (
const GmVector& ncoord,
const GmMatrix& X,
GmMatrix& dN,
bool transposed =
false)
const;
226 virtual double shapeCartesianPartialsFromJacobian(
const GmVector& ncoord,
const GmMatrix& J,
GmMatrix& dN,
bool transposed =
false)
const;
246 virtual double shapeCartesianPartialsFromPJ(
const GmMatrix& P,
const GmMatrix& J,
GmMatrix& dN,
bool transposed =
false)
const = 0;
278 jacobianAndPartials(ncoord, X, J, P, transposed);
288 virtual double scaledJacobianDet(
const GmMatrix& J)
const = 0;
321 virtual double borderScalingFactor(
int border,
const GmVector& borderCoord,
const GmVector& elementCoord,
322 const GmMatrix& X,
bool transposed =
false)
const = 0;
343 virtual double edgeScalingFactor(
int border,
const GmVector& borderCoord,
const GmVector& elementCoord,
344 const GmMatrix& X,
bool transposed =
false)
const = 0;
365 virtual double faceScalingFactor(
int border,
const GmVector& borderCoord,
const GmVector& elementCoord,
366 const GmMatrix& X,
bool transposed =
false)
const = 0;
369 virtual bool pointExtrapolationMatrix(
const GmMatrix& points,
GmMatrix& em)
const;
371 double interpolate(
const GmVector& nodeValues,
const GmVector& ncoord)
const;
373 static bool initShapeFunctions();
376 static const GmShape* linearShapeFromElementType(
GmCellType etype,
int P,
int Q);
384 #if defined ENABLE_TESTS 385 friend void shapeMappingTests();
387 virtual bool gradientBasedCartesianToNatural(
const GmVector& coord,
const GmMatrix& X,
double tol,
388 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const = 0;
390 bool gradientBasedCartesianToNaturalQuadHex(
const GmVector& coord,
const GmMatrix& X,
double tol,
391 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const;
392 bool gradientBasedCartesianToNaturalTriTet (
const GmVector& coord,
const GmMatrix& X,
double tol,
393 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const;
394 bool gradientBasedCartesianToNaturalBar (
const GmVector& coord,
const GmMatrix& X,
double tol,
395 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const;
396 bool gradientBasedCartesianToNaturalWedge (
const GmVector& coord,
const GmMatrix& X,
double tol,
397 int maxIter,
double natTol,
GmVector& ncoord,
bool* inside)
const;
411 GmVector& ncoord,
bool* inside)
const 413 Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
417 bool checkAndClampNaturalCoordinate(
double& ncoord,
double min,
double max,
double natTol)
const;
418 bool checkAndClampBarycentricNaturalCoordinates(
GmVector& ncoord,
double natTol)
const;
virtual bool cartesianToNatural(const GmVector &coord, const GmMatrix &X, bool preferGeometric, double gradTol, int maxGradIter, double outsideNatTol, GmVector &ncoord, bool *inside) const
Given a set of cartesian coordinates 'coord' and a matrix with node coordinates 'X',...
Definition: gmShape.h:169
A simple spin lock implementation based on a loop using test and set over an atomic int to change its...
Definition: gmSpinLock.h:36
static QMap< QPair< GmCellType, int >, GmMatrix > _emCache
A global cache, keyed by cell type + cacheKey(), storing extrapolation matrices.
Definition: gmShape.h:429
NOT a cell type. Stores the number of available types.
Definition: gmCellType.h:77
static int _maxGradIterConfig
The maximum number of iterations used by gradient based calls to cartesianToNatural.
Definition: gmShape.h:426
static QMutex _emCacheMutex
The mutex protecting _emCache.
Definition: gmShape.h:430
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmShape.h:404
Auxiliar class used to store the complete set of simulation data.
Definition: gmSimulationData.h:51
Integration rule base classe.
Definition: gmIntegrationRule.h:88
Shape function handling base classe.
Definition: gmShape.h:37
virtual void jacobian(const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, bool transposed=false) const
Calculates the jacobian matrix, relating cartesian coordinates to natural coordinates.
Definition: gmShape.h:275
#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
static double _gradTolConfig
The convergence tolerance used by gradient based calls to cartesianToNatural.
Definition: gmShape.h:425
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
virtual bool cartesianToNatural(const GmVector &coord, const GmMatrix &X, GmVector &ncoord, bool *inside) const
Given a set of cartesian coordinates 'coord' and a matrix with node coordinates 'X',...
Definition: gmShape.h:130
Declaration of the GmCell class.
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
static double _outsideNatTolConfig
The "snap" tolerance used by cartesianToNatural to consider a point outside the element.
Definition: gmShape.h:427
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
static bool _preferGeometricConfig
Should calls to cartesianToNatural prefer geometric algorithms (if available) over the gradient based...
Definition: gmShape.h:424
static GmSpinLock _hElementsLock
The spin-lock protecting the access/creation of shape pointers for hierarchical elements.
Definition: gmShape.h:422