MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <MonodomainDg0Assembler.hpp>

Inheritance diagram for MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

 MonodomainDg0Assembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainPde< ELEMENT_DIM, SPACE_DIM > *pPde, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBcc, unsigned numQuadPoints=2)
 ~MonodomainDg0Assembler ()

Static Public Attributes

static const unsigned E_DIM = ELEMENT_DIM
static const unsigned S_DIM = SPACE_DIM
static const unsigned P_DIM = 1u

Protected Types

typedef MonodomainDg0Assembler
< ELEMENT_DIM, SPACE_DIM > 
SelfType
typedef
SimpleDg0ParabolicAssembler
< ELEMENT_DIM, SPACE_DIM,
false, SelfType
BaseClassType

Protected Member Functions

virtual 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)
void ResetInterpolatedQuantities ()
void IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode)
virtual void PrepareForAssembleSystem (Vec existingSolutionOrGuess, double time)
void InitialiseForSolve (Vec initialSolution)

Protected Attributes

double mSourceTerm
MonodomainPde< ELEMENT_DIM,
SPACE_DIM > * 
mpMonodomainPde

Friends

class AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, 1u, false, BaseClassType >
 Allow the AbstractStaticAssembler to call our private/protected methods using static polymorphism.


Detailed Description

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

MonodomainDg0Assembler

This is essentially the same as the SimpleDg0ParabolicAssembler (which it inherits from), except that the source term (ie ionic current + stimulus) is interpolated from their nodal values, instead of computed at the gauss point, since they are only known at the nodes.

Also, the MonodomainAssembler automatically creates zero neumann boundary conditions when constructed and therefore does not need to take in a BoundaryConditionsContainer.

The user should call Solve() from the superclass AbstractDynamicAssemblerMixin.

Definition at line 52 of file MonodomainDg0Assembler.hpp.


Member Typedef Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
typedef MonodomainDg0Assembler<ELEMENT_DIM, SPACE_DIM> MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::SelfType [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
typedef SimpleDg0ParabolicAssembler<ELEMENT_DIM, SPACE_DIM, false, SelfType> MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::BaseClassType [protected]


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::MonodomainDg0Assembler ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
MonodomainPde< ELEMENT_DIM, SPACE_DIM > *  pPde,
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *  pBcc,
unsigned  numQuadPoints = 2 
) [inline]

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

Destructor.

Definition at line 116 of file MonodomainDg0Assembler.cpp.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 1 *(ELEMENT_DIM+1)> MonodomainDg0Assembler< 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, virtual]

ComputeVectorTerm()

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

Here, the SimpleDg0ParabolicAssembler version of this method is overloaded using the interpolated source term.

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 from SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, false, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

Definition at line 35 of file MonodomainDg0Assembler.cpp.

References MonodomainPde< ELEM_DIM, SPACE_DIM >::ComputeDuDtCoefficientFunction(), AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 >::mDtInverse, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainPde, and MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mSourceTerm.

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

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::PrepareForAssembleSystem ( Vec  existingSolutionOrGuess,
double  time 
) [inline, protected, virtual]

This method is called at the beginning of AssembleSystem() and should be overloaded in the concrete assembler class if there is any work to be done before assembling, for example integrating ODEs such as in the Monodomain assembler.

Parameters:
existingSolutionOrGuess 
time 

Reimplemented from AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 64 of file MonodomainDg0Assembler.cpp.

References AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 >::mDt, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainPde, and AbstractCardiacPde< ELEM_DIM, SPACE_DIM >::SolveCellSystems().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve ( Vec  initialSolution  )  [inline, protected, virtual]


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
const unsigned MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::E_DIM = ELEMENT_DIM [static]

The element dimension (to save typing).

Reimplemented from SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, false, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

Definition at line 56 of file MonodomainDg0Assembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
const unsigned MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::S_DIM = SPACE_DIM [static]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
const unsigned MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::P_DIM = 1u [static]

The problem dimension (to save typing).

Reimplemented from SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, false, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

Definition at line 58 of file MonodomainDg0Assembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mSourceTerm [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
MonodomainPde<ELEMENT_DIM,SPACE_DIM>* MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainPde [protected]


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

Generated on Tue Aug 4 16:11:27 2009 for Chaste by  doxygen 1.5.5