Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08: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
Modules | Namespaces | Classes | Functions
Collaboration diagram for PETScWrappers:

Modules

 SLEPcWrappers
 

Namespaces

namespace  PETScWrappers
 
namespace  PETScWrappers::MPI
 
namespace  internal
 

Classes

class  PETScWrappers::MPI::BlockSparseMatrix
 
class  PETScWrappers::MPI::BlockVector
 
class  PETScWrappers::FullMatrix
 
class  PETScWrappers::MatrixIterators::const_iterator
 
class  PETScWrappers::MatrixBase
 
class  PETScWrappers::MatrixFree
 
class  PETScWrappers::PreconditionBase
 
class  PETScWrappers::PreconditionJacobi
 
class  PETScWrappers::PreconditionBlockJacobi
 
class  PETScWrappers::PreconditionSOR
 
class  PETScWrappers::PreconditionSSOR
 
class  PETScWrappers::PreconditionICC
 
class  PETScWrappers::PreconditionILU
 
class  PETScWrappers::PreconditionLU
 
class  PETScWrappers::PreconditionBoomerAMG
 
class  PETScWrappers::PreconditionParaSails
 
class  PETScWrappers::PreconditionNone
 
class  PETScWrappers::PreconditionBDDC< dim >
 
class  PETScWrappers::SolverBase
 
class  PETScWrappers::SolverRichardson
 
class  PETScWrappers::SolverChebychev
 
class  PETScWrappers::SolverCG
 
class  PETScWrappers::SolverBiCG
 
class  PETScWrappers::SolverGMRES
 
class  PETScWrappers::SolverBicgstab
 
class  PETScWrappers::SolverCGS
 
class  PETScWrappers::SolverTFQMR
 
class  PETScWrappers::SolverTCQMR
 
class  PETScWrappers::SolverCR
 
class  PETScWrappers::SolverLSQR
 
class  PETScWrappers::SolverPreOnly
 
class  PETScWrappers::SparseDirectMUMPS
 
class  PETScWrappers::SparseMatrix
 
class  PETScWrappers::MPI::SparseMatrix
 
class  PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >
 
class  PETScWrappers::MPI::Vector
 
class  PETScWrappers::VectorBase
 

Functions

 PETScWrappers::NonlinearSolverData::NonlinearSolverData (const std::string &options_prefix="", const std::string &snes_type="", const std::string &snes_linesearch_type="", const real_type absolute_tolerance=0, const real_type relative_tolerance=0, const real_type step_tolerance=0, const int maximum_non_linear_iterations=-1, const int max_n_function_evaluations=-1)
 
template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType>
 PETScWrappers::DEAL_II_CXX20_REQUIRES ((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 >)) class NonlinearSolver
 
 PETScWrappers::TimeStepperData::TimeStepperData (const std::string &options_prefix="", const std::string &ts_type="", const real_type initial_time=0.0, const real_type final_time=0.0, const real_type initial_step_size=0.0, const int max_steps=-1, const bool match_step=false, const std::string &ts_adapt_type="none", const real_type minimum_step_size=-1.0, const real_type maximum_step_size=-1.0, const real_type absolute_tolerance=-1.0, const real_type relative_tolerance=-1.0, const bool ignore_algebraic_lte=true)
 

Detailed Description

The classes in this module are wrappers around functionality provided by the PETSc library. They provide a modern object-oriented interface that is compatible with the interfaces of the other linear algebra classes in deal.II. All classes and functions in this group reside in a namespace PETScWrappers.

These classes are only available if a PETSc installation was detected during configuration of deal.II. Refer to the README file for more details about this.

Function Documentation

◆ NonlinearSolverData()

PETScWrappers::NonlinearSolverData::NonlinearSolverData ( const std::string &  options_prefix = "",
const std::string &  snes_type = "",
const std::string &  snes_linesearch_type = "",
const real_type  absolute_tolerance = 0,
const real_type  relative_tolerance = 0,
const real_type  step_tolerance = 0,
const int  maximum_non_linear_iterations = -1,
const int  max_n_function_evaluations = -1 
)
inline

Initialization parameters for NonlinearSolverData.

Running parameters:

Parameters
options_prefixThe string indicating the options prefix for command line customization.
snes_typeThe string indicating the PETSc SNES solver type.
snes_linesearch_typeThe string indicating the PETSc linesearch type.
absolute_toleranceAbsolute error tolerance.
relative_toleranceRelative error tolerance.
step_toleranceStep tolerance.
maximum_non_linear_iterationsMaximum number of iterations allowed.
max_n_function_evaluationsMaximum number of function evaluations allowed.
Note
All parameters values specified here can be overridden by command line choices.

Definition at line 74 of file petsc_snes.h.

◆ DEAL_II_CXX20_REQUIRES()

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

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.

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

./myApp -snes_test_jacobian

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

The deal.II style approach still 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 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,

./myApp -ksp_type cg -pc_type gamg

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.

Type that holds real-valued numbers.

Used to represent norms.

Constructor.

Destructor.

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.

Return the PETSc SNES object.

Return the underlying MPI communicator.

Reset the solver, it does not change the customization.

Reset solver. Change customization according to data.

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

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

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.

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.

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

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

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.

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

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.

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.

The PETSc SNES object.

Pointers to the internal PETSc matrix objects.

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

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 260 of file petsc_snes.h.

◆ TimeStepperData()

PETScWrappers::TimeStepperData::TimeStepperData ( const std::string &  options_prefix = "",
const std::string &  ts_type = "",
const real_type  initial_time = 0.0,
const real_type  final_time = 0.0,
const real_type  initial_step_size = 0.0,
const int  max_steps = -1,
const bool  match_step = false,
const std::string &  ts_adapt_type = "none",
const real_type  minimum_step_size = -1.0,
const real_type  maximum_step_size = -1.0,
const real_type  absolute_tolerance = -1.0,
const real_type  relative_tolerance = -1.0,
const bool  ignore_algebraic_lte = true 
)
inline

Initialization parameters for TimeStepper.

Running parameters:

Parameters
options_prefixThe string indicating the options prefix for command line customization.
ts_typeThe string indicating the PETSc solver type.
initial_timeInitial simulation time.
final_timeFinal simulation time.
initial_step_sizeInitial step size.
max_stepsMaximum number of steps allowed.
match_stepWhether or not to exactly stop at final time or step over it.

Error parameters:

Parameters
ts_adapt_typeThe string indicating the PETSc time step adaptor type.
minimum_step_sizeMinimum step size allowed.
maximum_step_sizeMaximum step size allowed.
absolute_toleranceAbsolute error tolerance.
relative_toleranceRelative error tolerance.
ignore_algebraic_lteIgnore algebraic terms for error computations

Note that one between final_time or max_steps must be specified by the user, otherwise PETSc will complain. Adaptive time stepping is disabled by default. Negative values indicate using PETSc's default.

Note
All parameters values specified here can be overridden by command line choices.

Definition at line 87 of file petsc_ts.h.