GemaCoreLib
The GeMA Core library
Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
GmCellGeometry Class Reference

A class used to return static metadata information about a cell geometry, along with some methods for calculating geometric properties for a cell, given its coordinate matrix. More...

#include <gmCellGeometry.h>

Collaboration diagram for GmCellGeometry:
Collaboration graph
[legend]

Public Member Functions

 GmCellGeometry (GmCellType type, int P=0, int Q=0)
 Constructor. Receives the cell type for which geometry information will be returned. P and Q are needed only when working with hierarchical element types, where they are used together with the shape and integration rule factories. It is safe to create a geometry object with P & Q == 0, even for hierarchical types, if you are never going tocall the shape or integration rule factory functions.
 
GmCellType type () const
 Return the cell type associated to this object.
 
const char * typeName () const
 Returns the type name.
 
GmCellFamilyType family () const
 Returns the family to which this cell type belongs.
 
bool isInterface () const
 Returns true if this is an interface element, false otherwise.
 
bool isHierarchical () const
 Returns true if this is a hierarchical element, false otherwise.
 
int order () const
 Returns the cell interpolation order(1 = linear, 2 = quadratic, 3 = cubic, ...)
 
GmCellType linearElement () const
 Returns the type of the equivalent linear element (for a Quad9, returns a Quad4, for example).
 
int numNodes () const
 Returns the total number of nodes of this cell type.
 
int numVertices () const
 Returns the number of vertices of this cell type. This excludes center edge and face nodes - returns 4 for a quad8, for example. It also excludes additional dof nodes. More...
 
int numExtraDofNodes () const
 Returns the number of extra degrees of freedom (dof) nodes. This is usually 0, except for some kinds of interface elements. See also comments on numVertices().
 
int numCoord () const
 Returns the size of the cartesian coordinates for this cell type.
 
int numEdges () const
 Returns the number of edges of this cell type.
 
int numFaces () const
 Returns the number of faces of this cell type (0 for "bar" alike cells, 1 for surface cells)
 
bool isLine () const
 Is this element an 1D "line" ?
 
bool isSurface () const
 Is this element an 2D "surface".
 
bool isSolid () const
 Is this element an 3D "solid"?
 
int numEdgeNodes (int edgeIndex) const
 Returns the number of nodes in the cell's edge numbered edgeIndex (edgeIndex from 0 to numEdges()-1)
 
int edgeNode (int edgeIndex, int i) const
 Returns the ith node belonging to edge edgeIndex. More...
 
int firstEdgeNode (int edgeIndex) const
 Returns the first node in the cell's edge numbered edgeIndex (edgeIndex from 0 to numEdges()-1) Equivalent to edgeNode(edgeIndex, 0)
 
int lastEdgeNode (int edgeIndex) const
 Returns the last node in the cell's edge numbered edgeIndex (edgeIndex from 0 to numEdges()-1) Equivalent to edgeNode(edgeIndex, numEdgeNodes()-1)
 
bool nodeOnEdge (int nodeIndex, int edgeIndex) const
 Returns true if the given (local) node index belongs to the given edge definition.
 
int edgeFromNodes (int n1, int n2) const
 Returns the edge index of the edge starting with n1 and ending with n2 or vice-versa. Both n1 and n2 should be local node indices. Returns -1 if the two local nodes do not define an element edge.
 
int numFaceEdges (int faceIndex) const
 Returns the number of edges in the cell's face numbered faceIndex (faceIndex from 0 to numFaces()-1)
 
int faceEdge (int faceIndex, int i) const
 Returns the ith edge belonging to face faceIndex. More...
 
int numEdgeFaces (int edgeIndex) const
 Returns the number of faces in the cell's edge numbered edgeIndex (edgeIndex from 0 to numEdges()-1)
 
int edgeFace (int edgeIndex, int i) const
 Returns the ith face belonging to edge edgeIndex. More...
 
int numFaceBoundaryNodes (int faceIndex) const
 Returns the number of boundary nodes in the cell's face numbered faceIndex (faceIndex from 0 to numFaces()-1). This node count excludes internal face nodes.
 
