Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/time_stepping.h>
Classes | |
struct | Status |
Public Member Functions | |
EmbeddedExplicitRungeKutta ()=default | |
EmbeddedExplicitRungeKutta (const runge_kutta_method method, const double coarsen_param=1.2, const double refine_param=0.8, const double min_delta=1e-14, const double max_delta=1e100, const double refine_tol=1e-8, const double coarsen_tol=1e-12) | |
~EmbeddedExplicitRungeKutta () | |
void | free_memory () |
void | initialize (const runge_kutta_method method) |
double | evolve_one_time_step (const std::function< VectorType(const double, const VectorType &)> &f, const std::function< VectorType(const double, const double, const VectorType &)> &id_minus_tau_J_inverse, double t, double delta_t, VectorType &y) |
double | evolve_one_time_step (const std::function< VectorType(const double, const VectorType &)> &f, double t, double delta_t, VectorType &y) |
void | set_time_adaptation_parameters (const double coarsen_param, const double refine_param, const double min_delta, const double max_delta, const double refine_tol, const double coarsen_tol) |
const Status & | get_status () const |
Public Member Functions inherited from TimeStepping::RungeKutta< VectorType > | |
virtual | ~RungeKutta ()=default |
double | evolve_one_time_step (std::vector< std::function< VectorType(const double, const VectorType &)> > &F, std::vector< std::function< VectorType(const double, const double, const VectorType &)> > &J_inverse, double t, double delta_t, VectorType &y) |
Public Member Functions inherited from TimeStepping::TimeStepping< VectorType > | |
virtual | ~TimeStepping ()=default |
Private Member Functions | |
void | compute_stages (const std::function< VectorType(const double, const VectorType &)> &f, const double t, const double delta_t, const VectorType &y, std::vector< VectorType > &f_stages) |
Private Attributes | |
double | coarsen_param |
double | refine_param |
double | min_delta_t |
double | max_delta_t |
double | refine_tol |
double | coarsen_tol |
bool | last_same_as_first |
std::vector< double > | b1 |
std::vector< double > | b2 |
VectorType * | last_stage |
Status | status |
Additional Inherited Members | |
Protected Attributes inherited from TimeStepping::RungeKutta< VectorType > | |
unsigned int | n_stages |
std::vector< double > | b |
std::vector< double > | c |
std::vector< std::vector< double > > | a |
This is class is derived from RungeKutta and implement embedded explicit methods.
Definition at line 432 of file time_stepping.h.
|
default |
Default constructor. initialize(runge_kutta_method) and set_time_adaptation_parameters(double, double, double, double, double, double) need to be called before the object can be used.
TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::EmbeddedExplicitRungeKutta | ( | const runge_kutta_method | method, |
const double | coarsen_param = 1.2 , |
||
const double | refine_param = 0.8 , |
||
const double | min_delta = 1e-14 , |
||
const double | max_delta = 1e100 , |
||
const double | refine_tol = 1e-8 , |
||
const double | coarsen_tol = 1e-12 |
||
) |
Constructor. This function calls initialize(runge_kutta_method) and initialize the parameters needed for time adaptation.
|
inline |
Destructor.
Definition at line 459 of file time_stepping.h.
void TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::free_memory | ( | ) |
If necessary, deallocate memory allocated by the object.
|
virtual |
Initialize the embedded explicit Runge-Kutta method.
Implements TimeStepping::RungeKutta< VectorType >.
|
virtual |
This function is used to advance from time t
to t+ delta_t
. f
is the function \( f(t,y) \) that should be integrated, the input parameters are the time t and the vector y and the output is value of f at this point. id_minus_tau_J_inverse
is a function that computes \( inv(I-\tau J)\) where \( I \) is the identity matrix, \( \tau \) is given, and \( J \) is the Jacobian \( \frac{\partial J}{\partial y} \). The input parameters are the time, \( \tau \), and a vector. The output is the value of function at this point. evolve_one_time_step returns the time at the end of the time step.
Implements TimeStepping::RungeKutta< VectorType >.
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::evolve_one_time_step | ( | const std::function< VectorType(const double, const VectorType &)> & | f, |
double | t, | ||
double | delta_t, | ||
VectorType & | y | ||
) |
This function is used to advance from time t
to t+ delta_t
. This function is similar to the one derived from TimeStepping, but does not required id_minus_tau_J_inverse because it is not used for explicit methods. evolve_one_time_step returns the time at the end of the time step.
void TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::set_time_adaptation_parameters | ( | const double | coarsen_param, |
const double | refine_param, | ||
const double | min_delta, | ||
const double | max_delta, | ||
const double | refine_tol, | ||
const double | coarsen_tol | ||
) |
Set the parameters necessary for the time adaptation.
|
virtual |
Return the status of the current object.
Implements TimeStepping::TimeStepping< VectorType >.
|
private |
Compute the different stages needed.
|
private |
This parameter is the factor (>1) by which the time step is multiplied when the time stepping can be coarsen.
Definition at line 549 of file time_stepping.h.
|
private |
This parameter is the factor (<1) by which the time step is multiplied when the time stepping must be refined.
Definition at line 555 of file time_stepping.h.
|
private |
Smallest time step allowed.
Definition at line 560 of file time_stepping.h.
|
private |
Largest time step allowed.
Definition at line 565 of file time_stepping.h.
|
private |
Refinement tolerance: if the error estimate is larger than refine_tol, the time step is refined.
Definition at line 571 of file time_stepping.h.
|
private |
Coarsening tolerance: if the error estimate is smaller than coarse_tol, the time step is coarsen.
Definition at line 577 of file time_stepping.h.
|
private |
If the flag is true, the last stage is the same as the first stage and one evaluation of f can be saved.
Definition at line 583 of file time_stepping.h.
|
private |
Butcher tableau coefficients.
Definition at line 588 of file time_stepping.h.
|
private |
Butcher tableau coefficients.
Definition at line 593 of file time_stepping.h.
|
private |
If the last_same_as_first flag is set to true, the last stage is saved and reused as the first stage of the next time step.
Definition at line 599 of file time_stepping.h.
|
private |
Status structure of the object.
Definition at line 604 of file time_stepping.h.