MonodomainAssembler< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <MonodomainAssembler.hpp>

Inheritance diagram for MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

 MonodomainAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *pTissue, double dt, unsigned numQuadPoints=2)

Protected 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)
c_vector< double,
1 *(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, 1 > &rU, c_matrix< double, 1, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
c_vector< double, ELEMENT_DIM > ComputeVectorSurfaceTerm (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &rSurfaceElement, c_vector< double, ELEMENT_DIM > &rPhi, ChastePoint< SPACE_DIM > &rX)
void ResetInterpolatedQuantities (void)
void IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode)

Protected Attributes

MonodomainTissue< ELEMENT_DIM,
SPACE_DIM > * 
mpMonodomainTissue
HeartConfigmpConfig
double mIionic
double mIIntracellularStimulus
double mDt


Detailed Description

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

Assembler for assembling the LHS matrix and RHS vector of the linear systems solved in monodomain problems

Definition at line 41 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,
double  dt,
unsigned  numQuadPoints = 2 
) [inline]

Constructor

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

Definition at line 99 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, protected]

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

Definition at line 35 of file MonodomainAssembler.cpp.

References HeartConfig::GetCapacitance(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), HeartConfig::GetSurfaceAreaToVolumeRatio(), MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mDt, MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpConfig, and MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainTissue.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 1 *(ELEMENT_DIM+1)> MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorTerm ( 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, protected]

ComputeVectorTerm()

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

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

Definition at line 62 of file MonodomainAssembler.cpp.

References HeartConfig::GetCapacitance(), HeartConfig::GetSurfaceAreaToVolumeRatio(), MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mDt, MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mIIntracellularStimulus, MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mIionic, and MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpConfig.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, ELEMENT_DIM > MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorSurfaceTerm ( const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &  rSurfaceElement,
c_vector< double, ELEMENT_DIM > &  rPhi,
ChastePoint< SPACE_DIM > &  rX 
) [inline, protected, virtual]

ComputeVectorSurfaceTerm()

This method is called by AssembleOnSurfaceElement() and tells the assembler what to add to the element stiffness matrix arising from surface element contributions.

Parameters:
rSurfaceElement the element which is being considered.
rPhi The basis functions, rPhi(i) = phi_i, i=1..numBases
rX The point in space

Reimplemented from AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, 1, true, true, CARDIAC >.

Definition at line 79 of file MonodomainAssembler.cpp.

References AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, 1, true, true, CARDIAC >::mpBoundaryConditions.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::ResetInterpolatedQuantities ( void   )  [inline, protected, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::IncrementInterpolatedQuantities ( double  phiI,
const Node< SPACE_DIM > *  pNode 
) [inline, protected, virtual]


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainTissue [protected]

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mIionic [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mIIntracellularStimulus [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >::mDt [protected]


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

Generated on Mon Nov 1 12:37:03 2010 for Chaste by  doxygen 1.5.5