GemaCoreLib
The GeMA Core library
gmBarShape.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_BAR_SHAPE_H_
26 #define _GEMA_BAR_SHAPE_H_
27 
28 #include "gmShape.h"
29 
30 //--------------------------------------------------------------------------------------------------------------------------
31 // Bar Elements
32 //--------------------------------------------------------------------------------------------------------------------------
33 
42 {
43 public:
44  // See comments on the base class
45  virtual int numNaturalCoord() const { return 1; }
46 
47  // See comments on the base class
48  virtual void naturalCoordLimits(int coord, double* min, double* max) const { Q_UNUSED(coord); *min = -1.0; *max = 1.0; }
49 
50  // See comments on the base class
51  virtual void naturalCenter(GmVector& coord) const { coord.zeros(1); } // Center = 0.0 for all bars
52 
53  // See comments on the base class. Undefined for line elements.
54  virtual bool translateEdgePoint(int edge, const GmVector& srcEdgeCoord, GmVector& elementCoord) const
55  {
56  Q_UNUSED(edge); Q_UNUSED(srcEdgeCoord); Q_UNUSED(elementCoord); return false;
57  }
58 
59  // See comments on the base class. Undefined for line elements.
60  virtual bool translateFacePoint(int face, const GmVector& srcFaceCoord, GmVector& elementCoord) const
61  {
62  Q_UNUSED(face); Q_UNUSED(srcFaceCoord); Q_UNUSED(elementCoord); return false;
63  }
64 
65  // See comments on the base class
66  virtual double shapeCartesianPartialsFromPJ(const GmMatrix& P, const GmMatrix& J, GmMatrix& dN, bool transposed = false) const;
67 
68  // See comments on the base class
69  virtual double scaledJacobianDet(const GmMatrix& J) const {
70  assert(J.n_elem <= 3);
71  // return sqrt(dx*dx + dy*dy) for 2D
72  // sqrt(dx*dx + dy*dy + dz*dz) for 3D
73  return arma::norm(J);
74  }
75 
76  // See comments on the base class
77  virtual void jacobianAndPartials(const GmVector& ncoord, const GmMatrix& X, GmMatrix& J, GmMatrix& P, bool transposed = false) const;
78 
79  // See comments on the base class. Undefined for line elements.
80  virtual double borderScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
81  const GmMatrix& X, bool transposed = false) const
82  {
83  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
84  return 0.0;
85  }
86 
87  // See comments on the base class. Undefined for line elements.
88  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
89  const GmMatrix& X, bool transposed = false) const
90  {
91  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
92  return 0.0;
93  }
94 
95  // See comments on the base class. Undefined for line elements.
96  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
97  const GmMatrix& X, bool transposed = false) const
98  {
99  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
100  return 0.0;
101  }
102 
103 protected:
104  // See comments on the base class. Calls the default algorithm for Bar elements
105  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
106  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
107  {
108  return gradientBasedCartesianToNaturalBar(coord, X, tol, maxIter, natTol, ncoord, inside);
109  }
110 };
111 
120 {
121 public:
122  // See comments on the base class
123  virtual GmCellType elemType() const { return GM_BAR2; }
124 
125  // See comments on the base class
126  virtual int numFunctions() const { return 2; }
127 
128  virtual int numCartesianCoord() const { return 2; }
129 
130  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
131 
132  virtual void shapeValues (const GmVector& ncoord, GmVector& N) const;
133  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
134 
135 protected:
136  // See comments on the base class
137  virtual bool hasGeometryBasedCartesianToNatural() const { return true; }
138 
139  virtual bool geometryBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double natTol,
140  GmVector& ncoord, bool* inside) const;
141 };
142 
153 {
154 public:
155  // See comments on the base class
156  virtual GmCellType elemType() const { return GM_BAR3; }
157 
158  // See comments on the base class
159  virtual int numFunctions() const { return 3; }
160 
161  virtual int numCartesianCoord() const { return 2; }
162 
163  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
164 
165  virtual void shapeValues (const GmVector& ncoord, GmVector& N) const;
166  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
167 };
168 
177 {
178 public:
179  // See comments on the base class
180  virtual GmCellType elemType() const { return GM_BAR3D2; }
181  virtual int numCartesianCoord() const { return 3; }
182 
183 protected:
184  // See comments on the base class. At the moment we don't have a working version for 3D bars
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_BAR3D3; }
207  virtual int numCartesianCoord() const { return 3; }
208 
209 protected:
210  // See comments on the base class. At the moment we don't have a working version for 3D bars
211  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
212  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
213  {
214  Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
215  return false;
216  }
217 };
218 #endif
219 
GmShape specialization for a 3D Bar with 2 nodes object.
Definition: gmBarShape.h:176
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmBarShape.h:123
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmBarShape.h:128
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: gmBarShape.h:60
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmBarShape.h:51
GmShape specialization for a 2D Bar with 3 nodes object.
Definition: gmBarShape.h:152
virtual bool translateEdgePoint(int edge, const GmVector &srcEdgeCoord, GmVector &elementCoord) const
Given a "bar like" edge coordinate (from -1 to 1), at the given edge, fills elementCoord with the equ...
Definition: gmBarShape.h:54
virtual void jacobianAndPartials(const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed=false) const =0
This function does the same calculations as the jacobian() call, but also filling the extra parameter...
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmBarShape.h:180
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 int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmBarShape.h:159
A bar cell 2D with 3 nodes.
Definition: gmCellType.h:33
A bar cell 2D with 2 nodes.
Definition: gmCellType.h:32
Shape function handling base classe.
Definition: gmShape.h:37
GmShape specialization for a 3D Bar with 3 nodes object.
Definition: gmBarShape.h:202
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmBarShape.h:206
A bar cell 3D with 3 nodes.
Definition: gmCellType.h:35
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: gmBarShape.h:96
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 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: gmBarShape.h:48
bool gradientBasedCartesianToNaturalBar(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 Bar elements....
Definition: gmShape.cpp:297
GmShape specialization for Bars, containing common functions used by concrete Bar specializations,...
Definition: gmBarShape.h:41
#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 double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmBarShape.h:69
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmBarShape.h:137
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 edgeScalingFactor(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 edges. This has the same ef...
Definition: gmBarShape.h:88
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...
A bar cell 3D with 2 nodes.
Definition: gmCellType.h:34
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmBarShape.h:45
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmBarShape.h:156
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: gmBarShape.h:80
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: gmBarShape.h:126
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmBarShape.h:181
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmBarShape.h:161
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmBarShape.h:207
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.
GmShape specialization for a 2D Bar with 2 nodes object.
Definition: gmBarShape.h:119