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

#include <MonodomainSolver.hpp>

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

Public Member Functions

void PrepareForSetupLinearSystem (Vec currentSolution)
 
virtual void InitialiseForSolve (Vec initialSolution)
 
 MonodomainSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBoundaryConditions)
 
virtual ~MonodomainSolver ()
 
- 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)
 
virtual void FollowingSolveLinearSystem (Vec currentSolution)
 
LinearSystemGetLinearSystem ()
 

Private Member Functions

void SetupLinearSystem (Vec currentSolution, bool computeMatrix)
 

Private Attributes

MonodomainTissue< ELEMENT_DIM, SPACE_DIM > * mpMonodomainTissue
 
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > * mpBoundaryConditions
 
MonodomainAssembler< ELEMENT_DIM, SPACE_DIM > * mpMonodomainAssembler
 
NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, 1 > * mpNeumannSurfaceTermsAssembler
 
MonodomainCorrectionTermAssembler< ELEMENT_DIM, SPACE_DIM > * mpMonodomainCorrectionTermAssembler
 
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 MonodomainSolver< ELEMENT_DIM, SPACE_DIM >

A monodomain solver, which uses various assemblers to set up the monodomain linear system.

The discretised monodomain equation leads to the linear system (see FEM implementations document)

( (chi*C/dt) M + K ) V^{n+1} = (chi*C/dt) M V^{n} + M F^{n} + c_surf

where chi is the surface-area to volume ratio, C the capacitance, dt the timestep M the mass matrix, K the stiffness matrix, V^{n} the vector of voltages at time n, F^{n} the vector of (chi*Iionic + Istim) at each node, and c_surf a vector arising from any surface stimuli (usually zero).

This solver uses two assemblers, MonodomainAssembler to assemble the LHS matrix, [(chi*C/dt) M + K], and also to compute c_surf, and MassMatrixAssembler to assemble the mass matrix M used on the RHS. Note that MonodomainAssembler itself calls: MassMatrixAssembler for M (as this class does directly for RHS), and MonodomainStiffnessMatrixAssembler for K.

Also allows state variable interpolation (SVI) to be used on elements for which it will be needed, if the appropriate HeartConfig boolean is set. See https://chaste.github.io/docs/user-guides/state-variable-interpolation/ for more details on this. In this case the equation is ( (chi*C/dt) M + K ) V^{n+1} = (chi*C/dt) M V^{n} + M F^{n} + c_surf + c_correction and another assembler is used to create the c_correction.

Definition at line 74 of file MonodomainSolver.hpp.

Constructor & Destructor Documentation

◆ MonodomainSolver()

◆ ~MonodomainSolver()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::~MonodomainSolver ( )
virtual

Destructor

Definition at line 214 of file MonodomainSolver.cpp.

References PetscTools::Destroy().

Member Function Documentation

◆ InitialiseForSolve()

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

Overloaded InitialiseForSolve

Parameters
initialSolutioninitial solution

Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 134 of file MonodomainSolver.cpp.

References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve(), HeartConfig::Instance(), and PetscTools::SetupMat().

◆ PrepareForSetupLinearSystem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::PrepareForSetupLinearSystem ( Vec  currentSolution)
virtual

Overloaded PrepareForSetupLinearSystem() methods which gets the cell models to solve themselves

Parameters
currentSolutionsolution at current time

Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 173 of file MonodomainSolver.cpp.

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

◆ SetupLinearSystem()

Member Data Documentation

◆ mMassMatrix

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

The mass matrix, used to computing the RHS vector

Definition at line 98 of file MonodomainSolver.hpp.

◆ mpBoundaryConditions

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

Boundary conditions

Definition at line 83 of file MonodomainSolver.hpp.

◆ mpMonodomainAssembler

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

The monodomain assembler, used to set up the LHS matrix

Definition at line 86 of file MonodomainSolver.hpp.

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

◆ mpMonodomainCorrectionTermAssembler

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>* MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainCorrectionTermAssembler
private

If using state variable interpolation, points to an assembler to use in computing the correction term to apply to the RHS.

Definition at line 95 of file MonodomainSolver.hpp.

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

◆ mpMonodomainTissue

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

Monodomain tissue class (collection of cells, and conductivities)

Definition at line 80 of file MonodomainSolver.hpp.

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

◆ mpNeumannSurfaceTermsAssembler

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

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

Definition at line 89 of file MonodomainSolver.hpp.

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

◆ mVecForConstructingRhs

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

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

Definition at line 103 of file MonodomainSolver.hpp.

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


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