![]() |
Reference documentation for deal.II version GIT 7b2de2f2f9 2023-09-24 11:00:02+00:00
|
Modules | |
SLEPcWrappers | |
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) | |
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) | |
template<typename VectorType = PETScWrappers::VectorBase, typename PMatrixType = PETScWrappers::MatrixBase, typename AMatrixType = PMatrixType> | |
PETScWrappers::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 TimeStepper | |
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.
|
inline |
Initialization parameters for NonlinearSolverData.
Running parameters:
options_prefix | The string indicating the options prefix for command line customization. |
snes_type | The string indicating the PETSc SNES solver type. |
snes_linesearch_type | The string indicating the PETSc linesearch type. |
absolute_tolerance | Absolute error tolerance. |
relative_tolerance | Relative error tolerance. |
step_tolerance | Step tolerance. |
maximum_non_linear_iterations | Maximum number of iterations allowed. |
max_n_function_evaluations | Maximum number of function evaluations allowed. |
Definition at line 75 of file petsc_snes.h.
|
inline |
Initialization parameters for TimeStepper.
Running parameters:
options_prefix | The string indicating the options prefix for command line customization. |
ts_type | The string indicating the PETSc solver type. |
initial_time | Initial simulation time. |
final_time | Final simulation time. |
initial_step_size | Initial step size. |
max_steps | Maximum number of steps allowed. |
match_step | Whether or not to exactly stop at final time or step over it. |
Error parameters:
ts_adapt_type | The string indicating the PETSc time step adaptor type. |
minimum_step_size | Minimum step size allowed. |
maximum_step_size | Maximum step size allowed. |
absolute_tolerance | Absolute error tolerance. |
relative_tolerance | Relative error tolerance. |
ignore_algebraic_lte | Ignore 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.
Definition at line 88 of file petsc_ts.h.
PETScWrappers::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 the PETSc TS solver for Ordinary Differential Equations and Differential-Algebraic Equations. The TS solver is described in the PETSc manual.
This class supports two kinds of formulations. The explicit formulation:
\[ \begin{cases} \dot y = G(t,y)\, , \\ y(t_0) = y_0\, , \\ \end{cases} \]
and the implicit formulation:
\[ \begin{cases} F(t,y,\dot y) = 0\, , \\ y(t_0) = y_0\, . \\ \end{cases} \]
The interface to PETSc is realized by means of std::function callbacks like in the SUNDIALS::IDA and SUNDIALS::ARKode classes.
TimeStepper supports any vector and matrix type having constructors and methods:
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 explicit solvers (like for example explicit Runge-Kutta methods), the user only needs to provide the implementation of \(G\) via the TimeStepper::explicit_function. For implicit solvers, users have also the alternative of providing the \(F\) function via TimeStepper::implicit_function. IMEX methods are also supported by providing both callbacks.
The default linearization procedure of an implicit 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 that are ODE-solver specific. For details, consult the PETSc manual.
Users can also provide the implementations of the Jacobians. This can be accomplished in two ways:
TimeStepper::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 TimeStepper::set_matrix() can be checked using
See TimeStepper::set_matrix() and TimeStepper::set_matrices() for additional details.
The deal.II style approach still allows command line customization, like for example,
in case the user wants to change the default nonlinear solver to a trust region solver and iterate on the tangent system with CG, still using TimeStepper::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,
Type that holds real-valued numbers.
Used to represent time and norms tolerances.
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 TS 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 TimeStepper::setup_jacobian and TimeStepper::solve_with_jacobian, PETSc will approximate the linear system matrix-vector product using an internal matrix-free representation.
When used with TimeStepper::implicit_jacobian or TimeStepper::explicit_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.
Return current time.
Return current time step.
Return current step number.
Integrate the differential-algebraic equations starting from y
.
This function returns the final number of computed steps. Upon returning, the y
vector contains the solution of the DAE at the end time.
Integrate the differential-algebraic equations starting from y
.
This function returns the final number of computed steps. Upon returning, the y
vector contains the solution of the DAE at the end time.
Here we also set the matrix to precondition the tangent system.
Integrate the differential-algebraic equations starting from y
.
This function returns the final number of computed steps. Upon returning, the y
vector contains the solution of the DAE at the end time.
Here we also set the matrices to describe and precondition the tangent system.
Callback for the computation of the implicit residual \(F(t, y, \dot y)\).
Callback for the computation of the implicit Jacobian \(\dfrac{\partial F}{\partial y} + \alpha \dfrac{\partial F}{\partial \dot y}\).
All implicit solvers implementations are recast to use the above linearization. The \(\alpha\) parameter is time-step and solver-type specific.
Callback for the computation of the explicit residual \(G(t, y)\).
Callback for the computation of the explicit Jacobian \(\dfrac{\partial G}{\partial y}\).
Callback for monitoring the solution process.
This function is called by TimeStepper at the beginning of each time step.
Callback for the set up of the Jacobian system.
This callback gives full control to users to set up the linearized equations \(\dfrac{\partial F}{\partial y} + \alpha \dfrac{\partial F}{\partial \dot y}\).
All implicit solvers implementations are recast to use the above linearization. The \(\alpha\) parameter is time-step and solver-type specific.
Solvers must be provided via TimeStepper::solve_with_jacobian.
Callback for the solution of the tangent system set up with TimeStepper::setup_jacobian.
This is used as a preconditioner inside the Krylov process.
Callback to return an index set containing the algebraic components.
Implementation of this function is optional. If your equation is also algebraic (i.e., it contains algebraic constraints, or Lagrange multipliers), you should implement this function in order to return only these components of your system.
Callback to distribute solution to hanging nodes.
Implementation of this function is optional. It is called at the end of each successful stage. The same functionality can be equivalently implemented in TimeStepper::solve_with_jacobian.
Callback to set up mesh adaption.
Implementation of this function is optional. resize
must be set to true if mesh adaption is to be performed, false otherwise. The y
vector contains the current solution, t
the current time, @ step the step number. Solution transfer and actual mesh adaption must be performed in a separate callback, TimeStepper::interpolate
Callback to interpolate vectors and perform mesh adaption.
Implementation of this function is mandatory if TimeStepper::decide_for_coarsening_and_refinement is used. This function must perform mesh adaption and interpolate the discrete functions that are stored in all_in
onto the refined and/or coarsenend grid. Output vectors must be created inside the callback.
The PETSc object.
Pointers to the internal PETSc matrix objects.
Object to apply solve_with_jacobian.
This flag is set when changing the customization and used within solve.
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.
Internal data to handle recoverable errors.
Setup callbacks.
This function is called inside TimeStepper::solve routines and does not need to be called by the user. It is used to reinitialize the solver if mesh adaption has been performed.
Setup algebraic constraints.
This function is called inside TimeStepper::solve routines and does not need to be called by the user. It is used to reinitialize the solver if mesh adaption has been performed.
Definition at line 320 of file petsc_ts.h.