Chaste Commit::1fd4e48e3990e67db148bc1bc4cf6991a0049d0c
OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <OperatorSplittingMonodomainSolver.hpp>

+ Inheritance diagram for OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >:
+ Collaboration diagram for OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >:

Public Member Functions

void InitialiseForSolve (Vec initialSolution)
 
 OperatorSplittingMonodomainSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBoundaryConditions)
 
 ~OperatorSplittingMonodomainSolver ()
 
- Public Member Functions inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 1 >
 AbstractDynamicLinearPdeSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
void SetTimes (double tStart, double tEnd)
 
void SetTimeStep (double dt)
 
void SetInitialCondition (Vec initialCondition)
 
virtual Vec Solve ()
 
void SetMatrixIsNotAssembled ()
 
void SetTimeAdaptivityController (AbstractTimeAdaptivityController *pTimeAdaptivityController)
 
void SetOutputToVtk (bool output)
 
void SetOutputToParallelVtk (bool output)
 
void SetOutputToTxt (bool output)
 
void SetOutputDirectoryAndPrefix (std::string outputDirectory, std::string prefix)
 
void SetPrintingTimestepMultiple (unsigned multiple)
 
- Public Member Functions inherited from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
 AbstractLinearPdeSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
virtual ~AbstractLinearPdeSolver ()
 
virtual void FinaliseLinearSystem (Vec currentSolution)
 
LinearSystemGetLinearSystem ()
 

Private Member Functions

void SetupLinearSystem (Vec currentSolution, bool computeMatrix)
 
void PrepareForSetupLinearSystem (Vec currentSolution)
 
void FollowingSolveLinearSystem (Vec currentSolution)
 

Private Attributes

BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > * mpBoundaryConditions
 
MonodomainTissue< ELEMENT_DIM, SPACE_DIM > * mpMonodomainTissue
 
MonodomainAssembler< ELEMENT_DIM, SPACE_DIM > * mpMonodomainAssembler
 
NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, 1 > * mpNeumannSurfaceTermsAssembler
 
Mat mMassMatrix
 
Vec mVecForConstructingRhs
 

Additional Inherited Members

- Protected Member Functions inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 1 >
void InitialiseHdf5Writer ()
 
void WriteOneStep (double time, Vec solution)
 
- Protected Attributes inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 1 >
double mTstart
 
double mTend
 
bool mTimesSet
 
Vec mInitialCondition
 
bool mMatrixIsAssembled
 
bool mMatrixIsConstant
 
double mIdealTimeStep
 
double mLastWorkingTimeStep
 
AbstractTimeAdaptivityControllermpTimeAdaptivityController
 
bool mOutputToVtk
 
bool mOutputToParallelVtk
 
bool mOutputToTxt
 
std::string mOutputDirectory
 
std::string mFilenamePrefix
 
unsigned mPrintingTimestepMultiple
 
Hdf5DataWritermpHdf5Writer
 
std::vector< int > mVariableColumnIds
 
- Protected Attributes inherited from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
LinearSystemmpLinearSystem
 
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
 

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
class OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >

A monodomain solver that uses Strang operator splitting of the diffusion (conductivity) term and the reaction (ionic current) term, instead of solving the full reaction-diffusion PDE. This does NOT refer to operator splitting of the two PDEs in the bidomain equations. For details see for example Sundnes et al "Computing the Electrical Activity of the Heart".

The algorithm is, for solving from t=T to T+dt.

(i) Solve ODEs dV/dt = Iionic for t=T to T+dt/2 [giving updated V (internally, and in solution vector) and updated state variables] (ii) Solve PDE dV/dt = div sigma grad V for t=T to dt [using V from step i, –> updated V] (iii) Solve ODEs dV/dt = Iionic for t=T+dt/2 to T+dt [using V from step ii, –> final V]