int numFaceNodes (int faceIndex) const
 Returns the number of nodes in the cell's face numbered faceIndex (faceIndex from 0 to numFaces()-1)
 
int faceNode (int faceIndex, int i) const
 Returns the ith node belonging to the face faceIndex. Since boundary nodes are always returned first by this function, it can be used to retrieve all nodes or just the boundary ones. More...
 
int numIncidenceNodes (int nodeIndex) const
 Returns the number of incident nodes into the cell's node numbered nodeIndex (nodeIndex from 0 to numNodes()-1)
 
int incidenceNode (int nodeIndex, int i) const
 Returns the ith incidence node belonging to node nodeIndex. More...
 
int numVolumeInternalNodes () const
 Returns the number of internal nodes of this volume's cell type. Returns 0 for 1d/2d cells.
 
int volumeInternalNode (int i) const
 Returns the ith internal node belonging to the volume. More...
 
int numFaceTypes () const
 Returns the number of face types for a 3D element (generally 1), 1 for 2D elements and 0 for 1D elements. More...
 
GmCellFaceType faceType (int faceIndex) const
 Returns the face type for the cells face numbered faceIndex (faceIndex from 0 to numFaces()-1)
 
GmCellType edgeElement (int edgeIndex) const
 Returns the type of the equivalent Bar element (with either 2D or 3D coordinates) for an elements edge (edgeIndex from 0 to numEdges()-1). Returns GM_INV_CELL_TYPE for line elements.
 
GmCellType faceElement (int faceIndex) const
 Returns the type of the equivalent surface element (with 3D node coordinates) for a 3D element's face (faceIndex from 0 to numFaces()-1). Returns GM_INV_CELL_TYPE for line or surface elements.
 
int edgeElementNode (int edgeIndex, int i) const
 Returns the ith node belonging to the edge edgeIndex using the equivalent linear element (returned by edgeElement()) node ordering. More...
 
int faceElementNode (int faceIndex, int i) const
 Returns the ith node belonging to the face faceIndex using the equivalent surface element (returned by faceElement()) node ordering. More...
 
GmCellType edgeLinearElement (int edgeIndex) const
 Returns the equivalent linear element of the given element's edge.
 
GmCellType faceLinearElement (int faceIndex) const
 Returns the equivalent linear element of the given element's face.
 
double length (const GmMatrix &X) const
 Returns the length of a bar element with nodes defined by the X matrix (with node coordinates organized by column). Can also be used for calculating the length of a 2D interface element (the length is calculated at the average "plane" between the interface element edges) More...
 
double area (const GmMatrix &X) const
 Returns the area of a 2D element with nodes defined by the X matrix (with node coordinates organized by column). Can also be used for calculating the area of a 3D interface element (the area is calculated at the average "plane" between the interface element faces) More...
 
double volume (const GmMatrix &X) const
 Returns the volume of a 3D element with nodes defined by the X matrix (with node coordinates organized by column) More...
 
double characteristicLength (const GmMatrix &X) const
 Returns the cell characteristic length, defined as the length for 1D and 2D interface elements, the square root of the area for non iterface 2D and 3D interface elements and the cubic root of the volume for non interface 3D elements.
 
double characteristicDimension (const GmMatrix &X) const
 Returns the cell characteristic dimension, defined as the length for 1D and 2D interface elements, the area for non iterface 2D and 3D interface elements and the volume for non interface 3D elements.
 
void centroidCartesian (const GmMatrix &X, GmVector &coord) const
 Fills the coord vector with the cartesian coordinates of the cell centroid, with nodes defined by the X matrix (with node coordinates organized by column). The coord vector will have size equal to the number of rows in the matrix. More...
 
bool hasCapability (GmCellGeometryCapabilities capability) const
 Returns true if this cell type implements the geometric related capability described by the given paramter. Unlike mesh capabilities, the query possibilities are given by the GmCellGeometryCapabilities enum iinstead of by a string (this happends since this are core capabilities, not suited to be extended by plugin provided capabilities).
 
