GemaCoreLib
The GeMA Core library
gmShape.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
24 #ifndef _GEMA_SHAPE_H_
25 #define _GEMA_SHAPE_H_
26 
27 #include "gmCell.h"
28 
29 class GmSpinLock;
30 class GmSimulationData;
31 
38 {
39 public:
40 
42  virtual GmCellType elemType() const = 0;
43 
45  virtual int numFunctions() const = 0;
46 
48  virtual int numNaturalCoord() const = 0;
49 
51  virtual int numCartesianCoord() const = 0;
52 
54  virtual void naturalCoordLimits(int coord, double* min, double* max) const = 0;
55 
64  virtual void nodeNaturalCoord(int node, GmVector& coord) const = 0;
65 
71  virtual void naturalCenter(GmVector& coord) const = 0;
72 
73  void fillNaturalCoordinates(GmMatrix& coord, bool transposed = false) const;
74 
87  virtual bool translateEdgePoint(int edge, const GmVector& srcEdgeCoord, GmVector& elementCoord) const = 0;
88 
102  virtual bool translateFacePoint(int face, const GmVector& srcFaceCoord, GmVector& elementCoord) const = 0;
103 
130  virtual bool cartesianToNatural(const GmVector& coord, const GmMatrix& X, GmVector& ncoord, bool* inside) const
131  {
132  return cartesianToNatural(coord, X, _preferGeometricConfig, _gradTolConfig, _maxGradIterConfig,
133  _outsideNatTolConfig, ncoord, inside);
134  }
135 
169  virtual bool cartesianToNatural(const GmVector& coord, const GmMatrix& X, bool preferGeometric, double gradTol,
170  int maxGradIter, double outsideNatTol, GmVector& ncoord, bool* inside) const
171  {
172  if(preferGeometric && hasGeometryBasedCartesianToNatural())
173  return geometryBasedCartesianToNatural(coord, X, outsideNatTol, ncoord, inside);
174  else
175  return gradientBasedCartesianToNatural(coord, X, gradTol, maxGradIter, outsideNatTol, ncoord, inside);
176  }
177 
178  virtual void naturalToCartesian(const GmVector& ncoord, const GmMatrix& X, GmVector& coord) const;
179 
193  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const = 0;
194 
223  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const = 0;
224 
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;
227 
246  virtual double shapeCartesianPartialsFromPJ(const GmMatrix& P, const GmMatrix& J, GmMatrix& dN, bool transposed = false) const = 0;
247 
275  virtual void jacobian(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, bool transposed = false) const
276  {
277  GmMatrix P;
278  jacobianAndPartials(ncoord, X, J, P, transposed);
279  }
280 
288  virtual double scaledJacobianDet(const GmMatrix& J) const = 0;
289 
298  virtual void jacobianAndPartials(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed = false) const = 0;
299 
321  virtual double borderScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
322  const GmMatrix& X, bool transposed = false) const = 0;
323 
343  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
344  const GmMatrix& X, bool transposed = false) const = 0;
345 
365  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
366  const GmMatrix& X, bool transposed = false) const = 0;
367 
368  virtual const GmMatrix& gaussExtrapolationMatrix(const GmIntegrationRule* ir) const;
369  virtual bool pointExtrapolationMatrix(const GmMatrix& points, GmMatrix& em) const;
370 
371  double interpolate(const GmVector& nodeValues, const GmVector& ncoord) const;
372 
373  static bool initShapeFunctions();
374  static void setConfigOptions(const GmSimulationData* simData);
375  static const GmShape* shapeFromElementType (GmCellType etype, int P, int Q);
376  static const GmShape* linearShapeFromElementType(GmCellType etype, int P, int Q);
377 
378 protected:
379  virtual void jacIndependentNatCoord(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed) const;
380  virtual void jacDependentNatCoord (const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed) const;
381 
382  bool ginv(const GmMatrix& A, GmMatrix& inv) const;
383 
384 #if defined ENABLE_TESTS
385  friend void shapeMappingTests();
386 #endif
387  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
388  int maxIter, double natTol, GmVector& ncoord, bool* inside) const = 0;
389 
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;
398 
404  virtual bool hasGeometryBasedCartesianToNatural() const { return false; }
405 
410  virtual bool geometryBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double natTol,
411  GmVector& ncoord, bool* inside) const
412  {
413  Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
414  return false;
415  }
416 
417  bool checkAndClampNaturalCoordinate(double& ncoord, double min, double max, double natTol) const;
418  bool checkAndClampBarycentricNaturalCoordinates(GmVector& ncoord, double natTol) const;
419 
420 private:
421  static const void* _shapePointersList[GM_NUM_CELL_TYPES];
423 
425  static double _gradTolConfig;
426  static int _maxGradIterConfig;
427  static double _outsideNatTolConfig;
428 
431 };
432 
433 
434 #endif
435 
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