GemaCoreLib
The GeMA Core library
gmInt2dCellGeometryInfo.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_INT2D_CELL_GEOMETRY_INFO_H_
25 #define _GEMA_INT2D_CELL_GEOMETRY_INFO_H_
26 
27 #include "gmCellGeometryInfo.h"
29 
30 #include "gmGeometryUtils.h"
31 
32 #include "gmBarIntegrationRule.h"
33 #include "gmQuadIntegrationRule.h"
34 #include "gmInt2DShape.h"
35 
36 //---------------------------------------------------------
37 // Int2D integration rules
38 //---------------------------------------------------------
39 
43 
51 
52 //---------------------------------------------------------
53 // Int2DL4 & Int2DL6 rule defaults
54 //---------------------------------------------------------
55 
60 
65 
66 //---------------------------------------------------------
67 // Int2DQ6 & Int2DQ8 rule defaults
68 //---------------------------------------------------------
69 
74 
79 
80 
81 //---------------------------------------------------------
82 // GmInt2DL4CellGeometryInfo
83 //---------------------------------------------------------
84 
87  : public GmCellGeometryInfoSurfaceElement<GmInt2DIntegrationRuleSet, GmLinearInt2DIntegrationRuleSetDefaults,
88  GmInt2DEdgeIntegrationRuleSet, GmLinearInt2DEdgeIntegrationRuleSetDefaults>
89 {
90 public:
91  static const GmInt2DL4CellGeometryInfo* instance();
92 
93  // Shape factory. See comments on the base class.
94  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmInt2DL4Shape; }
95 
97  virtual double dimension(const GmMatrix& X) const
98  {
99  return GmGeometryUtils::edgeLength(0.5*(X.col(0)+X.col(3)), 0.5*(X.col(1)+X.col(2)));
100  }
101 
103  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const
104  {
105  GmGeometryUtils::edgeCentroid(0.5*(X.col(0)+X.col(3)), 0.5*(X.col(1)+X.col(2)), coord);
106  }
107 
108  // --------- Capabilities API ---------
109 
110  virtual bool isValid(const GmMatrix& X, double tol) const { Q_UNUSED(X); Q_UNUSED(tol); return GmGeometryUtils::polygon2DIsCCW(X); }
111 
112 protected:
117 };
118 
119 //---------------------------------------------------------
120 // GmInt2DL6CellGeometryInfo
121 //---------------------------------------------------------
122 
125 {
126 public:
127  static const GmInt2DL6CellGeometryInfo* instance();
128 
129  // Shape factory. See comments on the base class.
130  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmInt2DL6Shape; }
131 
132  // --------- Capabilities API ---------
133 
134  virtual bool isValid(const GmMatrix& X, double tol) const {
135  Q_UNUSED(tol);
136  //reorder nodes
137  constexpr int order[] = { 0,1,5,2,3,4 };
138  GmMatrix Xr = GmGeometryUtils::shuffle(X, order, 6);
140  }
141 
142 private:
145 };
146 
147 //---------------------------------------------------------
148 // GmInt2DQ6CellGeometryInfo
149 //---------------------------------------------------------
150 
153  : public GmCellGeometryInfoSurfaceElement<GmInt2DIntegrationRuleSet, GmQuadraticInt2DIntegrationRuleSetDefaults,
154  GmInt2DEdgeIntegrationRuleSet, GmQuadraticInt2DEdgeIntegrationRuleSetDefaults>
155 {
156 public:
157  static const GmInt2DQ6CellGeometryInfo* instance();
158 
159  // Shape factory. See comments on the base class.
160  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmInt2DQ6Shape; }
161 
163  virtual double dimension(const GmMatrix& X) const
164  {
165  GmMatrix Xm(2, 3);
166  Xm.col(0) = 0.5*(X.col(0)+X.col(3)); // First node
167  Xm.col(1) = 0.5*(X.col(1)+X.col(2)); // Last node
168  Xm.col(2) = 0.5*(X.col(4)+X.col(5)); // Mid node
170  }
171 
173  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const
174  {
175  GmMatrix Xm(2, 3);
176  Xm.col(0) = 0.5*(X.col(0)+X.col(3)); // First node
177  Xm.col(1) = 0.5*(X.col(1)+X.col(2)); // Last node
178  Xm.col(2) = 0.5*(X.col(4)+X.col(5)); // Mid node
180  }
181 
182  // --------- Capabilities API ---------
183 
184  virtual bool isValid(const GmMatrix& X, double tol) const {
185  Q_UNUSED(tol);
186  //reorder nodes
187  constexpr int order[] = {0,4,1,2,5,3};
188  GmMatrix Xr = GmGeometryUtils::shuffle(X, order, 6);
190  }
191 
192 protected:
197 };
198 
199 
200 //---------------------------------------------------------
201 // GmInt2DQ8CellGeometryInfo
202 //---------------------------------------------------------
203 
206 {
207 public:
208  static const GmInt2DQ8CellGeometryInfo* instance();
209 
210  // Shape factory. See comments on the base class.
211  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmInt2DQ8Shape; }
212 
213  // --------- Capabilities API ---------
214 
215  virtual bool isValid(const GmMatrix& X, double tol) const {
216  Q_UNUSED(tol);
217  //reorder nodes
218  constexpr int order[] = {0,4,1,7,2,5,3,6};
219  GmMatrix Xr = GmGeometryUtils::shuffle(X, order, 8);
221  }
222 
223 private:
226 };
227 
228 #endif
229 
GmCellGeometryIntegrationRuleSet< GmQuadGaussEdgeIntegrationRule, GmQuadLobattoEdgeIntegrationRule, GmCellGeometryNullIntegrationRule, GmCellGeometryNullIntegrationRule > GmInt2DEdgeIntegrationRuleSet
The set of possible EDGE integration rules, per integration rule type, for the family of 2D interface...
Definition: gmInt2dCellGeometryInfo.h:50
GmInt2DQ6CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Protected constructor. Only a single Int2DQ6 geometry info object is ever necessary,...
Definition: gmInt2dCellGeometryInfo.h:194
1D Lobatto Integration rules (for bar elements)
Definition: gmBarIntegrationRule.h:83
virtual bool isValid(const GmMatrix &X, double tol) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_VALID c...
Definition: gmInt2dCellGeometryInfo.h:134
Lobatto integration rule for Quad element borders.
Definition: gmQuadIntegrationRule.h:144
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 2, 2, 2, 2 > GmLinearInt2DIntegrationRuleSetDefaults
The set of default rules for a linear 2D interface ELEMENT. Template parameters 2 to 5 are the defaul...
Definition: gmInt2dCellGeometryInfo.h:59
void centroidByIntegration(const GmMatrix &X, GmVector &coord, const GmIntegrationRule *ir=NULL) const
Given a cell geometry defined by matrix X (with node coordinates organized by column) performs a nume...
Definition: gmCellGeometryInfo.cpp:114
Declaration of the 2D shapes for interface elements inheriting from GmShape class Previosuly part of ...
1D Newton cotes Integration rules (for bar elements)
Definition: gmBarIntegrationRule.h:59
Gauss rules.
Definition: gmIntegrationRule.h:75
virtual GmShape * shapeInstance(int P, int Q) const
Shape function factory. Should return a NEW instance of the shape function object for this type....
Definition: gmInt2dCellGeometryInfo.h:211
virtual bool isValid(const GmMatrix &X, double tol) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_VALID c...
Definition: gmInt2dCellGeometryInfo.h:184
GmMatrix shuffle(const GmMatrix &X, const int order[], int size)
Returns quad8 vertices in circular order.
Definition: gmGeometryUtils.cpp:612
A type used to represent a 'unexistant" integration rule. Should be used as template parameter for Gm...
Definition: gmCellGeometryIntegrationRuleSet.h:35
The Int2DL6 implementation.
Definition: gmInt2dCellGeometryInfo.h:124
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 3, 3, 3, 3 > GmQuadraticInt2DIntegrationRuleSetDefaults
The set of default rules for a quadratic 2D interface ELEMENT. Template parameters 2 to 5 are the def...
Definition: gmInt2dCellGeometryInfo.h:73
A helper class used to store the default rules for each rule type. Template parameters are the defaul...
Definition: gmCellGeometryIntegrationRuleSet.h:43
virtual GmShape * shapeInstance(int P, int Q) const
Shape function factory. Should return a NEW instance of the shape function object for this type....
Definition: gmInt2dCellGeometryInfo.h:94
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 2, 2, -1, -1 > GmLinearInt2DEdgeIntegrationRuleSetDefaults
The set of default rules for a linear 2D interface EDGE. Template parameters 2 to 5 are the default r...
Definition: gmInt2dCellGeometryInfo.h:64
virtual bool isValid(const GmMatrix &X, double tol) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_VALID c...
Definition: gmInt2dCellGeometryInfo.h:110
A bar cell 2D with 3 nodes.
Definition: gmCellType.h:33
GmShape specialization for a 2D quadratic interface element with 8 nodes object 2 nodes (x,...
Definition: gmInt2DShape.h:208
void edgeCentroid(const GmVector &p1, const GmVector &p2, GmVector &coord)
Returns the cartesian coordinate of the edge centroid.
Definition: gmGeometryUtils.cpp:283
Shape function handling base classe.
Definition: gmShape.h:37
virtual double dimension(const GmMatrix &X) const
Returns the Int2DQ6 length at the "mid edge".
Definition: gmInt2dCellGeometryInfo.h:163
GmCellGeometryIntegrationRuleSet< GmLineGaussIntegrationRule, GmLineLobattoIntegrationRule, GmLineNewtonIntegrationRule, GmLineNewtonIntegrationRule > GmInt2DIntegrationRuleSet
The set of possible ELEMENT integration rules, per integration rule type, for the family of 2D interf...
Definition: gmInt2dCellGeometryInfo.h:42
Declaration of the GmCellGeometryIntegrationRuleSet class & friends.
GmInt2DL6CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Private constructor. Only a single Int2DL6 geometry info object is ever necessary,...
Definition: gmInt2dCellGeometryInfo.h:144
static const GmCellGeometryInfo * geometryInfo(GmCellType type)
Returns the geometry info objct associated with type. NOT FOR GENERAL USE.
Definition: gmCellGeometry.h:522
A helper class used to store the class types for each rule type. Template parameters are the implemen...
Definition: gmCellGeometryIntegrationRuleSet.h:196
The Int2DL4 implementation.
Definition: gmInt2dCellGeometryInfo.h:86
GmShape specialization for a 2D linear interface element with 4 nodes object.
Definition: gmInt2DShape.h:133
bool polygon2DIsCCW(const GmMatrix &X)
Returns true if the polygon delimited by the cell vertices is oriented CCW.
Definition: gmGeometryUtils.cpp:527
virtual double dimension(const GmMatrix &X) const
Returns the Int2DL4 length at the "mid edge".
Definition: gmInt2dCellGeometryInfo.h:97
Utilitary functions for working with geometry.
An auxiliary class that can be used as base for surface elements. Implements the needed integration r...
Definition: gmCellGeometryInfo.h:348
double dimensionByIntegration(const GmMatrix &X, const GmIntegrationRule *ir=NULL) const
Given a cell geometry defined by matrix X (with node coordinates organized by column) performs a nume...
Definition: gmCellGeometryInfo.cpp:76
virtual void centroidCartesian(const GmMatrix &X, GmVector &coord) const
Returns the Int2DQ6 centroid as a Bar3 at the "mid edge".
Definition: gmInt2dCellGeometryInfo.h:173
virtual GmShape * shapeInstance(int P, int Q) const
Shape function factory. Should return a NEW instance of the shape function object for this type....
Definition: gmInt2dCellGeometryInfo.h:160
#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
Declaration of the GmQuadXxxxIntegrationRule family of classes, including border rules....
static const GmInt2DQ6CellGeometryInfo * instance()
Instance function for the Int2DQ6 single geometry object.
Definition: gmInt2dCellGeometryInfo.cpp:150
GmInt2DQ8CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Private constructor. Only a single Int2DQ8 geometry info object is ever necessary,...
Definition: gmInt2dCellGeometryInfo.h:225
double edgeLength(const GmVector &p1, const GmVector &p2)
Calculates the length of the edge determined by points p1 and p2.
Definition: gmGeometryUtils.cpp:157
virtual GmShape * shapeInstance(int P, int Q) const
Shape function factory. Should return a NEW instance of the shape function object for this type....
Definition: gmInt2dCellGeometryInfo.h:130
The Int2DQ8 implementation.
Definition: gmInt2dCellGeometryInfo.h:205
virtual void centroidCartesian(const GmMatrix &X, GmVector &coord) const
Returns the Int2DL4 centroid as a Bar2 at the "mid edge".
Definition: gmInt2dCellGeometryInfo.h:103
static const GmInt2DL4CellGeometryInfo * instance()
Instance function for the Int2DL4 single geometry object.
Definition: gmInt2dCellGeometryInfo.cpp:34
Declaration of the GmCellGeometryInfo base class.
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 3, 3, -1, -1 > GmQuadraticInt2DEdgeIntegrationRuleSetDefaults
The set of default rules for a quadratic 2D interface EDGE. Template parameters 2 to 5 are the defaul...
Definition: gmInt2dCellGeometryInfo.h:78
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
The Int2DQ6 implementation.
Definition: gmInt2dCellGeometryInfo.h:152
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
virtual bool isValid(const GmMatrix &X, double tol) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_VALID c...
Definition: gmInt2dCellGeometryInfo.h:215
GmInt2DL4CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Protected constructor. Only a single Int2DL4 geometry info object is ever necessary,...
Definition: gmInt2dCellGeometryInfo.h:114
1D Gauss Integration rules (for bar elements)
Definition: gmBarIntegrationRule.h:36
Declaration of the GmLineXxxxIntegrationRule family of classes Content previously in gmIntegrationRul...
Plane structure storing the full set of geometric metadata for a cell type.
Definition: gmCellGeometryInfo.h:57
GmShape specialization for a 2D coupled interface element with 6 nodes object with 4 nodes (x,...
Definition: gmInt2DShape.h:166
Gauss integration rule for Quad element borders.
Definition: gmQuadIntegrationRule.h:130