AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > Class Template Reference

#include <AbstractAssembler.hpp>

Inherited by AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > [virtual], AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 > [virtual], AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 2 > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainRhsMatrixAssembler< DIM > > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainWithBathRhsMatrixAssembler< DIM > > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM > > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE > > [virtual], AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE > [virtual], and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > > [virtual].

Collaboration diagram for AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:
Collaboration graph
[legend]

List of all members.

Public Member Functions

virtual LinearSystem ** GetLinearSystem ()=0
 AbstractAssembler ()
void SetBoundaryConditionsContainer (BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pBoundaryConditions)
virtual ~AbstractAssembler ()
virtual void ResetInterpolatedQuantities ()
virtual void IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode)

Protected Member Functions

virtual void SetMatrixIsConst (bool matrixIsConstant=true)
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)=0
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)=0
virtual c_vector< double,
PROBLEM_DIM *ELEMENT_DIM > 
ComputeVectorSurfaceTerm (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &rSurfaceElement, c_vector< double, ELEMENT_DIM > &rPhi, ChastePoint< SPACE_DIM > &rX)=0
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, bool assembleVector, bool assembleMatrix)=0
virtual void AssembleOnSurfaceElement (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &rSurfaceElement, c_vector< double, PROBLEM_DIM *ELEMENT_DIM > &rBSurfElem)=0
virtual void AssembleSystem (bool assembleVector, bool assembleMatrix, Vec currentSolutionOrGuess=NULL, double currentTime=0.0)=0
virtual void PrepareForSolve ()=0
virtual void PrepareForAssembleSystem (Vec currentSolutionOrGuess, double currentTime)
virtual void FinaliseAssembleSystem (Vec currentSolutionOrGuess, double currentTime)
virtual void FinaliseLinearSystem (Vec currentSolutionOrGuess, double currentTime, bool assembleVector, bool assembleMatrix)
virtual void ApplyDirichletConditions (Vec currentSolutionOrGuess, bool applyToMatrix)=0
virtual bool ProblemIsNonlinear ()=0
virtual Vec StaticSolve (Vec currentSolutionOrGuess=NULL, double currentTime=0.0, bool assembleMatrix=true)=0
virtual void InitialiseForSolve (Vec initialGuess)=0
virtual ReplicatableVectorrGetCurrentSolutionOrGuess ()=0
void ApplyNeummanBoundaryConditions ()

Protected Attributes

BoundaryConditionsContainer
< ELEMENT_DIM, SPACE_DIM,
PROBLEM_DIM > * 
mpBoundaryConditions

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
class AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >

AbstractAssembler

Base class from which all solvers for linear and nonlinear PDEs inherit. Templated over the PROBLEM_DIM so also handles problems with more than one unknown variable (ie those of the form u_xx + v = 0, v_xx + 2u = 1, where PROBLEM_DIM is equal to 2)

It defines a common interface for AssembleSystem, AssembleOnElement and AssembleOnSurfaceElement. Each of these work for any PROBLEM_DIM>=1. Each of these methods work in both the dynamic case (when there is a current solution available) and the static case. The same code is used for the nonlinear and linear cases. Default code is defined in AbstractStaticAssembler

user calls:

Solve(). In the linear case Solve() calls AssembleSystem() directly, in the nonlinear case Solve() calls the PETSc nonlinear solver which then calls AssembleResidual or AssembleJacobian, both of which call AssembleSystem():

AssembleSystem(). (implemented in AbstractStaticAssembler, loops over elements and adds to the linear system or residual vector or jacobian matrix) AssembleSystem() calls:

AssembleOnElement() and AssembleOnSurfaceElement(). (implemented in AbstractStaticAssembler. These loop over gauss points and create the element stiffness matrix and vector in the linear case ). They call:

ComputeMatrixTerm(), ComputeVectorTerm(), ComputeVectorSurfaceTerm() (implemented in the concrete assembler class (eg SimpleDg0ParabolicAssembler), which tells this assembler exactly what function of bases, position, pde constants etc to add to the element stiffness matrix/vector).

Definition at line 77 of file AbstractAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractAssembler (  )  [inline]

Default constructor.

Definition at line 412 of file AbstractAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractAssembler (  )  [inline, virtual]

Delete any memory allocated by this class.

