GemaCoreLib
The GeMA Core library
gmQuadShape.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_QUAD_SHAPE_H_
26 #define _GEMA_QUAD_SHAPE_H_
27 
28 #include "gmShape.h"
29 
30 //--------------------------------------------------------------------------------------------------------------------------
31 // Quad Elements
32 //--------------------------------------------------------------------------------------------------------------------------
33 
42 {
43 public:
44  // See comments on the base class
45  virtual int numNaturalCoord() const { return 2; }
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 = -1.0; *max = 1.0; }
52 
53  // See comments on the base class
54  virtual void naturalCenter(GmVector& coord) const { coord.zeros(2); } // Center = (0.0, 0.0) for all quads
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 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  jacIndependentNatCoord(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  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
83  const GmMatrix& X, bool transposed = false) const;
84 
85  // See comments on the base class. Undefined for 2D elements.
86  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
87  const GmMatrix& X, bool transposed = false) const
88  {
89  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
90  return 0.0;
91  }
92 
93  // fixedEdgeCoord() & Value() are made static since they are "reused" by Int3DQFace shape functions
94  static int fixedEdgeCoord(int border);
95  static double fixedEdgeValue(int border);
96 
97 protected:
98 
99  // See comments on the base class. Calls the default algorithm for Quad / Hex elements
100  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
101  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
102  {
103  return gradientBasedCartesianToNaturalQuadHex(coord, X, tol, maxIter, natTol, ncoord, inside);
104  }
105 };
106 
107 
118 {
119 public:
120  // See comments on the base class
121  virtual GmCellType elemType() const { return GM_QUAD4; }
122 
123  // See comments on the base class
124  virtual int numFunctions() const { return 4; }
125 
126  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
127 
128  virtual void shapeValues (const GmVector& ncoord, GmVector& N) const;
129  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
130 
131 protected:
132  // See comments on the base class
133  virtual bool hasGeometryBasedCartesianToNatural() const { return true; }
134 
135  virtual bool geometryBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double natTol,
136  GmVector& ncoord, bool* inside) const;
137 };
138 
149 {
150 public:
151  // See comments on the base class
152  virtual GmCellType elemType() const { return GM_QUAD8; }
153 
154  // See comments on the base class
155  virtual int numFunctions() const { return 8; }
156 
157  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
158 
159  // We previously had (by Erwan) a tentative implementation of a geometric method
160  // for cartesian to natural based on the paper "Distorsion measures ans inverse mapping
161  // for isoparametric 8-node plane finite elements with curved boundaries" - Lautersztajn (1998)
162  // but it didn't worked and was removed
163 
164  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
165  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
166 };
167 
179 {
180 public:
181  // See comments on the base class
182  virtual GmCellType elemType() const { return GM_QUAD9; }
183 
184  // See comments on the base class
185  virtual int numFunctions() const { return 9; }
186 
187  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
188 
189  virtual void shapeValues (const GmVector& ncoord, GmVector& N) const;
190  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
191 };
192 
203 {
204 public:
205  // See comments on the base class
206  virtual GmCellType elemType() const { return GM_QUAD3D4; }
207  virtual int numCartesianCoord() const { return 3; }
208 
209  // See comments on the base class
210  virtual void jacIndependentNatCoord(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed) const;
211 
212 protected:
213  // Implementation of the geometric cartesian to natural method for Quad4 only works with 2D coordinates
214  virtual bool hasGeometryBasedCartesianToNatural() const { return false; }
215 
216  // See comments on the base class. At the moment we don't have a working version for 3D quads
217  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
218  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
219  {
220  Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
221  return false;
222  }
223 
224 };
225 
236 {
237 public:
238  // See comments on the base class
239  virtual GmCellType elemType() const { return GM_QUAD3D8; }
240  virtual int numCartesianCoord() const { return 3; }
241 
242  virtual void jacIndependentNatCoord(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed) const;
243 
244 protected:
245  // See comments on the base class. At the moment we don't have a working version for 3D quads
246  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
247  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
248  {
249  Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
250  return false;
251  }
252 };
253 
254 
255 #endif
256 
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmQuadShape.h:207
A 3D quadratic quadrilateral cell with 8 nodes.
Definition: gmCellType.h:40
GmShape specialization for a 3D Quadrilateral with 4 nodes object.
Definition: gmQuadShape.h:202
GmShape specialization for a Quadrilateral with 4 nodes object.
Definition: gmQuadShape.h:117
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmQuadShape.h:54
A linear quadrilateral cell with only 4 nodes.
Definition: gmCellType.h:36
GmShape specialization for a Quadrilateral with 8 nodes object.
Definition: gmQuadShape.h:148
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: gmQuadShape.h:86
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: gmQuadShape.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...
A quadratic quadrilateral cell with 8 nodes.
Definition: gmCellType.h:37
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 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: gmQuadShape.h:70
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
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmQuadShape.h:45
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmQuadShape.h:214
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmQuadShape.h:48
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmQuadShape.h:121
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmQuadShape.h:155
GmShape specialization for a 3D Quadrilateral with 8 nodes object.
Definition: gmQuadShape.h:235
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: gmQuadShape.h:76
A 3D linear quadrilateral cell with only 4 nodes.
Definition: gmCellType.h:39
#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: gmQuadShape.h:182
A quadratic quadrilateral cell with 9 nodes.
Definition: gmCellType.h:38
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
GmShape specialization for Quadrilaterals, containing common functions used by concrete Quad speciali...
Definition: gmQuadShape.h:41
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: gmQuadShape.h:152
virtual double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmQuadShape.h:67
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmQuadShape.h:133
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: gmQuadShape.h:59
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 numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmQuadShape.h:240
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmQuadShape.h:124
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmQuadShape.h:239
GmShape specialization for a Quadrilateral with 9 nodes object.
Definition: gmQuadShape.h:178
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 int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmQuadShape.h:185
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmQuadShape.h:206
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...