CEPS  24.01
Cardiac ElectroPhysiology Simulator
AbstractTimedPdeSolver.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 
33 
34 // Forward declarations
36 
39 {
40 
41  public:
42 
45 
47  virtual ~AbstractTimedPdeSolver();
48 
50  void
51  solve() override;
52 
55  void
56  setupWithParameters(InputParameters* params) override;
57 
59  CepsUInt
60  getExpectedNumberOfOutputs() const override;
61 
63  virtual void
64  output(DHVecPtr solution, CepsBool immediateWriting = true);
65 
68  DHVecPtr
69  getSolution() const;
70 
72  void
74 
79  getErrors() const;
80 
81  protected:
82 
85  getTimedProblem() const;
86 
88  virtual void
89  displayProgress() const;
90 
92  virtual CepsBool
93  finished() const;
94 
96  void
98 
100  virtual void
101  allocateArrays();
102 
103  // Linear system building methods
104  void
105  updateAssemblers() override;
106 
110  virtual void
112 
114  virtual void
116 
117  protected:
125 
139 
144 
149 
152 
154  static constexpr CepsReal c_alphaSBDF[5] = {0., 1., 2. / 3., 6. / 11., 12. / 25.};
156  static constexpr CepsReal c_alphaCN[3] = {0., 1., 1. / 2.};
158  static constexpr CepsReal c_akSBDF[4][4] = {{1.}, {4. / 3., -1. / 3.}, {18. / 11., -9. / 11., 2. / 11.}, {48. / 25., -36. / 25., 16. / 25., -3. / 25.}};
159 
162 
164 };
std::basic_string< CepsChar > CepsString
C++ format string.
Definition: CepsTypes.hpp:128
CepsArray< _Type, 2U > CepsArray2
C++ array, 2 elements.
Definition: CepsTypes.hpp:162
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
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.
Base class for PDE solving.
Astract Problem which does depend on time.
Base class for solving PDE with time dependance.
virtual void allocateArrays()
Allocates the arrays with the correct number of multistep.
CepsBool m_cheatConvergence
If true, use exact solution for first steps of multi steps methods.
virtual void displayProgress() const
Tells time.
TimeStepper * m_timeStepper
Monitors time evolution.
DHVecPtr getSolution() const
Returns a copy of the distributed vector containing the current solution.
void determineSolverType(CepsString s)
Initialize the type from parameter string.
virtual void swapMultiStepPointers()
Swaps vectors of array pointers for next step (for multistep methods)
virtual CepsBool finished() const
Tells if simulation reached end time or, if enabled, all points have seen an AP.
void solve() override
Solves the all PDE (all iterations)
CepsVector< DHVecPtr > m_regSourcesNm
Sources at previous steps (if problem has source terms)
CepsString m_type
Numerical scheme descriptor.
DMatPtr m_dtMat
Matrix of time derivative of the system.
void enableErrorComputation()
Sets the comparison with analytic solution.
CepsVector< DHVecPtr > m_uNm
Solution vectors at different time steps. Several vectors are needed for multi-step methods....
CepsVector< DHVecPtr > m_rkF
RK methods: intermediate rates.
CepsUInt m_nbMultiSteps
Maximum number of steps needed to compute next solution (eg CN:2)
DVecPtr m_timeDerEqsVec
Vector of 1 or 0, depending on whether the dof equation has a time derivative or not.
CepsVector< DHVecPtr > m_lapSourcesNm
Sources at previous steps (if problem has laplace source terms)
CepsVector< DHVecPtr > m_rkU
RK methods: intermediate evaluations.
static constexpr CepsReal c_alphaCN[3]
LHS coeffcients of the CN scheme.
CepsUInt getExpectedNumberOfOutputs() const override
Number of files written.
virtual ~AbstractTimedPdeSolver()
Destructor.
AbstractTimedPdeSolver(AbstractTimedPdeProblem *pb)
Default constructor.
CepsReal m_scaleSystem
A factor for both matrix and vector (typically 1/dt)
void setupWithParameters(InputParameters *params) override
Set attributes from input file. Parameters are passed as arguments in case one wants to use other par...
CepsUInt m_nbIterSnapshot
Output perdiodicity.
AbstractTimedPdeProblem * getTimedProblem() const
Returns a pointer with the appropriate type of pb.
void updateAssemblers() override
Update assemblers and recompute everything is needed.
CepsArray2< CepsArray3< CepsArray3< CepsReal > > > getErrors() const
Gets the currently computed errors. First index selects absolute(0) or relative(1) second index selec...
static constexpr CepsReal c_alphaSBDF[5]
LHS coeffcients of the SBDF schemes.
CepsString m_rkType
RK Method if the solver is of type "RK".
virtual void output(DHVecPtr solution, CepsBool immediateWriting=true)
Prints the solution.
DVecPtr m_noTimeDerEqsVec
The opposite version of m_timeDerEqsVec, such as m_timeDerEqsVec + m_noTimeDerEqsVec = vec(1)
DHVecPtr m_uNp1
Vectors of unknowns at time t^n+1.
virtual void fillSourceTermsVector(CepsReal tn)
Adds contributions from all source terms into m_regSourcesNm[0].
static constexpr CepsReal c_akSBDF[4][4]
RHS coeffcients of the SBDF schemes.
Reads and stores simulation configuration.
A simple time stepper used in dynamic linear solvers.
Definition: TimeStepper.hpp:40