#include <AbstractNonlinearAssemblerSolverHybrid.hpp>
Inherits AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >.
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, unsigned numQuadPoints=2) | |
virtual | ~AbstractNonlinearAssemblerSolverHybrid () |
virtual Vec | Solve (Vec initialGuess, bool useAnalyticalJacobian=false) |
void | SetNonlinearSolver (AbstractNonlinearSolver *pNonlinearSolver) |
bool | VerifyJacobian (double tol) |
Protected Member Functions | |
void | ApplyDirichletConditions (Vec currentGuess, Vec residual, Mat *pMat) |
void | ComputeJacobianNumerically (const Vec currentGuess, Mat *pJacobian) |
Protected Attributes | |
AbstractTetrahedralMesh < ELEMENT_DIM, SPACE_DIM > * | mpMesh |
BoundaryConditionsContainer < ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * | mpBoundaryConditions |
AbstractNonlinearSolver * | mpSolver |
bool | mWeAllocatedSolverMemory |
bool | mUseAnalyticalJacobian |
NaturalNeumannSurfaceTermAssembler < ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * | mpNeumannSurfaceTermsAssembler |
The ASSEMBLER-SOLVER class.
An abstract solver for solving nonlinear elliptic PDEs.
Definition at line 98 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid | ( | AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * | pMesh, | |
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * | pBoundaryConditions, | |||
unsigned | numQuadPoints = 2 | |||
) | [inline] |
Constructor.
pMesh | The mesh | |
pBoundaryConditions | The boundary conditions to apply | |
numQuadPoints | number of quadrature points (defaults to 2) |
Definition at line 219 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannSurfaceTermsAssembler, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory.
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractNonlinearAssemblerSolverHybrid | ( | ) | [inline, virtual] |
Destructor.
Definition at line 245 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannSurfaceTermsAssembler, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory.
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions | ( | Vec | currentGuess, | |
Vec | residual, | |||
Mat * | pMat | |||
) | [inline, protected] |
Apply Dirichlet boundary conditions to either the residual or Jacobian.
currentGuess | the solution guess for the current nonlinear solver iteration | |
residual | the residual to apply boundary conditions to (can be NULL if bcs to be applied to Jacobian only) | |
pMat | the Jacobian to apply boundary conditions to (can be NULL if bcs to be applied to residual only) |
Definition at line 255 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual().
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian | ( | const Vec | currentGuess, | |
Mat * | pJacobian | |||
) | [inline] |
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()).
currentGuess | The solution guess for the current iteration. | |
pJacobian | Pointer 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 291 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleMatrix(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically(), PetscMatTools::Finalise(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian, AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SetCurrentSolution(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetMatrixToAssemble(), and PetscMatTools::SwitchWriteMode().
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::VerifyJacobian().
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically | ( | const Vec | currentGuess, | |
Mat * | pJacobian | |||
) | [inline, protected] |
Computes the Jacobian numerically i.e. an approximation, using numerical partial derivatives.
currentGuess | Independent variable, u in f(u), for example | |
pJacobian | A pointer to the Jacobian matrix |
Definition at line 312 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References PetscVecTools::AddToElement(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual(), PetscTools::CreateVec(), PetscMatTools::Finalise(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, PetscVecTools::Scale(), PetscMatTools::SetElement(), and PetscVecTools::WAXPY().
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian().
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual | ( | const Vec | currentGuess, | |
Vec | residualVector | |||
) | [inline] |
Compute the residual vector given the current solution guess.
currentGuess | The solution guess for the current nonlinear solve iteration. | |
residualVector | We fill this with the residual vector. |
NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.
Definition at line 270 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions(), AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleVector(), PetscVecTools::Finalise(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannSurfaceTermsAssembler, AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::SetCurrentSolution(), and AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::SetVectorToAssemble().
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically().
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetNonlinearSolver | ( | AbstractNonlinearSolver * | pNonlinearSolver | ) | [inline] |
SetNonlinearSolver - by default a SimplePetscNonlinearSolver is created and used in this class, this method can be called to use a different AbstractNonlinearSolver.
pNonlinearSolver | a nonlinear solver |
Definition at line 396 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory.
Vec AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve | ( | Vec | initialGuess, | |
bool | useAnalyticalJacobian = false | |||
) | [inline, virtual] |
Assemble and solve the system for a nonlinear elliptic PDE.
initialGuess | An initial guess for the iterative solver | |
useAnalyticalJacobian | Set to true to use an analytically calculated Jacobian matrix rather than a numerically approximated one. |
Definition at line 372 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References EXCEPTION, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian, and AbstractNonlinearSolver::Solve().
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::VerifyJacobian | ( | double | tol | ) | [inline] |
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.
tol | A tolerance which defaults to 1e-5 |
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 407 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References PetscMatTools::CheckEquality(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian(), PetscTools::CreateAndSetVec(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian, and PetscTools::SetupMat().
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions [protected] |
Boundary conditions container.
Definition at line 106 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions().
AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh [protected] |
The mesh.
Reimplemented from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >.
Definition at line 103 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::VerifyJacobian().
NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannSurfaceTermsAssembler [protected] |
An assembler for the surface integral terms in the residual
Definition at line 118 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractNonlinearAssemblerSolverHybrid().
AbstractNonlinearSolver* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver [protected] |
The nonlinear solver.
Definition at line 109 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetNonlinearSolver(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractNonlinearAssemblerSolverHybrid().
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian [protected] |
Whether to use an analytical expression for the Jacobian.
Definition at line 115 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::VerifyJacobian().
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory [protected] |
Whether memory has been allocated for the solver.
Definition at line 112 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetNonlinearSolver(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractNonlinearAssemblerSolverHybrid().