CEPS  24.01
Cardiac ElectroPhysiology Simulator
GeomCell.hpp
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 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~*/
30 #pragma once
31 
33 
40 class GeomCell :
42  public ceps::HoldsDimension,
43  public ceps::HoldsBoundary,
44  public ceps::HoldsAttributes,
45  public ceps::HoldsProcIds,
46  public ceps::HoldsNodes<GeomNode>,
47  public ceps::HoldsCells<GeomCell>
48 {
49 
50  public:
51 
53  using JacobianMatrixType = Eigen::Matrix<CepsReal, 3, Eigen::Dynamic>;
55  using InverseJacobianMatrixType = Eigen::Matrix<CepsReal, Eigen::Dynamic, 3>;
56 
58  GeomCell();
59 
61  explicit GeomCell(CepsVector<GeomNode*> nodes, const CepsCellGlobalIndex& gID=0);
62 
64  GeomCell(const GeomCell &) = default;
65 
67  GeomCell(GeomCell &&) = default;
68 
70  GeomCell&
71  operator=(const GeomCell &) = default;
72 
74  GeomCell&
75  operator=(GeomCell &&) = default;
76 
78  ~GeomCell() override = default;
79 
80  // -------------------------------------------------------------------------
81  // Nodes manipulation
82  // -------------------------------------------------------------------------
83 
85  GeomNode*
87 
89  const CepsCellType&
90  getCellType() const;
91 
93  virtual void
94  reset();
95 
96  // --------------------------------------------------------------------------------
97  // Jacobian
98  // --------------------------------------------------------------------------------
99 
101  CepsReal
102  getJacobianDeterminant() const;
103 
106  virtual const JacobianMatrixType&
108 
111  virtual const InverseJacobianMatrixType&
113 
115  virtual void
116  refresh();
117 
119  void
121 
124  getNormal() const;
125 
126  friend class TestGeomCells;
127 
129  void
131 
135 
137  CepsReal
138  getDiameter();
139 
140 
141  protected:
142 
144  void
145  initialize(CepsUInt dim, CepsCellType cellType, CepsUInt nNodes);
146 
148  virtual void
149  initializeNodes(GeomNode *const *nodes, CepsUInt nNodes);
150 
152  void
153  computeDiameter();
154 
155  // ===================================================================
156  // pure virtual funcs
157 
158  public:
159 
161  virtual CepsReal
162  getMeasure() const = 0;
163 
165  virtual CepsUInt
166  getNumberOfNodes() const = 0;
167 
168 
169  protected:
170 
172  virtual void
174 
176  virtual void
178 
180  virtual void
182 
184  virtual CepsBool
186 
187 
188  protected:
189 
198 };
199 
CepsCellType
Enum for different shapes of cells.
Definition: CepsEnums.hpp:128
CepsGlobalIndex CepsCellGlobalIndex
Indices of cells.
Definition: CepsTypes.hpp:222
std::set< _Type, _Compare, _Alloc > CepsSet
C++ set.
Definition: CepsTypes.hpp:209
std::vector< _Type, _Alloc > CepsVector
C++ vector.
Definition: CepsTypes.hpp:155
bool CepsBool
Booleans.
Definition: CepsTypes.hpp:124
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
Eigen::Matrix< CepsScalar, 3, 1 > CepsMathVertex
Vertex, eigen format.
Definition: CepsTypes.hpp:135
CepsStandardArgs getStandardArgsFrom(GeomCell *cell)
Returns a CepsStandardArgs structure with data from given cell.
Definition: GeomCell.cpp:203
Abstract class for geometrical cell. On top of index and attributes managament, the cell has informat...
Definition: GeomCell.hpp:48
JacobianMatrixType m_jacobian
Jacobian Matrix.
Definition: GeomCell.hpp:191
void setNodeSharingCellsIndices(const CepsSet< CepsCellGlobalIndex > &nCellIds)
Cells that share at least one node with self.
Definition: GeomCell.cpp:150
InverseJacobianMatrixType m_invJacobian
Inverse of Jacobian Matrix.
Definition: GeomCell.hpp:192
virtual CepsUInt getNumberOfNodes() const =0
Volume for 3D, area for 2D, length for 1D, 0 for points.
CepsBool m_invJcomputed
Flag to trigger computation of invJ.
Definition: GeomCell.hpp:194
void initialize(CepsUInt dim, CepsCellType cellType, CepsUInt nNodes)
internal init method
Definition: GeomCell.cpp:170
Eigen::Matrix< CepsReal, 3, Eigen::Dynamic > JacobianMatrixType
Typedef for jaccobian matrix.
Definition: GeomCell.hpp:53
virtual void computeJacobianMatrix()=0
Compute geometrical Jacobian matrix.
const CepsCellType & getCellType() const
Get the type of cell.
Definition: GeomCell.cpp:72
CepsMathVertex getNormal() const
Sets the normal vector, (useful for boundary cells only)
Definition: GeomCell.cpp:144
GeomCell(GeomCell &&)=default
copy constructor
virtual CepsReal getMeasure() const =0
Volume for 3D, area for 2D, length for 1D, 0 for points.
Eigen::Matrix< CepsReal, Eigen::Dynamic, 3 > InverseJacobianMatrixType
Typedef for inverse jaccobian matrix.
Definition: GeomCell.hpp:55
friend class TestGeomCells
Definition: GeomCell.hpp:126
virtual void reset()
Wipes content.
Definition: GeomCell.cpp:78
GeomCell()
default constructor
Definition: GeomCell.cpp:32
GeomCell & operator=(GeomCell &&)=default
assignment operator
CepsReal getJacobianDeterminant() const
Jacobian of geometrical deformation from reference cell.
Definition: GeomCell.cpp:95
~GeomCell() override=default
destructor
virtual CepsBool checkJacobianDeterminant()=0
Tells if determinant is negative, tries to swap points to make it positive.
const CepsSet< CepsCellGlobalIndex > & getNodeSharingCellsIndices() const
Cells that share at least one node with self.
Definition: GeomCell.cpp:156
GeomCell & operator=(const GeomCell &)=default
assignment operator
virtual void initializeNodes(GeomNode *const *nodes, CepsUInt nNodes)
internal init method
Definition: GeomCell.cpp:181
CepsReal m_detJacobian
determinant of Jacobian matrix
Definition: GeomCell.hpp:193
virtual const JacobianMatrixType & getJacobianMatrix()
Jacobian Matrix of deformation from reference cell (3 rows for R^3, element dimension columns)
Definition: GeomCell.cpp:101
GeomCell(const GeomCell &)=default
copy constructor
CepsReal getDiameter()
Size of cell.
Definition: GeomCell.cpp:162
CepsSet< CepsCellGlobalIndex > m_nodeSharingCells
IDs of cells that have a node in common.
Definition: GeomCell.hpp:197
virtual const InverseJacobianMatrixType & getInverseJacobianMatrix()
Jacobian Matrix of deformation from reference cell (3 rows for R^3, element dimension columns)
Definition: GeomCell.cpp:111
CepsReal m_diameter
Size of cell.
Definition: GeomCell.hpp:196
virtual void refresh()
Recomputes jacobian matrix and determinant.
Definition: GeomCell.cpp:124
virtual void computeInverseJacobianMatrix()=0
Compute inverse of geometrical Jacobian matrix.
CepsMathVertex m_normal
For boundary cells, outward normal vector.
Definition: GeomCell.hpp:195
void setNormal(CepsMathVertex n)
Sets the normal vector, (useful for boundary cells only)
Definition: GeomCell.cpp:138
virtual void computeJacobianDeterminant()=0
Compute geometrical Jacobian determinant.
void computeDiameter()
Computes size of cell (max edge length)
Definition: GeomCell.cpp:190
CepsCellType m_cellType
Cell type.
Definition: GeomCell.hpp:190
GeomNode * getNodeWithMinimalIndex() const
Returns node with smallest index (used eg to determine ownership)
Definition: GeomCell.cpp:51
Base class for nodes used in meshes.
Definition: GeomNode.hpp:53
An abstract class from which all objects that contain region attributes should derive.
Abstract class to describe if an entity is on a boundary or not Objects that can be located in a mesh...
An abstract class for objects that regroup pointers to cells (eg a mesh, a finite elements discretiza...
Definition: HoldsCells.hpp:46
Abstract class for objects that have a dimensionality (0D to 3D)
Abstract class for objects that have a global index.
Abstract class for objects that contain a CPU Id.
Structure used to pass arguments to SAFunc (see pde directory) The flags of the SAFunc allows extract...
Definition: CepsTypes.hpp:239