CEPS  24.01
Cardiac ElectroPhysiology Simulator
FiniteElements.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 
36 
51 {
52 
53  public:
54  //------------------------------------------------------------------------------
55  // Constructors
56 
58  FiniteElements() = delete;
59 
62 
64  FiniteElements(const FiniteElements&) = delete;
65 
67  ~FiniteElements() override;
68 
69  //------------------------------------------------------------------------------
70  // FiniteElements shortcuts
71 
75 
79 
83 
85  CepsUInt
86  getNumberOfHaloFENodes() const;
87 
89  CepsUInt
91 
93  const CepsVector<FENode*>&
95 
97  const CepsVector<FENode*>&
99 
101  FENode*
103 
105  CepsBool
106  isOwnedFENode(CepsGlobalIndex globalID) const;
107 
109  CepsBool
110  isHaloFENode(CepsGlobalIndex globalID) const;
111 
114  getDofsOnElement(FEBase* elem);
115 
122 
124  CepsIndex
125  getDofSpatialId(const DegreeOfFreedom* dof) const override;
126 
127  // ========== Access to matrices======================================================
128 
130  DMatPtr
131  getMassMatrix() const;
132 
134  DMatPtr
135  getStiffnessMatrix() const;
136 
137  // ========== METRICS ================================================================
138 
140  CepsReal
141  dotProduct(
142  DHVecPtr u,
143  DHVecPtr v,
144  CepsBool boundary = false,
145  const CepsSet<CepsAttribute>& attrs = {},
146  const CepsVector<Unknown*>& unknowns = {}
147  ) final;
148 
151  CepsReal
152  stiffProduct(
153  DHVecPtr u,
154  DHVecPtr v,
155  CepsBool boundary = false,
156  const CepsSet<CepsAttribute>& attrs = {},
157  const CepsVector<Unknown*>& unknowns = {}
158  );
159 
161  CepsReal
162  h1Norm(
163  DHVecPtr v,
164  CepsBool boundary = false,
165  const CepsSet<CepsAttribute>& attrs = {},
166  const CepsVector<Unknown*>& unknowns = {}
167  ) override;
168 
177  void
179  const CepsUnknownIndex& uId,
180  DVecPtr vec,
181  CepsVector<CepsIndex>& indices,
182  CepsVector<CepsReal>& values,
183  CepsBool sendToMaster = false,
184  CepsBool returnGeomIndices = false
185  ) override;
186 
187  protected:
189  void
191 
201  void
202  buildElements();
203 
206  void
207  buildElements(Mesh* mesh, CepsBool onBoundary);
208 
217  void
219  FEBase* fElem,
220  FEBase* nElem,
221  CepsUInt nID,
222  CepsBool boundary,
223  CepsBool& foundSameDof
224  );
225 
232  void
233  createNewNode(FEBase* el, CepsUInt indexDof, CepsBool isBoundary, Mesh* mesh);
234 
236  void
238 
240  CepsHash3
241  getNodeHash(const FENode* n) const;
242 
245  void
246  buildDofsTree() override;
247 
250  void
252 
253  protected:
254 
271 
272 
273 };
274 
279 CepsReal
281  const CepsVector<FEBase*>& cells,
282  const CepsSet<CepsAttribute>& attributes,
283  CepsBool onlyOfThisProc = false
284 );
285 
289 CepsReal
290 getDomainMeasure(const CepsVector<FEBase*>& cells, CepsBool onlyOfThisProc = false);
CepsCellType
Enum for different shapes of cells.
Definition: CepsEnums.hpp:128
CepsIndex CepsUnknownIndex
For unknowns.
Definition: CepsTypes.hpp:217
std::map< _Key, _Tp, _Compare, _Alloc > CepsMap
C++ map.
Definition: CepsTypes.hpp:196
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
CepsIndex CepsGlobalIndex
Many uses. Has to be signed for PETSc.
Definition: CepsTypes.hpp:218
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
CepsInt CepsIndex
Index rowid etc.
Definition: CepsTypes.hpp:111
std::shared_ptr< DistributedHaloVector > DHVecPtr
Typedef for pointer on Distributed Halo CepsVector.
std::shared_ptr< DistributedMatrix > DMatPtr
Short typedef for pointer on dist matrix.
std::shared_ptr< DistributedVector > DVecPtr
Short typedef for pointer on distributed vector.
CepsReal getDomainMeasure(const CepsVector< FEBase * > &cells, const CepsSet< CepsAttribute > &attributes, CepsBool onlyOfThisProc=false)
Compute the measure of a set of a cells.
Abstract Class for all numerical method (FE, FD, FV etc)
Base class for creating PDEs to solve.
A degree of freedom for any kind of problem The dof can be associated to a geometrical element or not...
Abstract class for finite elements.
Definition: FEBase.hpp:65
A nodal point on a finite element. It is different from a geom node as it may have different properti...
Definition: FENode.hpp:46
Holds all finite elements corresponding to each geometrical element.
FENode * getFENode(CepsGlobalIndex nID)
Link to a single node.
CepsReal dotProduct(DHVecPtr u, DHVecPtr v, CepsBool boundary=false, const CepsSet< CepsAttribute > &attrs={}, const CepsVector< Unknown * > &unknowns={}) final
returns where is the mass matrix
DMatPtr m_mass
Mass matrix, computed at instantiation.
void initializeRefElements()
Configure reference elements.
CepsHash3 getNodeHash(const FENode *n) const
Get the hash value from the coordinates of the node.
CepsVector< FEBase * > & getBoundaryFiniteElements()
Get vector containing all finite elements of the maximum valid dim.
ReferenceFE * getReferenceElementOfDim(CepsCellType type, CepsUInt dim) const
Access to reference elements.
CepsVector< DegreeOfFreedom * > getDofsOnElement(FEBase *elem)
ALl the dofs defined on an element.
CepsUInt getNumberOfOwnedFENodes() const
Local number of nodes belonging to this process, halo nodes NOT included.
void createNewNode(FEBase *el, CepsUInt indexDof, CepsBool isBoundary, Mesh *mesh)
Allocate new node that it is not a geometrical node.
void buildElements()
main routine of FiniteElements called in the constructor
~FiniteElements() override
Destructor.
FiniteElements()=delete
Default Constuctor deleted.
const CepsVector< FENode * > & getHaloFENodes()
vector of halo nodes
CepsVector< FEBase * > & getFiniteElements()
Get vector containing all finite elements of the maximum valid dim.
CepsBool isHaloFENode(CepsGlobalIndex globalID) const
Tells if node is halo by current CPU.
void setSpatialUnknownsDofsInteractions() override
Build dofs neighboring list for spatial unknowns. Different in FE (same element dofs) and VF,...
CepsIndex getDofSpatialId(const DegreeOfFreedom *dof) const override
Return the Geom Node Id for a node (not the FENode index returned by dof->getGeomId())
DistributedInfos< FEBase * > * m_bdrElems
Holds 1D-3D finite boundary elements.
DMatPtr m_stiff
Stiffness matrix, computed at instantiation.
const CepsVector< FENode * > & getOwnedFENodes()
vector of owned nodes
DistributedInfos< FENode *, CepsHash3 > * m_nodes
Holds nodes.
DistributedInfos< FEBase * > * m_volElems
Holds 1D-3D finite elements.
CepsVector< DegreeOfFreedom * > getDofsOnElementForUnknowns(FEBase *elem, const CepsVector< Unknown * > &us)
All dofs defined on an element, considering unknowns us.
void checkNeighborsBeforeAssign(FEBase *fElem, FEBase *nElem, CepsUInt nID, CepsBool boundary, CepsBool &foundSameDof)
Checks if a FE node already exists in a neighbor element.
void extractValuesForUnknown(const CepsUnknownIndex &uId, DVecPtr vec, CepsVector< CepsIndex > &indices, CepsVector< CepsReal > &values, CepsBool sendToMaster=false, CepsBool returnGeomIndices=false) override
Slicing method that extracts the data corresponding to a specific unknown.
FiniteElements(const FiniteElements &)=delete
Copy Constructor deleted.
CepsReal h1Norm(DHVecPtr v, CepsBool boundary=false, const CepsSet< CepsAttribute > &attrs={}, const CepsVector< Unknown * > &unknowns={}) override
sqrt(stiffproduct(v,v))
CepsUInt m_maxValidDim
Get the maximum valid dimension.
CepsReal stiffProduct(DHVecPtr u, DHVecPtr v, CepsBool boundary=false, const CepsSet< CepsAttribute > &attrs={}, const CepsVector< Unknown * > &unknowns={})
returns where is the stiffness matrix
void buildDofsTree() override
Builds the degrees of freedom tree from the PDE data, uses previously build elements and nodes....
CepsMap< CepsCellType, CepsArray4< ReferenceFE * > > m_refElems
Reference elements (point, segment, tri, tetra)
CepsUInt getNumberOfHaloFENodes() const
Global number of owned nodes.
void computeNeighbors()
Compute Neighbors.
CepsUInt m_order
Order used to build elements.
DMatPtr getStiffnessMatrix() const
Pointer on stiffness matrix.
DMatPtr getMassMatrix() const
Pointer on mass matrix.
CepsBool isOwnedFENode(CepsGlobalIndex globalID) const
Tells if node is owned by current CPU.
Geometrical information of 1,2 or 3D distributed cells.
Definition: Mesh.hpp:57
Base class for reference finite elements.
Definition: ReferenceFE.hpp:48
A triple hash to be used for coordinates (multiplied *10^12 then truncated)