GemaCoreLib
The GeMA Core library
gmCellGeometryInfo.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_CELL_GEOMETRY_INFO_H_
25 #define _GEMA_CELL_GEOMETRY_INFO_H_
26 
27 #include "gmCoreConfig.h"
28 #include "gmCellType.h"
29 #include "gmIntegrationRule.h"
30 #include "gmMatrix.h"
31 #include "gmVector.h"
32 
33 #include <QAtomicPointer>
34 #include <QMutex>
35 
36 class GmShape;
38 
43 {
44  // IMPORTANT: DO NOT CHANGE THE ORDER OF THE ENUM BELOW without changing the order
45  // used to fill the GmCellGeometryMetadata._capabilities vector for
46  // EVERY GmCellGeometryInfo derived class!!!!
47  // Also, when ADDING entries, rember to fill the said vector in EVERY derived class.
51  // ----------------------------
52  // No adding after this line
54 };
55 
58 {
60  struct EdgeNodeData
61  {
67 
72  };
73 
75  struct FaceNodeData
76  {
78 
82 
88 
95  };
96 
98  const char* _name;
100  bool _interface;
102  int _order;
104  int _nnodes;
107  int _ncoord;
108  int _nedges;
109  int _nfaces;
111 
112  // A vector storing true or false to specifiy if a cell type has the capabilities listed in GmCellGeometryCapabilities
113  bool _capabilities[GM_CELL_GEOMETRY_CAPABILITY_COUNT];
114 
117 
120 
123 
126 
129 
132 
134  GmCellGeometryMetadata(const GmCellGeometryMetadata& other) = default;
135 
138 };
139 
140 
141 /* \brief An interface that links a cell geometric metadata with the associated
142  shape function and integration rules through factory methods that can
143  create those associated objects. It also implements the geometric routines
144  specifica for a type (like length and area calculation routines).
145 
146  This interface should be instanced only by way of derived classes that make the
147  concrete link with the shape function / integration rules.
148 
149  If in the future we need to support an hypotetic cell type that has no associated
150  shape function or integration rule (a cell type that is not an element type), a
151  derived class that returns NULL for all factory methods can be implemented.
152 */
154 {
155 public:
158 
161 
167  virtual GmShape* shapeInstance(int P, int Q) const = 0;
168 
186  virtual GmIntegrationRule* integrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3, int P, int Q) const = 0;
187 
207  virtual GmBorderIntegrationRule* edgeIntegrationRule(GmIntegrationRuleType irType, int rule1, int P, int Q) const = 0;
208 
234  virtual GmBorderIntegrationRule* faceIntegrationRule(int faceType, GmIntegrationRuleType irType, int rule1, int rule2, int P, int Q) const = 0;
235 
236 
243  virtual double dimension(const GmMatrix& X) const = 0;
244 
249  virtual void centroidCartesian(const GmMatrix& X, GmVector& coord) const = 0;
250 
251  //
252  // Capability related methods. Must be implmented on every derived classes that
253  // do implement the capability and should be IN SYNC with the information
254  // stored in the associeted GmCellGeometryMetadata _capabilities field.
255  //------------------------------------------------------------------------------------
256 
264  virtual bool isValid(const GmMatrix& X, double tol) const { Q_UNUSED(tol); Q_UNUSED(X); return true; }
265 
273  virtual double quality(const GmMatrix& J, double tol) const {
274  assert(tol > 0);
275  if (arma::det(J) > tol) {
276  GmMatrix U, V;
277  GmVector s;
278  arma::svd_econ(U, s, V, J, "left");
279 
280  //Proportion between the smaller / largest eigen value.
281  return s[s.n_elem - 1] / s[0];
282  }
283  return 0;
284  }
285 
292  virtual bool contains(const GmMatrix& X, const GmVector& coord) const { Q_UNUSED(X); Q_UNUSED(coord); return false; }
293 
294  //------------------------------------------------------------------------------------
295  // Helper functions for dimension and centroid calculation by numeric integration
296  //------------------------------------------------------------------------------------
297  const GmIntegrationRule* cellIntegrationRule() const;
298 
299  double dimensionByIntegration(const GmMatrix& X, const GmIntegrationRule* ir = NULL) const;
300  void centroidByIntegration (const GmMatrix& X, GmVector& coord, const GmIntegrationRule* ir = NULL) const;
301 
302 private:
308 };
309 
310 
314 template <class RuleSet, class RuleSetDefaults>
316 {
317 public:
320 
322  virtual GmIntegrationRule* integrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3, int P, int Q) const
323  {
324  Q_UNUSED(P); Q_UNUSED(Q);
325  return RuleSet::instance<RuleSetDefaults>(this, irType, rule1, rule2, rule3);
326  }
327 
328  // Element edge integration rule factory. See comments on the base class.
329  virtual GmBorderIntegrationRule* edgeIntegrationRule(GmIntegrationRuleType irType, int rule1, int P, int Q) const
330  {
331  Q_UNUSED(irType); Q_UNUSED(rule1); Q_UNUSED(P); Q_UNUSED(Q);
332  return NULL;
333  }
334 
335  // Element face integration rule factory. See comments on the base class.
336  virtual GmBorderIntegrationRule* faceIntegrationRule(int faceType, GmIntegrationRuleType irType, int rule1, int rule2, int P, int Q) const
337  {
338  Q_UNUSED(faceType); Q_UNUSED(irType); Q_UNUSED(rule1); Q_UNUSED(rule2); Q_UNUSED(P); Q_UNUSED(Q);
339  return NULL;
340  }
341 };
342 
343 
347 template <class RuleSet, class RuleSetDefaults, class EdgeRuleSet, class EdgeRuleSetDefaults>
349 {
350 public:
353 
354  // Element integration rule factory. See comments on the base class.
355  virtual GmIntegrationRule* integrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3, int P, int Q) const
356  {
357  Q_UNUSED(P); Q_UNUSED(Q);
358  return RuleSet::instance<RuleSetDefaults>(this, irType, rule1, rule2, rule3);
359  }
360 
361  // Element edge integration rule factory. See comments on the base class.
362  virtual GmBorderIntegrationRule* edgeIntegrationRule(GmIntegrationRuleType irType, int rule1, int P, int Q) const
363  {
364  Q_UNUSED(P); Q_UNUSED(Q);
365  return EdgeRuleSet::edgeInstance<EdgeRuleSetDefaults>(this, irType, rule1);
366  }
367 
368  // Element face integration rule factory. See comments on the base class.
369  virtual GmBorderIntegrationRule* faceIntegrationRule(int faceType, GmIntegrationRuleType irType, int rule1, int rule2, int P, int Q) const
370  {
371  Q_UNUSED(faceType); Q_UNUSED(irType); Q_UNUSED(rule1); Q_UNUSED(rule2); Q_UNUSED(P); Q_UNUSED(Q);
372  return NULL;
373  }
374 };
375 
379 template <class RuleSet, class RuleSetDefaults, class EdgeRuleSet, class EdgeRuleSetDefaults, class FaceRuleSet, class FaceRuleSetDefaults>
381 {
382 public:
385 
386  // Element integration rule factory. See comments on the base class.
387  virtual GmIntegrationRule* integrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3, int P, int Q) const
388  {
389  Q_UNUSED(P); Q_UNUSED(Q);
390  return RuleSet::instance<RuleSetDefaults>(this, irType, rule1, rule2, rule3);
391  }
392 
393  // Element edge integration rule factory. See comments on the base class.
394  virtual GmBorderIntegrationRule* edgeIntegrationRule(GmIntegrationRuleType irType, int rule1, int P, int Q) const
395  {
396  Q_UNUSED(P); Q_UNUSED(Q);
397  return EdgeRuleSet::edgeInstance<EdgeRuleSetDefaults>(this, irType, rule1);
398  }
399 
400  // Element face integration rule factory. See comments on the base class.
401  virtual GmBorderIntegrationRule* faceIntegrationRule(int faceType, GmIntegrationRuleType irType, int rule1, int rule2, int P, int Q) const
402  {
403  Q_UNUSED(P); Q_UNUSED(Q);
404  return FaceRuleSet::faceInstance<FaceRuleSetDefaults>(this, faceType, irType, rule1, rule2);
405  }
406 };
407 
408 
409 #endif
410 
QVector< int > _edgeNodes
The set of local node numbers, ordered sequentialy in the edge. The extremes(first and last entries) ...
Definition: gmCellGeometryInfo.h:71
GmIntegrationRuleType
The type of desired integration rule (Gauss quadrature, Lobatto quadrature, etc)
Definition: gmIntegrationRule.h:67
Can the cell geometry info type check wheter a point is inside a cell?
Definition: gmCellGeometryInfo.h:50
GmCellType _eqFaceType
The 3D face equivalent surface element or GM_INV_CELL_TYPE if there is none. The invalid type is also...
Definition: gmCellGeometryInfo.h:81
const char * _name
The cell type name.
Definition: gmCellGeometryInfo.h:98
GmCellGeometryInfoSurfaceElement(GmCellGeometryMetadata &&metadata)
Constructor receiving as parameter a metadata object that will be MOVED to the new object.
Definition: gmCellGeometryInfo.h:352
int _nfaces
Number of cell faces.
Definition: gmCellGeometryInfo.h:109
An auxiliary class that can be used as base for line (bar) elements. Implements the needed integratio...
Definition: gmCellGeometryInfo.h:315
Aux structure storing information about one edge of a cell type.
Definition: gmCellGeometryInfo.h:60
An auxiliary class that can be used as base for solid elements. Implements the needed integration rul...
Definition: gmCellGeometryInfo.h:380
NOT a cell type. Stores the number of available types.
Definition: gmCellType.h:77
GmCellGeometryInfoSolidElement(GmCellGeometryMetadata &&metadata)
Constructor receiving as parameter a metadata object that will be MOVED to the new object.
Definition: gmCellGeometryInfo.h:384
QVector< int > _volumeNodeInfo
Mid volume node list.
Definition: gmCellGeometryInfo.h:131
Declaration of usefull configuration definitions for the Core library.
GmCellType _type
The cell type.
Definition: gmCellGeometryInfo.h:97
GmCellGeometryInfo(GmCellGeometryMetadata &&metadata)
Standard move constructor receiving the type metadata.
Definition: gmCellGeometryInfo.h:160
Declaration of the GmCellType enum.
int _order
The cell order (1: linear, 2: quadratic, 3: cubic, ...)
Definition: gmCellGeometryInfo.h:102
int _ncoord
The number of cartesian coordinates for nodes of this cell type.
Definition: gmCellGeometryInfo.h:107
virtual GmBorderIntegrationRule * edgeIntegrationRule(GmIntegrationRuleType irType, int rule1, int P, int Q) const
A factory function that returns a new border integration rule object suited for this kind of element....
Definition: gmCellGeometryInfo.h:362
Declaration of the GmMatrix class.
GmCellFamilyType _family
The family that this cell type belongs to.
Definition: gmCellGeometryInfo.h:99
GmCellGeometryInfoLineElement(GmCellGeometryMetadata &&metadata)
Constructor receiving as parameter a metadata object that will be MOVED to the new object.
Definition: gmCellGeometryInfo.h:319
QVector< QVector< int > > _edgeFaceInfo
Face information for each cell's edge.
Definition: gmCellGeometryInfo.h:122
Not a capability type. Stores the number of capabilities in the enum.
Definition: gmCellGeometryInfo.h:53
QVector< QVector< int > > _nodeIncidenceInfo
Node incidence information for each cell's node. Size equal to "nnodes".
Definition: gmCellGeometryInfo.h:128
Declaration of the GmIntegrationRule class and its helper decendants.
QVector< int > _faceNodes
The list of nodes in the face, ordered in a way that the first b entries are the boundary nodes order...
Definition: gmCellGeometryInfo.h:87
GmCellGeometryInfo(const GmCellGeometryInfo &other)
Standard copy constructor receiving the type metadata.
Definition: gmCellGeometryInfo.h:157
Integration rule base classe.
Definition: gmIntegrationRule.h:88
Shape function handling base classe.
Definition: gmShape.h:37
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: gmCellGeometryInfo.h:264
QVector< QVector< int > > _faceEdgeInfo
Edge information for each cell's face, ordered according to the cell's face description....
Definition: gmCellGeometryInfo.h:119
virtual GmIntegrationRule * integrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3, int P, int Q) const
A factory function that returns a new integration rule object suited for this kind of element.
Definition: gmCellGeometryInfo.h:387
int _nnodes
Number of cell nodes.
Definition: gmCellGeometryInfo.h:104
virtual GmIntegrationRule * integrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3, int P, int Q) const
Element integration rule factory. See comments on the base class.
Definition: gmCellGeometryInfo.h:322
bool _hierarchical
Is this cell an hierarchical element (needing P & Q order parameters?)
Definition: gmCellGeometryInfo.h:101
Can the cell geometry info type check if a cell is valid?
Definition: gmCellGeometryInfo.h:48
GmCellGeometryCapabilities
An enum storing the possible geometric capabilities for a cell geometry info type that can be queried...
Definition: gmCellGeometryInfo.h:42
Border integration rule base classe.
Definition: gmBorderIntegrationRule.h:37
int _nextraDofNodes
Number of extra dof nodes. Usually 0. See comments for GmCellGeometry::numVertices().
Definition: gmCellGeometryInfo.h:106
QVector< FaceNodeData > _faceNodeInfo
Node information for each cell's face. Size equal to "nfaces". Empty for 1D elements.
Definition: gmCellGeometryInfo.h:125
QVector< EdgeNodeData > _edgeNodeInfo
Node information for each cell's edge. Size equal to "nedges".
Definition: gmCellGeometryInfo.h:116
QVector< int > _faceTypeNodes
If faceType is not GM_INV_CELL_TYPE, stores the list of nodes in the equivalent surface.
Definition: gmCellGeometryInfo.h:94
int _nedges
Number of cell edges.
Definition: gmCellGeometryInfo.h:108
virtual double quality(const GmMatrix &J, double tol) const
Virtual method that should be implemented if this geometry info supports the GM_CELL_GEOMETRY_QUALITY...
Definition: gmCellGeometryInfo.h:273
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: gmCellGeometryInfo.h:292
virtual GmBorderIntegrationRule * edgeIntegrationRule(GmIntegrationRuleType irType, int rule1, int P, int Q) const
A factory function that returns a new border integration rule object suited for this kind of element....
Definition: gmCellGeometryInfo.h:394
An auxiliary class that can be used as base for surface elements. Implements the needed integration r...
Definition: gmCellGeometryInfo.h:348
Declaration of the GmVector class.
virtual GmBorderIntegrationRule * faceIntegrationRule(int faceType, GmIntegrationRuleType irType, int rule1, int rule2, int P, int Q) const
A factory function that returns a new border integration rule object suited for this kind of element....
Definition: gmCellGeometryInfo.h:369
#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 GmBorderIntegrationRule * faceIntegrationRule(int faceType, GmIntegrationRuleType irType, int rule1, int rule2, int P, int Q) const
A factory function that returns a new border integration rule object suited for this kind of element....
Definition: gmCellGeometryInfo.h:401
GmCellType
Mesh Cell types. Don't change type orders or add types without reading comments below.
Definition: gmCellType.h:30
bool _interface
Is this cell an interface element?
Definition: gmCellGeometryInfo.h:100
virtual GmBorderIntegrationRule * edgeIntegrationRule(GmIntegrationRuleType irType, int rule1, int P, int Q) const
A factory function that returns a new border integration rule object suited for this kind of element....
Definition: gmCellGeometryInfo.h:329
virtual GmIntegrationRule * integrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3, int P, int Q) const
A factory function that returns a new integration rule object suited for this kind of element.
Definition: gmCellGeometryInfo.h:355
Aux structure storing information about the set of face edges for a cell type.
Definition: gmCellGeometryInfo.h:75
GmCellType _eqEdgeType
The 2D/3D bar equivalent edge element or GM_INV_CELL_TYPE if there is none. Currently,...
Definition: gmCellGeometryInfo.h:66
Definition: gmCellGeometryInfo.h:153
static QMutex _cellIntegrationRulesMutex
The mutex protecting _cellIntegrationRules.
Definition: gmCellGeometryInfo.h:307
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
int _nboundaryNodes
The number "b" of boundary nodes for this cell face. Excludes internal face nodes.
Definition: gmCellGeometryInfo.h:77
virtual GmBorderIntegrationRule * faceIntegrationRule(int faceType, GmIntegrationRuleType irType, int rule1, int rule2, int P, int Q) const
A factory function that returns a new border integration rule object suited for this kind of element....
Definition: gmCellGeometryInfo.h:336
arma::mat GmMatrix
The basic type for a GeMA matrix object. Currently based on an Armadillo matrix.
Definition: gmMatrix.h:38
int _nvertices
Number of cell vertices. Excludes "quadratic" nodes and extra dof nodes. See comments for GmCellGeome...
Definition: gmCellGeometryInfo.h:105
GmCellFamilyType
Mesh cell type families. Use to group all quad, tri, etc elements in a "family". Its start value is a...
Definition: gmCellType.h:106
Plane structure storing the full set of geometric metadata for a cell type.
Definition: gmCellGeometryInfo.h:57
int _nfaceTypes
Number of face types for a 3D element. 1 For surface elements.
Definition: gmCellGeometryInfo.h:110
Can the cell geometry info type give a cell quality measure?
Definition: gmCellGeometryInfo.h:49
GmCellType _eqLinearType
The equivalent linear element.
Definition: gmCellGeometryInfo.h:103