Definition at line 305 of file AbstractAssembler.hpp.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions ( Vec  currentSolutionOrGuess,
bool  applyToMatrix 
) [protected, pure virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyNeummanBoundaryConditions (  )  [inline, protected]

Apply Neumann boundary conditions to the RHS vector by looping over surface elements (though actually looping over the boundary condition objects).

Note for PROBLEM_DIM>1. We assume that if an element has a boundary condition on any unknown there is a boundary condition on unknown 0. This can be so for any problem by adding zero constant conditions where required although this is a bit inefficient. Proper solution involves changing BCC to have a map of arrays boundary conditions rather than an array of maps.

Definition at line 385 of file AbstractAssembler.hpp.

References AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AssembleOnSurfaceElement(), GenericEventHandler< 13, HeartEventHandler >::BeginEvent(), GenericEventHandler< 13, HeartEventHandler >::EndEvent(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetLinearSystem(), AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::GetStiffnessMatrixGlobalIndices(), and AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions.

Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleSystem(), and AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DoMatrixBasedRhsAssembly().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::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,
bool  assembleVector,
bool  assembleMatrix 
) [protected, pure virtual]

Calculate the contribution of a single element to the linear system.

Parameters:
rElement The element to assemble on.
rAElem The element's contribution to the LHS matrix is returned in this n by n matrix, where n is the no. of nodes in this element. There is no need to zero this matrix before calling.
rBElem The element's contribution to the RHS vector is returned in this vector of length n, the no. of nodes in this element. There is no need to zero this vector before calling.
assembleVector a bool stating whether to assemble the load vector (in the linear case) or the residual vector (in the nonlinear case)
assembleMatrix a bool stating whether to assemble the stiffness matrix (in the linear case) or the Jacobian matrix (in the nonlinear case)

Called by AssembleSystem() Calls ComputeMatrixTerm() etc

Implemented in AbstractStaticAssembler

Implemented in AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainWithBathRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM > >, and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AssembleOnSurfaceElement ( const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &  rSurfaceElement,
c_vector< double, PROBLEM_DIM *ELEMENT_DIM > &  rBSurfElem 
) [protected, pure virtual]

Calculate the contribution of a single surface element with Neumann boundary condition to the linear system.

Parameters:
rSurfaceElement The element to assemble on.
rBSurfElem The element's contribution to the RHS vector is returned in this vector of length n, the no. of nodes in this element. There is no need to zero this vector before calling.

Implemented in AbstractStaticAssembler

Implemented in AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainWithBathRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM > >, and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

Referenced by AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyNeummanBoundaryConditions().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AssembleSystem ( bool  assembleVector,
bool  assembleMatrix,
Vec  currentSolutionOrGuess = NULL,
double  currentTime = 0.0 
) [protected, pure virtual]

AssembleSystem - the major method for all assemblers

Assemble the linear system for a linear PDE, or the residual or Jacobian for nonlinear PDEs. Loops over each element (and each each surface element if there are non-zero Neumann boundary conditions), calls AssembleOnElement() and adds the contribution to the linear system.

Takes in current solution and time if necessary but only used if the problem is a dynamic one. This method uses PROBLEM_DIM and can assemble linear systems for any number of unknown variables.

Parameters:
assembleVector Whether to assemble the RHS vector of the linear system (i.e. the residual vector for nonlinear problems).
assembleMatrix Whether to assemble the LHS matrix of the linear system (i.e. the jacobian matrix for nonlinear problems).
currentSolutionOrGuess The current solution in a linear dynamic problem, or the current guess in a nonlinear problem. Should be NULL for linear static problems. Defaults to NULL.
currentTime The current time for dynamic problems. Not used in static problems. Defaults to 0.0.

Implemented in AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >, AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainWithBathRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM > >, and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

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

This method returns the matrix to be added to element stiffness matrix for a given gauss point. The arguments are the bases, bases gradients, x and current solution computed at the Gauss point. The returned matrix will be multiplied by the gauss weight and jacobian determinent and added to the element stiffness matrix (see AssembleOnElement()).

--This method has to be implemented in the concrete class--

NOTE: for linear problems rGradU is NOT set up correctly because it should not be needed

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.

Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleOnElement().

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

This method returns the vector to be added to element stiffness vector for a given gauss point in BoundaryElement. The arguments are the bases, x and current solution computed at the Gauss point. The returned vector will be multiplied by the gauss weight and jacobian determinent and added to the element stiffness matrix (see AssembleOnElement()).

--This method has to be implemented in the concrete class--

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

Implemented in BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM >, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE >, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM >, and BidomainDg0Assembler< DIM, DIM >.

Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleOnSurfaceElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual c_vector<double,PROBLEM_DIM*(ELEMENT_DIM+1)> AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_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, PROBLEM_DIM > &  rU,
c_matrix< double, PROBLEM_DIM, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
) [protected, pure virtual]

This method returns the vector to be added to element stiffness vector for a given gauss point. The arguments are the bases, x and current solution computed at the Gauss point. The returned vector will be multiplied by the gauss weight and jacobian determinent and added to the element stiffness matrix (see AssembleOnElement()).

--This method has to be implemented in the concrete class--

NOTE: for linear problems rGradPhi and rGradU are NOT set up correctly because they should not be needed

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

Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleOnElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::FinaliseAssembleSystem ( Vec  currentSolutionOrGuess,
double  currentTime 
) [inline, protected, virtual]

This method is called at the end of AssembleSystem() and should be overloaded in the concrete assembler class if there is any further work to be done.

