BidomainSolver< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <BidomainSolver.hpp>

Inherits AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >.

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

List of all members.

Public Member Functions

 BidomainSolver (bool bathSimulation, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, BidomainTissue< SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 2 > *pBoundaryConditions, unsigned numQuadPoints=2)

Private Member Functions

void InitialiseForSolve (Vec initialSolution)
void SetupLinearSystem (Vec currentSolution, bool computeMatrix)

Private Attributes

Mat mMassMatrix
Vec mVecForConstructingRhs
BidomainAssembler< ELEMENT_DIM,
SPACE_DIM > * 
mpBidomainAssembler
BidomainNeumannSurfaceTermAssembler
< ELEMENT_DIM, SPACE_DIM > * 
mpBidomainNeumannSurfaceTermAssembler
BidomainCorrectionTermAssembler
< ELEMENT_DIM, SPACE_DIM > * 
mpBidomainCorrectionTermAssembler

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
class BidomainSolver< ELEMENT_DIM, SPACE_DIM >

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

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

[ (chi*C/dt) M + K1 K1 ] [ V^{n+1} ] = [ (chi*C/dt) M V^{n} + M F^{n} + c1_surf ] [ K1 K2 ] [ PhiE^{n+1}] [ c2_surf ]

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

This solver uses two assemblers, one to assemble the whole LHS matrix, and also to compute c1_surf and c2_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 vector [c_correction, 0] is added to the above, and another assembler is used to create the c_correction.

Definition at line 75 of file BidomainSolver.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BidomainSolver< ELEMENT_DIM, SPACE_DIM >::BidomainSolver ( bool  bathSimulation,
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
BidomainTissue< SPACE_DIM > *  pTissue,
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 2 > *  pBoundaryConditions,
unsigned  numQuadPoints = 2 
) [inline]

Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void BidomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve ( Vec  initialSolution  )  [inline, private, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void BidomainSolver< 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 59 of file BidomainSolver.cpp.

References BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::Assemble(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleVector(), DistributedVector::Begin(), GenericEventHandler< 16, HeartEventHandler >::BeginEvent(), DistributedVectorFactory::CreateDistributedVector(), DistributedVector::End(), GenericEventHandler< 16, HeartEventHandler >::EndEvent(), PetscMatTools::Finalise(), AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseForBath(), LinearSystem::FinaliseLhsMatrix(), LinearSystem::FinaliseRhsVector(), HeartConfig::GetCapacitance(), PdeSimulationTime::GetPdeTimeStepInverse(), HeartConfig::GetSurfaceAreaToVolumeRatio(), HeartConfig::Instance(), HeartRegionCode::IsRegionBath(), AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::mBathSimulation, BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mMassMatrix, BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBidomainAssembler, BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBidomainCorrectionTermAssembler, BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBidomainNeumannSurfaceTermAssembler, AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBidomainTissue, AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBoundaryConditions, AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpLinearSystem, AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mVecForConstructingRhs, DistributedVector::Restore(), LinearSystem::rGetLhsMatrix(), LinearSystem::rGetRhsVector(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetMatrixToAssemble(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetVectorToAssemble(), and LinearSystem::SwitchWriteModeLhsMatrix().


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Mat BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mMassMatrix [private]

Mass matrix, used to computing the RHS vector (actually: mass-matrix in voltage-voltage block, zero elsewhere)

Definition at line 81 of file BidomainSolver.hpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BidomainAssembler<ELEMENT_DIM,SPACE_DIM>* BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBidomainAssembler [private]

The bidomain assembler, used to set up the LHS matrix

Definition at line 90 of file BidomainSolver.hpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BidomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>* BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBidomainCorrectionTermAssembler [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 99 of file BidomainSolver.hpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BidomainNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM>* BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mpBidomainNeumannSurfaceTermAssembler [private]

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

Definition at line 93 of file BidomainSolver.hpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Vec BidomainSolver< ELEMENT_DIM, SPACE_DIM >::mVecForConstructingRhs [private]

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

Definition at line 87 of file BidomainSolver.hpp.

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


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