MonodomainAssembler< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <MonodomainAssembler.hpp>

Inherits AbstractCardiacFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, 1, false, true, CARDIAC >.

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

List of all members.

Public Member Functions

c_matrix< double,
1 *(ELEMENT_DIM+1),
1 *(ELEMENT_DIM+1)> 
ComputeMatrixTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, 1 > &rU, c_matrix< double, 1, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
 MonodomainAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *pTissue, unsigned numQuadPoints=2)

Protected Attributes

HeartConfigmpConfig
MassMatrixAssembler
< ELEMENT_DIM, SPACE_DIM > 
mMassMatrixAssembler
MonodomainStiffnessMatrixAssembler
< ELEMENT_DIM, SPACE_DIM > 
mStiffnessMatrixAssembler

Detailed Description

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

Assembler, mainly used for assembling the LHS matrix of the linear system that arises when the monodomain equations are discretised.

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 assembler is used for assembling the matrix A :=(chi*C/dt) M + K. Hence, this class inherits from AbstractCardiacFeVolumeIntegralAssembler and implements the method ComputeMatrixTerm().

Definition at line 57 of file MonodomainAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::MonodomainAssembler ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *  pTissue,
unsigned  numQuadPoints = 2 
) [inline]

Constructor

Parameters:
pMesh pointer to the mesh
pTissue pointer to the tissue
numQuadPoints number of quadrature points (defaults to 2)

Definition at line 51 of file MonodomainAssembler.cpp.

References HeartConfig::Instance(), and MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpConfig.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_matrix< double, 1 *(ELEMENT_DIM+1), 1 *(ELEMENT_DIM+1)> MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm ( c_vector< double, ELEMENT_DIM+1 > &  rPhi,
c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &  rGradPhi,
ChastePoint< SPACE_DIM > &  rX,
c_vector< double, 1 > &  rU,
c_matrix< double, 1, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
) [inline]

ComputeMatrixTerm()

This method is called by AssembleOnElement() and tells the assembler the contribution to add to the element stiffness matrix.

Parameters:
rPhi The basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhi Basis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rX The point in space
rU The unknown as a vector, u(i) = u_i
rGradU The gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElement Pointer to the element

Am and Cm are set as scaling factors for the mass matrix in its constructor.

Definition at line 36 of file MonodomainAssembler.cpp.

References PdeSimulationTime::GetPdeTimeStepInverse(), MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mMassMatrixAssembler, and MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mStiffnessMatrixAssembler.


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MassMatrixAssembler<ELEMENT_DIM, SPACE_DIM> MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mMassMatrixAssembler [protected]

This assembler uses another assembler, though just for calling the ComputeMatrixTerm() method.

Definition at line 66 of file MonodomainAssembler.hpp.

Referenced by MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
HeartConfig* MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpConfig [protected]

Local cache of the configuration singleton instance

Definition at line 62 of file MonodomainAssembler.hpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainStiffnessMatrixAssembler<ELEMENT_DIM, SPACE_DIM> MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mStiffnessMatrixAssembler [protected]

This assembler uses another assembler, though just for calling the ComputeMatrixTerm() method.

Definition at line 70 of file MonodomainAssembler.hpp.

Referenced by MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm().


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