GemaCoreLib
The GeMA Core library
gmInt2DShape.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_INT2D_SHAPE_H_
26 #define _GEMA_INT2D_SHAPE_H_
27 
28 #include "gmShape.h"
29 
30 //--------------------------------------------------------------------------------------------------------------------------
31 // 2D Interface 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 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(1); }
55 
56  // See comments on the base class.
57  virtual bool translateEdgePoint(int edge, const GmVector& srcEdgeCoord, GmVector& elementCoord) const
58  {
59  Q_UNUSED(edge);
60  elementCoord = srcEdgeCoord; // edge and element coordinates are the same!!!
61  return !(edge % 2); // Returns true for edges 0 & 2 and false for edges 1 & 3, since those should not be used for integration
62  }
63 
64  // See comments on the base class. Undefined for 2D interface elements.
65  virtual bool translateFacePoint(int face, const GmVector& srcFaceCoord, GmVector& elementCoord) const
66  {
67  Q_UNUSED(face); Q_UNUSED(srcFaceCoord); Q_UNUSED(elementCoord); return false;
68  }
69 
70  virtual double shapeCartesianPartialsFromPJ(const GmMatrix& P, const GmMatrix& J, GmMatrix& dN, bool transposed = false) const;
71 
72  // See comments on the base class
73  virtual double scaledJacobianDet(const GmMatrix& J) const
74  {
75  assert(J.n_elem == 2);
76  // Get dx/dxi and dy/dxi derivate. Remember that J might be transposed
77  double dx = J(0, 0);
78  double dy = (J.n_cols == 2) ? J(0, 1) : J(1, 0);
79  return sqrt(dx*dx + dy*dy);
80  }
81 
82  // See comments on the base class
83  virtual void jacobianAndPartials(const GmVector& ncoord, const GmMatrix& MX, GmMatrix& J,
84  GmMatrix& P, bool transposed = false) const;
85 
86  // Undefined
87  virtual double borderScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
88  const GmMatrix& X, bool transposed = false) const
89  {
90  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
91  return 0.0;
92  }
93 
94  // Undefined
95  virtual double edgeScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
96  const GmMatrix& X, bool transposed = false) const
97  {
98  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
99  return 0.0;
100  }
101 
102  // Undefined
103  virtual double faceScalingFactor(int border, const GmVector& borderCoord, const GmVector& elementCoord,
104  const GmMatrix& X, bool transposed = false) const
105  {
106  Q_UNUSED(border); Q_UNUSED(borderCoord); Q_UNUSED(elementCoord); Q_UNUSED(X); Q_UNUSED(transposed);
107  return 0.0;
108  }
109 
110  // See comments on the base class
111  void naturalToCartesian(const GmVector& ncoord, const GmMatrix& X, GmVector& coord) const;
112 
113 protected:
114  // See comments on the base class. At the moment we don't have a working version for interface elements
115  // Can we just reuse the Bar implementation for all elements?
116  virtual bool gradientBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double tol,
117  int maxIter, double natTol, GmVector& ncoord, bool* inside) const
118  {
119  Q_UNUSED(coord); Q_UNUSED(X); Q_UNUSED(tol); Q_UNUSED(maxIter); Q_UNUSED(natTol); Q_UNUSED(ncoord); Q_UNUSED(inside);
120  return false;
121  }
122 };
123 
134 {
135 public:
136  // See comments on the base class
137  virtual GmCellType elemType() const { return GM_INT2DL4; }
138 
139  // Returns the number of shape functions of this cohesive element (equal to the node number at middle plane)
140  virtual int numFunctions() const { return 2; }
141 
142  // See comments on the base class
143  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
144 
145  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
146  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
147 
148 protected:
149  // See comments on the base class
150  virtual bool hasGeometryBasedCartesianToNatural() const { return true; }
151 
152  virtual bool geometryBasedCartesianToNatural(const GmVector& coord, const GmMatrix& X, double natTol,
153  GmVector& ncoord, bool* inside) const;
154 };
155 
167 {
168 public:
169  // See comments on the base class
170  virtual GmCellType elemType() const { return GM_INT2DL6; }
171 };
172 
183 {
184 public:
185  // See comments on the base class
186  virtual GmCellType elemType() const { return GM_INT2DQ6; }
187 
188  // Returns the number of shape functions of this interface element (equal to the node number at middle plane)
189  virtual int numFunctions() const { return 3; }
190 
191  // See comments on the base class
192  virtual void nodeNaturalCoord(int node, GmVector& coord) const;
193 
194  virtual void shapeValues(const GmVector& ncoord, GmVector& N) const;
195  virtual void shapePartials(const GmVector& ncoord, GmMatrix& dN, bool transposed = false) const;
196 };
197 
209 {
210 public:
211  // See comments on the base class
212  virtual GmCellType elemType() const { return GM_INT2DQ8; }
213 };
214 
215 #endif
216 
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmInt2DShape.h:140
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this element type.
Definition: gmInt2DShape.h:45
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: gmInt2DShape.h:65
virtual double scaledJacobianDet(const GmMatrix &J) const
Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differe...
Definition: gmInt2DShape.h:73
virtual void naturalCenter(GmVector &coord) const
Fills the coord vector with the set of natural coordinates for the element center.
Definition: gmInt2DShape.h:54
virtual void naturalToCartesian(const GmVector &ncoord, const GmMatrix &X, GmVector &coord) const
Given a set of natural coordinates 'ncoord' and a matrix with node coordinates 'X',...
Definition: gmShape.cpp:106
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 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 GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:137
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:170
A 2D linear interface element with 4 nodes.
Definition: gmCellType.h:45
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: gmInt2DShape.h:57
GmShape specialization for a 2D quadratic interface element with 8 nodes object 2 nodes (x,...
Definition: gmInt2DShape.h:208
virtual GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:186
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 GmCellType elemType() const
Returns the type of this element.
Definition: gmInt2DShape.h:212
GmShape specialization for a 2D linear interface element with 4 nodes object.
Definition: gmInt2DShape.h:133
GmShape specialization for Cohesive elements, containing common functions used to model fracture proc...
Definition: gmInt2DShape.h:41
virtual int numFunctions() const
Returns the number of shape functions of this element type (equal to the number of nodes)
Definition: gmInt2DShape.h:189
virtual int numCartesianCoord() const
Returns the number of cartesian coordinates expected by this element type.
Definition: gmInt2DShape.h:48
#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
A 2D linear interface element with 6 nodes (2 nodes of pore pressure dof at the middle)
Definition: gmCellType.h:46
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 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: gmInt2DShape.h:87
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: gmInt2DShape.h:95
A 2D quadratic interface element with 6 nodes.
Definition: gmCellType.h:47
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: gmInt2DShape.h:51
virtual bool hasGeometryBasedCartesianToNatural() const
A virtual function that should be replaced to return true if the shape implements a geometry based al...
Definition: gmInt2DShape.h:150
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
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
GmShape specialization for a 2D quadratic interface element with 6 nodes object.
Definition: gmInt2DShape.h:182
Declaration of the GmShape base class.
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: gmInt2DShape.h:103
GmShape specialization for a 2D coupled interface element with 6 nodes object with 4 nodes (x,...
Definition: gmInt2DShape.h:166
A 2D quadratic interface element with 8 nodes (2 nodes of pore pressure dof at the middle)
Definition: gmCellType.h:48