Chaste  Release::2017.1
MonodomainAssembler< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <MonodomainAssembler.hpp>

+ Inheritance diagram for MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >:
+ Collaboration diagram for MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >:

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)
 
- Public Member Functions inherited from AbstractCardiacFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, 1, false, true, CARDIAC >
 AbstractCardiacFeVolumeIntegralAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM > *pTissue)
 
- Public Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
 AbstractFeVolumeIntegralAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
virtual ~AbstractFeVolumeIntegralAssembler ()
 
- Public Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
 AbstractFeAssemblerCommon ()
 
void SetCurrentSolution (Vec currentSolution)
 
virtual ~AbstractFeAssemblerCommon ()
 
- Public Member Functions inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
 AbstractFeAssemblerInterface ()
 
void SetMatrixToAssemble (Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
 
void SetVectorToAssemble (Vec &rVecToAssemble, bool zeroVectorBeforeAssembly)
 
void Assemble ()
 
void AssembleMatrix ()
 
void AssembleVector ()
 
virtual ~AbstractFeAssemblerInterface ()
 

Protected Attributes

MassMatrixAssembler< ELEMENT_DIM, SPACE_DIM > mMassMatrixAssembler
 
MonodomainStiffnessMatrixAssembler< ELEMENT_DIM, SPACE_DIM > mStiffnessMatrixAssembler
 
- Protected Attributes inherited from AbstractCardiacFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, 1, false, true, CARDIAC >
AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM > * mpCardiacTissue
 
HeartConfigmpConfig
 
- Protected Attributes inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
 
GaussianQuadratureRule< ELEMENT_DIM > * mpQuadRule
 
- Protected Attributes inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
ReplicatableVector mCurrentSolutionOrGuessReplicated
 
- Protected Attributes inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
Vec mVectorToAssemble
 
Mat mMatrixToAssemble
 
bool mAssembleMatrix
 
bool mAssembleVector
 
bool mZeroMatrixBeforeAssembly
 
bool mZeroVectorBeforeAssembly
 
PetscInt mOwnershipRangeLo
 
PetscInt mOwnershipRangeHi
 

Additional Inherited Members

- Protected Types inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
typedef LinearBasisFunction< ELEMENT_DIM > BasisFunction
 
- Protected Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
void ComputeTransformedBasisFunctionDerivatives (const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rReturnValue)
 
void DoAssemble ()
 
virtual c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(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, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
 
virtual c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> ComputeVectorTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
 
virtual void AssembleOnElement (Element< ELEMENT_DIM, SPACE_DIM > &rElement, c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1) > &rAElem, c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> &rBElem)
 
virtual bool ElementAssemblyCriterion (Element< ELEMENT_DIM, SPACE_DIM > &rElement)
 
- Protected Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
virtual double GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown)
 
virtual void ResetInterpolatedQuantities ()
 
virtual void IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode)
 
virtual void IncrementInterpolatedGradientQuantities (const c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, unsigned phiIndex, const Node< SPACE_DIM > *pNode)
 

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 64 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 
)

Constructor

Parameters
pMeshpointer to the mesh
pTissuepointer to the tissue

Definition at line 57 of file MonodomainAssembler.cpp.

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 
)

ComputeMatrixTerm()

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

Parameters
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rXThe point in space
rUThe unknown as a vector, u(i) = u_i
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElementPointer to the element
Returns
stencil matrix

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

Definition at line 43 of file MonodomainAssembler.cpp.

References PdeSimulationTime::GetPdeTimeStepInverse().

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 70 of file MonodomainAssembler.hpp.

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 74 of file MonodomainAssembler.hpp.


The documentation for this class was generated from the following files: