CEPS  24.01
Cardiac ElectroPhysiology Simulator
FEBase.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 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~*/
32 
33 FEBase::FEBase(GeomCell *cell, CepsBool owned, ReferenceFE *reference) :
34  ceps::HoldsGeomObject<GeomCell>(),
35  ceps::HoldsNodes<FENode>(),
36  ceps::HoldsCells<FEBase>(),
37  m_reference(reference)
38 {
40  "Cannot create finite element as ptr to geometrical element is null"
41  );
42  CEPS_ABORT_IF(ceps::isNullPtr(reference),
43  "Cannot create finite element as ptr to reference element is null"
44  );
45  // Check compatibility of reference element and geom cell
47  or m_reference->getDimension() != cell->getDimension(),
48  "Cannot create finite element as the passed geometric cell and reference element\n" <<
49  " do not share the same geometrical properties"
50  );
51  this->setGeomObject(cell,owned);
52 }
53 
55 {
56 }
57 
58 void
60 {
62  "Cannot create finite element as ptr to geometrical element is null"
63  );
65  this->setNodes({});
66  this->setCells({});
67  return;
68 }
69 
70 // --------------------------------------------------------------------------------
71 // Reference Element
72 // --------------------------------------------------------------------------------
73 
76 {
77  return m_reference;
78 }
79 
82 {
84 }
85 
88 {
90 }
91 
92 const CepsReal3D &
94 {
96 }
97 
100 {
102 }
103 
106 {
107  return getReferenceElement()->getCellType();
108 }
109 
110 CepsUInt
112 {
114 }
115 
116 // --------------------------------------------------------------------------------
117 // Geometrical
118 // --------------------------------------------------------------------------------
119 
120 CepsReal
122 {
123  return getGeomObject()->getMeasure();
124 }
125 
128 {
129  return getGeomObject()->getJacobianMatrix();
130 }
131 
132 CepsReal
134 {
136 }
137 
140 {
142 }
CepsCellType
Enum for different shapes of cells.
Definition: CepsEnums.hpp:128
#define CEPS_ABORT_IF(condition, message)
Stops the execution with a message if condition is true. If testing is enabled, only throws a runtime...
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
CepsMathVertex getCoordinatesForEigen() const
Get three coordinates.
Definition: CepsVertex.cpp:186
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
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
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
virtual CepsUInt getNumberOfNodes() const =0
Volume for 3D, area for 2D, length for 1D, 0 for points.
Eigen::Matrix< CepsReal, 3, Eigen::Dynamic > JacobianMatrixType
Typedef for jaccobian matrix.
Definition: GeomCell.hpp:53
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
CepsReal getJacobianDeterminant() const
Jacobian of geometrical deformation from reference cell.
Definition: GeomCell.cpp:95
virtual const JacobianMatrixType & getJacobianMatrix()
Jacobian Matrix of deformation from reference cell (3 rows for R^3, element dimension columns)
Definition: GeomCell.cpp:101
virtual const InverseJacobianMatrixType & getInverseJacobianMatrix()
Jacobian Matrix of deformation from reference cell (3 rows for R^3, element dimension columns)
Definition: GeomCell.cpp:111
Base class for reference finite elements.
Definition: ReferenceFE.hpp:48
CepsUInt getNumberOfDofs() const
Number of degrees of freedom.
Definition: ReferenceFE.cpp:60
CepsUInt getDimension() const
Get dimension of geom object.
CepsUInt getPolynomialOrder() const
Order.
Definition: ReferenceFE.cpp:54
CepsCellType getCellType() const
CellType.
Definition: ReferenceFE.cpp:48
const CepsVector< CepsReal3D > & getBasisVertices() const
Get vertices used to create basis functions.
Definition: ReferenceFE.cpp:92
CepsUInt getNumberOfGeomVertices() const
Get number of vertices of the geometric cell the element is based on.
const CepsReal3D & getBasisVertex(CepsUInt i) const
Get vertex 'i' used to create basis functions.
Definition: ReferenceFE.cpp:98
void setCells(const CepsVector< FEBase * > &cells)
Replaces currently held cell pointers.
const CepsUInt & getDimension() const
Get the dimension of the object.
_GeomObject * getProperties() const
Get the properties (in fact it's just the geom object)
virtual void setGeomObject(_GeomObject *object, CepsBool owned)
Set the geom object.
GeomCell * getGeomObject() const
Get the geom object.
FENode * getNodeAt(const CepsUInt &i)
Pointer to the i-th node.
void setNodes(const CepsVector< FENode * > &nodes)
Assign all nodes.
A namespace for all utility methods.
CepsBool isNullPtr(_Type *ptr)
Tells if passed pointer is null.
Definition: CepsMemory.hpp:45