Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07:20:02+00:00
\(\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
Classes | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
NonlinearSolverSelector< VectorType > Class Template Reference

#include <deal.II/numerics/nonlinear.h>

Classes

class  AdditionalData
 

Public Member Functions

 NonlinearSolverSelector ()
 
 NonlinearSolverSelector (const AdditionalData &additional_data)
 
 NonlinearSolverSelector (const AdditionalData &additional_data, const MPI_Comm &mpi_communicator)
 
void select (const typename AdditionalData::SolverType &type)
 
void set_data (const AdditionalData &additional_data)
 
void set_data (const typename TrilinosWrappers::NOXSolver< VectorType >::AdditionalData &additional_data, const Teuchos::RCP< Teuchos::ParameterList > &parameters=Teuchos::rcp(new Teuchos::ParameterList))
 
void set_data (const typename SUNDIALS::KINSOL< VectorType >::AdditionalData &additional_data)
 
void set_data (const typename PETScWrappers::NonlinearSolverData &additional_data)
 
void solve (VectorType &initial_guess_and_solution)
 

Public Attributes

std::function< void(VectorType &)> reinit_vector
 
std::function< void(const VectorType &src, VectorType &dst)> residual
 
std::function< void(const VectorType &current_u)> setup_jacobian
 
std::function< void(const VectorType &rhs, VectorType &dst, const double tolerance)> solve_with_jacobian
 

Private Member Functions

void solve_with_petsc (VectorType &initial_guess_and_solution)
 
void solve_with_petsc (PETScWrappers::MPI::Vector &initial_guess_and_solution)
 

Private Attributes

AdditionalData additional_data
 
MPI_Comm mpi_communicator
 
SUNDIALS::KINSOL< VectorType >::AdditionalData additional_data_kinsol
 
TrilinosWrappers::NOXSolver< VectorType >::AdditionalData additional_data_nox
 
Teuchos::RCP< Teuchos::ParameterList > parameters_nox
 
PETScWrappers::NonlinearSolverData additional_data_petsc_snes
 

Detailed Description

template<typename VectorType = Vector<double>>
class NonlinearSolverSelector< VectorType >

This class offers a unified interface to several nonlinear solver implementations like KINSOL, NOX, and SNES.

KINSOL nonlinear solvers are part of the SUNDIALS package, while NOX is part of Trilinos and SNES is part of PETSc, respectively. If no solver is manually specified, this class will automaticlaly choose one of the available solvers depending on the enabled dependencies.

By calling the solve function of this NonlinearSolverSelector, it selects the solve function of that Solver that was specified in the constructor of this class, similar to the SolverSelector class.

Usage

An example of code one would run with this class is the following:

// Generate a @p NonlinearSolverSelector that uses @p KINSOL
// Functions that are required for solving a nonlinear problem,
// which are utilized in both @p KINSOL and @p NOX.
nonlinear_solver.reinit_vector = [&](Vector<double> &x) {...};
nonlinear_solver.residual =
[&](const Vector<double> &evaluation_point,
nonlinear_solver.setup_jacobian =
[&](const Vector<double> &current_u) {...};
nonlinear_solver.solve_with_jacobian =
[&](const Vector<double> &rhs,
const double tolerance) {...};
// Calling the @p solve function with an initial guess.
nonlinear_solver.solve(current_solution);
std::function< void(const VectorType &src, VectorType &dst)> residual
Definition nonlinear.h:299
AdditionalData additional_data
Definition nonlinear.h:358

Definition at line 75 of file nonlinear.h.

Constructor & Destructor Documentation

◆ NonlinearSolverSelector() [1/3]

template<typename VectorType >
NonlinearSolverSelector< VectorType >::NonlinearSolverSelector ( )
default

Constructor, filling in default values

◆ NonlinearSolverSelector() [2/3]

template<typename VectorType >
NonlinearSolverSelector< VectorType >::NonlinearSolverSelector ( const AdditionalData additional_data)

Constructor, selecting the solver and other parametersspecified in additional_data.

Definition at line 490 of file nonlinear.h.

◆ NonlinearSolverSelector() [3/3]

template<typename VectorType >
NonlinearSolverSelector< VectorType >::NonlinearSolverSelector ( const AdditionalData additional_data,
const MPI_Comm mpi_communicator 
)

Constructor.

Parameters
additional_dataNonlinearSolverSelector configuration data
mpi_communicatorMPI communicator used by the nonlinear solver.

Definition at line 500 of file nonlinear.h.

Member Function Documentation

◆ select()

template<typename VectorType >
void NonlinearSolverSelector< VectorType >::select ( const typename AdditionalData::SolverType type)

Select a new nonlinear solver. All solver names used in this class are all lower case.

Definition at line 513 of file nonlinear.h.

◆ set_data() [1/4]

template<typename VectorType >
void NonlinearSolverSelector< VectorType >::set_data ( const AdditionalData additional_data)

Set the generic additional data for all nonlinear solvers.

Definition at line 404 of file nonlinear.h.

◆ set_data() [2/4]

template<typename VectorType >
void NonlinearSolverSelector< VectorType >::set_data ( const typename TrilinosWrappers::NOXSolver< VectorType >::AdditionalData &  additional_data,
const Teuchos::RCP< Teuchos::ParameterList > &  parameters = Teuchos::rcp(new Teuchos::ParameterList) 
)

Set the additional data for NOX. See TrilinosWrappers::NOXSolver for more information.

Definition at line 544 of file nonlinear.h.

◆ set_data() [3/4]

template<typename VectorType >
void NonlinearSolverSelector< VectorType >::set_data ( const typename SUNDIALS::KINSOL< VectorType >::AdditionalData &  additional_data)

Set the additional data for KINSOL. See SUNDIALS::KINSOL for more information.

Definition at line 559 of file nonlinear.h.

◆ set_data() [4/4]

template<typename VectorType = Vector<double>>
void NonlinearSolverSelector< VectorType >::set_data ( const typename PETScWrappers::NonlinearSolverData additional_data)

Set the additional data for SNES. See PETScWrappers::NonlinearSolverData for more information.

◆ solve()

template<typename VectorType >
void NonlinearSolverSelector< VectorType >::solve ( VectorType &  initial_guess_and_solution)

Solve the nonlinear system.

The content of initial_guess_and_solution is used as an initial guess and the final solution is stored in the same vector.

Definition at line 609 of file nonlinear.h.

◆ solve_with_petsc() [1/2]

template<typename VectorType >
void NonlinearSolverSelector< VectorType >::solve_with_petsc ( VectorType &  initial_guess_and_solution)
private

Solve with PETSc SNES. Internal functions specialized for PETSc Vectors. This is necessary to only instantiate the SNES solvers with PETSc Vectors, as this is the only vector type currently supported.

Definition at line 570 of file nonlinear.h.

◆ solve_with_petsc() [2/2]

void NonlinearSolverSelector< PETScWrappers::MPI::Vector >::solve_with_petsc ( PETScWrappers::MPI::Vector initial_guess_and_solution)
private

Definition at line 582 of file nonlinear.h.

Member Data Documentation

◆ reinit_vector

template<typename VectorType = Vector<double>>
std::function<void(VectorType &)> NonlinearSolverSelector< VectorType >::reinit_vector

A function object that users need to supply and that is intended to reinitize the given vector to its correct size, and block structure (if block vectors are used), along with any other properties necessary.

Definition at line 284 of file nonlinear.h.

◆ residual

template<typename VectorType = Vector<double>>
std::function<void(const VectorType &src, VectorType &dst)> NonlinearSolverSelector< VectorType >::residual

A function object that users should supply and that is intended to compute the residual dst = F(src).

Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions. Some of the underlying packages used by this class can deal with recoverable exceptions, whereas others cannot. As a consequence, if a callback throws an exception of type RecoverableUserCallbackError, then this exception may or may not be treated like any other exception.

Definition at line 299 of file nonlinear.h.

◆ setup_jacobian

template<typename VectorType = Vector<double>>
std::function<void(const VectorType &current_u)> NonlinearSolverSelector< VectorType >::setup_jacobian

A function object that users may supply and that is intended to prepare the linear solver for subsequent calls to solve_with_jacobian).

The job of setup_jacobian() is to prepare the linear solver for subsequent calls to solve_with_jacobian(), in the solution of linear systems \(Ax = b\). The exact nature of this system depends on the SolutionStrategy that has been selected.

In the cases strategy = SolutionStrategy::newton or SolutionStrategy::linesearch, \(A\) is the Jacobian \(J = \partial F/\partial u\). If strategy = SolutionStrategy::picard, \(A\) is the approximate Jacobian matrix \(L\).

Parameters
current_uCurrent value of \(u\)
Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions. Some of the underlying packages used by this class can deal with recoverable exceptions, whereas others cannot. As a consequence, if a callback throws an exception of type RecoverableUserCallbackError, then this exception may or may not be treated like any other exception.

Definition at line 327 of file nonlinear.h.

◆ solve_with_jacobian

template<typename VectorType = Vector<double>>
std::function< void(const VectorType &rhs, VectorType &dst, const double tolerance)> NonlinearSolverSelector< VectorType >::solve_with_jacobian

A function object that users may supply and that is intended to solve a linear system with the Jacobian matrix.

Note
PETSc SNES does not provide us with a linear tolerance to solve linear system with. We will provide a value of 1e-6 in that case.
Parameters
[in]rhsThe system right hand side to solve for.
[out]dstThe solution of \(J^{-1} * \texttt{src}\).
[in]toleranceThe tolerance with which to solve the linear system of equations.
Note
This variable represents a user provided callback. See there for a description of how to deal with errors and other requirements and conventions. Some of the underlying packages used by this class can deal with recoverable exceptions, whereas others cannot. As a consequence, if a callback throws an exception of type RecoverableUserCallbackError, then this exception may or may not be treated like any other exception.

Definition at line 352 of file nonlinear.h.

◆ additional_data

template<typename VectorType = Vector<double>>
AdditionalData NonlinearSolverSelector< VectorType >::additional_data
private

NonlinearSolverSelector configuration data.

Definition at line 358 of file nonlinear.h.

◆ mpi_communicator

template<typename VectorType = Vector<double>>
MPI_Comm NonlinearSolverSelector< VectorType >::mpi_communicator
private

The MPI communicator to be used by this solver, if any.

Definition at line 364 of file nonlinear.h.

◆ additional_data_kinsol

template<typename VectorType = Vector<double>>
SUNDIALS::KINSOL<VectorType>::AdditionalData NonlinearSolverSelector< VectorType >::additional_data_kinsol
private

KINSOL configuration data

Definition at line 370 of file nonlinear.h.

◆ additional_data_nox

template<typename VectorType = Vector<double>>
TrilinosWrappers::NOXSolver<VectorType>::AdditionalData NonlinearSolverSelector< VectorType >::additional_data_nox
private

NOX configuration data

Definition at line 378 of file nonlinear.h.

◆ parameters_nox

template<typename VectorType = Vector<double>>
Teuchos::RCP<Teuchos::ParameterList> NonlinearSolverSelector< VectorType >::parameters_nox
private
Initial value:
=
Teuchos::rcp(new Teuchos::ParameterList)

Definition at line 379 of file nonlinear.h.

◆ additional_data_petsc_snes

template<typename VectorType = Vector<double>>
PETScWrappers::NonlinearSolverData NonlinearSolverSelector< VectorType >::additional_data_petsc_snes
private

PETSc SNES configuration data

Definition at line 387 of file nonlinear.h.


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