![]() |
GemaCoreLib
The GeMA Core library
|
Base interface for mesh cells. More...
#include <gmCell.h>
Public Member Functions | |
virtual | ~GmCell () |
Virtual destructor. | |
virtual int | cellId () const =0 |
Returns a number identifying the cell. It should be possible to use this id as an index to GmCellMesh::cell() or to get values from a cell attribute accessor. | |
virtual GmCellMesh * | mesh () const =0 |
Returns the "father" mesh for this cell. | |
virtual bool | active () const =0 |
Is this cell active ? When looping over the cell list, this function should be probably called to ensure that only active cells are processed. Alternatively, use the GmForeachActiveCell() macro. | |
virtual void | pushProxy (lua_State *L, const GmLogCategory &logger)=0 |
Pushes a proxy object for the current cell onto the Lua stack. More... | |
virtual GmCellType | type () const =0 |
Returns the cell type. | |
const char * | typeStr () const |
Returns a name describing the cell type. | |
virtual GmCellGeometry | geometry () const |
Returns an object that contains detailed cell geometry information. | |
virtual int | numNodes () const =0 |
Returns the number of nodes of this cell. Equivalent to typeToNumNodes(type()) but usually much more efficient;. | |
virtual int | numGhostNodes () const |
Returns the number of ghost nodes of this cell, if any. | |
virtual int | totalNumNodes () const |
Returns the total number of nodes of the cell, including normal nodes & ghost nodes. | |
virtual int | nodeIndex (int localIndex) const =0 |
Returns the global node index used to locate node coordinates for each of the cell nodes. More... | |
virtual void | nodes (int *nodeList, bool ghost) const =0 |
Fills the vector nodeList with the set of cell nodes, including or not eventual ghost nodes depending on the boolean ghost parameter. The vector size must be equal to numNodes() if ghost is false or equal to totalNumNodes() if true. Returned indices follow the same rule detailed for nodeIndex() with respect to ghost nodes. | |
virtual int | propertyIndex (int propertySet) const =0 |
Given a property set number, returns the line of the property set that contains property values for this cell. The returned value can be used to get a property from a property accessor retrieved from this mesh. More... | |
virtual void | setActive (bool active)=0 |
Marks the cell as active or inactive. Inactive cells are not processed by several of the algorithms, such as analysis processes. | |
virtual bool | setNodes (const int *nodeList) |
Function allowing editable meshes to set the node list of newlly created cells. More... | |
virtual int | addGhostNode (int globalIndex) |
Adds a new ghost node to the cell, returning the new node local index or -1 if there was a problem adding the node. More... | |
virtual void | removeGhostNode (int localIndex) |
Removes a ghost node from the cell definition given its local index. More... | |
virtual bool | setProperties (const int *propList, int nprop) |
Function allowing editable meshes to set the property indices of newlly created cells. More... | |
GmCell * | adjacentCell (int sideIndex, int *adjSideIndex=NULL) const |
Returns the adjacent cell sharing the given side, represented by a cell edge index for surface meshes or by a cell face index for solid mehses. More... | |
double | length (const GmValueAccessor *coordAccessor, int edgeIndex=0, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Returns the length of a 1D element, the length of a 2D interface element, or the length of the given edge for other element types, with node coordinates given by the provided accessor. More... | |
double | length (const GmMatrix &X, int edgeIndex=0, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Returns the length of a 1D element, the length of a 2D interface element, or the length of the given edge for other element types, with node coordinates defined by the X matrix (with node coordinates organized by column) More... | |
double | area (const GmValueAccessor *coordAccessor, int faceIndex=0, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Returns the area of a 2D (non interface) element, the area of a 3D interface element, or the area of the given face for 3D (non interface) element types, with node coordinates given by the provided accessor. More... | |
double | area (const GmMatrix &X, int faceIndex=0, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Returns the area of a 2D (non interface) element, the area of a 3D interface element, or the area of the given face for 3D (non interface) element types, with node coordinates defined by the X matrix (with node coordinates organized by column) More... | |
double | volume (const GmValueAccessor *coordAccessor, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Returns the volume of a 3D element with node coordinates given by the provided accessor. More... | |
double | volume (const GmMatrix &X, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Returns the volume of a 3D element with node coordinates defined by the X matrix (with node coordinates organized by column) More... | |
double | characteristicLength (const GmValueAccessor *coordAccessor, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) 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. More... | |
double | characteristicLength (const GmMatrix &X, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) 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. More... | |
double | characteristicDimension (const GmValueAccessor *coordAccessor, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) 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. More... | |
double | characteristicDimension (const GmMatrix &X, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) 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. More... | |
void | centroidCartesian (const GmValueAccessor *coordAccessor, GmVector &coord, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Fills the coord vector with the cartesian coordinates for the cell centroid. More... | |
void | centroidCartesian (const GmMatrix &X, GmVector &coord, GmCellGeometryMode linearMode=GM_CELLGEOM_AUTO) const |
Fills the coord vector with the cartesian coordinates for the cell centroid. More... | |
virtual bool | isValid (const GmValueAccessor *nodeAccessor, double tol=1e-5) const |
Checks if cell is geometrically valid tol is the relative tolerance for face planar test. More... | |
virtual double | quality (const GmValueAccessor *nodeAccessor, double tol=1e-5) 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 cell geometry type does not implements the quality capability, this function will always return 0.0. | |
virtual bool | contains (const GmValueAccessor *nodeAccessor, const GmVector &coord) const |
Returns true if the cell contains the point specified by the given cartesian coordinates 'coord'. If the cell geometry type does not implements the contains capability, this function will always return false. | |
void | fillNodeMatrix (const GmValueAccessor *nodeAccessor, GmMatrix &m, bool transposed=false, GmCellFillMode mode=GM_CELL_NODES) const |
Given a node accessor, fills the matrix 'm' with the coordinates of the cell nodes. More... | |
void | fillFaceNodeMatrix (const GmValueAccessor *nodeAccessor, GmMatrix &m, int faceIndex, bool transposed=false, GmCellFillMode mode=GM_CELL_NODES) const |
Similar to fillNodeMatrix() but filling m with the node coordinates for the given 3D element face. Nodes are orderer according to the equivalent face element type, as returned by GmCellGeomtry::faceElement() or faceLinearElement() if mode is GM_CELL_VERTICES. NOT defined if the face type is GM_INV_CELL_TYPE (this includes all surface elements). More... | |
void | fillEdgeNodeMatrix (const GmValueAccessor *nodeAccessor, GmMatrix &m, int edgeIndex, bool transposed=false, GmCellFillMode mode=GM_CELL_NODES) const |
Similar to fillNodeMatrix() but filling m with the node coordinates for the given element edge. Nodes are orderer according to the equivalent edge element type, as returned by GmCellGeomtry::edgeElement() or edgeLinearElement() if mode is GM_CELL_VERTICES. NOT defined if the edge type is GM_INV_CELL_TYPE (this includes all bar elements). More... | |
void | fillDeformedNodeMatrix (const GmValueAccessor *nodeAccessor, const GmValueAccessor *uAccessor, GmMatrix &m, GmVector &u, bool transposed=false, GmCellFillMode mode=GM_CELL_NODES) const |
void | fillDeformedFaceNodeMatrix (const GmValueAccessor *nodeAccessor, const GmValueAccessor *uAccessor, GmMatrix &m, GmVector &u, int faceIndex, bool transposed=false, GmCellFillMode mode=GM_CELL_NODES) const |
Similar to fillDeformedNodeMatrix() but filling m & u with the node coordinates / deformation for the given element face. Nodes are orderer according to the. More... | |
void | fillDeformedEdgeNodeMatrix (const GmValueAccessor *nodeAccessor, const GmValueAccessor *uAccessor, GmMatrix &m, GmVector &u, int edgeIndex, bool transposed=false, GmCellFillMode mode=GM_CELL_NODES) const |
Similar to fillDeformedNodeMatrix() but filling m & u with the node coordinates / deformation for the given element edge. Nodes are orderer according to the equivalent edge element type, as returned by GmCellGeomtry::edgeElement() or edgeLinearElement() if mode is GM_CELL_VERTICES. More... | |
virtual int | meshId () const =0 |
Returns the meshId associated with mesh(). Ugly. Needed by the GeMA mesh implementation. | |
virtual void | replaceCellId (int id, bool keepActiveFlag) |
Updates the cell id. For those who really know what they are doing. This function is a very specific function that should not be called lightly. It exists on behalf of the state dumping code and for helping other mesh implementations wisshng to encode extra iformation on the cell id. More... | |
virtual void | replaceGhostNodes (int *ghostNodes, int numNodes) |
Replaces the full set of ghost nodes of the cell. This function is a very specific function that should not be called lightly. It exists on behalf of the state dumping code. More... | |
Static Public Member Functions | |
static int | typeToNumNodes (GmCellType type) |
Returns the number of nodes associated to the cell type. | |
static const char * | typeToStr (GmCellType type) |
Converts a cell type to a string. | |
static int | strToType (QString str) |
Converts a string to a cell type. Returns -1 if no match is found. | |
Base interface for mesh cells.
Mesh cells represent a "patch" of the domain covered by a mesh. The GmMeshCell interface provides methods to discover the list of nodes that define the cell and methods to query the set of properties associated to the cell.
Cell properties are divided into property sets. Each property set is like a table, containing multiple properties (columns). Each row is a set of property values that can be associated with mesh cells.
Each cell stores the index of a row in a property set. Multiple cells can point to the same row in a property set table.
Each mesh can be associated with more then one property set, so each cell contains a value (row) index for each property set that is associated with its parent mesh. The list of property sets associated with a cell is indexed from 0 to GmCellMesh::propertySets().size() - 1. The mapping between this local index and a property set object is done by the mesh. The cell function propertyIndex() returns the "row" index of the cell for the specified property set, definied by it's "local" index.
To use this class, a mesh plugin should inherit from GmCell and implement the functions needed to recover node information and property values, as can be seen in the example class from the standard mesh plugin.
|
inlinevirtual |
Adds a new ghost node to the cell, returning the new node local index or -1 if there was a problem adding the node.
The added ghost node should already belong to the mesh. The given global index is a number from 0 to mesh->numGhostNodes()-1. It accepts both values WITH or WITHOUT the highest bit set. See comments on nodeIndex().
This function should be reimplemented by mesh implementations that support ghost nodes (hasCapability("ghostNodes") returns true) and should probably call GmCellMesh::ghostNodesUpdated().
IMPORTANT: Keep in mind that no checking is done to see if the given node geometrically belongs to the cell.
GmCell * GmCell::adjacentCell | ( | int | sideIndex, |
int * | adjSideIndex = NULL |
||
) | const |
Returns the adjacent cell sharing the given side, represented by a cell edge index for surface meshes or by a cell face index for solid mehses.
The sideIndex is a value between 0 and the number of element edges/face - 1, defining the position of the edge/face of interest in the cell. Returns NULL if the edge/face belongs to the mesh boundary, if the mesh does NOT support the topology capability or if the topological data structure was not yet built and could not be created. If you need to distinguish between those cases, consider using the GmCellMeshTopology::adjacentCell() method directly.
If adjSideIndex is different from NULL, it will be filled with the side index of the adjacent edge/face in the adjacent cell.
double GmCell::area | ( | const GmValueAccessor * | coordAccessor, |
int | faceIndex = 0 , |
||
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Returns the area of a 2D (non interface) element, the area of a 3D interface element, or the area of the given face for 3D (non interface) element types, with node coordinates given by the provided accessor.
Depending on the element type, the area is either calculated by geometric relations or by numeric integration using the default integration rule for the element / face Returns 0.0 for 1D, and 2D interface elements.
coordAccessor | The mesh coordinate accessor |
faceIndex | The face index when calculating the area of the face of a 3D (non interface) element. Should be 0 otherwise. |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. |
double GmCell::area | ( | const GmMatrix & | X, |
int | faceIndex = 0 , |
||
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Returns the area of a 2D (non interface) element, the area of a 3D interface element, or the area of the given face for 3D (non interface) element types, with node coordinates defined by the X matrix (with node coordinates organized by column)
Depending on the element type, the area is either calculated by geometric relations or by numeric integration using the default integration rule for the element / face Returns 0.0 for 1D, and 2D interface elements.
For interface elements, the area is calculated at the mid plane (which only really makes a difference if the X matrix is filled taking into account deformations).
X | The cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true) IMPORTANT: X expects to be the full element node matrix, even if we are calculating the area of a single 3D element face. |
faceIndex | The face index when calculating the area of the face of a 3D (non interface) element. Should be 0 otherwise. |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. Do not use a linear mode if deformations where added to the X matrix! |
void GmCell::centroidCartesian | ( | const GmValueAccessor * | coordAccessor, |
GmVector & | coord, | ||
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Fills the coord vector with the cartesian coordinates for the cell centroid.
If coord holds a vector created by DECLARE_REF_VECTOR(), its size MUST already be equal to the expected result size.
coordAccessor | The mesh coordinate accessor |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. |
void GmCell::centroidCartesian | ( | const GmMatrix & | X, |
GmVector & | coord, | ||
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Fills the coord vector with the cartesian coordinates for the cell centroid.
If coord holds a vector created by DECLARE_REF_VECTOR(), its size MUST already be equal to the expected result size.
X | The cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true) |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. Do not use a linear mode if deformations where added to the X matrix! |
double GmCell::characteristicDimension | ( | const GmValueAccessor * | coordAccessor, |
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | 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.
coordAccessor | The mesh coordinate accessor |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. |
double GmCell::characteristicDimension | ( | const GmMatrix & | X, |
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | 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.
X | The cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true) |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. Do not use a linear mode if deformations where added to the X matrix! |
double GmCell::characteristicLength | ( | const GmValueAccessor * | coordAccessor, |
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | 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.
coordAccessor | The mesh coordinate accessor |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. |
double GmCell::characteristicLength | ( | const GmMatrix & | X, |
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | 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.
X | The cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true) |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. Do not use a linear mode if deformations where added to the X matrix! |
void GmCell::fillDeformedEdgeNodeMatrix | ( | const GmValueAccessor * | nodeAccessor, |
const GmValueAccessor * | uAccessor, | ||
GmMatrix & | m, | ||
GmVector & | u, | ||
int | edgeIndex, | ||
bool | transposed = false , |
||
GmCellFillMode | mode = GM_CELL_NODES |
||
) | const |
Similar to fillDeformedNodeMatrix() but filling m & u with the node coordinates / deformation for the given element edge. Nodes are orderer according to the equivalent edge element type, as returned by GmCellGeomtry::edgeElement() or edgeLinearElement() if mode is GM_CELL_VERTICES.
See comments for fillEdgeNodeMatrix() & fillDeformedNodeMatrix() for the remaining parameters.
void GmCell::fillDeformedFaceNodeMatrix | ( | const GmValueAccessor * | nodeAccessor, |
const GmValueAccessor * | uAccessor, | ||
GmMatrix & | m, | ||
GmVector & | u, | ||
int | faceIndex, | ||
bool | transposed = false , |
||
GmCellFillMode | mode = GM_CELL_NODES |
||
) | const |
Similar to fillDeformedNodeMatrix() but filling m & u with the node coordinates / deformation for the given element face. Nodes are orderer according to the.
equivalent face element type, as returned by GmCellGeomtry::faceElement() or faceLinearElement() if mode is GM_CELL_VERTICES.
See comments for fillFaceNodeMatrix() & fillDeformedNodeMatrix() for the remaining parameters.
void GmCell::fillEdgeNodeMatrix | ( | const GmValueAccessor * | nodeAccessor, |
GmMatrix & | m, | ||
int | edgeIndex, | ||
bool | transposed = false , |
||
GmCellFillMode | mode = GM_CELL_NODES |
||
) | const |
Similar to fillNodeMatrix() but filling m with the node coordinates for the given element edge. Nodes are orderer according to the equivalent edge element type, as returned by GmCellGeomtry::edgeElement() or edgeLinearElement() if mode is GM_CELL_VERTICES. NOT defined if the edge type is GM_INV_CELL_TYPE (this includes all bar elements).
See comments for fillNodeMatrix() for the remaining parameters. IMPORTANT: This function does NOT handle ghost nodes, so "mode" must be either GM_CELL_NODES or GM_CELL_VERTICES. Also, this function is NOT defined if the face has an unknown edge element type.
void GmCell::fillFaceNodeMatrix | ( | const GmValueAccessor * | nodeAccessor, |
GmMatrix & | m, | ||
int | faceIndex, | ||
bool | transposed = false , |
||
GmCellFillMode | mode = GM_CELL_NODES |
||
) | const |
Similar to fillNodeMatrix() but filling m with the node coordinates for the given 3D element face. Nodes are orderer according to the equivalent face element type, as returned by GmCellGeomtry::faceElement() or faceLinearElement() if mode is GM_CELL_VERTICES. NOT defined if the face type is GM_INV_CELL_TYPE (this includes all surface elements).
See comments for fillNodeMatrix() for the remaining parameters. IMPORTANT: This function does NOT handle ghost nodes, so "mode" must be either
GM_CELL_NODES or GM_CELL_VERTICES. Also, this function is NOT defined if the face has an unknown face element type.
void GmCell::fillNodeMatrix | ( | const GmValueAccessor * | nodeAccessor, |
GmMatrix & | m, | ||
bool | transposed = false , |
||
GmCellFillMode | mode = GM_CELL_NODES |
||
) | const |
Given a node accessor, fills the matrix 'm' with the coordinates of the cell nodes.
The result 'm' wil be filled with an "n x d" matrix where 'n' is the number of nodes in the cell and 'd' is the number of dimensions of each coordinate, thus line 'i' holds the coordinates of node 'i'.
If transposed == true, the matrix will be transposed with size "d x n", with column 'i' holding the coordinates of a node 'i'.
The mode parameter controls which nodes will be added to the matrix. The default GM_CELL_NODES fills it with common geometry nodes. The GM_CELL_VERTICES mode fills the matrix with coordinates for the element vertices + any additional extra dof nodes, GM_CELL_GHOST fills it with ghost nodes only and the GM_CELL_ALL mode fills the matrix with both geometry and ghost nodes.
IMPORTANT: If m holds a matrix created by DECLARE_REF_MATRIX(), its size MUST already be equal to the expected result size.
|
virtual |
Checks if cell is geometrically valid tol is the relative tolerance for face planar test.
Tries to check whether the cell has a valid geometry or not. For every cell type, checks that all node indices are valid and do not contain duplicate nodes. Depending on the cell geomety implementing the valid capability, other checks, like orientation, convexity or self intersection checks, might also be performed.
double GmCell::length | ( | const GmValueAccessor * | coordAccessor, |
int | edgeIndex = 0 , |
||
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Returns the length of a 1D element, the length of a 2D interface element, or the length of the given edge for other element types, with node coordinates given by the provided accessor.
Depending on the element type, the length is either calculated by geometric relations or by numeric integration using the default integration rule for the element / edge
coordAccessor | The mesh coordinate accessor |
edgeIndex | The edge index when calculating the length of an edge. Should be 0 otherwise. |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. |
double GmCell::length | ( | const GmMatrix & | X, |
int | edgeIndex = 0 , |
||
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Returns the length of a 1D element, the length of a 2D interface element, or the length of the given edge for other element types, with node coordinates defined by the X matrix (with node coordinates organized by column)
Depending on the element type, the length is either calculated by geometric relations or by numeric integration using the default integration rule for the element / edge
For interface elements, the length is calculated at the "mid plane" (which only really makes a difference if the X matrix is filled taking into account deformations).
X | The cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true) IMPORTANT: X expects to be the full element node matrix, even if we are calculating the length of a single element edge. |
edgeIndex | The edge index when calculating the length of an edge. Should be 0 otherwise. |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. Do not use a linear mode if deformations where added to the X matrix! |
|
pure virtual |
Returns the global node index used to locate node coordinates for each of the cell nodes.
Local indices are numbered from 0 to numNodes()-1 for common cell nodes, i.e., for those nodes that take part in the cell geometry definition. Ghost nodes, on the other hand are numbered from numNodes() to totalNumNodes()-1.
For the set of 'geometry' nodes, their precise local order in the cell depends on cell type. For the set of 'ghost' nodes, their order follows the order in which ghost nodes where added to the cell.
The returned global index is a number from 0 to mesh->numNodes()-1 if localIndex is a 'geometry' index and a value from 0 to mesh->numGhostNodes()-1 WITH the highest bit set if localIndex is a 'ghost' index (that means that the first two ghost nodes have, respectively, an index equal to 0x800000 and 0x80000001 if node indexes are stored in 32 bit integers)
|
pure virtual |
Given a property set number, returns the line of the property set that contains property values for this cell. The returned value can be used to get a property from a property accessor retrieved from this mesh.
The propertySet number is a local to the mesh identification of the property sets that are bound to the mesh. This id varies between 0 and GmCellMesh::propertySets().size() - 1. See class comments for more details.
|
pure virtual |
Pushes a proxy object for the current cell onto the Lua stack.
The proxy object should be a sub class of GmLuaCell and can include specific functions not in the standard cell interfaces. This interface is based on pushing the object on the Lua stack and not by returning a GmLuaCell object since the later option would create a circular dependency between the GemaCore library and the GemaLuaCore library.
|
inlinevirtual |
Removes a ghost node from the cell definition given its local index.
This operation does NOT remove the ghost node from the mesh, only from the cell defintion.
This function should be reimplemented by mesh implementations that support ghost nodes (hasCapability("ghostNodes") returns true).
|
inlinevirtual |
Updates the cell id. For those who really know what they are doing. This function is a very specific function that should not be called lightly. It exists on behalf of the state dumping code and for helping other mesh implementations wisshng to encode extra iformation on the cell id.
This function can preserve the current active bit depending on the keepActiveFlag parameter. Should be used to code other information on the id high bits. The "base id", coinciding with the cell vector location, should be the same.
|
inlinevirtual |
Replaces the full set of ghost nodes of the cell. This function is a very specific function that should not be called lightly. It exists on behalf of the state dumping code.
The indices in the list should be values from 0 to mesh->numGhostNodes()-1 WITH the highest bit set. Its call does NOT calls GmCellMesh::ghostNodesUpdated().
This function should be reimplemented by mesh implementations that support both ghost nodes and state dumping (hasCapability("ghostNodes") returns true and hasCapability("stateDump") also) and are implemented based on the standard Gema mesh.
|
inlinevirtual |
Function allowing editable meshes to set the node list of newlly created cells.
The nodeList parameter should contain a list with size equal to numNodes() and referencing the nodes of the element in the correct order.
Returns true on success, false on failure (invalid data, for example).
|
inlinevirtual |
Function allowing editable meshes to set the property indices of newlly created cells.
The propList parameter should contain a list with size equal to GmCellMesh::propertySets().size() and referencing valid lines in the appropriate property sets. The parameter nprop should be equal to the list size. Property values are copied.
Returns true on success, false on failure (invalid data, for example).
double GmCell::volume | ( | const GmValueAccessor * | coordAccessor, |
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Returns the volume of a 3D element with node coordinates given by the provided accessor.
Depending on the element type, the volume is either calculated by geometric relations or by numeric integration using the default integration rule for the element Returns 0.0 for 1D, 2D and 3D interface elements.
coordAccessor | The mesh coordinate accessor |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. |
double GmCell::volume | ( | const GmMatrix & | X, |
GmCellGeometryMode | linearMode = GM_CELLGEOM_AUTO |
||
) | const |
Returns the volume of a 3D element with node coordinates defined by the X matrix (with node coordinates organized by column)
Depending on the element type, the volume is either calculated by geometric relations or by numeric integration using the default integration rule for the element Returns 0.0 for 1D, 2D and 3D interface elements.
X | The cell's node matrix. Can be created by calling GmCell::fillNodeMatrix(coordAc, X, true) |
linearMode | When equal to GM_CELLGEOM_LINEAR, assumes that the geometry of the element is linear, even if the cell itself is quadratic. When set to GM_CELLGEOM_AUTO, gets that information from the cell's mesh. Do not use a linear mode if deformations where added to the X matrix! |