AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE > Class Template Reference

#include <AbstractStaticAssembler.hpp>

Inheritance diagram for AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

LinearSystem ** GetLinearSystem ()
 AbstractStaticAssembler (unsigned numQuadPoints=2)
void SetNumberOfQuadraturePointsPerDimension (unsigned numQuadPoints)
void SetMesh (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
virtual ~AbstractStaticAssembler ()

Protected Types

typedef LinearBasisFunction
< ELEMENT_DIM > 
BasisFunction
typedef LinearBasisFunction
< ELEMENT_DIM-1 > 
SurfaceBasisFunction

Protected Member Functions

void ComputeTransformedBasisFunctionDerivatives (const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rReturnValue)
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)
virtual void AssembleOnSurfaceElement (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &rSurfaceElement, c_vector< double, PROBLEM_DIM *ELEMENT_DIM > &rBSurfElem)
virtual void AssembleSystem (bool assembleVector, bool assembleMatrix, Vec currentSolutionOrGuess=NULL, double currentTime=0.0)
virtual void PrepareForSolve ()
ReplicatableVectorrGetCurrentSolutionOrGuess ()
virtual double GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown)

Protected Attributes

AbstractTetrahedralMesh
< ELEMENT_DIM, SPACE_DIM > * 
mpMesh
GaussianQuadratureRule
< ELEMENT_DIM > * 
mpQuadRule
GaussianQuadratureRule
< ELEMENT_DIM-1 > * 
mpSurfaceQuadRule
ReplicatableVector mCurrentSolutionOrGuessReplicated
LinearSystemmpLinearSystem


Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
class AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >

AbstractStaticAssembler

Implmentation of main assembler methods so that the virtual base class does not need to contain data (which should improve performance).

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 default code 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

See also documentation for AbstractAssembler

Definition at line 104 of file AbstractStaticAssembler.hpp.


Member Typedef Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
typedef LinearBasisFunction<ELEMENT_DIM> AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::BasisFunction [protected]

Basis function for use with normal elements

Definition at line 118 of file AbstractStaticAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
typedef LinearBasisFunction<ELEMENT_DIM-1> AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::SurfaceBasisFunction [protected]

Basis function for use with boundary elements

Definition at line 120 of file AbstractStaticAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AbstractStaticAssembler ( unsigned  numQuadPoints = 2  )  [inline]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::~AbstractStaticAssembler (  )  [inline, virtual]


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::ComputeTransformedBasisFunctionDerivatives ( const ChastePoint< ELEMENT_DIM > &  rPoint,
const c_matrix< double, ELEMENT_DIM, SPACE_DIM > &  rInverseJacobian,
c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &  rReturnValue 
) [inline, protected]

Compute the derivatives of all basis functions at a point within an element. This method will transform the results, for use within gaussian quadrature for example.

This is almost identical to LinearBasisFunction::ComputeTransformedBasisFunctionDerivatives, except that it is also templated over SPACE_DIM and can handle cases such as 1d in 3d space.

Todo:
Template LinearBasisFunction over SPACE_DIM?
Parameters:
rPoint The point at which to compute the basis functions. The results are undefined if this is not within the canonical element.
rInverseJacobian The inverse of the Jacobian matrix mapping the real element into the canonical element.
rReturnValue A reference to a vector, to be filled in
Returns:
The derivatives of the basis functions, in local index order. Each entry is a vector (c_vector<double, SPACE_DIM> instance) giving the derivative along each axis.

Definition at line 291 of file AbstractStaticAssembler.hpp.

