SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > Class Template Reference

#include <SimpleDg0ParabolicAssembler.hpp>

Inheritance diagram for SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >:

Inheritance graph
[legend]
Collaboration diagram for SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >:

Collaboration graph
[legend]

List of all members.

Public Member Functions

 SimpleDg0ParabolicAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, AbstractLinearParabolicPde< ELEMENT_DIM, SPACE_DIM > *pPde, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBoundaryConditions, unsigned numQuadPoints=2)
virtual void PrepareForSolve ()
Vec Solve (Vec currentSolutionOrGuess=NULL, double currentTime=0.0)

Static Public Attributes

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

Protected Member Functions

virtual 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)
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)
virtual c_vector< double,
ELEMENT_DIM > 
ComputeVectorSurfaceTerm (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &rSurfaceElement, c_vector< double, ELEMENT_DIM > &rPhi, ChastePoint< SPACE_DIM > &rX)

Private Types

typedef
SimpleDg0ParabolicAssembler
< ELEMENT_DIM, SPACE_DIM,
NON_HEART, CONCRETE > 
SelfType
typedef
AbstractLinearAssembler
< ELEMENT_DIM, SPACE_DIM,
1, NON_HEART, SelfType
BaseClassType

Private Attributes

AbstractLinearParabolicPde
< ELEMENT_DIM, SPACE_DIM > * 
mpParabolicPde

Friends

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


Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE = boost::mpl::void_>
class SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >

SimpleDg0ParabolicAssembler

Assembler for solving AbstractLinearParabolicPdes

Definition at line 52 of file SimpleDg0ParabolicAssembler.hpp.


Member Typedef Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE = boost::mpl::void_>
typedef SimpleDg0ParabolicAssembler<ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE> SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::SelfType [private]

This type (to save typing).

Reimplemented in MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 65 of file SimpleDg0ParabolicAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE = boost::mpl::void_>
typedef AbstractLinearAssembler<ELEMENT_DIM, SPACE_DIM, 1, NON_HEART, SelfType> SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::BaseClassType [private]

Base class type (to save typing).

Reimplemented in MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 66 of file SimpleDg0ParabolicAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE>
SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::SimpleDg0ParabolicAssembler ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
AbstractLinearParabolicPde< ELEMENT_DIM, SPACE_DIM > *  pPde,
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *  pBoundaryConditions,
unsigned  numQuadPoints = 2 
) [inline]

Constructor stores the mesh, pde and boundary conditons, and calls base constructor.

Parameters:
pMesh pointer to the mesh
pPde pointer to the PDE
pBoundaryConditions pointer to the boundary conditions
numQuadPoints number of quadrature points (defaults to 2)

Definition at line 39 of file SimpleDg0ParabolicAssemblerImplementation.hpp.

References SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::mpParabolicPde, AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetBoundaryConditionsContainer(), AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 >::SetMatrixIsConstant(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE > >::SetMesh().


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE>
c_matrix< double, 1 *(ELEMENT_DIM+1), 1 *(ELEMENT_DIM+1)> SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::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, virtual]

The term to be added to the element stiffness matrix:

grad_phi[row] . ( pde_diffusion_term * grad_phi[col]) + (1.0/mDt) * pPde->ComputeDuDtCoefficientFunction(rX) * rPhi[row] * rPhi[col]

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
Todo:
should this be rU?
Parameters:
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 77 of file SimpleDg0ParabolicAssemblerImplementation.hpp.

References AbstractLinearParabolicPde< ELEM_DIM, SPACE_DIM >::ComputeDiffusionTerm(), AbstractLinearParabolicPde< ELEM_DIM, SPACE_DIM >::ComputeDuDtCoefficientFunction(), AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 >::mDtInverse, and SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::mpParabolicPde.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE>
c_vector< double, 1 *(ELEMENT_DIM+1)> SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::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]

The term to be added 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

Reimplemented in MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 93 of file SimpleDg0ParabolicAssemblerImplementation.hpp.

References AbstractLinearParabolicPde< ELEM_DIM, SPACE_DIM >::ComputeDuDtCoefficientFunction(), AbstractLinearParabolicPde< ELEM_DIM, SPACE_DIM >::ComputeLinearSourceTerm(), AbstractLinearParabolicPde< ELEM_DIM, SPACE_DIM >::ComputeNonlinearSourceTerm(), AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 >::mDtInverse, and SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::mpParabolicPde.

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

The term arising from boundary conditions to be added to the element stiffness vector.

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

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 108 of file SimpleDg0ParabolicAssemblerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::GetNeumannBCValue(), and AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE>
void SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::PrepareForSolve (  )  [inline, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE>
Vec SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::Solve ( Vec  currentSolutionOrGuess = NULL,
double  currentTime = 0.0 
) [inline]

Solve the static pde.

The mesh, pde and boundary conditions container must be set before Solve() is called.

Parameters:
currentSolutionOrGuess either the current solution or initial guess (defaults to NULL)
currentTime the current time (defaults to 0.0)

Reimplemented from AbstractDynamicAssemblerMixin< ELEMENT_DIM, SPACE_DIM, 1 >.

Definition at line 68 of file SimpleDg0ParabolicAssemblerImplementation.hpp.

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


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE = boost::mpl::void_>
const unsigned SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::E_DIM = ELEMENT_DIM [static]

The element dimension (to save typing).

Reimplemented in MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 57 of file SimpleDg0ParabolicAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE = boost::mpl::void_>
const unsigned SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::S_DIM = SPACE_DIM [static]

The space dimension (to save typing).

Reimplemented in MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 58 of file SimpleDg0ParabolicAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE = boost::mpl::void_>
const unsigned SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::P_DIM = 1u [static]

The problem dimension (to save typing).

Reimplemented in MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 59 of file SimpleDg0ParabolicAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, bool NON_HEART, class CONCRETE = boost::mpl::void_>
AbstractLinearParabolicPde<ELEMENT_DIM, SPACE_DIM>* SimpleDg0ParabolicAssembler< ELEMENT_DIM, SPACE_DIM, NON_HEART, CONCRETE >::mpParabolicPde [private]


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

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