GemaCoreLib
The GeMA Core library
gmPyraShape.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_PYRA_SHAPE_H_
25 #define _GEMA_PYRA_SHAPE_H_
26 
27 #include "gmShape.h"
28 
29 //--------------------------------------------------------------------------------------------------------------------------
30 // Pyra Elements
31 //--------------------------------------------------------------------------------------------------------------------------
32 
41 {
42 public:
43  // See comments on the base class
44  virtual int numNaturalCoord() const { return 3; }
45 
46  // See comments on the base class
47  virtual int numCartesianCoord() const { return 3; }
48 
49  // See comments on the base class
50  virtual void naturalCoordLimits(int coord, double* min, double* max) const { Q_UNUSED(coord); *min = -1.0; *max = 1.0; }
51 
52  // See comments on the base class
53  virtual void naturalCenter(GmVector& coord) const { coord.zeros(3); } // Center = (0.0, 0.0, 0.0) for all pyramids
54 
55  virtual bool translateEdgePoint(int edge, const GmVector& srcEdgeCoord, GmVector& elementCoord) const;
56  virtual bool translateFacePoint(int face, const GmVector& srcFaceCoord, GmVector& elementCoord) const;
57 
58  virtual double shapeCartesianPartialsFromPJ(const GmMatrix& P, const GmMatrix& J, GmMatrix& dN, bool transposed = false) const;
59 
60  // See comments on the base class
61  virtual double scaledJacobianDet(const GmMatrix& J) const { return arma::det(J); }
62 
63  // See comments on the base class
64  virtual void jacobianAndPartials(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed = false) const
65  {
66  jacIndependentNatCoord(ncoord, X, J, P, transposed);
67  }
68 
69  // See comments on the base class
70  virtual double borderScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
71  const GmMatrix& X, bool transposed = false) const
72  {
73  return faceScalingFactor(border, borderCoord, elementCoord, X, transposed);
74  }
75 
76  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
77  const GmMatrix& X, bool transposed = false) const;
78 
79  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
80  const GmMatrix& X, bool transposed = false) const;
81 
82 protected:
83  const int* fixedEdgeCoord(int border) const;
84  const double* fixedEdgeValue(int border) const;
85  int fixedFaceCoord(int border) const;
86  double fixedFaceValue(int border) const;
87 
88  // See comments on the base class. Calls the default algorithm for Tri / Tet elements
89  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
90  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
91  {
92  return gradientBasedCartesianToNaturalTriTet(coord, X, tol, maxIter, natTol, ncoord, inside);
93  }
94 };
95 
105 {
106 public:
107  // See comments on the base class
108  virtual GmCellType elemType() const { return GM_PYRA5; }
109 
110  // See comments on the base class
111  virtual int numFunctions() const { return 5; }
112 
113  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
114 
115  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
116  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
117 
118 };
119 
120 
130 {
131 public:
132  // See comments on the base class
133  virtual GmCellType elemType() const { return GM_PYRA13; }
134 
135  // See comments on the base class
136  virtual int numFunctions() const { return 13; }
137 
138  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
139 
140  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
141  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
142 };
143 
144 
145 #endif
146 
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...
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: gmPyraShape.h:50
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: gmPyraShape.h:70
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmPyraShape.h:111
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 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...
GmShape specialization for a Tetrahedron with 4 nodes object.
Definition: gmPyraShape.h:129
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 double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmPyraShape.h:61
Shape function handling base classe.
Definition: gmShape.h:37
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmPyraShape.h:53
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmPyraShape.h:47
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
GmShape specialization for Pyramid, containing common functions used by concrete Pyra specializations...
Definition: gmPyraShape.h:40
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...
A pyramid with 13 nodes - quadratic.
Definition: gmCellType.h:65
#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
GmCellType
Mesh Cell types. Don't change type orders or add types without reading comments below.
Definition: gmCellType.h:30
GmShape specialization for a Pyramid with 5 nodes object.
Definition: gmPyraShape.h:104
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 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: gmPyraShape.h:64
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmPyraShape.h:108
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmPyraShape.h:44
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 numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmPyraShape.h:136
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.
A pyramid with 5 nodes - linear.
Definition: gmCellType.h:64
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmPyraShape.h:133
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...