GemaCoreLib
The GeMA Core library
gmTetShape.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 
25 #ifndef _GEMA_TET_SHAPE_H_
26 #define _GEMA_TET_SHAPE_H_
27 
28 #include "gmShape.h"
29 
30 //--------------------------------------------------------------------------------------------------------------------------
31 // Tet Elements
32 //--------------------------------------------------------------------------------------------------------------------------
33 
42 {
43 public:
44  // See comments on the base class
45  virtual int numNaturalCoord() const { return 4; }
46 
47  // See comments on the base class
48  virtual int numCartesianCoord() const { return 3; }
49 
50  // See comments on the base class
51  virtual void naturalCoordLimits(int coord, double* min, double* max) const { Q_UNUSED(coord); *min = 0.0; *max = 1.0; }
52 
53  // See comments on the base class
54  virtual void naturalCenter(GmVector& coord) const { coord.set_size(4); coord.fill(1/4.0); } // Center = (1/4, 1/4, 1/4, 1/4) for all tetrahedrons
55 
56  virtual bool translateEdgePoint(int edge, const GmVector& srcEdgeCoord, GmVector& elementCoord) const;
57  virtual bool translateFacePoint(int face, const GmVector& srcFaceCoord, GmVector& elementCoord) const;
58 
59  virtual double shapeCartesianPartialsFromPJ(const GmMatrix& P, const GmMatrix& J, GmMatrix& dN, bool transposed = false) const;
60 
61  // See comments on the base class
62  virtual double scaledJacobianDet(const GmMatrix& J) const { return (1.0/6.0) * arma::det(J); }
63 
64  // See comments on the base class
65  virtual void jacobianAndPartials(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed = false) const
66  {
67  jacDependentNatCoord(ncoord, X, J, P, transposed);
68  }
69 
70  // See comments on the base class
71  virtual double borderScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
72  const GmMatrix& X, bool transposed = false) const
73  {
74  return faceScalingFactor(border, borderCoord, elementCoord, X, transposed);
75  }
76 
77  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
78  const GmMatrix& X, bool transposed = false) const;
79 
80  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
81  const GmMatrix& X, bool transposed = false) const;
82 
83 protected:
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;
88 
89  // See comments on the base class. Calls the default algorithm for Tri / Tet elements
90  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
91  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
92  {
93  return gradientBasedCartesianToNaturalTriTet(coord, X, tol, maxIter, natTol, ncoord, inside);
94  }
95 };
96 
106 {
107 public:
108  // See comments on the base class
109  virtual GmCellType elemType() const { return GM_TET4; }
110 
111  // See comments on the base class
112  virtual int numFunctions() const { return 4; }
113 
114  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
115 
116  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
117  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
118 
119 protected:
120  // See comments on the base class
121  virtual bool hasGeometryBasedCartesianToNatural() const { return true; }
122 
123  virtual bool geometryBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double natTol,
124  GmVector& ncoord, bool* inside) const;
125 };
126 
127 
137 {
138 public:
139  // See comments on the base class
140  virtual GmCellType elemType() const { return GM_TET10; }
141 
142  // See comments on the base class
143  virtual int numFunctions() const { return 10; }
144 
145  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
146 
147  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
148  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
149 };
150 
151 
152 #endif
153 
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...
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 double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmTetShape.h:62
A tetrahedron with 4 nodes.
Definition: gmCellType.h:52
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: gmTetShape.h:71
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 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 void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmTetShape.h:54
Shape function handling base classe.
Definition: gmShape.h:37
GmShape specialization for Tetrahedrons, containing common functions used by concrete Hex specializat...
Definition: gmTetShape.h:41
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmTetShape.h:143
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 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: gmTetShape.h:65
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 tetrahedron with 10 nodes.
Definition: gmCellType.h:61
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmTetShape.h:48
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...
#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 GmCellType elemType() const
Returns the type of this element.
Definition: gmTetShape.h:140
GmShape specialization for a Tetrahedron with 4 nodes object.
Definition: gmTetShape.h:105
GmShape specialization for a Tetrahedron with 4 nodes object.
Definition: gmTetShape.h:136
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmTetShape.h:112
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 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 hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmTetShape.h:121
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: gmTetShape.h:51
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 numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmTetShape.h:45
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 GmCellType elemType() const
Returns the type of this element.
Definition: gmTetShape.h:109
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...