Abstract class for PETSc's ODE solvers.
More...
#include <petsc.hpp>
Abstract class for PETSc's ODE solvers.
Definition at line 777 of file petsc.hpp.
The type of the ODE. Use ODE_SOLVER_LINEAR if the jacobians are linear and independent of time.
Enumerator |
---|
ODE_SOLVER_LINEAR |
|
ODE_SOLVER_GENERAL |
|
Definition at line 782 of file petsc.hpp.
mfem::PetscODESolver::PetscODESolver |
( |
MPI_Comm |
comm, |
|
|
const std::string & |
prefix = std::string() |
|
) |
| |
mfem::PetscODESolver::~PetscODESolver |
( |
| ) |
|
|
virtual |
Initialize the ODE solver.
Definition at line 3415 of file petsc.cpp.
mfem::PetscODESolver::operator TS |
( |
| ) |
const |
|
inline |
Conversion function to PETSc's TS type.
Definition at line 807 of file petsc.hpp.
void mfem::PetscODESolver::Run |
( |
Vector & |
x, |
|
|
double & |
t, |
|
|
double & |
dt, |
|
|
double |
tf |
|
) |
| |
|
virtual |
Perform time integration from time t [in] to time tf [in].
- Parameters
-
[in,out] | x | Approximate solution. |
[in,out] | t | Time associated with the approximate solution x. |
[in,out] | dt | Time step size. |
[in] | tf | Requested final time. |
The default implementation makes consecutive calls to Step() until reaching tf. The following rules describe the common behavior of the method:
- The input x [in] is the approximate solution for the input time t [in].
- The input dt [in] is the initial time step size.
- The output dt [out] is the last time step taken by the method which may be smaller or larger than the input dt [in] value, e.g. because of time step control.
- The output value of t [out] is not smaller than tf [in].
Reimplemented from mfem::ODESolver.
Definition at line 3537 of file petsc.cpp.
Specifies the desired format of the Jacobian in case a PetscParMatrix is not returned by the GetGradient methods
Definition at line 3482 of file petsc.cpp.
void mfem::PetscODESolver::Step |
( |
Vector & |
x, |
|
|
double & |
t, |
|
|
double & |
dt |
|
) |
| |
|
virtual |
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
- Parameters
-
[in,out] | x | Approximate solution. |
[in,out] | t | Time associated with the approximate solution x. |
[in,out] | dt | Time step size. |
The following rules describe the common behavior of the method:
- The input x [in] is the approximate solution for the input time t [in].
- The input dt [in] is the desired time step size, defining the desired target time: t [target] = t [in] + dt [in].
- The output x [out] is the approximate solution for the output time t [out].
- The output dt [out] is the last time step taken by the method which may be smaller or larger than the input dt [in] value, e.g. because of time step control.
- The method may perform more than one time step internally; in this case dt [out] is the last internal time step size.
- The output value of t [out] may be smaller or larger than t [target], however, it is not smaller than t [in] + dt [out], if at least one internal time step was performed.
- The value x [out] may be obtained by interpolation using internally stored data.
- In some cases, the contents of x [in] may not be used, e.g. when x [out] from a previous Step() call was obtained by interpolation.
- In consecutive calls to this method, the output t [out] of one Step() call has to be the same as the input t [in] to the next Step() call.
- If the previous rule has to be broken, e.g. to restart a time stepping sequence, then the ODE solver must be re-initialized by calling Init() between the two Step() calls.
Implements mfem::ODESolver.
Definition at line 3512 of file petsc.cpp.
The documentation for this class was generated from the following files: