CEPS  24.01
Cardiac ElectroPhysiology Simulator
FECardiacBCAssembler.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 
34  CardiacProblem* pb,
35  FiniteElements* fe,
36  const CepsVector<Unknown*>& unknowns
37 ) : FECardiacBCAssembler(pb, pb->getTissueAttributes(), fe, unknowns)
38 {
39 }
40 
42  CardiacProblem* pb,
43  const CepsSet<CepsAttribute> &attributes,
44  FiniteElements* fe,
45  const CepsVector<Unknown*>& unknowns
46 ) :
47  FEDivKGradBCAssembler(fe,attributes,unknowns),
49 {
50 }
51 
52 void
54  FEBase* elem,
55  CepsReal3D xQ,
56  CepsReal t,
57  const CepsMathDynamic1D& phi,
58  const CepsMathDynamic2D& gradPhi
59 )
60 {
61  // a * sigma * grad_n u + b * u = g
62  // => sigma * grad_n u = g / a - (b / a) * u
63 
64  // Fetch the unknown for which we assemble
66  Unknown* u = m_problem->getUnknown(uid);
67 
68  // Check that we can assemble on that element for this unknown.
69  // If not, the vectors following structures are empty
70  // (we still need to call setupBlocks though !)
71  auto us = extractUnknownsFrom({u},elem);
73  if (us.empty())
74  return;
75 
76  // Evaluate the BC coefficients at the quadrature point
77  CepsStandardArgs args;
78  args.t = t;
79  args.x = xQ;
80  args.cellId = elem->getGeomObject()->getGlobalIndex();
81  args.attr = elem->getGeomObject()->getAttributes();
82  args.unknownId = uid;
83  CepsReal alpha = m_bc->getAlpha(args);
84  CepsReal beta = m_bc->getBeta(args);
85  CepsReal gamma = m_bc->eval(args);
86 
87  // Get Am and Cm at position xQ
88  CepsReal Am = m_Am->eval(args);
89  CepsReal Cm = m_Cm->eval(args);
90 
91  // Now computing the submatrix/subvector
92  m_bMat[0][0] += 1. / (Am * Cm) * beta / alpha * phi * phi.transpose();
93  m_bVec[0] += 1. / (Am * Cm) * gamma / alpha * phi;
94 
95  return;
96 }
Eigen::Matrix< CepsScalar, Eigen::Dynamic, 1 > CepsMathDynamic1D
Dynamic 1D array, eigen format.
Definition: CepsTypes.hpp:139
CepsIndex CepsUnknownIndex
For unknowns.
Definition: CepsTypes.hpp:217
std::set< _Type, _Compare, _Alloc > CepsSet
C++ set.
Definition: CepsTypes.hpp:209
std::vector< _Type, _Alloc > CepsVector
C++ vector.
Definition: CepsTypes.hpp:155
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
Eigen::Matrix< CepsScalar, Eigen::Dynamic, Eigen::Dynamic > CepsMathDynamic2D
Dynamic 2D array, eigen format.
Definition: CepsTypes.hpp:140
AbstractPdeProblem * m_problem
Direct link to problem, you may write getXXXProblem methods in derived classes that return the right ...
Unknown * getUnknown(const CepsString &label) const
Get an unknown by its name.
CepsUnknownIndex getUnknownId() const
Get the id of the unknown that this BC is for.
_Result eval(CepsStandardArgs args) override
Evaluate a boundary condition, essentially call the functor inside.
CoeffFunc getBeta() const
Get beta in .
CoeffFunc getAlpha() const
Get alpha in .
A abstract class that regroups common parameters of cardiac problems.
CepsVector< CepsMathDynamic1D > m_bVec
Block vector,.
CepsVector< Unknown * > extractUnknownsFrom(const CepsVector< Unknown * > &us, FEBase *elem) const
From a vector of unknowns, get the ones defined on that element.
void setupBlocksOnElementForUnknown(FEBase *elem, Unknown *unknown)
Get all the dofs on an element related to the given unknowns. This routine MUST be called before sumB...
CepsVector< CepsVector< CepsMathDynamic2D > > m_bMat
Block submatrices. The first two vectors account for unknowns: if assembly on element if for unknowns...
Abstract class for finite elements.
Definition: FEBase.hpp:65
Class to fill linear systems with coefficients linked to boundary conditions.
FECardiacBCAssembler()=delete
No default constructor.
virtual void computeBlocksOnElementAtQuadPoint(FEBase *element, CepsReal3D xQ, CepsReal t, const CepsMathDynamic1D &phi, const CepsMathDynamic2D &gradPhi) override
The function that is called to get the coefficients of the submatrix on a given finite element....
Common elements of cardiac assemblers.
ScalarSAFunc * m_Cm
Shortcut to Cm coefficient.
ScalarSAFunc * m_Am
Shortcut to Am coefficient.
Class to fill linear systems with coefficients linked to boundary conditions.
ScalarBoundaryCondition * m_bc
Current boundary condition.
Holds all finite elements corresponding to each geometrical element.
A class used to defined an unknown of a PDE problem The unknown can be defined on a specific region,...
Definition: Unknown.hpp:45
CepsSet< CepsAttribute > & getAttributes()
Returns the attributes of the entity.
_GeomObject * getGeomObject() const
Get the geom object.
const CepsGlobalIndex & getGlobalIndex() const
Get the index
Structure used to pass arguments to SAFunc (see pde directory) The flags of the SAFunc allows extract...
Definition: CepsTypes.hpp:239
CepsSet< CepsAttribute > attr
attributes
Definition: CepsTypes.hpp:243
CepsCellGlobalIndex cellId
cell index
Definition: CepsTypes.hpp:245
CepsReal3D x
space
Definition: CepsTypes.hpp:241
CepsReal t
time
Definition: CepsTypes.hpp:240
CepsUnknownIndex unknownId
unknown index
Definition: CepsTypes.hpp:247