GemaCoreLib
The GeMA Core library
gmTriShape.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_TRI_SHAPE_H_
26 #define _GEMA_TRI_SHAPE_H_
27 
28 #include "gmShape.h"
29 
30 
31 //--------------------------------------------------------------------------------------------------------------------------
32 // Tri Elements
33 //--------------------------------------------------------------------------------------------------------------------------
42 {
43 public:
44  // See comments on the base class
45  virtual int numNaturalCoord() const { return 3; }
46 
47  // See comments on the base class
48  virtual int numCartesianCoord() const { return 2; }
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(3); coord.fill(1/3.0); } // Center = (1/3, 1/3, 1/3) for all triangles
55 
56  virtual bool translateEdgePoint(int edge, const GmVector& srcEdgeCoord, GmVector& elementCoord) const;
57 
58  // See comments on the base class. Undefined for surface elements.
59  virtual bool translateFacePoint(int face, const GmVector& srcFaceCoord, GmVector& elementCoord) const
60  {
61  Q_UNUSED(face); Q_UNUSED(srcFaceCoord); Q_UNUSED(elementCoord); return false;
62  }
63 
64  virtual double shapeCartesianPartialsFromPJ(const GmMatrix& P, const GmMatrix& J, GmMatrix& dN, bool transposed = false) const;
65 
66  // See comments on the base class
67  virtual double scaledJacobianDet(const GmMatrix& J) const { return 0.5 * arma::det(J); }
68 
69  // See comments on the base class
70  virtual void jacobianAndPartials(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed = false) const
71  {
72  jacDependentNatCoord(ncoord, X, J, P, transposed);
73  }
74 
75  // See comments on the base class
76  virtual double borderScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
77  const GmMatrix& X, bool transposed = false) const
78  {
79  return edgeScalingFactor(border, borderCoord, elementCoord, X, transposed);
80  }
81 
82  // See comments on the base class. Undefined for 2D elements.
83  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
84  const GmMatrix& X, bool transposed = false) const
85  {
86  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
87  return 0.0;
88  }
89 
90  // fixedEdgeCoord() is made static since it is "reused" by Int3DTFace shape functions
91  static int fixedEdgeCoord(int border);
92 
93 protected:
94 
95  // See comments on the base class. Calls the default algorithm for Tri / Tet elements
96  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
97  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
98  {
99  return gradientBasedCartesianToNaturalTriTet(coord, X, tol, maxIter, natTol, ncoord, inside);
100  }
101 };
102 
112 {
113 public:
114  // See comments on the base class
115  virtual GmCellType elemType() const { return GM_TRI3; }
116 
117  // See comments on the base class
118  virtual int numFunctions() const { return 3; }
119 
120  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
121 
122  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
123  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
124 
125  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
126  const GmMatrix& X, bool transposed = false) const;
127 
128 protected:
129  // See comments on the base class
130  virtual bool hasGeometryBasedCartesianToNatural() const { return true; }
131 
132  virtual bool geometryBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double natTol,
133  GmVector& ncoord, bool* inside) const;
134 };
135 
146 {
147 public:
148  // See comments on the base class
149  virtual GmCellType elemType() const { return GM_TRI6; }
150 
151  // See comments on the base class
152  virtual int numFunctions() const { return 6; }
153 
154  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
155 
156  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
157  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
158 
159  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
160  const GmMatrix& X, bool transposed = false) const;
161 };
162 
172 {
173 public:
174  // See comments on the base class
175  virtual GmCellType elemType() const { return GM_TRI3D3; }
176  virtual int numCartesianCoord() const { return 3; }
177 
178  virtual void jacDependentNatCoord(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed) const;
179 
180 protected:
181  // Implementation of the geometric cartesian to natural method for Tri3 only works with 2D coordinates
182  virtual bool hasGeometryBasedCartesianToNatural() const { return false; }
183 
184  // See comments on the base class. At the moment we don't have a working version for 3D triangles
185  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
186  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
187  {
188  Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
189  return false;
190  }
191 };
192 
203 {
204 public:
205  // See comments on the base class
206  virtual GmCellType elemType() const { return GM_TRI3D6; }
207  virtual int numCartesianCoord() const { return 3; }
208 
209  virtual void jacDependentNatCoord(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed) const;
210 
211 protected:
212  // See comments on the base class. At the moment we don't have a working version for 3D triangles
213  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
214  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
215  {
216  Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
217  return false;
218  }
219 };
220 
221 #endif
222 
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:206
virtual double faceScalingFactor(int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
Returns the scaling factor needed when calculating integrals over element faces. This has the same ef...
Definition: gmTriShape.h:83
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmTriShape.h:182
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 int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmTriShape.h:118
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: gmTriShape.h:70
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmTriShape.h:176
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: gmTriShape.h:76
GmShape specialization for a Triangle with 6 nodes object.
Definition: gmTriShape.h:145
A 3D quadratic triangle with 6 nodes.
Definition: gmCellType.h:44
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 double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmTriShape.h:67
A 3D linear triangle with only 3 nodes.
Definition: gmCellType.h:43
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...
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 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 quadratic triangle with 6 nodes.
Definition: gmCellType.h:42
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmTriShape.h:48
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmTriShape.h:45
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmTriShape.h:152
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmTriShape.h:130
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:175
A linear triangle with only 3 nodes.
Definition: gmCellType.h:41
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmTriShape.h:207
#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 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
GmShape specialization for a 3D Triangle with 6 nodes object.
Definition: gmTriShape.h:202
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 translateFacePoint(int face, const GmVector &srcFaceCoord, GmVector &elementCoord) const
Given a "quad like" or "tri like" face coordinate (-1 to 1 pair for quad faces and 0 to 1 barycentric...
Definition: gmTriShape.h:59
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:115
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmTriShape.h:149
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
GmShape specialization for a Triangle with 3 nodes object.
Definition: gmTriShape.h:111
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 void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmTriShape.h:54
GmShape specialization for Triangles, containing common functions used by concrete Tri specialization...
Definition: gmTriShape.h:41
GmShape specialization for a 3D Triangle with 3 nodes object.
Definition: gmTriShape.h:171
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 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: gmTriShape.h:51