MonodomainSolver< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <MonodomainSolver.hpp>

Inherits AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 1 >.

Collaboration diagram for MonodomainSolver< ELEMENT_DIM, SPACE_DIM >:
Collaboration graph
[legend]

List of all members.

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

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, 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 64 of file MonodomainSolver.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
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 
) [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::~MonodomainSolver (  )  [inline]

Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve ( Vec  initialSolution  )  [inline, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::PrepareForSetupLinearSystem ( Vec  currentSolution  )  [inline, virtual]

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

Parameters:
currentSolution solution at current time

Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 169 of file MonodomainSolver.cpp.

References PdeSimulationTime::GetPdeTimeStep(), PdeSimulationTime::GetTime(), and MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainTissue.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem ( Vec  currentSolution,
bool  computeMatrix 
) [inline, 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

Parameters:
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 35 of file MonodomainSolver.cpp.

References AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::Assemble(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleMatrix(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleVector(), DistributedVector::Begin(), GenericEventHandler< 16, HeartEventHandler >::BeginEvent(), DistributedVectorFactory::CreateDistributedVector(), DistributedVector::End(), GenericEventHandler< 16, HeartEventHandler >::EndEvent(), PetscMatTools::Finalise(), LinearSystem::FinaliseLhsMatrix(), LinearSystem::FinalisePrecondMatrix(), LinearSystem::FinaliseRhsVector(), HeartConfig::GetCapacitance(), PdeSimulationTime::GetPdeTimeStepInverse(), HeartConfig::GetSurfaceAreaToVolumeRatio(), HeartConfig::Instance(), MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mMassMatrix, MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mNumQuadPoints, AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpLinearSystem, 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, DistributedVector::Restore(), LinearSystem::rGetLhsMatrix(), LinearSystem::rGetPrecondMatrix(), LinearSystem::rGetRhsVector(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetMatrixToAssemble(), LinearSystem::SetPrecondMatrixIsDifferentFromLhs(), HeartConfig::SetUseMassLumping(), and AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetVectorToAssemble().


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Mat MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mMassMatrix [private]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
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 76 of file MonodomainSolver.hpp.

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

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 79 of file MonodomainSolver.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>* MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainAssembler [private]
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 91 of file MonodomainSolver.hpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::mpMonodomainTissue [private]
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 85 of file MonodomainSolver.hpp.

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

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 99 of file MonodomainSolver.hpp.

Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::MonodomainSolver(), MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem(), and MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::~MonodomainSolver().


The documentation for this class was generated from the following files:
Generated on Thu Dec 22 13:06:32 2011 for Chaste by  doxygen 1.6.3