#include <AbstractStaticAssembler.hpp>
Public Member Functions | |
AbstractStaticAssembler (unsigned numQuadPoints=2) | |
void | SetNumberOfQuadraturePointsPerDimension (unsigned numQuadPoints) |
void | SetMesh (AbstractMesh< ELEMENT_DIM, SPACE_DIM > *pMesh) |
virtual | ~AbstractStaticAssembler () |
Protected Types | |
typedef LinearBasisFunction < ELEMENT_DIM > | BasisFunction |
typedef LinearBasisFunction < ELEMENT_DIM-1 > | SurfaceBasisFunction |
Protected Member Functions | |
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 () |
LinearSystem ** | GetLinearSystem () |
ReplicatableVector & | rGetCurrentSolutionOrGuess () |
virtual double | GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown) |
Protected Attributes | |
AbstractMesh< ELEMENT_DIM, SPACE_DIM > * | mpMesh |
GaussianQuadratureRule < ELEMENT_DIM > * | mpQuadRule |
GaussianQuadratureRule < ELEMENT_DIM-1 > * | mpSurfaceQuadRule |
ReplicatableVector | mCurrentSolutionOrGuessReplicated |
LinearSystem * | mpLinearSystem |
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 101 of file AbstractStaticAssembler.hpp.
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 114 of file AbstractStaticAssembler.hpp.
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 116 of file AbstractStaticAssembler.hpp.
AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AbstractStaticAssembler | ( | unsigned | numQuadPoints = 2 |
) | [inline] |
Default constructor. Uses linear basis functions.
numQuadPoints | Number of quadrature points to use per dimension. |
Definition at line 598 of file AbstractStaticAssembler.hpp.
virtual AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::~AbstractStaticAssembler | ( | ) | [inline, virtual] |
Delete any memory allocated by this class.
Definition at line 642 of file AbstractStaticAssembler.hpp.
virtual 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.
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. | |
currentSolutionOrGuess | For the parabolic linear case, the solution at the current timestep. NULL for the static linear case. In the nonlinear case, the current guess. | |
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) |
Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 152 of file AbstractStaticAssembler.hpp.
Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleSystem().
virtual 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.
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. |
Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 295 of file AbstractStaticAssembler.hpp.
virtual 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()
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 378 of file AbstractStaticAssembler.hpp.
Referenced by AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::StaticSolve().
virtual void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::PrepareForSolve | ( | ) | [inline, protected, virtual] |
This method is called at the beginning of Solve(). Subclass assemblers can use it to check everything has been set up correctly
Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 533 of file AbstractStaticAssembler.hpp.
Referenced by AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::Solve().
LinearSystem** AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::GetLinearSystem | ( | ) | [inline, protected, 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 570 of file AbstractStaticAssembler.hpp.
ReplicatableVector& AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::rGetCurrentSolutionOrGuess | ( | ) | [inline, protected, 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 579 of file AbstractStaticAssembler.hpp.
virtual 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
Definition at line 587 of file AbstractStaticAssembler.hpp.
Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleOnElement().
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. (///
numQuadPoints | Number of quadrature points to use per dimension. |
Definition at line 621 of file AbstractStaticAssembler.hpp.
void AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::SetMesh | ( | AbstractMesh< ELEMENT_DIM, SPACE_DIM > * | pMesh | ) | [inline] |
Set the mesh.
Definition at line 633 of file AbstractStaticAssembler.hpp.
AbstractMesh<ELEMENT_DIM, SPACE_DIM>* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpMesh [protected] |
Mesh to be solved on
Definition at line 105 of file AbstractStaticAssembler.hpp.
Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AbstractStaticAssembler(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleOnElement(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleOnSurfaceElement(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleSystem(), AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::InitialiseForSolve(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::PrepareForSolve(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::SetMesh(), and AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::Solve().
GaussianQuadratureRule<ELEMENT_DIM>* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpQuadRule [protected] |
Quadrature rule for use on normal elements
Definition at line 108 of file AbstractStaticAssembler.hpp.
Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AbstractStaticAssembler(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::SetNumberOfQuadraturePointsPerDimension(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::~AbstractStaticAssembler().
GaussianQuadratureRule<ELEMENT_DIM-1>* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpSurfaceQuadRule [protected] |
Quadrature rule for use on boundary elements
Definition at line 110 of file AbstractStaticAssembler.hpp.
Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AbstractStaticAssembler(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::SetNumberOfQuadraturePointsPerDimension(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::~AbstractStaticAssembler().
ReplicatableVector AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mCurrentSolutionOrGuessReplicated [protected] |
The CURRENT SOLUTION as a replicated vector for linear dynamic problems. (Empty for a static problem). The CURRENT GUESS for nonlinear problems.
Definition at line 122 of file AbstractStaticAssembler.hpp.
Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleOnElement(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleSystem(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::GetCurrentSolutionOrGuessValue(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::rGetCurrentSolutionOrGuess().
LinearSystem* AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::mpLinearSystem [protected] |
The linear system that is assembled in linear pde problems. Not used in nonlinear problems
Definition at line 128 of file AbstractStaticAssembler.hpp.
Referenced by AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AbstractStaticAssembler(), AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::ApplyDirichletConditions(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::AssembleSystem(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::GetLinearSystem(), AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::InitialiseForSolve(), AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::StaticSolve(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, MonodomainRhsMatrixAssembler< DIM > >::~AbstractStaticAssembler().