GemaCoreLib
The GeMA Core library
gmTriCellGeometryInfo.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_TRI_CELL_GEOMETRY_INFO_H_
25 #define _GEMA_TRI_CELL_GEOMETRY_INFO_H_
26 
27 #include "gmCellGeometryInfo.h"
29 
30 #include "gmGeometryUtils.h"
31 
32 #include "gmTriIntegrationRule.h"
33 #include "gmTriShape.h"
34 
35 //---------------------------------------------------------
36 // Tri integration rules
37 //---------------------------------------------------------
38 
42 
46 
47 
48 //---------------------------------------------------------
49 // Tri3 & Tri3D3 rule defaults
50 //---------------------------------------------------------
51 
56 
61 
62 //---------------------------------------------------------
63 // Tri6 & Tri3D6 rule defaults
64 //---------------------------------------------------------
65 
70 
75 
76 
77 //---------------------------------------------------------
78 // GmTri3CellGeometryInfo
79 //---------------------------------------------------------
80 
83  : public GmCellGeometryInfoSurfaceElement<GmTriIntegrationRuleSet, GmLinearTriIntegrationRuleSetDefaults,
84  GmTriEdgeIntegrationRuleSet, GmLinearTriEdgeIntegrationRuleSetDefaults>
85 {
86 public:
87  static const GmTri3CellGeometryInfo* 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 GmTri3Shape; }
91 
93  virtual double dimension(const GmMatrix& X) const
94  {
95  return GmGeometryUtils::triangleArea(X.col(0), X.col(1), X.col(2));
96  }
97 
99  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const
100  {
101  GmGeometryUtils::triangleCentroid(X.col(0), X.col(1), X.col(2), coord);
102  }
103 
104  // --------- Capabilities API ---------
105 
106  virtual bool isValid(const GmMatrix& X, double tol) const {
108  }
109 
110  // See comments on the base class
111  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { return GmGeometryUtils::pointInConvexPolygon(X, coord); }
112 
113 protected:
118 };
119 
120 //---------------------------------------------------------
121 // GmTri3D3CellGeometryInfo
122 //---------------------------------------------------------
123 
126 {
127 public:
128  static const GmTri3D3CellGeometryInfo* instance();
129 
130  // Shape factory. See comments on the base class.
131  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmTri3D3Shape; }
132 
133  // --------- Capabilities API ---------
134 
135  virtual bool isValid(const GmMatrix& X, double tol) const { Q_UNUSED(tol); return !GmGeometryUtils::tri3DIsColinear(X); }
136  // contains() capability is currently not supported in 3D, so the implementation
137  // from GmTri3CellGeometryInfo must be overridden. Quality is ok.
138  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { Q_UNUSED(X); Q_UNUSED(coord); return false; }
139 
140 private:
143 };
144 
145 //---------------------------------------------------------
146 // GmTri6CellGeometryInfo
147 //---------------------------------------------------------
148 
149 //tri6 nodes in circular order
150 static int tri6order[] = {0,3,1,4,2,5};
151 
154  : public GmCellGeometryInfoSurfaceElement<GmTriIntegrationRuleSet, GmQuadraticTriIntegrationRuleSetDefaults,
155  GmTriEdgeIntegrationRuleSet, GmQuadraticTriEdgeIntegrationRuleSetDefaults>
156 {
157 public:
158  static const GmTri6CellGeometryInfo* instance();
159 
160  // Shape factory. See comments on the base class.
161  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmTri6Shape; }
162 
164  virtual double dimension(const GmMatrix& X) const { return dimensionByIntegration(X); }
165 
167  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const { centroidByIntegration(X, coord); }
168 
169  // --------- Capabilities API ---------
170 
171  virtual bool isValid(const GmMatrix& X, double tol) const {
172  Q_UNUSED(tol);
173  GmMatrix Xr = GmGeometryUtils::shuffle(X, tri6order, 6);
175  }
176 
177  // See comments on the base class
178  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { return GmGeometryUtils::pointInTesselatedPolygon_T6(X, coord); }
179 
180 protected:
185 };
186 
187 //---------------------------------------------------------
188 // GmTri3D6CellGeometryInfo
189 //---------------------------------------------------------
190 
193 {
194 public:
195  static const GmTri3D6CellGeometryInfo* instance();
196 
197  // Shape factory. See comments on the base class.
198  virtual GmShape* shapeInstance(int P, int Q) const { Q_UNUSED(P); Q_UNUSED(Q); return new GmTri3D6Shape; }
199 
200  // --------- Capabilities API ---------
201 
202  virtual bool isValid(const GmMatrix& X, double tol) const {
204  && !GmGeometryUtils::tri3DIsColinear(X.cols(0,2));;
205  }
206 
207  // contains() capability is currently not supported in 3D, so the implementation
208  // from GmTri6CellGeometryInfo must be overridden.
209  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { Q_UNUSED(X); Q_UNUSED(coord); return false; }
210 
211 private:
214 };
215 
216 
217 #endif
218 
Declaration of the 2D tri shapes inheriting from GmShape class Previosuly part of the gmShape2D....
Gauss Integration rule for triangle elements. Base on Felippa, IFEM, section 24.2.
Definition: gmTriIntegrationRule.h:42
GmTri3D6CellGeometryInfo(const GmTri6CellGeometryInfo &tri6)
Private constructor. Only a single Tri3D6 geometry info object is ever necessary, created by GmTriCel...
Definition: gmTriCellGeometryInfo.h:213
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 3, -1, 1, 1 > GmLinearTriIntegrationRuleSetDefaults
The set of default rules for a linear Tri ELEMENT. Template parameters 2 to 5 are the default rule nu...
Definition: gmTriCellGeometryInfo.h:55
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
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 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: gmTriCellGeometryInfo.h:209
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: gmTriCellGeometryInfo.h:111
Gauss rules.
Definition: gmIntegrationRule.h:75
void triangleCentroid(const GmVector &p1, const GmVector &p2, const GmVector &p3, GmVector &coord)
Returns the cartesian coordinate of the triangle centroid.
Definition: gmGeometryUtils.cpp:293
GmCellGeometryIntegrationRuleSet< GmTriGaussIntegrationRule, GmCellGeometryNullIntegrationRule, GmTriNewtonIntegrationRule, GmTriNewtonIntegrationRule > GmTriIntegrationRuleSet
The set of possible ELEMENT integration rules, per integration rule type, for the family of Tri eleme...
Definition: gmTriCellGeometryInfo.h:41
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
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: gmTriCellGeometryInfo.h:138
GmShape specialization for a Triangle with 6 nodes object.
Definition: gmTriShape.h:145
A helper class used to store the default rules for each rule type. Template parameters are the defaul...
Definition: gmCellGeometryIntegrationRuleSet.h:43
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: gmTriCellGeometryInfo.h:135
The Tri3D3 implementation.
Definition: gmTriCellGeometryInfo.h:125
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: gmTriCellGeometryInfo.h:202
Newton Integration rule for triangle elements. Accepts both open and closed rules.
Definition: gmTriIntegrationRule.h:65
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: gmTriCellGeometryInfo.h:198
GmTri6CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Protected constructor. Only a single Tri6 geometry info object is ever necessary, created by GmTriCel...
Definition: gmTriCellGeometryInfo.h:182
Shape function handling base classe.
Definition: gmShape.h:37
GmTri3D3CellGeometryInfo(const GmTri3CellGeometryInfo &tri3)
Private constructor. Only a single Tri3D3 geometry info object is ever necessary, created by GmTriCel...
Definition: gmTriCellGeometryInfo.h:142
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: gmTriCellGeometryInfo.h:178
Declaration of the GmCellGeometryIntegrationRuleSet class & friends.
The Tri6 implementation.
Definition: gmTriCellGeometryInfo.h:153
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: gmTriCellGeometryInfo.h:131
virtual double dimension(const GmMatrix &X) const
Returns the Tri3 area.
Definition: gmTriCellGeometryInfo.h:93
A helper class used to store the class types for each rule type. Template parameters are the implemen...
Definition: gmCellGeometryIntegrationRuleSet.h:196
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 2, -1, -1, -1 > GmLinearTriEdgeIntegrationRuleSetDefaults
The set of default rules for a linear Tri EDGE. Template parameters 2 to 5 are the default rule numbe...
Definition: gmTriCellGeometryInfo.h:60
bool polygon2DIsCCW(const GmMatrix &X)
Returns true if the polygon delimited by the cell vertices is oriented CCW.
Definition: gmGeometryUtils.cpp:527
GmCellGeometryIntegrationRuleSet< GmTriGaussEdgeIntegrationRule, GmCellGeometryNullIntegrationRule, GmCellGeometryNullIntegrationRule, GmCellGeometryNullIntegrationRule > GmTriEdgeIntegrationRuleSet
The set of possible EDGE integration rules, per integration rule type, for the family of Tri elements...
Definition: gmTriCellGeometryInfo.h:45
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 GmShape * shapeInstance(int P, int Q) const
Shape function factory. Should return a NEW instance of the shape function object for this type....
Definition: gmTriCellGeometryInfo.h:161
Gauss integration rule for Tri element borders.
Definition: gmTriIntegrationRule.h:95
#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
double triangleArea(const GmVector &p1, const GmVector &p2, const GmVector &p3)
Calculates the area of the triangle determined by points p1, p2 and p3 (2D or 3D),...
Definition: gmGeometryUtils.cpp:165
virtual double dimension(const GmMatrix &X) const
Returns the Tri6 area by numeric integration.
Definition: gmTriCellGeometryInfo.h:164
GmShape specialization for a 3D Triangle with 6 nodes object.
Definition: gmTriShape.h:202
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 6, -1, 2, 2 > GmQuadraticTriIntegrationRuleSetDefaults
The set of default rules for a quadratic Tri ELEMENT. Template parameters 2 to 5 are the default rule...
Definition: gmTriCellGeometryInfo.h:69
static const GmTri3CellGeometryInfo * instance()
Instance function for the Tri3 single geometry object.
Definition: gmTriCellGeometryInfo.cpp:33
bool pointInTesselatedPolygon_T6(const GmMatrix &X, const GmVector &p)
Returns true if the point p is inside the tesselated T6 polygon or false otherwise.
Definition: gmGeometryUtils.cpp:666
The Tri3 implementation.
Definition: gmTriCellGeometryInfo.h:82
The Tri3D6 implementation.
Definition: gmTriCellGeometryInfo.h:192
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
Declaration of the GmTriXxxxIntegrationRule family of classes, including border rules....
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: gmTriCellGeometryInfo.h:171
Declaration of the GmCellGeometryInfo base class.
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: gmTriCellGeometryInfo.h:90
static const GmTri6CellGeometryInfo * instance()
Instance function for the Tri6 single geometry object.
Definition: gmTriCellGeometryInfo.cpp:111
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
GmCellGeometryIntegrationRuleSetDefaultRules< GM_GAUSS_RULE_TYPE, 3, -1, -1, -1 > GmQuadraticTriEdgeIntegrationRuleSetDefaults
The set of default rules for a quadratic Tri EDGE. Template parameters 2 to 5 are the default rule nu...
Definition: gmTriCellGeometryInfo.h:74
GmShape specialization for a Triangle with 3 nodes object.
Definition: gmTriShape.h:111
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
GmTri3CellGeometryInfo(GmCellGeometryMetadata &&metadata)
Protected constructor. Only a single Tri3 geometry info object is ever necessary, created by GmTriCel...
Definition: gmTriCellGeometryInfo.h:115
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: gmTriCellGeometryInfo.h:106
GmShape specialization for a 3D Triangle with 3 nodes object.
Definition: gmTriShape.h:171
Plane structure storing the full set of geometric metadata for a cell type.
Definition: gmCellGeometryInfo.h:57
virtual void centroidCartesian(const GmMatrix &X, GmVector &coord) const
Returns the Tri6 centroid by numeric integration.
Definition: gmTriCellGeometryInfo.h:167
virtual void centroidCartesian(const GmMatrix &X, GmVector &coord) const
Returns the Tri3 centroid by calculating the mean of the vertices (simplex centroid)
Definition: gmTriCellGeometryInfo.h:99