Parameters:
currentSolutionOrGuess 
currentTime 

Reimplemented in BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >, and BidomainDg0Assembler< DIM, DIM >.

Definition at line 212 of file AbstractAssembler.hpp.

Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleSystem(), and AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DoMatrixBasedRhsAssembly().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::FinaliseLinearSystem ( Vec  currentSolutionOrGuess,
double  currentTime,
bool  assembleVector,
bool  assembleMatrix 
) [inline, protected, virtual]

Can be overloaded if the user needs to edit the linear system after the boundary conditions have been added but before it is solved.

Parameters:
currentSolutionOrGuess 
currentTime 
assembleVector 
assembleMatrix 

Reimplemented in BidomainWithBathAssembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 224 of file AbstractAssembler.hpp.

Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleSystem().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual LinearSystem** AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetLinearSystem (  )  [pure virtual]

Accessor method that subclasses can use to get to useful data.

Implemented in AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainWithBathRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM > >, and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

Referenced by AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyNeummanBoundaryConditions(), AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DoMatrixBasedRhsAssembly(), BidomainWithBathAssembler< ELEMENT_DIM, SPACE_DIM >::FinaliseLinearSystem(), and AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::IncrementInterpolatedQuantities ( double  phiI,
const Node< SPACE_DIM > *  pNode 
) [inline, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve ( Vec  initialGuess  )  [protected, pure virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::PrepareForAssembleSystem ( Vec  currentSolutionOrGuess,
double  currentTime 
) [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:
currentSolutionOrGuess 
currentTime 

Reimplemented in BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >, and BidomainDg0Assembler< DIM, DIM >.

Definition at line 202 of file AbstractAssembler.hpp.

Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleSystem(), and AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DoMatrixBasedRhsAssembly().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::PrepareForSolve (  )  [protected, pure virtual]

This method is called at the beginning of Solve(). Subclass assemblers can use it to check everything has been set up correctly

Implemented in AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainWithBathRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, BidomainRhsMatrixAssembler< DIM > >, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM > >, and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >.

Referenced by AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual bool AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ProblemIsNonlinear (  )  [protected, pure virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ResetInterpolatedQuantities ( void   )  [inline, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual ReplicatableVector& AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::rGetCurrentSolutionOrGuess (  )  [protected, pure virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetBoundaryConditionsContainer ( BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *  pBoundaryConditions  )  [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual void AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetMatrixIsConst ( bool  matrixIsConstant = true  )  [inline, protected, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
virtual Vec AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::StaticSolve ( Vec  currentSolutionOrGuess = NULL,
double  currentTime = 0.0,
bool  assembleMatrix = true 
) [protected, pure virtual]

Perform the work of a single solve, but without any initialisation. Static assemblers must implement this method.

Parameters:
currentSolutionOrGuess either the current solution (dynamic problem) or initial guess (static problem); optional in some cases
currentTime for a dynamic problem, the current time
assembleMatrix whether to assemble the matrix (it may have been done by a previous call)
Returns:
the solution vector

Implemented in AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >, AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >, AbstractLinearAssembler< DIM, DIM, 2, false, BidomainWithBathRhsMatrixAssembler< DIM > >, AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, false, MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >, AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, true, SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE > >, AbstractLinearAssembler< DIM, DIM, 2, false, BidomainRhsMatrixAssembler< DIM > >, AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, 2, false, BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM > >, AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >, AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, false, MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM > >, and AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >.

Referenced by AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve().


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions [protected]

Boundary conditions to be applied

Definition at line 37 of file AbstractAssembler.hpp.

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::ApplyDirichletConditions(), AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::ApplyDirichletConditions(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyNeummanBoundaryConditions(), BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::BidomainDg0Assembler(), BidomainRhsMatrixAssembler< DIM >::BidomainRhsMatrixAssembler(), BidomainWithBathRhsMatrixAssembler< DIM >::BidomainWithBathRhsMatrixAssembler(), BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::CheckCompatibilityCondition(), SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorSurfaceTerm(), SimpleLinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM, CONCRETE >::ComputeVectorSurfaceTerm(), SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::ComputeVectorSurfaceTerm(), MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorSurfaceTerm(), BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorSurfaceTerm(), BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::FinaliseAssembleSystem(), MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::MonodomainDg0Assembler(), MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM >::MonodomainRhsMatrixAssembler(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::PrepareForSolve(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetBoundaryConditionsContainer(), BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::SetFixedExtracellularPotentialNodes(), AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::StaticSolve(), BidomainRhsMatrixAssembler< DIM >::~BidomainRhsMatrixAssembler(), BidomainWithBathRhsMatrixAssembler< DIM >::~BidomainWithBathRhsMatrixAssembler(), and MonodomainRhsMatrixAssembler< ELEMENT_DIM, SPACE_DIM >::~MonodomainRhsMatrixAssembler().


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

Generated by  doxygen 1.6.2