Reference documentation for deal.II version 9.5.0
\(\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
Classes | Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
TimeStepping::EmbeddedExplicitRungeKutta< VectorType > Class Template Reference

#include <deal.II/base/time_stepping.h>

Inheritance diagram for TimeStepping::EmbeddedExplicitRungeKutta< VectorType >:
[legend]

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 () override
 
void free_memory ()
 
void initialize (const runge_kutta_method method) override
 
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) override
 
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 Statusget_status () const override
 
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) override
 

Protected Attributes

unsigned int n_stages
 
std::vector< double > b
 
std::vector< double > c
 
std::vector< std::vector< double > > a
 

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 = false
 
std::vector< double > b1
 
std::vector< double > b2
 
VectorType * last_stage = nullptr
 
Status status
 

Detailed Description

template<typename VectorType>
class TimeStepping::EmbeddedExplicitRungeKutta< VectorType >

This class is derived from RungeKutta and implements embedded explicit methods.

Definition at line 670 of file time_stepping.h.

Constructor & Destructor Documentation

◆ EmbeddedExplicitRungeKutta() [1/2]

template<typename VectorType >
TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::EmbeddedExplicitRungeKutta ( )
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.

◆ EmbeddedExplicitRungeKutta() [2/2]

template<typename VectorType >
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.

◆ ~EmbeddedExplicitRungeKutta()

template<typename VectorType >
TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::~EmbeddedExplicitRungeKutta ( )
inlineoverride

Destructor.

Definition at line 697 of file time_stepping.h.

Member Function Documentation

◆ free_memory()

template<typename VectorType >
void TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::free_memory ( )

If necessary, deallocate memory allocated by the object.

◆ initialize()

template<typename VectorType >
void TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::initialize ( const runge_kutta_method  method)
overridevirtual

Initialize the embedded explicit Runge-Kutta method.

Implements TimeStepping::RungeKutta< VectorType >.

◆ evolve_one_time_step() [1/3]

template<typename VectorType >
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::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 
)
overridevirtual

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 f}{\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.

Note
id_minus_tau_J_inverse is ignored since the method is explicit.

Implements TimeStepping::RungeKutta< VectorType >.

◆ evolve_one_time_step() [2/3]

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

◆ set_time_adaptation_parameters()

template<typename VectorType >
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.

◆ get_status()

template<typename VectorType >
const Status & TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::get_status ( ) const
overridevirtual

Return the status of the current object.

Implements TimeStepping::TimeStepping< VectorType >.

◆ compute_stages()

template<typename VectorType >
void TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::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

Compute the different stages needed.

◆ evolve_one_time_step() [3/3]

template<typename VectorType >
double TimeStepping::RungeKutta< VectorType >::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 
)
overridevirtualinherited

This function is used to advance from time t to t+ delta_t. F is a vector of functions \( 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. J_inverse is a vector functions that compute the inverse of the Jacobians associated to the implicit problems. The input parameters are the time, \( \tau \), and a vector. The output is the value of function at this point. This function returns the time at the end of the time step. When using Runge-Kutta methods, F and J_inverse can only contain one element.

Implements TimeStepping::TimeStepping< VectorType >.

Member Data Documentation

◆ coarsen_param

template<typename VectorType >
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::coarsen_param
private

This parameter is the factor (>1) by which the time step is multiplied when the time stepping can be coarsen.

Definition at line 803 of file time_stepping.h.

◆ refine_param

template<typename VectorType >
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::refine_param
private

This parameter is the factor (<1) by which the time step is multiplied when the time stepping must be refined.

Definition at line 809 of file time_stepping.h.

◆ min_delta_t

template<typename VectorType >
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::min_delta_t
private

Smallest time step allowed.

Definition at line 814 of file time_stepping.h.

◆ max_delta_t

template<typename VectorType >
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::max_delta_t
private

Largest time step allowed.

Definition at line 819 of file time_stepping.h.

◆ refine_tol

template<typename VectorType >
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::refine_tol
private

Refinement tolerance: if the error estimate is larger than refine_tol, the time step is refined.

Definition at line 825 of file time_stepping.h.

◆ coarsen_tol

template<typename VectorType >
double TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::coarsen_tol
private

Coarsening tolerance: if the error estimate is smaller than coarse_tol, the time step is coarsen.

Definition at line 831 of file time_stepping.h.

◆ last_same_as_first

template<typename VectorType >
bool TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::last_same_as_first = false
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 837 of file time_stepping.h.

◆ b1

template<typename VectorType >
std::vector<double> TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::b1
private

Butcher tableau coefficients.

Definition at line 842 of file time_stepping.h.

◆ b2

template<typename VectorType >
std::vector<double> TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::b2
private

Butcher tableau coefficients.

Definition at line 847 of file time_stepping.h.

◆ last_stage

template<typename VectorType >
VectorType* TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::last_stage = nullptr
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 853 of file time_stepping.h.

◆ status

template<typename VectorType >
Status TimeStepping::EmbeddedExplicitRungeKutta< VectorType >::status
private

Status structure of the object.

Definition at line 858 of file time_stepping.h.

◆ n_stages

template<typename VectorType >
unsigned int TimeStepping::RungeKutta< VectorType >::n_stages
protectedinherited

Number of stages of the Runge-Kutta method.

Definition at line 282 of file time_stepping.h.

◆ b

template<typename VectorType >
std::vector<double> TimeStepping::RungeKutta< VectorType >::b
protectedinherited

Butcher tableau coefficients.

Definition at line 287 of file time_stepping.h.

◆ c

template<typename VectorType >
std::vector<double> TimeStepping::RungeKutta< VectorType >::c
protectedinherited

Butcher tableau coefficients.

Definition at line 292 of file time_stepping.h.

◆ a

template<typename VectorType >
std::vector<std::vector<double> > TimeStepping::RungeKutta< VectorType >::a
protectedinherited

Butcher tableau coefficients.

Definition at line 297 of file time_stepping.h.


The documentation for this class was generated from the following file: