Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <SimpleLinearEllipticSolver.hpp>

+ Inheritance diagram for SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >:
+ Collaboration diagram for SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >:

Public Member Functions

 SimpleLinearEllipticSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, AbstractLinearEllipticPde< ELEMENT_DIM, SPACE_DIM > *pPde, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBoundaryConditions)
 
void InitialiseForSolve (Vec initialSolution=nullptr)
 
- Public Member Functions inherited from AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, 1, NORMAL >
 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 ()
 
- Public Member Functions inherited from AbstractStaticLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, 1 >
 AbstractStaticLinearPdeSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
Vec Solve (Vec initialGuess=nullptr)
 
- Public Member Functions inherited from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
 AbstractLinearPdeSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
virtual ~AbstractLinearPdeSolver ()
 
virtual void PrepareForSetupLinearSystem (Vec currentSolution)
 
virtual void FinaliseLinearSystem (Vec currentSolution)
 
virtual void FollowingSolveLinearSystem (Vec currentSolution)
 
LinearSystemGetLinearSystem ()
 

Protected Member Functions

virtual c_matrix< double, 1 *(ELEMENT_DIM+1), 1 *(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, 1 > &rU, c_matrix< double, 1, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
 
virtual c_vector< double, 1 *(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, 1 > &rU, c_matrix< double, 1, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
 
void SetupLinearSystem (Vec currentSolution, bool computeMatrix)
 
- 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)
 

Protected Attributes

AbstractLinearEllipticPde< ELEMENT_DIM, SPACE_DIM > * mpEllipticPde
 
- Protected Attributes inherited from AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, 1, NORMAL >
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
 
- Protected Attributes inherited from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
LinearSystemmpLinearSystem
 
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
 

Additional Inherited Members

- Protected Types inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
typedef LinearBasisFunction< ELEMENT_DIM > BasisFunction
 

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
class SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >

SimpleLinearEllipticSolver.

Solver for solving AbstractLinearEllipticPdes.

Definition at line 49 of file SimpleLinearEllipticSolver.hpp.

Constructor & Destructor Documentation

◆ SimpleLinearEllipticSolver()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::SimpleLinearEllipticSolver ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
AbstractLinearEllipticPde< ELEMENT_DIM, SPACE_DIM > *  pPde,
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *  pBoundaryConditions 
)

Constructor.

Parameters
pMeshpointer to the mesh
pPdepointer to the PDE
pBoundaryConditionspointer to the boundary conditions

Definition at line 74 of file SimpleLinearEllipticSolver.cpp.

References SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde.

Member Function Documentation

◆ ComputeMatrixTerm()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_matrix< double, 1 *(ELEMENT_DIM+1), 1 *(ELEMENT_DIM+1)> SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm ( c_vector< double, ELEMENT_DIM+1 > &  rPhi,
c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &  rGradPhi,
ChastePoint< SPACE_DIM > &  rX,
c_vector< double, 1 > &  rU,
c_matrix< double, 1, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
)
protectedvirtual
Returns
the term to be added to the element stiffness matrix - see AbstractFeVolumeIntegralAssembler

grad_phi[row] . ( pde_diffusion_term * grad_phi[col])

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

Definition at line 39 of file SimpleLinearEllipticSolver.cpp.

◆ ComputeVectorTerm()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 1 *(ELEMENT_DIM+1)> SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::ComputeVectorTerm ( c_vector< double, ELEMENT_DIM+1 > &  rPhi,
c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &  rGradPhi,
ChastePoint< SPACE_DIM > &  rX,
c_vector< double, 1 > &  rU,
c_matrix< double, 1, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
)
protectedvirtual
Returns
the term arising from boundary conditions to be added to the element stiffness vector - see AbstractFeVolumeIntegralAssembler
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

Definition at line 62 of file SimpleLinearEllipticSolver.cpp.

◆ InitialiseForSolve()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve ( Vec  initialSolution = nullptr)
virtual

Overloaded InitaliseForSolve() which just calls the base class but also sets the matrix as symmetric and sets Conjugate Gradients as the solver

Parameters
initialSolutioninitialSolution (used in base class version of this method)

Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 85 of file SimpleLinearEllipticSolver.cpp.

References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve().

Referenced by CellBasedEllipticPdeSolver< DIM >::InitialiseForSolve().

◆ SetupLinearSystem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem ( Vec  currentSolution,
bool  computeMatrix 
)
inlineprotectedvirtual

Delegate to AbstractAssemblerSolverHybrid::SetupGivenLinearSystem.

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).

Implements AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 110 of file SimpleLinearEllipticSolver.hpp.

References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpLinearSystem, and AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, 1, NORMAL >::SetupGivenLinearSystem().

Member Data Documentation

◆ mpEllipticPde

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractLinearEllipticPde<ELEMENT_DIM,SPACE_DIM>* SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde
protected

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