Chaste  Release::2017.1
AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL > Class Template Reference

#include <AbstractAssemblerSolverHybrid.hpp>

+ Inheritance diagram for AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >:
+ Collaboration diagram for AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >:

Public Member Functions

 AbstractAssemblerSolverHybrid (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pBoundaryConditions)
 
virtual ~AbstractAssemblerSolverHybrid ()
 
void SetupGivenLinearSystem (Vec currentSolution, bool computeMatrix, LinearSystem *pLinearSystem)
 
- Public Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
 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 Attributes

NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > mNaturalNeumannSurfaceTermAssembler
 
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * mpBoundaryConditions
 
- Protected Attributes inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
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
 

Additional Inherited Members

- Protected Types inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
typedef LinearBasisFunction< ELEMENT_DIM > BasisFunction
 
- Protected Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
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)
 

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, InterpolationLevel INTERPOLATION_LEVEL>
class AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >

A class which inherits from AbstractFeVolumeIntegralAssembler and implements a method SetupGivenLinearSystem(), which sets up the given linear system using the assembler part of this class, which can be called by SetUpLinearSystem() on a concrete solver.

It assumes natural Neumann boundary conditions are needed and uses a NaturalNeumannSurfaceTermAssembler for this part of the vector.

See SimpleLinearEllipticSolver for an example of a concrete class

Definition at line 55 of file AbstractAssemblerSolverHybrid.hpp.

Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, InterpolationLevel INTERPOLATION_LEVEL>
AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::AbstractAssemblerSolverHybrid ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *  pBoundaryConditions 
)
inline

Constructor.

Parameters
pMeshpointer to the mesh
pBoundaryConditionspointer to the boundary conditions.

Definition at line 78 of file AbstractAssemblerSolverHybrid.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, InterpolationLevel INTERPOLATION_LEVEL>
virtual AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::~AbstractAssemblerSolverHybrid ( )
inlinevirtual

Destructor.

Definition at line 91 of file AbstractAssemblerSolverHybrid.hpp.

Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::SetupGivenLinearSystem ( Vec  currentSolution,
bool  computeMatrix,
LinearSystem pLinearSystem 
)

Implementation of AbstractLinearPdeSolver::SetupLinearSystem, using the assembler that this class also inherits from. Concrete classes inheriting from both this class and AbstractLinearPdeSolver can then have a one-line implementation of AbstractLinearPdeSolver::SetupLinearSystem which calls this method.

Parameters
currentSolutionThe current solution which can be used in setting up the linear system if needed (NULL if there isn't a current solution)
computeMatrixWhether to compute the LHS matrix of the linear system (mainly for dynamic solves)
pLinearSystemThe linear system to set up.

Definition at line 111 of file AbstractAssemblerSolverHybrid.hpp.

References AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::Assemble(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleVector(), LinearSystem::FinaliseLhsMatrix(), LinearSystem::FinaliseRhsVector(), AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::mNaturalNeumannSurfaceTermAssembler, AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::mpBoundaryConditions, LinearSystem::rGetLhsMatrix(), LinearSystem::rGetRhsVector(), AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SetCurrentSolution(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetMatrixToAssemble(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetVectorToAssemble(), and LinearSystem::SwitchWriteModeLhsMatrix().

Referenced by AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, 1, NORMAL >::~AbstractAssemblerSolverHybrid().

Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, InterpolationLevel INTERPOLATION_LEVEL>
NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM> AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::mNaturalNeumannSurfaceTermAssembler
protected

An assembler for Neumann surface integrals, which are assumed to arise from natural Neumann boundary conditions, ie such that this surface integral is (for a 1-unknown problem) integral(g phi_i dS), where g is the Neumann boundary condition function

Definition at line 64 of file AbstractAssemblerSolverHybrid.hpp.

Referenced by AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::SetupGivenLinearSystem().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, InterpolationLevel INTERPOLATION_LEVEL>
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::mpBoundaryConditions
protected

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