Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Private Attributes | List of all members
PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType > Class Template Reference

#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, NonlinearSolverA
 
SmartPointer< PMatrixType, NonlinearSolverP
 
bool need_dummy_assemble
 
std::exception_ptr pending_exception
 

Detailed Description

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
class PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >

Interface to PETSc SNES solver for nonlinear equations. The SNES solver is described in the PETSc manual.

This class solves the nonlinear system of algebraic equations \(F(x) = 0\).

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:

class VectorType : public Subscriptor
...
explicit VectorType(Vec);
...
Vec & petsc_vector();
...
class MatrixType : public Subscriptor
...
explicit MatrixType(Mat);
...
Mat & petsc_matrix();
...

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

In alternative, users can also provide the implementation of the Jacobian. This can be accomplished in two ways:

In case both approaches are implemented, the deal.II style will be used.

The second approach is more in style with the deal.II philosophy and it also allows command line customization, like for example,

./myApp -snes_type newtontr -ksp_type cg

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 first approach has instead the advantage that only the matrix assembly procedure has to be provided, thus allowing quicker implementations and faster turnaround for experimenting with linear solver preconditioning configurations via command line customizations, like for example,

./myApp -ksp_type cg -pc_type gamg

See NonlinearSolver::set_matrix and NonlinearSolver::set_matrices for additional details.

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: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 258 of file petsc_snes.h.

Member Typedef Documentation

◆ real_type

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
using PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::real_type = PetscReal

Type that holds real-valued numbers.

Used to represent norms.

Definition at line 266 of file petsc_snes.h.

Constructor & Destructor Documentation

◆ NonlinearSolver()

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::NonlinearSolver ( const NonlinearSolverData data = NonlinearSolverData(),
const MPI_Comm  mpi_comm = PETSC_COMM_WORLD 
)

Constructor.

◆ ~NonlinearSolver()

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::~NonlinearSolver ( )

Destructor.

Member Function Documentation

◆ operator SNES()

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

◆ petsc_snes()

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
SNES PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::petsc_snes ( )

Return the PETSc SNES object.

◆ get_mpi_communicator()

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
MPI_Comm PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::get_mpi_communicator ( ) const

Return the underlying MPI communicator.

◆ reinit() [1/2]

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
void PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::reinit ( )

Reset the solver. It does not change the customization.

◆ reinit() [2/2]

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
void PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::reinit ( const NonlinearSolverData data)

Reset solver. Change customization according to data.

◆ set_matrix()

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

◆ set_matrices()

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

◆ solve() [1/3]

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

◆ solve() [2/3]

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

◆ solve() [3/3]

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

Member Data Documentation

◆ residual

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
std::function<void(const VectorType &x, VectorType &res)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::residual

Callback for the computation of the nonlinear residual \(F(x)\).

Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions.

Definition at line 374 of file petsc_snes.h.

◆ jacobian

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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}\).

Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions.

Definition at line 386 of file petsc_snes.h.

◆ monitor

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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 time step. Input arguments are the current step number and the current value for ||F(x)||.

Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions.

Definition at line 403 of file petsc_snes.h.

◆ setup_jacobian

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
std::function<void(const VectorType &x)> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::setup_jacobian

Callback to set up 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.

Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions.

Definition at line 418 of file petsc_snes.h.

◆ solve_with_jacobian

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions.

Definition at line 432 of file petsc_snes.h.

◆ energy

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
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.

Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions.

Definition at line 451 of file petsc_snes.h.

◆ snes

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
SNES PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::snes
private

The PETSc SNES object.

Definition at line 457 of file petsc_snes.h.

◆ A

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
SmartPointer<AMatrixType, NonlinearSolver> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::A
private

Pointers to the internal PETSc matrix objects.

Definition at line 462 of file petsc_snes.h.

◆ P

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
SmartPointer<PMatrixType, NonlinearSolver> PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::P
private

Definition at line 463 of file petsc_snes.h.

◆ need_dummy_assemble

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
bool PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::need_dummy_assemble
private

This flag is used to support versions of PETSc older than 3.13.

Definition at line 468 of file petsc_snes.h.

◆ pending_exception

template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
std::exception_ptr PETScWrappers::NonlinearSolver< VectorType, PMatrixType, AMatrixType >::pending_exception
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 475 of file petsc_snes.h.


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