CEPS  24.01
Cardiac ElectroPhysiology Simulator
Mesh.cpp
Go to the documentation of this file.
1 /*~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2  This file is part of CEPS.
3 
4  CEPS is free software: you can redistribute it and/or modify
5  it under the terms of the GNU General Public License as published by
6  the Free Software Foundation, either version 3 of the License, or
7  (at your option) any later version.
8 
9  CEPS is distributed in the hope that it will be useful,
10  but WITHOUT ANY WARRANTY; without even the implied warranty of
11  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  GNU General Public License for more details.
13 
14  You should have received a copy of the GNU General Public License
15  along with CEPS (see file LICENSE at root of project).
16  If not, see <https://www.gnu.org/licenses/>.
17 
18 
19  Copyright 2019-2024 Inria, Universite de Bordeaux
20 
21  Authors, in alphabetical order:
22 
23  Pierre-Elliott BECUE, Florian CARO, Yves COUDIERE(*), Andjela DAVIDOVIC,
24  Charlie DOUANLA-LONTSI, Marc FUENTES, Mehdi JUHOOR, Michael LEGUEBE(*),
25  Pauline MIGERDITICHAN, Valentin PANNETIER(*), Nejib ZEMZEMI.
26  * : currently active authors
27 
28 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~*/
31 #include "geometry/Mesh.hpp"
32 
34  m_dim (dim),
35  m_geom (nullptr),
36  m_nbCells (0),
37  m_nbBoundaryCells (0),
38  m_nbNodes (0),
39  m_connectivityComputed (false),
40  m_maxCellDiameter (0.e0),
41  m_maxNodeConnectivity (0)
42 {
43 }
44 
46 {
51 
56 }
57 
60 {
61  return m_dim;
62 }
63 
64 void
66 {
67  CEPS_ABORT_IF(not c,"adding a nullptr cell to the mesh");
69  "trying to add a cell of dimension " << c->getDimension()
70  << "\nin a mesh of dimension " << m_dim
71  );
72 
75  "Cell #" << cID << " was already added to the mesh"
76  );
77 
78  // Register the cell
79  m_cells.append(c,getCellHash(c),cID);
80 
82 
83  return;
84 }
85 
86 void
88 {
89  CEPS_ABORT_IF(not c, "adding a nullptr cell to the mesh");
90  CEPS_ABORT_IF(c->getDimension () != m_dim-1,
91  "trying to add a cell of dimension " << c->getDimension ()
92  << "\nin a mesh of dimension " << m_dim
93  );
94 
97  "Cell #" << cID << " was already added to the mesh"
98  );
99 
100  // Register the cell
101  m_bCells.append(c,getCellHash(c),cID);
103  return;
104 }
105 
107 {
108  CEPS_ABORT_IF(not n, "adding a nullptr node to the mesh");
111  "trying to add node " << nID << " that is already added to the mesh"
112  );
113 
114  // Register the node
115  m_ownedNodes.append(n,getNodeHash(n),nID);
116  return;
117 }
118 
119 void
121 {
122  CEPS_ABORT_IF(not n,
123  "adding a nullptr node to the mesh"
124  );
127  "trying to add node " << nID << " that is already added to the mesh"
128  );
129  // Register the node
130  m_haloNodes.append(n,getNodeHash(n),nID);
131  return;
132 }
133 
134 void
136 {
137  for (GeomCell *&cell : m_cells.all())
138  cell->refresh ();
139  for (GeomCell *&boundaryCell : m_bCells.all())
140  boundaryCell->refresh();
141  return;
142 }
143 
144 void
146 {
147  m_nbCells = n;
148 }
149 
150 void
152 {
153  m_nbBoundaryCells = n;
154 }
155 
157 {
158  m_nbNodes = n;
159 }
160 
161 GeomNode*
163 {
164  return m_ownedNodes.atGlobal(geomId,nullptr);
165 }
166 
167 GeomNode*
169 {
170  return m_haloNodes.atGlobal(geomId,nullptr);
171 }
172 
173 GeomNode*
175 {
176  return m_ownedNodes.atGlobal(geomId,getHaloNode(geomId));
177 }
178 
179 CepsInt
181 {
182  return m_cells.hasGlobal(geomId) ? m_cells.globalToLocal(geomId) : -1;
183 }
184 
185 CepsInt
187 {
188  return m_bCells.hasGlobal(geomId) ? m_bCells.globalToLocal(geomId) : -1;
189 }
190 
191 GeomCell*
193 {
194  return m_cells.atGlobal(geomId,nullptr);
195 }
196 
197 GeomCell*
199 {
200  return m_bCells.atGlobal(geomId,nullptr);
201 }
202 
205 {
206  return m_cells.all();
207 }
208 
211 {
212  return m_bCells.all();
213 }
214 
217 {
218  return m_ownedNodes.all();
219 }
220 
223 {
224  return m_haloNodes.all();
225 }
226 
227 CepsUInt
229 {
230  return m_nbCells;
231 }
232 
233 CepsUInt
235 {
236  return m_nbBoundaryCells;
237 }
238 
239 CepsUInt
241 {
242  return m_nbNodes;
243 }
244 
245 CepsUInt
247 {
248  return m_cells.size();
249 }
250 
251 CepsUInt
253 {
254  return m_bCells.size();
255 }
256 
257 CepsUInt
259 {
260  return m_ownedNodes.size ();
261 }
262 
263 CepsUInt
265 {
266  return m_haloNodes.size();
267 }
268 
269 CepsUInt
271 {
273  return m_maxNodeConnectivity;
274  return this->computeMaxNodeConnectivity();
275 }
276 
277 CepsReal
279 {
280  return m_maxCellDiameter;
281 }
282 
283 void
285 {
286  m_geom = geom;
287 }
288 
289 Geometry*
291 {
292  return m_geom;
293 }
294 
295 void
296 Mesh::scale(CepsReal scaleFactor)
297 {
298  for (GeomNode* node: m_ownedNodes.all())
299  node->scale(scaleFactor);
300  for (GeomNode* node: m_haloNodes.all())
301  node->scale(scaleFactor);
302  // Recompute all jacobians
303  refreshCells ();
304 }
305 
306 CepsHash
308 {
309  return n->getGlobalIndex();
310 }
311 
312 CepsHash
314 {
315  return c->getGlobalIndex();
316 }
317 
318 CepsUInt
320 {
322 
323  // iterate through our nodes and halo node
324  for (GeomNode *node : m_ownedNodes.all())
325  m_maxNodeConnectivity = std::max (m_maxNodeConnectivity, node->getNumberOfNodes());
326  for (GeomNode *node : m_haloNodes.all())
327  m_maxNodeConnectivity = std::max (m_maxNodeConnectivity, node->getNumberOfNodes());
328 
329  m_connectivityComputed = true;
330  return m_maxNodeConnectivity;
331 }
332 
333 std::ostream&
334 operator<<(std::ostream& os, const Mesh& m)
335 {
336  os << "Cells: #" << m.m_cells.size() << std::endl;
337 
338  for (GeomCell *el : m.m_cells.all())
339  {
340  os << el->getGlobalIndex () << " - nodes : ";
341  for (GeomNode *node : el->getNodes())
342  os << node->getGlobalIndex () << " ";
343  os << std::endl;
344  }
345  os << CEPS_STRING_SEPARATOR << std::endl;
346 
347  os << "Boundary Cells: #" << m.m_bCells.size() << std::endl;
348  for (GeomCell *el : m.m_bCells.all())
349  {
350  os << el->getGlobalIndex() << " - nodes : ";
351  for (GeomNode *node : el->getNodes ())
352  os << node->getGlobalIndex() << " ";
353  os << std::endl;
354  }
355 
356  os << CEPS_STRING_SEPARATOR << std::endl;
357 
358  os << "Owned Nodes: #" << m.m_ownedNodes.size() << std::endl;
359  for (GeomNode *node : m.m_ownedNodes.all())
360  os << node->getGlobalIndex() << " - coordinates : " << *node << std::endl;
361 
362  os << "Halo Nodes: #" << m.m_haloNodes.size() << std::endl;
363  for (GeomNode *node : m.m_haloNodes.all())
364  os << node->getGlobalIndex() << " - coordinates : " << *node << std::endl;
365 
366  return os;
367 }
#define CEPS_ABORT_IF(condition, message)
Stops the execution with a message if condition is true. If testing is enabled, only throws a runtime...
#define CEPS_STRING_SEPARATOR
CepsString used as separator.
CepsGlobalIndex CepsCellGlobalIndex
Indices of cells.
Definition: CepsTypes.hpp:222
CepsSize CepsHash
Hashes for distributed data.
Definition: CepsTypes.hpp:220
std::vector< _Type, _Alloc > CepsVector
C++ vector.
Definition: CepsTypes.hpp:155
std::make_unsigned_t< CepsInt > CepsUInt
Unsigned version on CepsInt.
Definition: CepsTypes.hpp:109
float CepsReal
Need single precision floating point.
Definition: CepsTypes.hpp:100
CepsGlobalIndex CepsNodeGlobalIndex
Indices of nodes.
Definition: CepsTypes.hpp:224
int32_t CepsInt
Need 32 bit integer.
Definition: CepsTypes.hpp:106
std::ostream & operator<<(std::ostream &os, const Mesh &m)
Displays some info on the mesh.
Definition: Mesh.cpp:334
Abstract class for geometrical cell. On top of index and attributes managament, the cell has informat...
Definition: GeomCell.hpp:48
CepsReal getDiameter()
Size of cell.
Definition: GeomCell.cpp:162
Base class for nodes used in meshes.
Definition: GeomNode.hpp:53
Encapsulates all the geometrical data.
Definition: Geometry.hpp:50
CepsUInt size() const
Number of elements.
const CepsVector< _Type > & all() const
access to data
const CepsLocalIndex & globalToLocal(const CepsGlobalIndex &gId) const
Global index to local index conversion.
_Type & atGlobal(const CepsGlobalIndex &gId, const _Type &def)
Access value with global ID, with default value.
void destroyObjects()
Wipes content, keeps structure.
CepsBool hasGlobal(const CepsGlobalIndex &gId) const
Tells if item with global index gID can be found on this CPU.
void append(const _Type &x, const _Hash &hashValue)
Emplace an element in the map.
Geometrical information of 1,2 or 3D distributed cells.
Definition: Mesh.hpp:57
CepsReal m_maxCellDiameter
Largest cell size.
Definition: Mesh.hpp:250
void addBoundaryCell(GeomCell *boundaryCell)
Add a boundary cell to the mesh.
Definition: Mesh.cpp:87
CepsInt getLocalIndexOfBoundaryCell(CepsCellGlobalIndex geomIndex) const
Geom global to local mapping for boundary cells.
Definition: Mesh.cpp:186
const CepsVector< GeomNode * > & getOwnedNodes()
CepsVector of local nodes stored on this process.
Definition: Mesh.cpp:216
CepsUInt getNbBoundaryCells() const
Total number of boundary cells of the mesh.
Definition: Mesh.cpp:234
CepsUInt getMaxNodeConnectivity()
Local maximum node connectivity (nb of adjacent nodes)
Definition: Mesh.cpp:270
CepsUInt m_nbNodes
whole mesh nodes count
Definition: Mesh.hpp:248
CepsUInt computeMaxNodeConnectivity()
Locally computes the maximum node connectivity (nb of adjacent nodes)
Definition: Mesh.cpp:319
LocalGlobalMapping< GeomNode * > m_haloNodes
Nodes in halo.
Definition: Mesh.hpp:260
CepsUInt getLocalNbCells() const
Number of non-boundary cells stored on this process.
Definition: Mesh.cpp:246
void addHaloNode(GeomNode *node)
Add a halo node to the mesh.
Definition: Mesh.cpp:120
void refreshCells()
Forces re-evaluation of each cell structure (boundary elmts included).
Definition: Mesh.cpp:135
CepsHash getNodeHash(GeomNode *n)
gets an identifier for the node (here its global index)
Definition: Mesh.cpp:307
CepsUInt getNbNodes() const
Total number of nodes of the mesh.
Definition: Mesh.cpp:240
CepsBool m_connectivityComputed
computation flag
Definition: Mesh.hpp:249
CepsUInt getDimension() const
Geometrical dimension.
Definition: Mesh.cpp:59
void setNbCells(CepsUInt n)
Set the total number of cells of this mesh.
Definition: Mesh.cpp:145
void setNbNodes(CepsUInt n)
Set the total number of nodes of this mesh.
Definition: Mesh.cpp:156
GeomNode * getHaloNode(CepsNodeGlobalIndex geomIndex) const
Get node referenced by geom index.
Definition: Mesh.cpp:168
Geometry * getGeometry() const
Geometry this mesh is tied too.
Definition: Mesh.cpp:290
LocalGlobalMapping< GeomCell * > m_bCells
Boundary cells (dim-1)
Definition: Mesh.hpp:258
void scale(CepsReal scaleFactor)
Each node coordinate is scaled by scaleFactor.
Definition: Mesh.cpp:296
const CepsVector< GeomCell * > & getCells()
CepsVector of cells stored on this process.
Definition: Mesh.cpp:204
const CepsVector< GeomNode * > & getHaloNodes()
CepsVector of halo nodes stored on this process.
Definition: Mesh.cpp:222
~Mesh()
Destructor.
Definition: Mesh.cpp:45
Mesh(CepsUInt dim)
Default constructor.
Definition: Mesh.cpp:33
CepsUInt m_maxNodeConnectivity
maximum nb adjacent nodes.
Definition: Mesh.hpp:262
CepsUInt getNbCells() const
Total number of non-boundary cells of the mesh.
Definition: Mesh.cpp:228
void addNode(GeomNode *node)
Add a node to the mesh.
Definition: Mesh.cpp:106
void addCell(GeomCell *cell)
Add an cell to the mesh.
Definition: Mesh.cpp:65
CepsHash getCellHash(GeomCell *c)
gets an identifier for the node, (here its global index)
Definition: Mesh.cpp:313
CepsReal getMaxCellDiameter() const
Maximum size of cells.
Definition: Mesh.cpp:278
CepsUInt getLocalNbNodes() const
Number of nodes stored on this process.
Definition: Mesh.cpp:258
Geometry * m_geom
Link to geometry.
Definition: Mesh.hpp:245
CepsUInt getLocalNbBoundaryCells() const
Number of boundary cells stored on this process.
Definition: Mesh.cpp:252
void setGeometry(Geometry *geom)
Link geometry.
Definition: Mesh.cpp:284
GeomCell * getCell(CepsCellGlobalIndex geomIndex) const
pointer on requested cell, nullptr if not owned
Definition: Mesh.cpp:192
LocalGlobalMapping< GeomNode * > m_ownedNodes
Nodes of that proc.
Definition: Mesh.hpp:259
GeomCell * getBoundaryCell(CepsCellGlobalIndex geomIndex) const
pointer on requested boundary cell, nullptr if not owned
Definition: Mesh.cpp:198
CepsUInt getLocalNbHaloNodes() const
Number of halo nodes stored on this process.
Definition: Mesh.cpp:264
GeomNode * getNodeOrHaloNode(CepsNodeGlobalIndex geomIndex) const
Get node referenced by geom index.
Definition: Mesh.cpp:174
CepsUInt m_dim
Dimension of the cells.
Definition: Mesh.hpp:243
void setNbBoundaryCells(CepsUInt n)
Set the total number of boundary cells of this mesh.
Definition: Mesh.cpp:151
CepsUInt m_nbCells
whole mesh cells count
Definition: Mesh.hpp:246
CepsUInt m_nbBoundaryCells
whole mesh boundary cells count
Definition: Mesh.hpp:247
GeomNode * getNode(CepsNodeGlobalIndex globalIndex) const
Get node referenced by global index.
Definition: Mesh.cpp:162
const CepsVector< GeomCell * > & getBoundaryCells()
CepsVector of boundary cells stored on this process.
Definition: Mesh.cpp:210
LocalGlobalMapping< GeomCell * > m_cells
Cells of main dimension.
Definition: Mesh.hpp:257
CepsInt getLocalIndexOfCell(CepsCellGlobalIndex geomIndex) const
Geom global to local mapping for cells.
Definition: Mesh.cpp:180
const CepsUInt & getDimension() const
Get the dimension of the object.
const CepsGlobalIndex & getGlobalIndex() const
Get the index
CepsVector< NodeType * > & getNodes()
Reference to the node pointers array.
const _Type & max(const CepsVector< _Type, _Alloc > &vec)
Returns the maximum of the vector, const version.
Definition: CepsVector.hpp:464
void destroyObject(_Type &)
Destroy[delete] any type.
Definition: CepsMemory.hpp:116