#include <MonodomainDg0Assembler.hpp>
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. |
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.
typedef MonodomainDg0Assembler<ELEMENT_DIM, SPACE_DIM> MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::SelfType [protected] |
This type (to save typing).
Reimplemented from SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, false, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.
Definition at line 69 of file MonodomainDg0Assembler.hpp.
typedef SimpleDg0ParabolicAssembler<ELEMENT_DIM, SPACE_DIM, false, SelfType> MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::BaseClassType [protected] |
Base class type (to save typing).
Reimplemented from SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, false, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.
Definition at line 70 of file MonodomainDg0Assembler.hpp.
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] |
Constructor stores the mesh and pde and sets up boundary conditions.
pMesh | pointer to the mesh | |
pPde | pointer to the PDE | |
pBcc | pointer to the boundary conditions | |
numQuadPoints | number of quadrature points (defaults to 2) |
Definition at line 98 of file MonodomainDg0Assembler.cpp.
References AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainPde, AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 >::SetMatrixIsConstant(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > > >::SetMesh().
MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::~MonodomainDg0Assembler | ( | ) | [inline] |
Destructor.
Definition at line 116 of file MonodomainDg0Assembler.cpp.
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] |
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.
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.
void MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::ResetInterpolatedQuantities | ( | void | ) | [inline, protected, virtual] |
Overridden ResetInterpolatedQuantities() method.
Reimplemented from AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 49 of file MonodomainDg0Assembler.cpp.
References MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mSourceTerm.
void MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::IncrementInterpolatedQuantities | ( | double | phiI, | |
const Node< SPACE_DIM > * | pNode | |||
) | [inline, protected, virtual] |
Overridden IncrementInterpolatedQuantities() method.
phiI | ||
pNode |
Reimplemented from AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 56 of file MonodomainDg0Assembler.cpp.
References MonodomainPde< ELEM_DIM, SPACE_DIM >::ComputeNonlinearSourceTermAtNode(), Node< SPACE_DIM >::GetIndex(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > > >::mCurrentSolutionOrGuessReplicated, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainPde, and MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mSourceTerm.
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.
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().
void MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve | ( | Vec | initialSolution | ) | [inline, protected, virtual] |
Create the linear system object if it hasn't been already. Can use an initial solution as PETSc template, or base it on the mesh size.
initialSolution | an initial guess |
Reimplemented from AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > > >.
Definition at line 72 of file MonodomainDg0Assembler.cpp.
References AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >::InitialiseForSolve(), HeartConfig::Instance(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > > >::mpLinearSystem, LinearSystem::SetAbsoluteTolerance(), LinearSystem::SetKspType(), LinearSystem::SetMatrixIsSymmetric(), LinearSystem::SetPcType(), and LinearSystem::SetRelativeTolerance().
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.
const unsigned MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::S_DIM = SPACE_DIM [static] |
The space dimension (to save typing).
Reimplemented from SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, false, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.
Definition at line 57 of file MonodomainDg0Assembler.hpp.
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.
double MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mSourceTerm [protected] |
The source term.
Definition at line 63 of file MonodomainDg0Assembler.hpp.
Referenced by MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorTerm(), MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::IncrementInterpolatedQuantities(), and MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::ResetInterpolatedQuantities().
MonodomainPde<ELEMENT_DIM,SPACE_DIM>* MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::mpMonodomainPde [protected] |
The PDE to be solved.
Definition at line 66 of file MonodomainDg0Assembler.hpp.
Referenced by MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorTerm(), MonodomainMatrixBasedAssembler< ELEMENT_DIM, SPACE_DIM >::ConstructVectorForMatrixBasedRhsAssembly(), MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::IncrementInterpolatedQuantities(), MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::MonodomainDg0Assembler(), and MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::PrepareForAssembleSystem().