Chaste Release::3.1
|
#include <MonodomainSolver.hpp>
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, unsigned numQuadPoints=2) | |
~MonodomainSolver () | |
Private Member Functions | |
void | SetupLinearSystem (Vec currentSolution, bool computeMatrix) |
Private Attributes | |
MonodomainTissue< ELEMENT_DIM, SPACE_DIM > * | mpMonodomainTissue |
unsigned | mNumQuadPoints |
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 |
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, one to assemble the LHS matrix, (chi*C/dt) M + K, and also to compute c_surf, and one to assemble the mass matrix M.
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 wiki page ChasteGuides/StateVariableInterpolation 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 71 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, | ||
unsigned | numQuadPoints = 2 |
||
) |
Constructor
pMesh | pointer to the mesh |
pTissue | pointer to the tissue |
pBoundaryConditions | pointer to the boundary conditions |
numQuadPoints | number of quadrature points (defaults to 2) |
Definition at line 185 of file MonodomainSolver.cpp.
References HeartConfig::Instance(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 1 >::mMatrixIsConstant, MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mNumQuadPoints, 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().
MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::~MonodomainSolver | ( | ) |
void MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve | ( | Vec | initialSolution | ) | [virtual] |
Overloaded InitialiseForSolve
initialSolution | initial solution |
Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 137 of file MonodomainSolver.cpp.
References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve(), HeartConfig::Instance(), and PetscTools::SetupMat().
void MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::PrepareForSetupLinearSystem | ( | Vec | currentSolution | ) | [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 176 of file MonodomainSolver.cpp.
References PdeSimulationTime::GetPdeTimeStep(), and PdeSimulationTime::GetTime().
void MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem | ( | Vec | currentSolution, |
bool | computeMatrix | ||
) | [private, virtual] |
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().
Mat MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mMassMatrix [private] |
The mass matrix, used to computing the RHS vector
Definition at line 101 of file MonodomainSolver.hpp.
unsigned MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mNumQuadPoints [private] |
Number of quadrature points per dimension (only saved so it can be passed to the assembler)
Definition at line 83 of file MonodomainSolver.hpp.
Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver().
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBoundaryConditions [private] |
Boundary conditions
Definition at line 86 of file MonodomainSolver.hpp.
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 89 of file MonodomainSolver.hpp.
Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver().
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 98 of file MonodomainSolver.hpp.
Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver().
MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainTissue [private] |
Monodomain tissue class (collection of cells, and conductivities)
Definition at line 77 of file MonodomainSolver.hpp.
Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver().
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 92 of file MonodomainSolver.hpp.
Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver().
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 106 of file MonodomainSolver.hpp.
Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver().