Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL > Class Template Reference

#include <AbstractFeVolumeIntegralAssembler.hpp>

+ Inheritance diagram for AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >:
+ Collaboration diagram for AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >:

Public Member Functions

 AbstractFeVolumeIntegralAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
virtual ~AbstractFeVolumeIntegralAssembler ()
 
- Public Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
 AbstractFeAssemblerCommon ()
 
void SetCurrentSolution (Vec currentSolution)
 
virtual ~AbstractFeAssemblerCommon ()
 
- Public Member Functions inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
 AbstractFeAssemblerInterface ()
 
void SetMatrixToAssemble (Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
 
void SetVectorToAssemble (Vec &rVecToAssemble, bool zeroVectorBeforeAssembly)
 
void Assemble ()
 
void AssembleMatrix ()
 
void AssembleVector ()
 
virtual ~AbstractFeAssemblerInterface ()
 

Protected Types

typedef LinearBasisFunction< ELEMENT_DIM > BasisFunction
 

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 ()
 
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 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 bool ElementAssemblyCriterion (Element< ELEMENT_DIM, SPACE_DIM > &rElement)
 
- Protected Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
virtual double GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown)
 
virtual void ResetInterpolatedQuantities ()
 
virtual void IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode)
 
virtual void IncrementInterpolatedGradientQuantities (const c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, unsigned phiIndex, const Node< SPACE_DIM > *pNode)
 

Protected Attributes

AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
 
GaussianQuadratureRule< ELEMENT_DIM > * mpQuadRule
 
- Protected Attributes inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
ReplicatableVector mCurrentSolutionOrGuessReplicated
 
- Protected Attributes inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
Vec mVectorToAssemble
 
Mat mMatrixToAssemble
 
bool mAssembleMatrix
 
bool mAssembleVector
 
bool mZeroMatrixBeforeAssembly
 
bool mZeroVectorBeforeAssembly
 
PetscInt mOwnershipRangeLo
 
PetscInt mOwnershipRangeHi
 

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
class AbstractFeVolumeIntegralAssembler< 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 element-wise integrals and adding it to the full matrix or vector.

This class is used for VOLUME integrals. For surface integrals there is a similar class, AbstractFeSurfaceIntegralAssembler.

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 INTEGRAND, as a function of the basis functions.

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.

This class inherits from AbstractFeAssemblerCommon which is where some member variables (the matrix/vector to be created, for example) are defined.

Definition at line 78 of file AbstractFeVolumeIntegralAssembler.hpp.

Member Typedef Documentation

◆ BasisFunction

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> AbstractFeVolumeIntegralAssembler< 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 89 of file AbstractFeVolumeIntegralAssembler.hpp.

Constructor & Destructor Documentation

◆ AbstractFeVolumeIntegralAssembler()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AbstractFeVolumeIntegralAssembler ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh)

◆ ~AbstractFeVolumeIntegralAssembler()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::~AbstractFeVolumeIntegralAssembler ( )
inlinevirtual

Member Function Documentation

◆ AssembleOnElement()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeVolumeIntegralAssembler< 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 
)
protectedvirtual

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

Parameters
rElementThe element to assemble on.
rAElemThe 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.
rBElemThe 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 345 of file AbstractFeVolumeIntegralAssembler.hpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEMENT_DIM >::GetWeight(), Node< SPACE_DIM >::rGetLocation(), and GaussianQuadratureRule< ELEMENT_DIM >::rGetQuadPoint().

◆ ComputeMatrixTerm()

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)> AbstractFeVolumeIntegralAssembler< 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 
)
inlineprotectedvirtual
Returns
the matrix to be added to element stiffness matrix for a given Gauss point, ie, essentially the INTEGRAND 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 determinant 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: When INTERPOLATION_LEVEL==NORMAL, rGradU does not get set up and should not be used.

Parameters
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases.
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i).
rXThe point in space.
rUThe unknown as a vector, u(i) = u_i.
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j).
pElementPointer to the element.

Reimplemented in LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 143 of file AbstractFeVolumeIntegralAssembler.hpp.

References NEVER_REACHED.

◆ ComputeTransformedBasisFunctionDerivatives()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeVolumeIntegralAssembler< 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 
)
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
rPointThe point at which to compute the basis functions. The results are undefined if this is not within the canonical element.
rInverseJacobianThe inverse of the Jacobian matrix mapping the real element into the canonical element.
rReturnValueA reference to a vector, to be filled in (returned) Returns (via rReturnValue) 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 332 of file AbstractFeVolumeIntegralAssembler.hpp.

References LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctionDerivatives().

◆ ComputeVectorTerm()

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)> AbstractFeVolumeIntegralAssembler< 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 
)
inlineprotectedvirtual
Returns
the vector to be added to element stiffness vector for a given Gauss point, ie, essentially the INTEGRAND 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 determinant 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: When INTERPOLATION_LEVEL==NORMAL, rGradPhi and rGradU do not get set up and should not be used.

Parameters
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rXThe point in space
rUThe unknown as a vector, u(i) = u_i
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElementPointer to the element

Reimplemented in LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 178 of file AbstractFeVolumeIntegralAssembler.hpp.

References NEVER_REACHED.

◆ DoAssemble()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble ( )
protectedvirtual

◆ ElementAssemblyCriterion()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual bool AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ElementAssemblyCriterion ( Element< ELEMENT_DIM, SPACE_DIM > &  rElement)
inlineprotectedvirtual
Returns
true if we should 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
rElementthe element

Definition at line 219 of file AbstractFeVolumeIntegralAssembler.hpp.

Member Data Documentation

◆ mpMesh

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>* AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpMesh
protected

Mesh to be solved on.

Definition at line 83 of file AbstractFeVolumeIntegralAssembler.hpp.

◆ mpQuadRule

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

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