AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL > Class Template Reference

#include <AbstractFeObjectAssembler.hpp>

Inheritance diagram for AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >:

Inheritance graph
[legend]
Collaboration diagram for AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >:

Collaboration graph
[legend]

List of all members.

Public Member Functions

 AbstractFeObjectAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, unsigned numQuadPoints=2)
void SetApplyNeummanBoundaryConditionsToVector (BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pBcc)
void OnlyAssembleOnSurfaceElements (bool onlyAssembleOnSurfaceElements=true)
void SetMatrixToAssemble (Mat &rMatToAssemble)
void SetVectorToAssemble (Vec &rVecToAssemble, bool zeroVectorBeforeAssembly)
void SetCurrentSolution (Vec currentSolution)
void Assemble ()
void AssembleMatrix ()
void AssembleVector ()
virtual ~AbstractFeObjectAssembler ()

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)
void DoAssemble ()
void ApplyNeummanBoundaryConditionsToVector ()
virtual double GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown)
template<size_t SMALL_MATRIX_SIZE>
void AddMultipleValuesToMatrix (unsigned *matrixRowAndColIndices, c_matrix< double, SMALL_MATRIX_SIZE, SMALL_MATRIX_SIZE > &rSmallMatrix)
template<size_t SMALL_VECTOR_SIZE>
void AddMultipleValuesToVector (unsigned *vectorIndices, c_vector< double, SMALL_VECTOR_SIZE > &smallVector)
virtual void ResetInterpolatedQuantities ()
virtual void IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode)
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)
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)
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)
virtual bool ElementAssemblyCriterion (Element< ELEMENT_DIM, SPACE_DIM > &rElement)
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)
virtual void AssembleOnSurfaceElement (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &rSurfaceElement, c_vector< double, PROBLEM_DIM *ELEMENT_DIM > &rBSurfElem)

Protected Attributes

Vec mVectorToAssemble
Mat mMatrixToAssemble
bool mAssembleMatrix
bool mAssembleVector
bool mZeroVectorBeforeAssembly
bool mApplyNeummanBoundaryConditionsToVector
bool mOnlyAssembleOnSurfaceElements
PetscInt mOwnershipRangeLo
PetscInt mOwnershipRangeHi
GaussianQuadratureRule
< ELEMENT_DIM > * 
mpQuadRule
GaussianQuadratureRule
< ELEMENT_DIM-1 > * 
mpSurfaceQuadRule
ReplicatableVector mCurrentSolutionOrGuessReplicated
AbstractTetrahedralMesh
< ELEMENT_DIM, SPACE_DIM > * 
mpMesh
BoundaryConditionsContainer
< ELEMENT_DIM, SPACE_DIM,
PROBLEM_DIM > * 
mpBoundaryConditions


Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
class AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >

An abstract class for creating finite element vectors or matrices that are defined by integrals over the computational domain of functions of basis functions (for example, stiffness or mass matrices), that require assembly by looping over each element in the mesh and computing the elementwise integrals and adding it to the full matrix or vector.

This class can be used to assemble a matrix OR a vector OR one of each. The template booleans CAN_ASSEMBLE_VECTOR and CAN_ASSEMBLE_MATRIX should be chosen accordingly.

The class provides the functionality to loop over elements, perform element-wise integration (using Gaussian quadrature and linear basis functions), and add the results to the final matrix or vector. The concrete class which inherits from this must implement either COMPUTE_MATRIX_TERM or COMPUTE_VECTOR_TERM or both, which should return the INTERGRAND, as a function of the basis functions.

Adding integrals of surface elements can also be done (but only to vectors).

The final template parameter defines how much interpolation (onto quadrature points) is required by the concrete class.

CARDIAC: only interpolates the first component of the unknown (ie the voltage) NORMAL: interpolates the position X and all components of the unknown u NONLINEAR: interpolates X, u and grad(u). Also computes the gradient of the basis functions when assembling vectors.

Definition at line 89 of file AbstractFeObjectAssembler.hpp.


Member Typedef Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
typedef LinearBasisFunction<ELEMENT_DIM> AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::BasisFunction [protected]

Basis function for use with normal elements

Definition at line 129 of file AbstractFeObjectAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
typedef LinearBasisFunction<ELEMENT_DIM-1> AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SurfaceBasisFunction [protected]

Basis function for use with boundary elements

Definition at line 131 of file AbstractFeObjectAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AbstractFeObjectAssembler ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
unsigned  numQuadPoints = 2 
) [inline]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::~AbstractFeObjectAssembler (  )  [inline, virtual]

Destructor

Definition at line 465 of file AbstractFeObjectAssembler.hpp.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::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:
#1319 Template LinearBasisFunction over SPACE_DIM and remove this method?
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 690 of file AbstractFeObjectAssembler.hpp.

References LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctionDerivatives().

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble (  )  [inline, protected]

The main assembly method. Private, should only be called through Assemble(), AssembleMatrix() or AssembleVector() which set mAssembleMatrix, mAssembleVector accordingly

Definition at line 546 of file AbstractFeObjectAssembler.hpp.

References AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ApplyNeummanBoundaryConditionsToVector(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement(), GenericEventHandler< 13, HeartEventHandler >::BeginEvent(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ElementAssemblyCriterion(), GenericEventHandler< 13, HeartEventHandler >::EndEvent(), EXCEPTION, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::GetStiffnessMatrixGlobalIndices(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mApplyNeummanBoundaryConditionsToVector, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mAssembleMatrix, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mAssembleVector, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mMatrixToAssemble, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOnlyAssembleOnSurfaceElements, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpMesh, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mVectorToAssemble, and AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mZeroVectorBeforeAssembly.

Referenced by AbstractFeObjectAssembler< DIM, DIM, 2, false, true, CARDIAC >::Assemble(), AbstractFeObjectAssembler< DIM, DIM, 2, false, true, CARDIAC >::AssembleMatrix(), and AbstractFeObjectAssembler< DIM, DIM, 2, false, true, CARDIAC >::AssembleVector().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ApplyNeummanBoundaryConditionsToVector (  )  [inline, protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual double AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::GetCurrentSolutionOrGuessValue ( unsigned  nodeIndex,
unsigned  indexOfUnknown 
) [inline, protected, virtual]

Useful inline function for getting an entry from the current solution vector

Parameters:
nodeIndex node index
indexOfUnknown index of unknown

Definition at line 185 of file AbstractFeObjectAssembler.hpp.

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
template<size_t SMALL_MATRIX_SIZE>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AddMultipleValuesToMatrix ( unsigned *  matrixRowAndColIndices,
c_matrix< double, SMALL_MATRIX_SIZE, SMALL_MATRIX_SIZE > &  rSmallMatrix 
) [inline, protected]

Add multiple values to the matrix. Templated over the size of the small matrix.

Parameters:
matrixRowAndColIndices mapping from index of the ublas matrix (see param below) to index of the matrix
rSmallMatrix Ublas matrix containing the values to be added
N.B. Values which are not local (ie the row is not owned) will be skipped.

Definition at line 886 of file AbstractFeObjectAssembler.hpp.

References AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mMatrixToAssemble, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOwnershipRangeHi, and AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOwnershipRangeLo.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
template<size_t SMALL_VECTOR_SIZE>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AddMultipleValuesToVector ( unsigned *  vectorIndices,
c_vector< double, SMALL_VECTOR_SIZE > &  smallVector 
) [inline, protected]

Add multiple values to the vector. Templated over the size of the small vector.

Parameters:
vectorIndices mapping from index of the ublas vector (see param below) to index of the vector
smallVector Ublas vector containing the values to be added
N.B. Values which are not local (ie the row is not owned) will be skipped.

Definition at line 944 of file AbstractFeObjectAssembler.hpp.

References AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOwnershipRangeHi, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOwnershipRangeLo, and AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mVectorToAssemble.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ResetInterpolatedQuantities ( void   )  [inline, protected, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::IncrementInterpolatedQuantities ( double  phiI,
const Node< SPACE_DIM > *  pNode 
) [inline, protected, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual c_matrix<double,PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1)> AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::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 
) [inline, protected, virtual]

This method returns the matrix to be added to element stiffness matrix for a given gauss point, ie, essentially the INTERGRAND in the integral definition of the matrix. 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 if CAN_ASSEMBLE_MATRIX is true **

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.

Definition at line 265 of file AbstractFeObjectAssembler.hpp.

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual c_vector<double,PROBLEM_DIM*(ELEMENT_DIM+1)> AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::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 
) [inline, protected, virtual]

This method returns the vector to be added to element stiffness vector for a given gauss point, ie, essentially the INTERGRAND in the integral definition of the vector. 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 if CAN_ASSEMBLE_VECTOR is true **

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

Definition at line 298 of file AbstractFeObjectAssembler.hpp.

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual c_vector<double, PROBLEM_DIM*ELEMENT_DIM> AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeVectorSurfaceTerm ( const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &  rSurfaceElement,
c_vector< double, ELEMENT_DIM > &  rPhi,
ChastePoint< SPACE_DIM > &  rX 
) [inline, protected, virtual]

This method returns the vector to be added to element stiffness vector for a given gauss point in BoundaryElement, ie, essentially the INTERGRAND in the boundary integral part of the definition of the vector. 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 if CAN_ASSEMBLE_VECTOR is true and SetApplyNeumannBoundaryConditionsToVector will be called.**

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

Reimplemented in BidomainAssembler< ELEMENT_DIM, SPACE_DIM >, MonodomainAssembler< ELEMENT_DIM, SPACE_DIM >, SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >, SimpleLinearParabolicSolver< ELEMENT_DIM, SPACE_DIM >, SimpleNonlinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >, and SimpleLinearEllipticSolver< DIM, DIM >.

Definition at line 327 of file AbstractFeObjectAssembler.hpp.

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnSurfaceElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual bool AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ElementAssemblyCriterion ( Element< ELEMENT_DIM, SPACE_DIM > &  rElement  )  [inline, protected, virtual]

Whether to include this (volume) element when assembling. Returns true here but can be overridden by the concrete assembler if not all elements should be included

Parameters:
rElement the element

Definition at line 342 of file AbstractFeObjectAssembler.hpp.

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::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 
) [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.
Called by AssembleSystem() Calls ComputeMatrixTerm() etc

Todo:
#1320 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.

Definition at line 704 of file AbstractFeObjectAssembler.hpp.

References LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeMatrixTerm(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeTransformedBasisFunctionDerivatives(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeVectorTerm(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::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(), ReplicatableVector::GetSize(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::IncrementInterpolatedQuantities(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mAssembleMatrix, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mAssembleVector, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mCurrentSolutionOrGuessReplicated, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpMesh, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpQuadRule, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ResetInterpolatedQuantities(), ChastePoint< DIM >::rGetLocation(), and Node< SPACE_DIM >::rGetLocation().

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::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:
#1321 Improve efficiency of Neumann BC implementation.

Definition at line 829 of file AbstractFeObjectAssembler.hpp.

References LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeVectorSurfaceTerm(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEMENT_DIM >::GetWeight(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::IncrementInterpolatedQuantities(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpMesh, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpSurfaceQuadRule, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ResetInterpolatedQuantities(), ChastePoint< DIM >::rGetLocation(), and GaussianQuadratureRule< ELEMENT_DIM >::rGetQuadPoint().

Referenced by AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ApplyNeummanBoundaryConditionsToVector().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SetApplyNeummanBoundaryConditionsToVector ( BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *  pBcc  )  [inline]

Whether to apply surface integrals to the vector.

Parameters:
pBcc A boundary conditions container (boundary elements containing a Neumann bc will be looped over). Requires CAN_ASSEMBLE_VECTOR==true.

Definition at line 643 of file AbstractFeObjectAssembler.hpp.

References AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mApplyNeummanBoundaryConditionsToVector, and AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpBoundaryConditions.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::OnlyAssembleOnSurfaceElements ( bool  onlyAssembleOnSurfaceElements = true  )  [inline]

If SetApplyNeummanBoundaryConditionsToVector() has been called, can call this to only apply the boundary conditions, not to assemble over the volume elements. Obviously requires CAN_ASSEMBLE_VECTOR==true and ComputeVectorSurfaceTerm to be implemented

Parameters:
onlyAssembleOnSurfaceElements whether to only assemble on the surface elements (defaults to true)

Definition at line 402 of file AbstractFeObjectAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SetMatrixToAssemble ( Mat &  rMatToAssemble  )  [inline]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SetVectorToAssemble ( Vec &  rVecToAssemble,
bool  zeroVectorBeforeAssembly 
) [inline]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SetCurrentSolution ( Vec  currentSolution  )  [inline]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::Assemble (  )  [inline]

Assemble everything that the class can assemble

Definition at line 433 of file AbstractFeObjectAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleMatrix (  )  [inline]

Assemble the matrix. Requires CAN_ASSEMBLE_MATRIX==true and ComputeMatrixTerm() to be implemented

Definition at line 443 of file AbstractFeObjectAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleVector (  )  [inline]

Assemble the vector. Requires CAN_ASSEMBLE_VECTOR==true and ComputeVectorTerm() to be implemented

Definition at line 454 of file AbstractFeObjectAssembler.hpp.


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
Vec AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mVectorToAssemble [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
Mat AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mMatrixToAssemble [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
bool AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mAssembleMatrix [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
bool AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mAssembleVector [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
bool AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mZeroVectorBeforeAssembly [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
bool AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mApplyNeummanBoundaryConditionsToVector [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
bool AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOnlyAssembleOnSurfaceElements [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
PetscInt AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOwnershipRangeLo [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
PetscInt AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mOwnershipRangeHi [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
GaussianQuadratureRule<ELEMENT_DIM>* AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpQuadRule [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
GaussianQuadratureRule<ELEMENT_DIM-1>* AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpSurfaceQuadRule [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
ReplicatableVector AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mCurrentSolutionOrGuessReplicated [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpMesh [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpBoundaryConditions [protected]


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

Generated on Mon Nov 1 12:35:43 2010 for Chaste by  doxygen 1.5.5