Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <MonodomainPurkinjeSolver.hpp>

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

Public Member Functions

void PrepareForSetupLinearSystem (Vec currentSolution)
 
virtual void InitialiseForSolve (Vec initialSolution)
 
 MonodomainPurkinjeSolver (MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 2 > *pBoundaryConditions)
 
virtual ~MonodomainPurkinjeSolver ()
 
- Public Member Functions inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 2 >
 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)
 
void SetIdentityBlockToLhsMatrix ()
 

Private Attributes

MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > * mpMixedMesh
 
MonodomainTissue< ELEMENT_DIM, SPACE_DIM > * mpMonodomainTissue
 
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 2 > * mpBoundaryConditions
 
MonodomainPurkinjeVolumeAssembler< ELEMENT_DIM, SPACE_DIM > * mpVolumeAssembler
 
MonodomainPurkinjeCableAssembler< ELEMENT_DIM, SPACE_DIM > * mpCableAssembler
 
NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, 2 > * mpNeumannSurfaceTermsAssembler
 
Mat mMassMatrix
 
Vec mVecForConstructingRhs
 

Additional Inherited Members

- Protected Member Functions inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 2 >
void InitialiseHdf5Writer ()
 
void WriteOneStep (double time, Vec solution)
 
- Protected Attributes inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 2 >
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 MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM >

Solver class for Monodomain problems on tissues containing Purkinje fibres.

In such problems, there are a subset of nodes of the mesh which are labelled as Purkinje nodes (and 1D elements connecting them). There are two variables to be computed, the (normal, myocardium) transmembrane voltage (V), defined at ALL nodes of the mesh, and the Purkinje voltage (Vp), defined at the purkinje nodes. Hence, the Purkinje nodes have two variables defined on them.

For parallelisation/implementation reasons, we choose to have two variables to be defined at ALL nodes, introducing dummy variables Vp for non-Purkinje nodes. We set Vp = 0 at non-Purkinje nodes. Hence we solve for {V,Vp} at every node in the mesh, and PROBLEM_DIM=2. The linear system to be assembled, written as usual in block form but actually in striped form in the code, is:

[ A1 0 ][V ] = [b1] [ 0 A2 ][Vp] = [b2]

where each block is of size num_nodes.

A1 and b1 are obtained by integrating over 3D myocardium elements and are therefore exactly the same as in a normal monodomain problem.

Suppose the nodes are ordered such that all the Purkinje nodes are last. Then the matrix A2 needs to be of the form A2 = [I 0 ] [0 Ap] where the first block is of size num_non_purkinje_nodes, and the second is of size num_purkinje_nodes. Ap is obtained by assembling over 1D purkinje elements. After assembly we add in the identity block, which is just represents the equations Vp=0 for the dummy variables (non-purkinje nodes). Finally, we similarly have b2 = [0 ] [bp] where bp involves a loop over 1D purkinje elements.

This class implements the above, and is based on (but doesn't inherit from, as the PROBLEM_DIMENSION is different) MonodomainSolver.

Definition at line 89 of file MonodomainPurkinjeSolver.hpp.

Constructor & Destructor Documentation

◆ MonodomainPurkinjeSolver()

◆ ~MonodomainPurkinjeSolver()

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

Destructor

Definition at line 238 of file MonodomainPurkinjeSolver.cpp.

References PetscTools::Destroy().

Member Function Documentation

◆ InitialiseForSolve()

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

◆ PrepareForSetupLinearSystem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainPurkinjeSolver< 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 203 of file MonodomainPurkinjeSolver.cpp.

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

◆ SetIdentityBlockToLhsMatrix()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM >::SetIdentityBlockToLhsMatrix ( )
private

For the block in the LHS matrix corresponding to the Purkinje voltage at myocardium nodes: this block is zero after all elements are assembled so, so we set it to be the identity block. This is done by just checking which rows of the matrix has zero diagonal values.

Definition at line 140 of file MonodomainPurkinjeSolver.cpp.

References PetscTools::Destroy(), PetscVecTools::GetElement(), PetscMatTools::SetElement(), and PetscMatTools::SwitchWriteMode().

◆ SetupLinearSystem()

Member Data Documentation

◆ mMassMatrix

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

The mass matrix, used to computing the RHS vector

Definition at line 123 of file MonodomainPurkinjeSolver.hpp.

◆ mpBoundaryConditions

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,2>* MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM >::mpBoundaryConditions
private

Boundary conditions

Definition at line 103 of file MonodomainPurkinjeSolver.hpp.

◆ mpCableAssembler

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainPurkinjeCableAssembler<ELEMENT_DIM,SPACE_DIM>* MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM >::mpCableAssembler
private

The cable element assembler, used to set up cable integral parts of the LHS matrix

Definition at line 114 of file MonodomainPurkinjeSolver.hpp.

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

◆ mpMixedMesh

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM >::mpMixedMesh
private

Saved pointer to the mesh in this class, as the pointer saved in the parent class (AbstractDynamicLinearPdeSolver::mpMesh) is not declared to be a pointer to a mixed mesh

Definition at line 97 of file MonodomainPurkinjeSolver.hpp.

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

◆ mpMonodomainTissue

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

Monodomain tissue class (collection of cells, and conductivities)

Definition at line 100 of file MonodomainPurkinjeSolver.hpp.

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

◆ mpNeumannSurfaceTermsAssembler

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,2>* MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM >::mpNeumannSurfaceTermsAssembler
private

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

Definition at line 117 of file MonodomainPurkinjeSolver.hpp.

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

◆ mpVolumeAssembler

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainPurkinjeVolumeAssembler<ELEMENT_DIM,SPACE_DIM>* MonodomainPurkinjeSolver< ELEMENT_DIM, SPACE_DIM >::mpVolumeAssembler
private

The volume assembler, used to set up volume integral parts of the LHS matrix

Definition at line 109 of file MonodomainPurkinjeSolver.hpp.

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

◆ mVecForConstructingRhs

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Vec MonodomainPurkinjeSolver< 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 128 of file MonodomainPurkinjeSolver.hpp.

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


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