bool isValid (const GmMatrix &X, double tol) const
 
double quality (const GmMatrix &J, double tol) const
 Returns a normalized quality measure of the cell geometry, from 0 to 1, where 0.0 means very bad and 1.0 very good. If the GM_CELL_GEOMETRY_QUALITY capability is not implemented, this function simply returns 0.0. The J matrix is the Jacobian.
 
bool contains (const GmMatrix &X, const GmVector &coord) const
 Returns true if the cell contains the point specified by the given cartesian coordinates 'coord'. If the GM_CELL_GEOMETRY_CONTAINS capability is not implemented, this function simply returns false. The X matrix holds the cell node coordinates organized by column.
 

Static Public Member Functions

static int maxNumNodes ()
 A static function that returns the maximum number of nodes among all known cell types.
 
static const char * familyTypeToStr (GmCellFamilyType family)
 Converts a cell family type to a string.
 
static int strToFamilyType (QString str)
 Converts a string to a cell family type. Returns -1 if no match is found.
 
static const char * faceTypeToStr (GmCellFaceType type)
 Converts a cell face type to a string.
 
static int strToFaceType (QString str)
 Converts a string to a cell face type. Returns -1 if no match is found.
 
static void registerCellType (const GmCellGeometryInfo *cellInfo)
 Registers a new cell type in the global geometry registry.
 
static const GmCellGeometryInfogeometryInfo (GmCellType type)
 Returns the geometry info objct associated with type. NOT FOR GENERAL USE.
 

Private Member Functions

GmShapeshapeInstance () const
 Shape function factory. Returns a NEW instance of the shape function object for this type. Returns NULL if the type has no associated shape function.
 
GmIntegrationRuleintegrationRule (GmIntegrationRuleType irType, int rule1=-1, int rule2=-1, int rule3=-1) const
 A factory function that returns a NEW integration rule object suited for this kind of element. More...
 
GmBorderIntegrationRuleedgeIntegrationRule (GmIntegrationRuleType irType, int rule1=-1) const
 A factory function that returns a NEW border integration rule object suited for this kind of element. Edge rules are used to integrate functions over an element edge. This function is undefined for line elements and will return NULL in those cases. More...
 
GmBorderIntegrationRulefaceIntegrationRule (int faceType, GmIntegrationRuleType irType, int rule1=-1, int rule2=-1) const
 A factory function that returns a NEW border integration rule object suited for this kind of element. Face rules are used to integrate functions over an element face. This function is undefined for line or 2D elements and will return NULL in those cases. More...
 
GmBorderIntegrationRuleborderIntegrationRule (int faceType, GmIntegrationRuleType irType, int rule1=-1, int rule2=-1) const
 Returns a NEW border integration rule object suited for this kind of element. Border rules are used to integrate functions over an element border. That means that for a solid element we want a rule for integrating a face and for a surface element we want a rule to integrate a line. More...
 
void checkX (const GmMatrix &X) const
 

Static Private Member Functions

static void checkMetadata (const GmCellGeometryMetadata *data, bool registryFullyInitialised)
 Performs basic consistency checks for the given metadata object Additional checks for Incidence, Edge and Face compatibility are done by checkGeometryIEFCompatibility()
 
static void checkGeometryIEFCompatibility (GmCellType type)
 Checks that the edges, faces and node incidence informations are compatible among themselves.
 
static void checkElementRegistry ()
 Perform a full consistency check for all registered elements.
 

Private Attributes

const GmCellGeometryInfo_info
 The pointer to the cell geometry info object storing the cell type metadata.
 
int _P
 
int _Q
 The P & Q parameters for hierarchical element types. 0 otherwise.
 

Static Private Attributes

static const GmCellGeometryInfo_infoRegistry [GM_NUM_CELL_TYPES]
 The global registry with each cell type geometry info. More...
 
static int _maxNumNodes = 0
 The maximum number of nodes among all registered cell types. More...
 

Friends

class GmShape
 
class GmIntegrationRule
 
class GmBorderIntegrationRule
 

Detailed Description

A class used to return static metadata information about a cell geometry, along with some methods for calculating geometric properties for a cell, given its coordinate matrix.

Creating a cell geometry object is fast, so it can be created freely on stack from the cell type. For historical reasons, it's implementation is a wrapper storing a pointer to a (unique) GmCellGeometryInfo object that contains the cell type metadata + specific geometric methods. In user code, one should always use a GmCellGeometry object instead of a GmCellGeometryInfo one.

The cell geometry object also serves as a connecting bridge between the cell type and its associated shape and integration rule classes. This is done by a set of factory functions that are used by GmShape and GmIntegration rule factory functions to associate a type with a specific shape function class or integratio rule class. This factory functions are not for general use though.

When creating new cell types, the added cell should be registered with a call to GmCellGeometry::registerCellType().

Cell Geometric descriptions:

As a general rule, while decribing the cyclic order of nodes inside edges or faces, GeMA always uses the CCW (counter clock-wise) convention and references "linear nodes" before "quadratic nodes", "quadratic nodes" before "cubic nodes", etc... This particular ordering is used to facilitate several algorithms in the topological mesh implementation and also while implementing super-parametric physics elements. Also, for 2D elements, edge numbering follows the number of its first node.

IMPORTANT: Please keep in mind that the following documentation is a C++ reference, and so node, edge and face numbering starts with 0 and not 1. When using this documentation in Lua for building simulation models, please consider that while referencing an edge or face number, you should add one to the values seen in here (i.e. the edge between the first and second nodes of a quad4 is edge 1 in Lua and not edge 0 as in C++).

Coordinate systems:

2D:

|y | | |________x

3D:

|z /y | / | / |/________x

Member Function Documentation

◆ area()

double GmCellGeometry::area ( const GmMatrix X) const
inline

Returns the area of a 2D element with nodes defined by the X matrix (with node coordinates organized by column). Can also be used for calculating the area of a 3D interface element (the area is calculated at the average "plane" between the interface element faces)

Returns 0.0 for 1D, 2D interface elements and 3D (non interface) elements.

Parameters
XThe cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true)

◆ borderIntegrationRule()

GmBorderIntegrationRule* GmCellGeometry::borderIntegrationRule ( int  faceType,
GmIntegrationRuleType  irType,
int  rule1 = -1,
int  rule2 = -1 
) const
inlineprivate

Returns a NEW border integration rule object suited for this kind of element. Border rules are used to integrate functions over an element border. That means that for a solid element we want a rule for integrating a face and for a surface element we want a rule to integrate a line.

For 2D elements it is equivalent to edgeIntegrationRule() and for 3D elements to faceIntegrationRule(). This function is undefined for line elements. See the description of the aforementioned functions for details on the parameters.

◆ centroidCartesian()

void GmCellGeometry::centroidCartesian ( const GmMatrix X,
GmVector coord 
) const
inline

Fills the coord vector with the cartesian coordinates of the cell centroid, with nodes defined by the X matrix (with node coordinates organized by column). The coord vector will have size equal to the number of rows in the matrix.

By definition, the centroid is the point defined by the arithmetic mean position of all the points of the element. If coord holds a vector created by DECLARE_REF_VECTOR(), its size MUST already be equal to the expected result size. Whenever possible, the centroid is calculated geometrically. For quadratic elements, the calculation is done by numeric integration. If the cell type is an interface element, the centroid will be calculated at the average "plane" between the interface element edges/faces.

Parameters
XThe cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true)
coordThe vector filled with the centroid cartesian coordinates

◆ edgeElementNode()

int GmCellGeometry::edgeElementNode ( int  edgeIndex,
int  i 
) const
inline

Returns the ith node belonging to the edge edgeIndex using the equivalent linear element (returned by edgeElement()) node ordering.

IMPORTANT: The result of this function is undefined if edgeElement(edgeIndex) returns GM_INV_CELL_TYPE.

Returned values are local indices. The parameter edgeIndex should be in the range 0..numEdges()-1 The parameter i should be in the range 0..GmCellGeometry(edgeElement(edgeIndex)).numNodes()-1

◆ edgeFace()

int GmCellGeometry::edgeFace ( int  edgeIndex,
int  i 
) const
inline

Returns the ith face belonging to edge edgeIndex.

Returned values are local indices according to the element definition. The parameter edgeIndex should be in the range 0..numEdges()-1 The parameter i should be in the range 0..numEdgeFaces(edgeIndex)-1

◆ edgeIntegrationRule()

GmBorderIntegrationRule* GmCellGeometry::edgeIntegrationRule ( GmIntegrationRuleType  irType,
int  rule1 = -1 
) const
inlineprivate

A factory function that returns a NEW border integration rule object suited for this kind of element. Edge rules are used to integrate functions over an element edge. This function is undefined for line elements and will return NULL in those cases.

Parameter irType defines the kind on integration rule that should be used (Gauss or Lobatto, for example). Parameter rule1 is used to pass an specific rule type for derived classes and will be ignored if irrelevant. A value of -1 means that the rule parameter is unused.

If the given parameters are different from -1 and do not express a valid rule, or irType is not valid or unimplemented for this kind of element, returns NULL. In particular, this routine will return a default (valid) rule for the element when called with irType == auto and all rule values equal to -1.

◆ edgeNode()

int GmCellGeometry::edgeNode ( int  edgeIndex,
int  i 
) const
inline

Returns the ith node belonging to edge edgeIndex.

Returned values are local indices. Extremes (i = 0 and i == numEdgeNodes()-1) are the edge vertices. The parameter edgeIndex should be in the range 0..numEdges()-1 The parameter i should be in the range 0..numEdgeNodes(edgeIndex)-1

◆ faceEdge()

int GmCellGeometry::faceEdge ( int  faceIndex,
int  i 
) const
inline

Returns the ith edge belonging to face faceIndex.

Returned values are local indices according to the element definition. The parameter faceIndex should be in the range 0..numFaces()-1 The parameter i should be in the range 0..numFaceEdges(faceIndex)-1

◆ faceElementNode()

int GmCellGeometry::faceElementNode ( int  faceIndex,
int  i 
) const
inline

Returns the ith node belonging to the face faceIndex using the equivalent surface element (returned by faceElement()) node ordering.

IMPORTANT: The result of this function is undefined if faceElement(faceIndex) returns GM_INV_CELL_TYPE.

Returned values are local indices. The parameter faceIndex should be in the range 0..numFaces()-1 The parameter i should be in the range 0..GmCellGeometry(faceElement(faceIndex)).numNodes()-1

◆ faceIntegrationRule()

GmBorderIntegrationRule* GmCellGeometry::faceIntegrationRule ( int  faceType,
GmIntegrationRuleType  irType,
int  rule1 = -1,
int  rule2 = -1 
) const
inlineprivate

A factory function that returns a NEW border integration rule object suited for this kind of element. Face rules are used to integrate functions over an element face. This function is undefined for line or 2D elements and will return NULL in those cases.

