GemaCoreLib
The GeMA Core library
gmQuadCellGeometryInfo.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_QUAD_CELL_GEOMETRY_INFO_H_
25 #define _GEMA_QUAD_CELL_GEOMETRY_INFO_H_
26 
27 #include "gmCellGeometryInfo.h"
29 
30 #include "gmGeometryUtils.h"
31 
32 #include "gmQuadIntegrationRule.h"
33 #include "gmQuadShape.h"
34 
35 //---------------------------------------------------------
36 // Quad integration rules
37 //---------------------------------------------------------
38 
42 
46 
47 
48 //---------------------------------------------------------
49 // Quad4 & Quad3D4 rule defaults
50 //---------------------------------------------------------
51 
56 
61 
62 //---------------------------------------------------------
63 // Quad8, Quad3D8 & Quad9 rule defaults
64 //---------------------------------------------------------
65 
70 
75 
76 
77 //---------------------------------------------------------
78 // GmQuad4CellGeometryInfo
79 //---------------------------------------------------------
80 
83  : public GmCellGeometryInfoSurfaceElement<GmQuadIntegrationRuleSet, GmLinearQuadIntegrationRuleSetDefaults,
84  GmQuadEdgeIntegrationRuleSet, GmLinearQuadEdgeIntegrationRuleSetDefaults>
85 {
86 public:
87  static const GmQuad4CellGeometryInfo* instance();
88 
89  // Shape factory. See comments on the base class.
90  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmQuad4Shape; }
91 
93  virtual double dimension(const GmMatrix& X) const
94  {
95  return GmGeometryUtils::quadArea(X.col(0), X.col(1), X.col(2), X.col(3));
96  }
97 
99  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const
100  {
101  GmGeometryUtils::quadCentroid(X.col(0), X.col(1), X.col(2), X.col(3), coord);
102  }
103 
104  // --------- Capabilities API ---------
105 
106  // See comments on the base class
107  virtual bool isValid(const GmMatrix& X, double tol) const
108  {
109  Q_UNUSED(tol);
111  }
112 
113  // See comments on the base class
114  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { return GmGeometryUtils::pointInConvexPolygon(X, coord); }
115 
116 protected:
121 };
122 
123 //---------------------------------------------------------
124 // GmQuad3D4CellGeometryInfo
125 //---------------------------------------------------------
126 
129 {
130 public:
131  static const GmQuad3D4CellGeometryInfo* instance();
132 
133  // Shape factory. See comments on the base class.
134  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmQuad3D4Shape; }
135 
136  // --------- Capabilities API ---------
137 
138  virtual bool isValid(const GmMatrix& X, double tol) const {return GmGeometryUtils::polygon3DIsPlanar(X, tol); }
139  // contains() capabilities are currently not supported in 3D, so the implementation
140  // from Quad4CellGeometryInfo must be overridden. Quality is ok.
141  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { Q_UNUSED(X); Q_UNUSED(coord); return false; }
142 
143 private:
146 };
147 
148 
149 //---------------------------------------------------------
150 // GmQuad8CellGeometryInfo
151 //---------------------------------------------------------
152 
153 //vertices in circular order
154 static int quad8order[] = { 0,4,1,5,2,6,3,7 };
155 
158  : public GmCellGeometryInfoSurfaceElement<GmQuadIntegrationRuleSet, GmQuadraticQuadIntegrationRuleSetDefaults,
159  GmQuadEdgeIntegrationRuleSet, GmQuadraticQuadEdgeIntegrationRuleSetDefaults>
160 {
161 public:
162  static const GmQuad8CellGeometryInfo* instance();
163 
164  // Shape factory. See comments on the base class.
165  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmQuad8Shape; }
166 
168  virtual double dimension(const GmMatrix& X) const { return dimensionByIntegration(X); }
169 
171  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const { centroidByIntegration(X, coord); }
172 
173  // --------- Capabilities API ---------
174 
175  virtual bool isValid(const GmMatrix & X, double tol) const {
176  Q_UNUSED(tol);
177  //reorder nodes
178  GmMatrix Xr = GmGeometryUtils::shuffle(X, quad8order, 8);
180  }
181 
182  // See comments on the base class
183  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { return GmGeometryUtils::pointInTesselatedPolygon_Q8(X, coord); }
184 
185 protected:
190 };
191 
192 //---------------------------------------------------------
193 // GmQuad3D8CellGeometryInfo
194 //---------------------------------------------------------
195 
198 {
199 public:
200  static const GmQuad3D8CellGeometryInfo* instance();
201 
202  // Shape factory. See comments on the base class.
203  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmQuad3D8Shape; }
204 
205  // --------- Capabilities API ---------
206 
207  virtual bool isValid(const GmMatrix& X, double tol) const {
208  return GmGeometryUtils::polygon3DIsPlanar(GmGeometryUtils::shuffle(X, quad8order, 8), tol);
209  }
210  // No capabilities are currently supported in 3D, so the contains() implementation
211  // from Quad8CellGeometryInfo must be overridden.
212  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { Q_UNUSED(X); Q_UNUSED(coord); return false; }
213 
214 private:
217 };
218 
219 //---------------------------------------------------------
220 // GmQuad9CellGeometryInfo
221 //---------------------------------------------------------
222 
225  : public GmCellGeometryInfoSurfaceElement<GmQuadIntegrationRuleSet, GmQuadraticQuadIntegrationRuleSetDefaults,
226  GmQuadEdgeIntegrationRuleSet, GmQuadraticQuadEdgeIntegrationRuleSetDefaults>
227 {
228 public:
229  static const GmQuad9CellGeometryInfo* instance();
230 
231  // Shape factory. See comments on the base class.
232  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmQuad9Shape; }
233 
235  virtual double dimension(const GmMatrix& X) const { return dimensionByIntegration(X); }
236 
238  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const { centroidByIntegration(X, coord); }
239 
240  // --------- Capabilities API ---------
241  virtual bool isValid(const GmMatrix & X, double tol) const {
242  Q_UNUSED(tol);
243  //reorder nodes, ignores node 9 (center)
244  GmMatrix Xr = GmGeometryUtils::shuffle(X, quad8order, 8);
246  }
247 private:
252 };
253 
254 
255 #endif
256 
bool pointInConvexPolygon(const GmMatrix &X, const GmVector &p)
Returns true if the point p is inside the polygon delimited by the cell vertices (convex hull),...
Definition: gmGeometryUtils.cpp:645
bool polygon3DIsPlanar(const GmMatrix &X, double tol)
Returns true if the polygon delimited by the cell vertices is planar.
Definition: gmGeometryUtils.cpp:587
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: gmQuadCellGeometryInfo.h:165
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
Lobatto integration rule for Quad element borders.
Definition: gmQuadIntegrationRule.h:144
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
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: gmQuadCellGeometryInfo.h:90
The Quad8 implementation.
Definition: gmQuadCellGeometryInfo.h:157
Gauss rules.
Definition: gmIntegrationRule.h:75
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: gmQuadCellGeometryInfo.h:241
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: gmQuadCellGeometryInfo.h:134
virtual void centroidCartesian(const GmMatrix &X, GmVector &coord) const
Returns the Quad8 centroid by numeric integration.
Definition: gmQuadCellGeometryInfo.h:171
GmCellGeometryIntegrationRuleSet< GmQuadGaussEdgeIntegrationRule, GmQuadLobattoEdgeIntegrationRule, GmCellGeometryNullIntegrationRule, GmCellGeometryNullIntegrationRule > GmQuadEdgeIntegrationRuleSet
The set of possible EDGE integration rules, per integration rule type, for the family of Quad element...
Definition: gmQuadCellGeometryInfo.h:45
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
double quadArea(const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4)
Calculates the area of a quad defined by points p1 to p4 (2D or 3D), defined in CCW order.
Definition: gmGeometryUtils.cpp:196
GmShape specialization for a Quadrilateral with 8 nodes object.
Definition: gmQuadShape.h:148
The Quad3D4 implementation.
Definition: gmQuadCellGeometryInfo.h:128
GmQuad3D4CellGeometryInfo(const GmQuad4CellGeometryInfo &quad4)
Private constructor. Only a single Quad3D4 geometry info object is ever necessary,...
Definition: gmQuadCellGeometryInfo.h:145
virtual void centroidCartesian(const GmMatrix &X, GmVector &coord) const
Returns the Quad9 centroid by numeric integration.
Definition: gmQuadCellGeometryInfo.h:238
2D Isotropic / anisotropic Lobatto integration rules for Quad elements
Definition: gmQuadIntegrationRule.h:88
A helper class used to store the default rules for each rule type. Template parameters are the defaul...
Definition: gmCellGeometryIntegrationRuleSet.h:43
virtual void centroidCartesian(const GmMatrix &X, GmVector &coord) const
Returns the Quad4 centroid by decomposing it into two triangles and weighting each coordinate by the ...
Definition: gmQuadCellGeometryInfo.h:99
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: gmQuadCellGeometryInfo.h:175
2D Isotropic / anisotropic Gauss integration rules for Quad elements
Definition: gmQuadIntegrationRule.h:42
virtual bool contains(const GmMatrix &X, const GmVector &coord) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_CONTAIN...
Definition: gmQuadCellGeometryInfo.h:183
GmQuad9CellGeometryInfo(const GmQuad8CellGeometryInfo &quad8)
Private constructor. Only a single Quad9 geometry info object is ever necessary, created by GmQuadCel...
Definition: gmQuadCellGeometryInfo.h:249
Shape function handling base classe.
Definition: gmShape.h:37
Declaration of the GmCellGeometryIntegrationRuleSet class & friends.
virtual double dimension(const GmMatrix &X) const
Returns the Quad9 area by numeric integration.
Definition: gmQuadCellGeometryInfo.h:235
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 3, 3, 3, 3 > GmQuadraticQuadIntegrationRuleSetDefaults
The set of default rules for a quadratic Quad ELEMENT. Template parameters 2 to 5 are the default rul...
Definition: gmQuadCellGeometryInfo.h:69
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: gmQuadCellGeometryInfo.h:207
A helper class used to store the class types for each rule type. Template parameters are the implemen...
Definition: gmCellGeometryIntegrationRuleSet.h:196
GmShape specialization for a 3D Quadrilateral with 8 nodes object.
Definition: gmQuadShape.h:235
bool pointInTesselatedPolygon_Q8(const GmMatrix &X, const GmVector &p)
Returns true if the point p is inside the tesselated Q8 polygon or false otherwise.
Definition: gmGeometryUtils.cpp:700
bool polygon2DIsCCW(const GmMatrix &X)
Returns true if the polygon delimited by the cell vertices is oriented CCW.
Definition: gmGeometryUtils.cpp:527
Utilitary functions for working with geometry.
The Quad3D8 implementation.
Definition: gmQuadCellGeometryInfo.h:197
2D Isotropic / anisotropic Newton cotes integration rules for Quad elements
Definition: gmQuadIntegrationRule.h:65
An auxiliary class that can be used as base for surface elements. Implements the needed integration r...
Definition: gmCellGeometryInfo.h:348
virtual bool contains(const GmMatrix &X, const GmVector &coord) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_CONTAIN...
Definition: gmQuadCellGeometryInfo.h:141
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
GmQuad3D8CellGeometryInfo(const GmQuad8CellGeometryInfo &quad8)
Private constructor. Only a single Quad3D8 geometry info object is ever necessary,...
Definition: gmQuadCellGeometryInfo.h:216
void quadCentroid(const GmVector &p1, const GmVector &p2, const GmVector &p3, const GmVector &p4, GmVector &coord)
Returns the cartesian coordinate of the quad centroid.
Definition: gmGeometryUtils.cpp:304
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: gmQuadCellGeometryInfo.h:138
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 3, 3, -1, -1 > GmQuadraticQuadEdgeIntegrationRuleSetDefaults
The set of default rules for a quadratic Quad EDGE. Template parameters 2 to 5 are the default rule n...
Definition: gmQuadCellGeometryInfo.h:74
The Quad4 implementation.
Definition: gmQuadCellGeometryInfo.h:82
GmCellGeometryIntegrationRuleSet< GmQuadGaussIntegrationRule, GmQuadLobattoIntegrationRule, GmQuadNewtonIntegrationRule, GmQuadNewtonIntegrationRule > GmQuadIntegrationRuleSet
The set of possible ELEMENT integration rules, per integration rule type, for the family of Quad elem...
Definition: gmQuadCellGeometryInfo.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 bool contains(const GmMatrix &X, const GmVector &coord) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_CONTAIN...
Definition: gmQuadCellGeometryInfo.h:212
virtual double dimension(const GmMatrix &X) const
Returns the Quad8 area by numeric integration.
Definition: gmQuadCellGeometryInfo.h:168
Declaration of the GmQuadXxxxIntegrationRule family of classes, including border rules....
virtual double dimension(const GmMatrix &X) const
Returns the Quad4 area.
Definition: gmQuadCellGeometryInfo.h:93
virtual bool contains(const GmMatrix &X, const GmVector &coord) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_CONTAIN...
Definition: gmQuadCellGeometryInfo.h:114
GmQuad8CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Protected constructor. Only a single Quad8 geometry info object is ever necessary,...
Definition: gmQuadCellGeometryInfo.h:187
Declaration of the 2D Quad shapes inheriting from GmShape class Previosuly part of the gmShape2D....
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: gmQuadCellGeometryInfo.h:107
bool polygon2DIsSelfIntersecting(const GmMatrix &X)
Returns true if the polygon delimited by the cell vertices is self-intersecting, or false otherwise.
Definition: gmGeometryUtils.cpp:498
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 2, 2, 2, 2 > GmLinearQuadIntegrationRuleSetDefaults
The set of default rules for a linear Quad ELEMENT. Template parameters 2 to 5 are the default rule n...
Definition: gmQuadCellGeometryInfo.h:55
static const GmQuad8CellGeometryInfo * instance()
Instance function for the Quad8 single geometry object.
Definition: gmQuadCellGeometryInfo.cpp:114
static const GmQuad4CellGeometryInfo * instance()
Instance function for the Quad4 single geometry object.
Definition: gmQuadCellGeometryInfo.cpp:34
Declaration of the GmCellGeometryInfo base class.
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
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
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: gmQuadCellGeometryInfo.h:203
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: gmQuadCellGeometryInfo.h:232
The Quad9 implementation.
Definition: gmQuadCellGeometryInfo.h:224
GmQuad4CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Protected constructor. Only a single Quad4 geometry info object is ever necessary,...
Definition: gmQuadCellGeometryInfo.h:118
Plane structure storing the full set of geometric metadata for a cell type.
Definition: gmCellGeometryInfo.h:57
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 2, 2, -1, -1 > GmLinearQuadEdgeIntegrationRuleSetDefaults
The set of default rules for a linear Quad EDGE. Template parameters 2 to 5 are the default rule numb...
Definition: gmQuadCellGeometryInfo.h:60
Gauss integration rule for Quad element borders.
Definition: gmQuadIntegrationRule.h:130