Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/petsc_snes.h>
Public Types | |
using | real_type = PetscReal |
Public Member Functions | |
NonlinearSolver (const NonlinearSolverData &data=NonlinearSolverData(), const MPI_Comm mpi_comm=PETSC_COMM_WORLD) | |
~NonlinearSolver () | |
operator SNES () const | |
SNES | petsc_snes () |
MPI_Comm | get_mpi_communicator () const |
void | reinit () |
void | reinit (const NonlinearSolverData &data) |
void | set_matrix (PMatrixType &P) |
void | set_matrices (AMatrixType &A, PMatrixType &P) |
unsigned int | solve (VectorType &x) |
unsigned int | solve (VectorType &x, PMatrixType &P) |
unsigned int | solve (VectorType &x, AMatrixType &A, PMatrixType &P) |
Public Attributes | |
std::function< void(const VectorType &x, VectorType &res)> | residual |
std::function< void(const VectorType &x, AMatrixType &A, PMatrixType &P)> | jacobian |
std::function< void(const VectorType &x, const unsigned int step_number, const real_type f_norm)> | monitor |
std::function< void(const VectorType &x)> | setup_jacobian |
std::function< void(const VectorType &src, VectorType &dst)> | solve_with_jacobian |
std::function< void(const VectorType &x, real_type &energy_value)> | energy |
Private Attributes | |
SNES | snes |
SmartPointer< AMatrixType, NonlinearSolver > | A |
SmartPointer< PMatrixType, NonlinearSolver > | P |
bool | need_dummy_assemble |
std::exception_ptr | pending_exception |
Interface to PETSc SNES solver for nonlinear equations. The SNES solver is described in the [PETSc manual](https://petsc.org/release/manual/snes/). This class solves the nonlinear system of algebraic equations @f$F(x) = 0@f$. The interface to PETSc is realized by means of std::function callbacks like in the TrilinosWrappers::NOXSolver and SUNDIALS::KINSOL classes. NonlinearSolver supports any vector and matrix type having constructors and methods: @code class VectorType : public Subscriptor ... explicit VectorType(Vec); ... Vec & petsc_vector(); ... @endcode @code class MatrixType : public Subscriptor ... explicit MatrixType(Mat); ... Mat & petsc_matrix(); ... @endcode In particular, the supported types are the ones that can *wrap* PETSc's native vector and matrix classes, that are able to modify them in place, and that can return PETSc native types when requested. To use the solvers the user needs to provide the implementation of @f$F@f$ via the NonlinearSolver::residual callback. The default linearization procedure of a solver instantiated with this class consists in using Jacobian-Free-Newton-Krylov; the action of tangent matrices inside a linear solver process are approximated via matrix-free finite-differencing of the nonlinear residual equations. For details, consult the [PETSc manual](https://petsc.org/release/manual/snes/#sec-nlmatrixfree). Users can also provide the implementations of the Jacobian. This can be accomplished in two ways: - PETSc style using NonlinearSolver::jacobian - deal.II style using NonlinearSolver::setup_jacobian and NonlinearSolver::solve_with_jacobian. The preconditioning matrix can be specified using NonlinearSolver::set_matrix(). In case both approaches are implemented, the deal.II style will be used. NonlinearSolver::set_matrices() must be used in case the user wants to provide the iteration matrix of the tangent system in the deal.II style approach, thus replacing the matrix-free linearization. The correctness of the constructed Jacobians passed using NonlinearSolver::set_matrix() can be checked using @code ./myApp -snes_test_jacobian @endcode See NonlinearSolver::set_matrix() and NonlinearSolver::set_matrices() for additional details. The deal.II style approach still allows command line customization, like for example, @code ./myApp -snes_type newtontr -ksp_type cg @endcode in case the user wants to change the default nonlinear solver to a trust region solver and iterate on the matrix-free tangent system with CG, still using NonlinearSolver::solve_with_jacobian as a preconditioner. The PETSc style approach has instead the advantage that only the matrix assembly procedure has to be implemented, thus allowing quicker implementations and faster turnaround for experimenting with linear solver preconditioning configurations via command line customizations, like for example, @code ./myApp -ksp_type cg -pc_type gamg @endcode In case the nonlinear equations are derived from energy minimization arguments, it may be beneficial to perform linesearch or test trust-region model reductions using the energy functional. In such cases, users can set an optional NonlinearSolver::energy callback. @note This class, function, or variable is a template, and it can only be instantiated if the following condition is true: @code (concepts::is_dealii_petsc_vector_type<VectorType>
|| std::constructible_from<VectorType,Vec>) && (concepts::is_dealii_petsc_matrix_type<PMatrixType> || std::constructible_from<PMatrixType,Mat>) && (concepts::is_dealii_petsc_matrix_type<AMatrixType> || std::constructible_from<AMatrixType, Mat>) If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.
Definition at line 269 of file petsc_snes.h.
using PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::real_type = PetscReal |
Type that holds real-valued numbers.
Used to represent norms.
Definition at line 277 of file petsc_snes.h.
PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::NonlinearSolver | ( | const NonlinearSolverData & | data = NonlinearSolverData(), |
const MPI_Comm | mpi_comm = PETSC_COMM_WORLD ) |
Constructor.
PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::~NonlinearSolver | ( | ) |
Destructor.
PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::operator SNES | ( | ) | const |
Conversion operator to gain access to the underlying PETSc type. If you do this, you cut this class off some information it may need, so this conversion operator should only be used if you know what you do.
SNES PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::petsc_snes | ( | ) |
Return the PETSc SNES object.
MPI_Comm PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::get_mpi_communicator | ( | ) | const |
Return the underlying MPI communicator.
void PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::reinit | ( | ) |
Reset the solver, it does not change the customization.
void PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::reinit | ( | const NonlinearSolverData & | data | ) |
Reset solver. Change customization according to data
.
void PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::set_matrix | ( | PMatrixType & | P | ) |
Set the preconditioning matrix only.
When used with NonlinearSolver::setup_jacobian and NonlinearSolver::solve_with_jacobian, PETSc will approximate the linear system matrix-vector product using an internal matrix-free representation.
When used with NonlinearSolver::jacobian PETSc will use the same matrix for both preconditioning and matrix-vector products.
void PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::set_matrices | ( | AMatrixType & | A, |
PMatrixType & | P ) |
Set both the linear system matrix and the preconditioning matrix that PETSc will use (can be the same matrix). In this case, the Jacobian-Free-Newton-Krylov approach will not be used.
unsigned int PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::solve | ( | VectorType & | x | ) |
Solve the nonlinear system of equations \(F(x) = 0\).
This function returns the number of iterations. The vector x
must contain the initial guess. Upon returning, the x
vector contains the solution.
unsigned int PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::solve | ( | VectorType & | x, |
PMatrixType & | P ) |
Solve the nonlinear system of equations \(F(x) = 0\).
This function returns the number of iterations. The vector x
must contain the initial guess. Upon returning, the x
vector contains the solution.
Here we also set the matrix to precondition the tangent system.
unsigned int PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::solve | ( | VectorType & | x, |
AMatrixType & | A, | ||
PMatrixType & | P ) |
Solve the nonlinear system of equations \(F(x) = 0\).
This function returns the number of iterations. The vector x
must contain the initial guess. Upon returning, the x
vector contains the solution.
Here we also set the matrices to describe and precondition the tangent system.
std::function<void(const VectorType &x, VectorType &res)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::residual |
Callback for the computation of the nonlinear residual \(F(x)\).
Definition at line 387 of file petsc_snes.h.
std::function<void(const VectorType &x, AMatrixType &A, PMatrixType &P)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::jacobian |
Callback for the computation of the Jacobian \(\dfrac{\partial F}{\partial x}\).
Definition at line 399 of file petsc_snes.h.
std::function<void(const VectorType &x, const unsigned int step_number, const real_type f_norm)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::monitor |
Callback for monitoring the solution process.
This function is called by NonlinearSolver at the beginning of each step. Input arguments are the current step number and the current value for \(||F(x)||\).
Definition at line 416 of file petsc_snes.h.
std::function<void(const VectorType &x)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::setup_jacobian |
Callback for the set up of the Jacobian system.
This callback gives full control to users to set up the tangent operator \(\dfrac{\partial F}{\partial x}\).
Solvers must be provided via NonlinearSolver::solve_with_jacobian.
Definition at line 431 of file petsc_snes.h.
std::function<void(const VectorType &src, VectorType &dst)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::solve_with_jacobian |
Callback for the solution of the tangent system set up with NonlinearSolver::setup_jacobian.
This is used as a preconditioner inside the Krylov process.
Definition at line 445 of file petsc_snes.h.
std::function<void(const VectorType &x, real_type &energy_value)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::energy |
Callback for the computation of the energy function.
This is usually not needed, since by default SNES assumes that the objective function to be minimized is \(\frac{1}{2} || F(x) ||^2 \).
However, if the nonlinear equations are derived from energy arguments, it may be useful to use this callback to perform linesearch or to test for the reduction in a trust region step.
The value of the energy function must be returned in energy_value
.
Definition at line 464 of file petsc_snes.h.
|
private |
The PETSc SNES object.
Definition at line 470 of file petsc_snes.h.
|
private |
Pointers to the internal PETSc matrix objects.
Definition at line 475 of file petsc_snes.h.
|
private |
Definition at line 476 of file petsc_snes.h.
|
private |
This flag is used to support versions of PETSc older than 3.13.
Definition at line 481 of file petsc_snes.h.
|
mutableprivate |
A pointer to any exception that may have been thrown in user-defined call-backs and that we have to deal after the KINSOL function we call has returned.
Definition at line 488 of file petsc_snes.h.