CEPS  24.01
Cardiac ElectroPhysiology Simulator
ReferenceFE.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 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~*/
31 #pragma once
32 
35 
48 {
49 
50 
51  public:
53  ReferenceFE() = delete;
54 
56  ReferenceFE(CepsUInt dimension, CepsFeLagrangeType type);
57 
59  ReferenceFE(const ReferenceFE& that) = default;
60 
62  ReferenceFE &
63  operator=(const ReferenceFE& that) = default;
64 
66  ~ReferenceFE();
67 
70  getCellType() const;
71 
73  CepsUInt
74  getPolynomialOrder() const;
75 
77  CepsUInt
78  getNumberOfDofs() const;
79 
82  computeBasisFunctions(const CepsReal3D& point);
83 
87 
90  getBasisVertices() const;
91 
93  const CepsReal3D&
94  getBasisVertex(CepsUInt i) const;
95 
98  CepsInt
100 
105  GeomNode*
107 
116  CepsReal3D
118 
120  CepsUInt
121  getNumberOfGeomVertices() const;
122 
124  CepsUInt
125  getDimension() const;
126 
128  CepsUInt
129  getNumberOfBasisVertices() const;
130 
132  CepsUInt
134 
137  getFEType() const;
138 
140  CepsReal
141  computeBasisFunction(CepsUInt iPhi, const CepsReal3D& point);
142 
144  CepsReal3D
146 
147 
148  protected:
149 
151  void
153 
154 
155  protected :
156 
166 
168 
172 };
CepsCellType
Enum for different shapes of cells.
Definition: CepsEnums.hpp:128
Eigen::Matrix< CepsScalar, Eigen::Dynamic, 1 > CepsMathDynamic1D
Dynamic 1D array, eigen format.
Definition: CepsTypes.hpp:139
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
CepsArray3< CepsReal > CepsReal3D
Three real scalars, used like this for compatibility in polynomials.
Definition: CepsTypes.hpp:178
int32_t CepsInt
Need 32 bit integer.
Definition: CepsTypes.hpp:106
Eigen::Matrix< CepsScalar, Eigen::Dynamic, Eigen::Dynamic > CepsMathDynamic2D
Dynamic 2D array, eigen format.
Definition: CepsTypes.hpp:140
std::pair< CepsCellType, CepsUInt > CepsFeLagrangeType
Typedef for descriptor of Lagrange finite elements type.
Definition: FEType.hpp:35
Abstract class for geometrical cell. On top of index and attributes managament, the cell has informat...
Definition: GeomCell.hpp:48
Base class for nodes used in meshes.
Definition: GeomNode.hpp:53
Base class for reference finite elements.
Definition: ReferenceFE.hpp:48
CepsReal3D computeBasisFunctionDerivative(CepsUInt iPhi, const CepsReal3D &point)
Get values of single basis function derivative at given point.
CepsUInt m_order
Order of Lagrangian base polynomial.
CepsUInt getNumberOfBasisFunctions() const
Get number of basis functions.
CepsUInt getNumberOfDofs() const
Number of degrees of freedom.
Definition: ReferenceFE.cpp:60
CepsUInt m_nGeomVertices
Nb of geom vertices the elem is made of.
CepsUInt m_dim
Geometrical dimension.
CepsUInt getDimension() const
Get dimension of geom object.
CepsReal3D getTransformedBasisVertex(CepsUInt i, GeomCell *cell) const
Get vertex 'i' used to create basis functions and transform it on GeomCell.
CepsUInt getPolynomialOrder() const
Order.
Definition: ReferenceFE.cpp:54
CepsInt getBasisVertexGeomCellIndex(CepsUInt i) const
Returns the index in geom cell of a basis vertex. -1 if basis vertex does not match a geom vertex.
ReferenceFE & operator=(const ReferenceFE &that)=default
Copy assignement.
void buildBasisFunction()
Build all basis function and Vertex used.
CepsVector< CepsInt > m_geomCellIndex
CepsVector< Polynomial< 3U > > m_basisFunctions
List of basis functions.
CepsMathDynamic1D computeBasisFunctions(const CepsReal3D &point)
Get values of basis functions at given point.
Definition: ReferenceFE.cpp:66
CepsVector< CepsSet< CepsUInt > > m_ownerDeductionRule
Dirty: map that contains rules used to attribute ownership of nodes that do not coincide with geometr...
CepsUInt getNumberOfBasisVertices() const
Get number of vertices used to create basis functions.
CepsVector< CepsReal3D > m_vertices
List of vertex used to build basis functions.
ReferenceFE()=delete
Deleted default constructor.
CepsCellType getCellType() const
CellType.
Definition: ReferenceFE.cpp:48
CepsFeLagrangeType m_type
CepsFeLagrangeType.
const CepsVector< CepsReal3D > & getBasisVertices() const
Get vertices used to create basis functions.
Definition: ReferenceFE.cpp:92
~ReferenceFE()
Destructor.
Definition: ReferenceFE.cpp:44
CepsFeLagrangeType getFEType() const
The type of Lagragian Finite Element return. May change in the future ?
CepsUInt getNumberOfGeomVertices() const
Get number of vertices of the geometric cell the element is based on.
GeomNode * getNodeOfOwnerForBasisVertex(CepsUInt i, GeomCell *cell) const
Returns a node that determines ownership of a basis vertex, especially for those that do not match a ...
CepsReal computeBasisFunction(CepsUInt iPhi, const CepsReal3D &point)
Get values of single basis function at given point.
ReferenceFE(const ReferenceFE &that)=default
Copy constructor.
const CepsReal3D & getBasisVertex(CepsUInt i) const
Get vertex 'i' used to create basis functions.
Definition: ReferenceFE.cpp:98
CepsMathDynamic2D computeBasisFunctionDerivatives(const CepsReal3D &point)
Get values of basis function derivatives at given point.
Definition: ReferenceFE.cpp:77
CepsUInt m_nbDofs
Number of degrees of freedom.