Reference documentation for deal.II version 9.5.0
|
#include <deal.II/lac/petsc_ts.h>
Public Types | |
using | real_type = PetscReal |
Public Member Functions | |
TimeStepper (const TimeStepperData &data=TimeStepperData(), const MPI_Comm mpi_comm=PETSC_COMM_WORLD) | |
~TimeStepper () | |
operator TS () const | |
TS | petsc_ts () |
MPI_Comm | get_mpi_communicator () const |
void | reinit () |
void | reinit (const TimeStepperData &data) |
void | set_matrix (PMatrixType &P) |
void | set_matrices (AMatrixType &A, PMatrixType &P) |
real_type | get_time () |
real_type | get_time_step () |
unsigned int | solve (VectorType &y) |
unsigned int | solve (VectorType &y, PMatrixType &P) |
unsigned int | solve (VectorType &y, AMatrixType &A, PMatrixType &P) |
Public Attributes | |
std::function< void(const real_type t, const VectorType &y, const VectorType &y_dot, VectorType &res)> | implicit_function |
std::function< void(const real_type t, const VectorType &y, const VectorType &y_dot, const real_type alpha, AMatrixType &A, PMatrixType &P)> | implicit_jacobian |
std::function< void(const real_type t, const VectorType &y, VectorType &res)> | explicit_function |
std::function< void(const real_type t, const VectorType &y, AMatrixType &A, PMatrixType &P)> | explicit_jacobian |
std::function< void(const real_type t, const VectorType &y, const unsigned int step_number)> | monitor |
std::function< void(const real_type t, const VectorType &y, const VectorType &ydot, const real_type alpha)> | setup_jacobian |
std::function< void(const VectorType &src, VectorType &dst)> | solve_with_jacobian |
std::function< IndexSet()> | algebraic_components |
Private Attributes | |
TS | ts |
SmartPointer< AMatrixType, TimeStepper > | A |
SmartPointer< PMatrixType, TimeStepper > | P |
bool | need_dae_tolerances |
bool | need_dummy_assemble |
std::exception_ptr | pending_exception |
bool | error_in_function |
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.
In alternative, users can also provide the implementations of the Jacobians. This can be accomplished in two ways:
In case both approaches are coded, 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,
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 first approach has instead the advantage that only the matrix assembly procedure has to be coded, thus allowing quicker implementations and faster turnaround for experimenting with linear solver preconditioning configurations via command line customizations, like for example,
See TimeStepper::set_matrix and TimeStepper::set_matrices for additional details.
Definition at line 319 of file petsc_ts.h.
using PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::real_type = PetscReal |
Type that holds real-valued numbers.
Used to represent time and norms tolerances.
Definition at line 327 of file petsc_ts.h.
PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::TimeStepper | ( | const TimeStepperData & | data = TimeStepperData() , |
const MPI_Comm | mpi_comm = PETSC_COMM_WORLD |
||
) |
Constructor.
PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::~TimeStepper | ( | ) |
Destructor.
PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::operator TS | ( | ) | 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.
TS PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::petsc_ts | ( | ) |
Return the PETSc TS object.
MPI_Comm PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::get_mpi_communicator | ( | ) | const |
Return the underlying MPI communicator.
void PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::reinit | ( | ) |
Reset the solver, it does not change the customization.
void PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::reinit | ( | const TimeStepperData & | data | ) |
Reset solver. Change customization according to data
.
void PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::set_matrix | ( | PMatrixType & | P | ) |
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.
void PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::set_matrices | ( | AMatrixType & | A, |
PMatrixType & | P | ||
) |
Set both the linear system matrix and the preconditioning matrix that PETSc will use.
real_type PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::get_time | ( | ) |
Return current time.
real_type PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::get_time_step | ( | ) |
Return current time step.
unsigned int PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::solve | ( | VectorType & | y | ) |
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.
unsigned int PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::solve | ( | VectorType & | y, |
PMatrixType & | P | ||
) |
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.
unsigned int PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::solve | ( | VectorType & | y, |
AMatrixType & | A, | ||
PMatrixType & | P | ||
) |
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.
std::function<void(const real_type t, const VectorType &y, const VectorType &y_dot, VectorType & res)> PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::implicit_function |
Callback for the computation of the implicit residual \(F(t, y, \dot y)\).
Definition at line 453 of file petsc_ts.h.
std::function<void(const real_type t, const VectorType &y, const VectorType &y_dot, const real_type alpha, AMatrixType & A, PMatrixType & P)> PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::implicit_jacobian |
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.
Definition at line 475 of file petsc_ts.h.
std::function<void(const real_type t, const VectorType &y, VectorType &res)> PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::explicit_function |
Callback for the computation of the explicit residual \(G(t, y)\).
Definition at line 486 of file petsc_ts.h.
std::function<void(const real_type t, const VectorType &y, AMatrixType & A, PMatrixType & P)> PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::explicit_jacobian |
Callback for the computation of the explicit Jacobian \(\dfrac{\partial G}{\partial y}\).
Definition at line 501 of file petsc_ts.h.
std::function<void(const real_type t, const VectorType & y, const unsigned int step_number)> PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::monitor |
Callback for monitoring the solution process.
This function is called by TimeStepper at the beginning of each time step.
Definition at line 517 of file petsc_ts.h.
std::function<void(const real_type t, const VectorType &y, const VectorType &ydot, const real_type alpha)> PETScWrappers::TimeStepper< 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 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.
Definition at line 542 of file petsc_ts.h.
std::function<void(const VectorType &src, VectorType &dst)> PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::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.
Definition at line 556 of file petsc_ts.h.
std::function<IndexSet()> PETScWrappers::TimeStepper< VectorType, PMatrixType, AMatrixType >::algebraic_components |
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.
Definition at line 566 of file petsc_ts.h.
|
private |
The PETSc object.
Definition at line 572 of file petsc_ts.h.
|
private |
Pointers to the internal PETSc matrix objects.
Definition at line 577 of file petsc_ts.h.
|
private |
Definition at line 578 of file petsc_ts.h.
|
private |
This flag is set when changing the customization and used within solve.
Definition at line 583 of file petsc_ts.h.
|
private |
This flag is used to support versions of PETSc older than 3.13.
Definition at line 588 of file petsc_ts.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 595 of file petsc_ts.h.
|
private |
Internal data to handle recoverable errors.
Definition at line 600 of file petsc_ts.h.