GemaMesh
The GeMA Mesh Plugin
uibhmGeometry.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 _UIBHM_GEOMETRY_H_
25 #define _UIBHM_GEOMETRY_H_
26 
27 #include "uibhmTemplate.h"
28 
29 #include <gmpGemaMeshData.h>
30 #include <gmpGemaCellMeshData.h>
31 #include <gmpGemaCellMesh.h>
32 
33 
43 template <bool Hybrid, template <class> class Vector>
44 class GMP_GEMAMESH_API_EXPORT UibhmGemaCellMeshDataGeometry
45 {
46 public:
47 
49  template <typename U> using VectorType = Vector<U>;
50 
51  // Casting a cell from GmCell* to RealCell* works based on two (ugly but true) assumptions:
52  // 1) That our cells are implemented by GmpGemaCell or by GmpGemaGhostCell. This is already implicit by
53  // working with a GmpGemaCellMeshData object
54  // 2) That the size and memory layout of a cell object class independs on the base cell class (GmCell or
55  // GmElement), celltype and number of nodes template parameters as further detailed and tested in
56  // GmpGemaCellMeshData::cellMemory().
57  // We might make this a little less ugly by removing the cell type and number of nodes temaplate
58  // parameters from the base cell class and adding them to a last level in the inheritance chain for a
59  // cell!
60  // IMPORTANT: obviously, any type or number of nodes dependent method can NOT be called on a cell casted
61  // to a RealCell.
63 
66 
68  enum { IsHybrid = Hybrid };
69 
70  // The basic constructor as expected by UibhmTopology
72  : UibhmGemaCellMeshDataGeometry(dynamic_cast<GmElementMesh*>(mesh) ? ((ElemMeshBase*)mesh)->theNodeData() : ((CellMeshBase*)mesh)->theNodeData(),
73  dynamic_cast<GmElementMesh*>(mesh) ? ((ElemMeshBase*)mesh)->theCellData() : ((CellMeshBase*)mesh)->theCellData())
74  {}
75 
77 
78  bool updateCellTypes(int numAddedCells);
79 
81  int numNodes() const { return _nodeData->_numNodes; }
82 
84  int numCells() const { return _cellData->_numCells; }
85 
87  GmCell* cell(int cellId) { assert(cellId >= 0 && cellId < _cellData->_numCells); return _cellData->_cells[cellId]; }
88 
90  bool solid() const { return _solid; }
91 
96  bool hybrid() const { return Hybrid; }
97 
103  bool strange() const { return _strange; }
104 
109  const UibhmGemaTemplate* const* typeTemplates() const { return _typeTemplateList; }
110 
111  void faceVertices(int cellId, int face, QVarLengthArray<int, 10>& nodeList) const;
112  void faceNodes (int cellId, int face, bool mid, QVarLengthArray<int, 10>& nodeList) const;
113 
120  void faceEdgeVertices(int cellId, int face, int edge, int* first, int* last) const
121  {
122  const UibhmGemaTemplate* t = tt(cellId);
123  assert(t);
124  t->faceEdgeVertices(face, edge, first, last, cn(cellId));
125  }
126 
133  void faceEdgeNodes(int cellId, int face, int edge, QVarLengthArray<int, 4>& nodeList) const
134  {
135  const UibhmGemaTemplate* t = tt(cellId);
136  assert(t);
137  t->faceEdgeNodes(face, edge, nodeList, cn(cellId));
138  }
139 
140  int findNodeIndex(int cellId, int v) const;
141  int findFaceEdge (int cellId, int face, int v1, int v2) const;
142  int findFaceEdge (int cellId, int face, int* v1, int* v2) const;
143 
149  void edgeVertices(int cellId, int edge, int* first, int* last) const
150  {
151  const UibhmGemaTemplate* t = tt(cellId);
152  assert(t);
153  t->edgeVertices(edge, first, last, cn(cellId));
154  }
155 
161  void edgeNodes(int cellId, int edge, QVarLengthArray<int, 4>& nodeList) const
162  {
163  const UibhmGemaTemplate* t = tt(cellId);
164  assert(t);
165  t->edgeNodes(edge, nodeList, cn(cellId));
166  }
167 
169  int numSides(int cellId) const { const UibhmGemaTemplate* t = tt(cellId); return t ? t->numSides() : 0; }
170 
171  void sideVertices(int cellId, int side, QVarLengthArray<int, 10>& nodeList) const;
172 
176  int midCellNode(int cellId) const
177  {
178  const UibhmGemaTemplate* t = tt(cellId);
179  if(t && t->midNode() != -1)
180  return cn(cellId)[t->midNode()];
181  return -1;
182  }
183 
189  bool quadraticMidEdgeNodeAdjacentVertices(int cellId, int localNode, int* v1, int* v2) const
190  {
191  const UibhmGemaTemplate* t = tt(cellId);
192  assert(t);
193  return t->quadraticMidEdgeNodeAdjacentVertices(localNode, v1, v2, cn(cellId));
194  }
195 
199  void vertexAdjacentVertices(int cellId, int localNode, QVarLengthArray<int, 4>& nodeList) const
200  {
201  const UibhmGemaTemplate* t = tt(cellId);
202  assert(t);
203  t->vertexAdjacentVertices(localNode, nodeList, cn(cellId));
204  }
205 
211  void vertexAdjacentNodes(int cellId, int localNode, QVarLengthArray<int, 4>& nodeList) const
212  {
213  const UibhmGemaTemplate* t = tt(cellId);
214  assert(t);
215  t->vertexAdjacentNodes(localNode, nodeList, cn(cellId));
216  }
217 
219  static size_t encodeEdge(int v1, int v2)
220  {
221  if(v1 > v2) qSwap(v1, v2);
222  return (((size_t)v1) << 32) | ((size_t)v2);
223  }
224 
229  static void decodeEdge(size_t edge, int* v1, int* v2)
230  {
231  *v1 = (int)(edge >> 32);
232  *v2 = (int)(edge & 0xFFFFFFFF);
233  }
234 
236  const UibhmGemaTemplate* tt(int cellId) const
237  {
238  assert(cellId >= 0 && cellId < numCells());
239  if constexpr(Hybrid)
240  return _typeTemplateList[_cellData->_cells[cellId]->type()];
241  else
242  {
243  Q_UNUSED(cellId);
244  return _typeTemplate;
245  }
246  }
247 
249  int co(int cellId) const { return ((RealCell*)_cellData->_cells[cellId])->offset(); }
250 
252  const int* cn(int cellId) const { return _cellData->cellNodesPtr(co(cellId)); }
253 
256  bool _solid;
257  bool _strange;
259 
264  const UibhmGemaTemplate* _typeTemplateList[GM_NUM_CELL_TYPES];
265 };
266 
267 
268 #endif
void vertexAdjacentNodes(int cellId, int localNode, QVarLengthArray< int, 4 > &nodeList) const
Given a local vertex number (a geometry node, NOT a quadratic node), returns the set of global nodes ...
Definition: uibhmGeometry.h:211
void vertexAdjacentVertices(int cellId, int localNode, QVarLengthArray< int, 4 > &nodeList) const
Given a local vertex number (a geometry node, NOT a quadratic node), returns the set of global vertic...
Definition: uibhmGeometry.h:199
GmpGemaMeshData< Vector > * _nodeData
The mesh node data object.
Definition: uibhmGeometry.h:254
int numSides(int cellId) const
Returns the number of sides(edges or faces) for the given element. Returns 0 for strange elements.
Definition: uibhmGeometry.h:169
GM_NUM_CELL_TYPES
GmpGemaCellMeshData< Vector > * _cellData
The mesh cell data object.
Definition: uibhmGeometry.h:255
int numNodes() const
Return the number of mesh nodes.
Definition: uibhmGeometry.h:81
bool _strange
Flag set to true if the mesh includes "strange" cells.
Definition: uibhmGeometry.h:257
const UibhmGemaTemplate * _typeTemplate
The type template if the mesh has a single cell type that is not a strange type (not hybrid,...
Definition: uibhmGeometry.h:258
void edgeNodes(int e, QVarLengthArray< int, 4 > &list) const
Returns the list of edge nodes for the given edge definied by its local element number.
Definition: uibhmTemplate.h:177
static size_t encodeEdge(int v1, int v2)
Encodes the edge definde by v1, v2 into a 64 bits integer number.
Definition: uibhmGeometry.h:219
int midCellNode(int cellId) const
If the element has a mid surface or a mid volume quadratic node, returns the node index....
Definition: uibhmGeometry.h:176
Vector< U > VectorType
The vector type used by the mesh.
Definition: uibhmGeometry.h:49
bool hybrid() const
Returns true if this mesh contains (or supports) hybrid cell types. It will always return true if the...
Definition: uibhmGeometry.h:96
int numSides() const
Returns numEdges() for surface elements and numFaces() for solid elements.
Definition: uibhmTemplate.h:69
GM_BAR2
Declaration of the GmpGemaCellMeshData structure.
const UibhmGemaTemplate * tt(int cellId) const
Returns the template type for the given cell. Will return NULL for "strange" cells.
Definition: uibhmGeometry.h:236
GmCell * cell(int cellId)
Returns the cell object from its id.
Definition: uibhmGeometry.h:87
Declaration of the GmpGemaCellMesh class.
Basic implementation for a "Gema cell", adding the information about the offset of its nodes in the g...
Definition: gmpGemaCell.h:37
A class used to store the type information needed by the Uibhm structure. The data is created based o...
Definition: uibhmTemplate.h:45
Declaration of the GmpGemaMeshData structure.
int co(int cellId) const
Returns the offset of the first node stored by the given cell.
Definition: uibhmGeometry.h:249
void edgeVertices(int cellId, int edge, int *first, int *last) const
Returns the pair of vertices for the given element's edge, where the edge is defined by its local ind...
Definition: uibhmGeometry.h:149
void vertexAdjacentNodes(int v, QVarLengthArray< int, 4 > &list) const
Given a local vertex number, returns the set of local nodes adjacent to the given one....
Definition: uibhmTemplate.h:165
int numCells() const
Return the number of mesh cells (including any 'strange' elements)
Definition: uibhmGeometry.h:84
const UibhmGemaTemplate *const * typeTemplates() const
Returns a vector with GM_NUM_CELL_TYPES entries storing the templates for the used cell types....
Definition: uibhmGeometry.h:109
Declaration of the UibhmGemaTemplate class.
const int * cn(int cellId) const
Returns the cell node list.
Definition: uibhmGeometry.h:252
void edgeNodes(int cellId, int edge, QVarLengthArray< int, 4 > &nodeList) const
Returns the list of nodes for the given element's edge, where the edge is defined by its local index ...
Definition: uibhmGeometry.h:161
Auxiliar structure used to share data between GmpGemaCellMesh and GmpMeshLoader.
Definition: gmpGemaCellMeshData.h:36
bool strange() const
Returns true if the mesh includes "strange" cells. A surface mesh has strange cells if it includes ba...
Definition: uibhmGeometry.h:103
void vertexAdjacentVertices(int v, QVarLengthArray< int, 4 > &list) const
Given a local vertex number, returns the set of local vertices adjacent to the given one....
Definition: uibhmTemplate.h:153
static void decodeEdge(size_t edge, int *v1, int *v2)
Decodes the result of an encodeEdge() operation returning both nodes. Since encoded values are always...
Definition: uibhmGeometry.h:229
void faceEdgeNodes(int cellId, int face, int edge, QVarLengthArray< int, 4 > &nodeList) const
Returns the list of nodes for the given element's face edge, where the face is defined by its local i...
Definition: uibhmGeometry.h:133
void edgeVertices(int e, int *first, int *last) const
Returns the pair of edge vertices for the given edge definied by its local element number.
Definition: uibhmTemplate.h:171
A class used to provide all the geometric information needed by the topology and query Uibhm classes....
Definition: uibhmGeometry.h:44
bool quadraticMidEdgeNodeAdjacentVertices(int cellId, int localNode, int *v1, int *v2) const
Given a mid-edge quadratic node defined by its local cell index, returns the pair of (global) vertice...
Definition: uibhmGeometry.h:189
bool solid() const
Returns true if the mesh includes any solid cell (a cell with more than one face)
Definition: uibhmGeometry.h:90
int midNode() const
For quadratic elements with a mid face (surface element) or mid volume (solid element) node,...
Definition: uibhmTemplate.h:109
Auxiliar structure used to share data between GmpGemaMesh and GmpMeshLoader.
Definition: gmpGemaMeshData.h:31
bool quadraticMidEdgeNodeAdjacentVertices(int node, int *v1, int *v2) const
Given a mid-edge quadratic node defined by its local cell index v, returns the pair of vertices of th...
Definition: uibhmTemplate.h:140
void faceEdgeVertices(int f, int e, int *first, int *last) const
Returns the pair of edge vertices for the given face, local (face) edge pair.
Definition: uibhmTemplate.h:183
void faceEdgeNodes(int f, int e, QVarLengthArray< int, 4 > &list) const
Returns the set of edge nodes for the given face, local (face) edge pair. The list endpoints are the ...
Definition: uibhmTemplate.h:193
bool _solid
Flag set to true if the mesh includes solid cells.
Definition: uibhmGeometry.h:256
void faceEdgeVertices(int cellId, int face, int edge, int *first, int *last) const
Returns the pair of vertices for the given element's face edge, where the face is defined by its loca...
Definition: uibhmGeometry.h:120