The faceType parameter can be used to identify the desired face type for elements with more than one face type (like a wedge that has 3 quadrilateral faces and 2 triangular faces). In that case, faceType should be a GmCellFaceType value. If the element has only one face type, -1 can be passed to the faceType parameter (if the value is not -1, it MUST be equal to one of the type's face type).

Parameter irType defines the kind on integration rule that should be used (Gauss or Lobatto, for example). Parameters rule1 and rule2 are used to pass an specific rule type for derived classes and will be ignored if irrelevant. A value of -1 means that the rule parameter is unused.

If the given parameters are different from -1 and do not express a valid rule, or irType is not valid or unimplemented for this kind of element, returns NULL. In particular, this routine will return a default (valid) rule for the element when called with irType == auto and all rule values equal to -1.

◆ faceNode()

int GmCellGeometry::faceNode ( int  faceIndex,
int  i 
) const
inline

Returns the ith node belonging to the face faceIndex. Since boundary nodes are always returned first by this function, it can be used to retrieve all nodes or just the boundary ones.

Returned values are local indices. The parameter faceIndex should be in the range 0..numFaces()-1 The parameter i should be in the range 0..numFaceNodes(faceIndex)-1

◆ incidenceNode()

int GmCellGeometry::incidenceNode ( int  nodeIndex,
int  i 
) const
inline

Returns the ith incidence node belonging to node nodeIndex.

Returned values are local indices. The parameter nodeIndex should be in the range 0..numNodes()-1 The parameter i should be in the range 0..numIncidenceNodes(nodeIndex)-1

◆ integrationRule()

GmIntegrationRule* GmCellGeometry::integrationRule ( GmIntegrationRuleType  irType,
int  rule1 = -1,
int  rule2 = -1,
int  rule3 = -1 
) const
inlineprivate

A factory function that returns a NEW integration rule object suited for this kind of element.

Parameter irType defines the kind on integration rule that should be used (Gauss or Lobatto, for example). Parameters rule1, rule2 and rule3 are used to pass an specific rule type for derived classes and will be ignored if irrelevant. A value of -1 means that the rule parameter is unused.

If the given parameters are different from -1 and do not express a valid rule, or irType is not valid or unimplemented for this kind of element, returns NULL. In particular, this routine will return a default (valid) rule for the element when called with irType == auto and all rule values equal to -1.

◆ isValid()

bool GmCellGeometry::isValid ( const GmMatrix X,
double  tol 
) const
inline

Performs some geometry checks, like orientation, convexity or self intersection checks to check a cell validity. If the GM_CELL_GEOMETRY_VALID capability is not implemented, this function simply returns true. The X matrix holds the cell node coordinates organized by column.

◆ length()

double GmCellGeometry::length ( const GmMatrix X) const
inline

Returns the length of a bar element with nodes defined by the X matrix (with node coordinates organized by column). Can also be used for calculating the length of a 2D interface element (the length is calculated at the average "plane" between the interface element edges)

Returns 0.0 for 2D (non interface) and 3D elements.

Parameters
XThe cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true)

◆ numFaceTypes()

int GmCellGeometry::numFaceTypes ( ) const
inline

Returns the number of face types for a 3D element (generally 1), 1 for 2D elements and 0 for 1D elements.

In general, 3D elements have a single face type (QUAD for hexahedrons, TRI for tetrahedrons), but some special elements might have heterogeneous face types, like a wedge, combining quadrilateral and triangular faces.

◆ numVertices()

int GmCellGeometry::numVertices ( ) const
inline

Returns the number of vertices of this cell type. This excludes center edge and face nodes - returns 4 for a quad8, for example. It also excludes additional dof nodes.

In practice, for linear elements, numNodes() == numVertices() + numExtraDofNodes(). Also, if the element is NOT linear, its number of vertices should be the same as the number of vertices for its equivalent linear element. At present, this is also expected for the number of extra dof nodes.

◆ volume()

double GmCellGeometry::volume ( const GmMatrix X) const
inline

Returns the volume of a 3D element with nodes defined by the X matrix (with node coordinates organized by column)

Returns 0.0 for 1D, 2D and 3D interface elements.

Parameters
XThe cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true)

◆ volumeInternalNode()

int GmCellGeometry::volumeInternalNode ( int  i) const
inline

Returns the ith internal node belonging to the volume.

Returned values are local indices. The parameter i should be in the range 0..numVolumeInternalNodes()-1

Member Data Documentation

◆ _infoRegistry

const GmCellGeometryInfo * GmCellGeometry::_infoRegistry
staticprivate

The global registry with each cell type geometry info.

The global cell geometry info registry vector. Static initialized to NULL and so already available during dynamic initialization steps where registerCellType() can be called.

◆ _maxNumNodes

int GmCellGeometry::_maxNumNodes = 0
staticprivate

The maximum number of nodes among all registered cell types.

The maximum number of nodes among all known cell types. Initialized by the cell registering process.


The documentation for this class was generated from the following files: