Chaste Commit::1fd4e48e3990e67db148bc1bc4cf6991a0049d0c
|
#include <MonodomainSolver.hpp>
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 |
AbstractTimeAdaptivityController * | mpTimeAdaptivityController |
bool | mOutputToVtk |
bool | mOutputToParallelVtk |
bool | mOutputToTxt |
std::string | mOutputDirectory |
std::string | mFilenamePrefix |
unsigned | mPrintingTimestepMultiple |
Hdf5DataWriter * | mpHdf5Writer |
std::vector< int > | mVariableColumnIds |
Protected Attributes inherited from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > | |
LinearSystem * | mpLinearSystem |
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * | mpMesh |
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.
MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver | ( | AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * | pMesh, |
MonodomainTissue< ELEMENT_DIM, SPACE_DIM > * | pTissue, | ||
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > * | pBoundaryConditions | ||
) |
Constructor
pMesh | pointer to the mesh |
pTissue | pointer to the tissue |
pBoundaryConditions | pointer to the boundary conditions |
Definition at line 180 of file MonodomainSolver.cpp.
References HeartConfig::Instance(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 1 >::mMatrixIsConstant, AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainAssembler, MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainCorrectionTermAssembler, MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainTissue, MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpNeumannSurfaceTermsAssembler, MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mVecForConstructingRhs, and AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM >::SetCacheReplication().
|
virtual |
|
virtual |
Overloaded InitialiseForSolve
initialSolution | initial 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().
|
virtual |
Overloaded PrepareForSetupLinearSystem() methods which gets the cell models to solve themselves
currentSolution | solution 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().
|
privatevirtual |
Implementation of SetupLinearSystem() which uses the assembler to compute the LHS matrix, but sets up the RHS vector using the mass-matrix (constructed using a separate assembler) multiplied by a vector
currentSolution | Solution at current time |
computeMatrix | Whether to compute the matrix of the linear system |
Implements AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 42 of file MonodomainSolver.cpp.
References AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::Assemble(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleMatrix(), DistributedVector::Begin(), GenericEventHandler< 16, HeartEventHandler >::BeginEvent(), DistributedVectorFactory::CreateDistributedVector(), DistributedVector::End(), GenericEventHandler< 16, HeartEventHandler >::EndEvent(), PetscMatTools::Finalise(), HeartConfig::GetCapacitance(), PdeSimulationTime::GetPdeTimeStepInverse(), HeartConfig::GetSurfaceAreaToVolumeRatio(), HeartConfig::Instance(), DistributedVector::Restore(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetMatrixToAssemble(), and HeartConfig::SetUseMassLumping().
|
private |
The mass matrix, used to computing the RHS vector
Definition at line 98 of file MonodomainSolver.hpp.
|
private |
Boundary conditions
Definition at line 83 of file MonodomainSolver.hpp.
|
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().
|
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().
|
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().
|
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().
|
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().