Notes (a) Stages (iii) and (i) can normally be solved together in one go, except just before/after printing the voltage to file. However for simplicity of code this has not been implemented (b) Therefore, the effective ODE timestep will be: min(ode_dt, pde_dt/2), where ode_dt and pde_dt are those given via HeartConfig. (c) This solver is FOR COMPARING ACCURACY, NOT PERFORMANCE. It has not been optimised and may or may not perform well in parallel. (d) We don't implement the simpler form of operator splitting, Godunov splitting, where the ODEs are solved for one timestep and the PDEs are solved for one timestep, since this is formally equivalent to the default implementation where the ionic current is interpolated from the nodal values (ie ICI - see ICI/SVI discussion in documentation)

Definition at line 72 of file OperatorSplittingMonodomainSolver.hpp.

Constructor & Destructor Documentation

◆ OperatorSplittingMonodomainSolver()

◆ ~OperatorSplittingMonodomainSolver()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::~OperatorSplittingMonodomainSolver ( )

Destructor

Definition at line 188 of file OperatorSplittingMonodomainSolver.cpp.

References PetscTools::Destroy().

Member Function Documentation

◆ FollowingSolveLinearSystem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::FollowingSolveLinearSystem ( Vec  currentSolution)
privatevirtual

Called after solving the linear system, used to solve the cell models for second half timestep (step (iii) above)

Parameters
currentSolutionthe latest solution vector (ie the solution of the linear system).

Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 117 of file OperatorSplittingMonodomainSolver.cpp.

References PdeSimulationTime::GetNextTime(), PdeSimulationTime::GetPdeTimeStep(), and PdeSimulationTime::GetTime().

◆ InitialiseForSolve()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve ( Vec  initialSolution)
virtual

◆ PrepareForSetupLinearSystem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::PrepareForSetupLinearSystem ( Vec  currentSolution)
privatevirtual

Called before setting up the linear system, used to solve the cell models for first half timestep (step (i) above)

Parameters
currentSolutionthe latest solution vector

Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 109 of file OperatorSplittingMonodomainSolver.cpp.

References PdeSimulationTime::GetPdeTimeStep(), and PdeSimulationTime::GetTime().

◆ SetupLinearSystem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem ( Vec  currentSolution,
bool  computeMatrix 
)
privatevirtual

Member Data Documentation

◆ mMassMatrix

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Mat OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mMassMatrix
private

The mass matrix, used to computing the RHS vector

Definition at line 89 of file OperatorSplittingMonodomainSolver.hpp.

◆ mpBoundaryConditions

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBoundaryConditions
private

Boundary conditions

Definition at line 77 of file OperatorSplittingMonodomainSolver.hpp.

◆ mpMonodomainAssembler

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>* OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainAssembler
private

The monodomain assembler, used to set up the LHS matrix

Definition at line 83 of file OperatorSplittingMonodomainSolver.hpp.

Referenced by OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::OperatorSplittingMonodomainSolver().

◆ mpMonodomainTissue

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainTissue
private

Monodomain tissue class (collection of cells, and conductivities)

Definition at line 80 of file OperatorSplittingMonodomainSolver.hpp.

Referenced by OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::OperatorSplittingMonodomainSolver().

◆ mpNeumannSurfaceTermsAssembler

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,1>* OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpNeumannSurfaceTermsAssembler
private

Assembler for surface integrals coming from any non-zero Neumann boundary conditions

Definition at line 86 of file OperatorSplittingMonodomainSolver.hpp.

Referenced by OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::OperatorSplittingMonodomainSolver().

◆ mVecForConstructingRhs

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Vec OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mVecForConstructingRhs
private

The vector multiplied by the mass matrix. Ie, if the linear system to be solved is Ax=b, this vector is z where b=Mz.

In the normal solver this is equal to alpha*V + beta*Iionic + gamma*Istim, here there is no Iionic term.

Definition at line 98 of file OperatorSplittingMonodomainSolver.hpp.

Referenced by OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::OperatorSplittingMonodomainSolver().


The documentation for this class was generated from the following files: