CEPS  24.01
Cardiac ElectroPhysiology Simulator
FEBase.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 
33 #include "common/CepsCommon.hpp"
36 
37 class FENode;
38 
61 class FEBase :
62  public ceps::HoldsGeomObject<GeomCell>,
63  public ceps::HoldsNodes<FENode>,
64  public ceps::HoldsCells<FEBase>
65 {
66 
67  public:
69  FEBase() = delete;
70 
72  FEBase(GeomCell *geom, CepsBool owned, ReferenceFE *reference);
73 
75  FEBase(const FEBase &that) = default;
76 
78  FEBase&
79  operator=(const FEBase &that);
80 
82  ~FEBase() override;
83 
84  void
85  setGeomObject(GeomCell *object, CepsBool owned = false) override;
86 
87  // --------------------------------------------------------------------------------
88  // Reference Element
89  // --------------------------------------------------------------------------------
92  getReferenceElement() const;
93 
95  CepsUInt
96  getNumberOfDofs() const;
97 
100  getBasisVertices() const;
101 
103  const CepsReal3D&
104  getBasisVertex(CepsUInt i) const;
105 
108  getOrigin() const;
109 
112  getCellType() const;
113 
115  CepsUInt
116  getPolynomialOrder() const;
117 
118  // --------------------------------------------------------------------------------
119  // Geometrical
120  // --------------------------------------------------------------------------------
121 
123  CepsReal
124  getMeasure() const;
125 
128  getJacobianMatrix() const;
129 
131  CepsReal
132  getJacobianDeterminant() const;
133 
136  getInverseJacobian() const;
137 
138  protected:
140 };
CepsCellType
Enum for different shapes of cells.
Definition: CepsEnums.hpp:128
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
CepsArray3< CepsReal > CepsReal3D
Three real scalars, used like this for compatibility in polynomials.
Definition: CepsTypes.hpp:178
Abstract class for finite elements.
Definition: FEBase.hpp:65
FEBase()=delete
Deleted default construtcor.
CepsUInt getNumberOfDofs() const
Number of DoFs.
Definition: FEBase.cpp:81
void setGeomObject(GeomCell *object, CepsBool owned=false) override
Set the geom object.
Definition: FEBase.cpp:59
const CepsVector< CepsReal3D > & getBasisVertices() const
Direct access to positionning coefficients.
Definition: FEBase.cpp:87
FEBase & operator=(const FEBase &that)
Copy assignement.
CepsMathVertex getOrigin() const
Access to coordinate of the first node (the origin node) as eigen object.
Definition: FEBase.cpp:99
const GeomCell::JacobianMatrixType & getJacobianMatrix() const
Jacobian matrix from the geometrical element.
Definition: FEBase.cpp:127
~FEBase() override
Destructor.
Definition: FEBase.cpp:54
CepsUInt getPolynomialOrder() const
Polynomial order.
Definition: FEBase.cpp:111
FEBase(const FEBase &that)=default
Copy constructor.
CepsReal getJacobianDeterminant() const
Jacobian determinant from the geometrical element.
Definition: FEBase.cpp:133
CepsReal getMeasure() const
Size of the element.
Definition: FEBase.cpp:121
const CepsReal3D & getBasisVertex(CepsUInt i) const
Direct access to positionning coefficients.
Definition: FEBase.cpp:93
const GeomCell::InverseJacobianMatrixType & getInverseJacobian() const
Inverse of the jacobian matrix from the geometrical element.
Definition: FEBase.cpp:139
ReferenceFE * getReferenceElement() const
Pointer on Reference element.
Definition: FEBase.cpp:75
ReferenceFE * m_reference
Link to reference element.
Definition: FEBase.hpp:139
CepsCellType getCellType() const
Cell type.
Definition: FEBase.cpp:105
A nodal point on a finite element. It is different from a geom node as it may have different properti...
Definition: FENode.hpp:46
Abstract class for geometrical cell. On top of index and attributes managament, the cell has informat...
Definition: GeomCell.hpp:48
Eigen::Matrix< CepsReal, 3, Eigen::Dynamic > JacobianMatrixType
Typedef for jaccobian matrix.
Definition: GeomCell.hpp:53
Eigen::Matrix< CepsReal, Eigen::Dynamic, 3 > InverseJacobianMatrixType
Typedef for inverse jaccobian matrix.
Definition: GeomCell.hpp:55
Base class for reference finite elements.
Definition: ReferenceFE.hpp:48
An abstract class for objects that regroup pointers to cells (eg a mesh, a finite elements discretiza...
Definition: HoldsCells.hpp:46
For objects that have pointers to either a node or a cell.