BidomainAssembler< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <BidomainAssembler.hpp>

Inheritance diagram for BidomainAssembler< ELEMENT_DIM, SPACE_DIM >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

 BidomainAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, BidomainTissue< SPACE_DIM > *pTissue, unsigned numQuadPoints=2)
 ~BidomainAssembler ()

Protected Member Functions

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

Protected Attributes

HeartConfigmpConfig
double mIionic
double mIIntracellularStimulus
double mIExtracellularStimulus


Detailed Description

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

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

Definition at line 43 of file BidomainAssembler.hpp.


Constructor & Destructor Documentation

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

Constructor stores the mesh and pde and sets up boundary conditions.

Parameters:
pMesh pointer to the mesh
pTissue pointer to the tissue
numQuadPoints number of quadrature points in each dimension

Definition at line 161 of file BidomainAssembler.cpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BidomainAssembler< ELEMENT_DIM, SPACE_DIM >::~BidomainAssembler (  )  [inline]

Destructor.

Definition at line 143 of file BidomainAssembler.hpp.


Member Function Documentation

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

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_matrix< double, 2 *(ELEMENT_DIM+1), 2 *(ELEMENT_DIM+1)> BidomainAssembler< 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, 2 > &  rU,
c_matrix< double, 2, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
) [inline, protected, virtual]

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

Reimplemented in BidomainWithBathAssembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 56 of file BidomainAssembler.cpp.

References HeartConfig::GetCapacitance(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), PdeSimulationTime::GetPdeTimeStepInverse(), HeartConfig::GetSurfaceAreaToVolumeRatio(), AbstractCardiacFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, 2, true, true, CARDIAC >::mpCardiacTissue, and BidomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpConfig.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 2 *(ELEMENT_DIM+1)> BidomainAssembler< 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, 2 > &  u,
c_matrix< double, 2, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
) [inline, protected, virtual]

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
u 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

Reimplemented in BidomainWithBathAssembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 112 of file BidomainAssembler.cpp.

References HeartConfig::GetCapacitance(), PdeSimulationTime::GetPdeTimeStepInverse(), HeartConfig::GetSurfaceAreaToVolumeRatio(), BidomainAssembler< ELEMENT_DIM, SPACE_DIM >::mIIntracellularStimulus, BidomainAssembler< ELEMENT_DIM, SPACE_DIM >::mIionic, and BidomainAssembler< ELEMENT_DIM, SPACE_DIM >::mpConfig.

Referenced by BidomainWithBathAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorTerm().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 2 *ELEMENT_DIM > BidomainAssembler< 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, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >.

Definition at line 139 of file BidomainAssembler.cpp.

References AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpBoundaryConditions.


Member Data Documentation

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

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

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double BidomainAssembler< ELEMENT_DIM, SPACE_DIM >::mIExtracellularStimulus [protected]

Extracellular stimulus to be interpolated from cache

Definition at line 53 of file BidomainAssembler.hpp.


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

Generated on Mon Apr 18 11:36:23 2011 for Chaste by  doxygen 1.5.5