#include <AbstractNonlinearAssemblerSolverHybrid.hpp>
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 |
Definition at line 125 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 247 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, 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 272 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References 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 281 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 >::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 335 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual(), PetscTools::CreateVec(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh.
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. |
Definition at line 297 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >::AssembleVector(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >::SetApplyNeummanBoundaryConditionsToVector(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >::SetCurrentSolution(), and AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >::SetVectorToAssemble().
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically().
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. |
Definition at line 311 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions(), AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >::AssembleMatrix(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian, AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >::SetCurrentSolution(), and AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >::SetMatrixToAssemble().
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::VerifyJacobian().
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 411 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().
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 441 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory.
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 |
Definition at line 453 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References 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().
AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh [protected] |
Mesh class
Reimplemented from AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >.
Definition at line 129 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().
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions [protected] |
Boundary conditions container
Reimplemented from AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >.
Definition at line 132 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual().
AbstractNonlinearSolver* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver [protected] |
The nonlinear solver.
Definition at line 135 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 >::mWeAllocatedSolverMemory [protected] |
Whether memory has been allocated for the solver.
Definition at line 138 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().
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian [protected] |
Whether to use an analytical expression for the Jacobian.
Definition at line 141 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().