References LinearBasisFunction< ELEM_DIM >::ComputeBasisFunctionDerivatives().

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::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 
) [inline, protected, 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

Todo:
This assumes that the Jacobian is constant on an element. This is true for linear basis functions, but not for any other type of basis function.

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 305 of file AbstractStaticAssembler.hpp.

References LinearBasisFunction< ELEM_DIM >::ComputeBasisFunctions(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeMatrixTerm(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::ComputeTransformedBasisFunctionDerivatives(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeVectorTerm(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::GetCurrentSolutionOrGuessValue(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), GaussianQuadratureRule< ELEM_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEM_DIM >::GetWeight(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::IncrementInterpolatedQuantities(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mCurrentSolutionOrGuessReplicated, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpMesh, AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ProblemIsNonlinear(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ResetInterpolatedQuantities(), ChastePoint< DIM >::rGetLocation(), Node< SPACE_DIM >::rGetLocation(), GaussianQuadratureRule< ELEM_DIM >::rGetQuadPoint(), and ReplicatableVector::size().

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleOnSurfaceElement ( const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &  rSurfaceElement,
c_vector< double, PROBLEM_DIM *ELEMENT_DIM > &  rBSurfElem 
) [inline, protected, 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.

Todo:
: add interpolation of u as well

Todo:
Improve efficiency of Neumann BC implementation.

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 438 of file AbstractStaticAssembler.hpp.

References LinearBasisFunction< ELEM_DIM >::ComputeBasisFunctions(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeVectorSurfaceTerm(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::IncrementInterpolatedQuantities(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpMesh, AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ResetInterpolatedQuantities(), and ChastePoint< DIM >::rGetLocation().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleSystem ( bool  assembleVector,
bool  assembleMatrix,
Vec  currentSolutionOrGuess = NULL,
double  currentTime = 0.0 
) [inline, protected, 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.

Called by Solve() Calls AssembleOnElement()

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.
currentTime The current time for dynamic problems. Not used in static problems.

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Reimplemented in AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >, and AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >.

Definition at line 493 of file AbstractStaticAssembler.hpp.

References LinearSystem::AddLhsMultipleValues(), LinearSystem::AddRhsMultipleValues(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyNeummanBoundaryConditions(), LinearSystem::AssembleFinalLhsMatrix(), LinearSystem::AssembleIntermediateLhsMatrix(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleOnElement(), LinearSystem::AssembleRhsVector(), GenericEventHandler< 11, HeartEventHandler >::BeginEvent(), GenericEventHandler< 11, HeartEventHandler >::EndEvent(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::FinaliseAssembleSystem(), AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::FinaliseLinearSystem(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), LinearSystem::GetSize(), AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::GetStiffnessMatrixGlobalIndices(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mCurrentSolutionOrGuessReplicated, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpLinearSystem, AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpMesh, AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::PrepareForAssembleSystem(), ReplicatableVector::ReplicatePetscVector(), LinearSystem::rGetLhsMatrix(), LinearSystem::rGetRhsVector(), ReplicatableVector::size(), LinearSystem::ZeroLhsMatrix(), and LinearSystem::ZeroRhsVector().

Referenced by AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::StaticSolve().

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
LinearSystem ** AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::GetLinearSystem (  )  [inline, virtual]

Accessor method that subclasses of AbstractAssembler (but not us) can use to get to useful data.

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 646 of file AbstractStaticAssembler.hpp.

References AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpLinearSystem.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
ReplicatableVector & AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::rGetCurrentSolutionOrGuess (  )  [inline, protected, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
double AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::GetCurrentSolutionOrGuessValue ( unsigned  nodeIndex,
unsigned  indexOfUnknown 
) [inline, protected, virtual]

Get the value of the current solution (or guess) vector at the given node

Parameters:
nodeIndex 
indexOfUnknown 

Definition at line 660 of file AbstractStaticAssembler.hpp.

References AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mCurrentSolutionOrGuessReplicated.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::SetNumberOfQuadraturePointsPerDimension ( unsigned  numQuadPoints  )  [inline]

Set the number of quadrature points to use, per dimension.

This method will throw an exception if the requested number of quadrature points is not supported. (///

Todo:
: There may be a small memory leak if this occurs.)
Parameters:
numQuadPoints Number of quadrature points to use per dimension (defaults to 2)

Definition at line 683 of file AbstractStaticAssembler.hpp.

References AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpQuadRule, and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpSurfaceQuadRule.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::SetMesh ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh  )  [inline]

Set the mesh.

Parameters:
pMesh Pointer to a mesh

Definition at line 693 of file AbstractStaticAssembler.hpp.

References AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpMesh.


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpMesh [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
GaussianQuadratureRule<ELEMENT_DIM>* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpQuadRule [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
GaussianQuadratureRule<ELEMENT_DIM-1>* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpSurfaceQuadRule [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
ReplicatableVector AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mCurrentSolutionOrGuessReplicated [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
LinearSystem* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpLinearSystem [protected]


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

Generated on Tue Aug 4 16:10:46 2009 for Chaste by  doxygen 1.5.5