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

#include <AbstractNonlinearAssemblerSolverHybrid.hpp>

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

Public Member Functions

void ComputeResidual (const Vec currentGuess, Vec residualVector)
 
void ComputeJacobian (const Vec currentGuess, Mat *pJacobian)
 
 AbstractNonlinearAssemblerSolverHybrid (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pBoundaryConditions)
 
virtual ~AbstractNonlinearAssemblerSolverHybrid ()
 
virtual Vec Solve (Vec initialGuess, bool useAnalyticalJacobian=false)
 
void SetNonlinearSolver (AbstractNonlinearSolver *pNonlinearSolver)
 
bool VerifyJacobian (double tol)
 
- Public Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >
 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 Member Functions

void ApplyDirichletConditions (Vec currentGuess, Vec residual, Mat *pMat)
 
void ComputeJacobianNumerically (const Vec currentGuess, Mat *pJacobian)
 
- Protected Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >
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
 
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * mpBoundaryConditions
 
AbstractNonlinearSolvermpSolver
 
bool mWeAllocatedSolverMemory
 
bool mUseAnalyticalJacobian
 
NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * mpNeumannSurfaceTermsAssembler
 
- Protected Attributes inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >
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, NONLINEAR >
typedef LinearBasisFunction< ELEMENT_DIM > BasisFunction
 

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
class AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >

The ASSEMBLER-SOLVER class.

An abstract solver for solving nonlinear elliptic PDEs.

Definition at line 132 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

Constructor & Destructor Documentation

◆ AbstractNonlinearAssemblerSolverHybrid()

◆ ~AbstractNonlinearAssemblerSolverHybrid()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractNonlinearAssemblerSolverHybrid ( )
virtual

Destructor.

Definition at line 275 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

Member Function Documentation

◆ ApplyDirichletConditions()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions ( Vec  currentGuess,
Vec  residual,
Mat pMat 
)
protected

Apply Dirichlet boundary conditions to either the residual or Jacobian.

Parameters
currentGuessthe solution guess for the current nonlinear solver iteration
residualthe residual to apply boundary conditions to (can be NULL if bcs to be applied to Jacobian only)
pMatthe Jacobian to apply boundary conditions to (can be NULL if bcs to be applied to residual only)

Definition at line 285 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

◆ ComputeJacobian()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian ( const Vec  currentGuess,
Mat pJacobian 
)

Compute the Jacobian matrix given a current guess at the solution. Choose whether to use a numerical or analytical method based on a flag provided by the user (in Solve()).

Parameters
currentGuessThe solution guess for the current iteration.
pJacobianPointer to object to fill with the Jacobian matrix.

NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.

Definition at line 321 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscMatTools::Finalise(), and PetscMatTools::SwitchWriteMode().

◆ ComputeJacobianNumerically()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically ( const Vec  currentGuess,
Mat pJacobian 
)
protected

Computes the Jacobian numerically i.e. an approximation, using numerical partial derivatives.

Parameters
currentGuessIndependent variable, u in f(u), for example
pJacobianA pointer to the Jacobian matrix
Todo:
This loop is setting the column "global_index_outer" of the Jacobian matrix to the result vector in a non-sparse way.

Definition at line 342 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscVecTools::AddToElement(), PetscTools::CreateVec(), PetscTools::Destroy(), PetscMatTools::Finalise(), CompareDoubles::IsNearZero(), PetscVecTools::Scale(), PetscMatTools::SetElement(), and PetscVecTools::WAXPY().

◆ ComputeResidual()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual ( const Vec  currentGuess,
Vec  residualVector 
)

Compute the residual vector given the current solution guess.

Parameters
currentGuessThe solution guess for the current nonlinear solve iteration.
residualVectorWe fill this with the residual vector.

NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.

Definition at line 300 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscVecTools::Finalise().

◆ SetNonlinearSolver()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetNonlinearSolver ( AbstractNonlinearSolver pNonlinearSolver)

SetNonlinearSolver - by default a SimplePetscNonlinearSolver is created and used in this class, this method can be called to use a different AbstractNonlinearSolver.

Parameters
pNonlinearSolvera nonlinear solver

Definition at line 432 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

◆ Solve()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
Vec AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve ( Vec  initialGuess,
bool  useAnalyticalJacobian = false 
)
virtual

Assemble and solve the system for a nonlinear elliptic PDE.

Parameters
initialGuessAn initial guess for the iterative solver
useAnalyticalJacobianSet to true to use an analytically calculated Jacobian matrix rather than a numerically approximated one.
Returns
A PETSc vector giving the solution at each mesh node.

Definition at line 408 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References EXCEPTION.

◆ VerifyJacobian()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::VerifyJacobian ( double  tol)

A helper method for use when writing concrete assemblers. Once the user has calculated (on paper) the weak form and the form of the ComputeMatrixTerm method, they can check whether the analytic Jacobian matches the numerical Jacobian to verify the correctness of the code.

Parameters
tolA tolerance which defaults to 1e-5
Returns
true if the componentwise difference between the matrices is less than the tolerance, false otherwise.

This method should NOT be run in simulations - it is only to verify the correctness of the concrete assembler code. Allocates dense matrices.

Definition at line 443 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscMatTools::CheckEquality(), PetscTools::CreateAndSetVec(), PetscTools::Destroy(), and PetscTools::SetupMat().

Member Data Documentation

◆ mpBoundaryConditions

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions
protected

Boundary conditions container.

Definition at line 140 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

◆ mpMesh

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh
protected

◆ mpNeumannSurfaceTermsAssembler

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannSurfaceTermsAssembler
protected

◆ mpSolver

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractNonlinearSolver* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver
protected

◆ mUseAnalyticalJacobian

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian
protected

Whether to use an analytical expression for the Jacobian.

Definition at line 149 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

◆ mWeAllocatedSolverMemory

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory
protected

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