15#ifndef dealii_time_stepping_h
16#define dealii_time_stepping_h
173 template <
typename VectorType>
194 std::vector<std::function<VectorType(
const double,
const VectorType &)>>
196 std::vector<std::function<
197 VectorType(
const double,
const double,
const VectorType &)>> &J_inverse,
220 template <
typename VectorType>
248 std::vector<std::function<VectorType(
const double,
const VectorType &)>>
250 std::vector<std::function<
251 VectorType(
const double,
const double,
const VectorType &)>> &J_inverse,
254 VectorType &y)
override;
269 const std::function<VectorType(
const double,
const VectorType &)> &f,
271 VectorType(
const double,
const double,
const VectorType &)>
272 &id_minus_tau_J_inverse,
286 std::vector<double>
b;
291 std::vector<double>
c;
296 std::vector<std::vector<double>>
a;
305 template <
typename VectorType>
344 const std::function<VectorType(
const double,
const VectorType &)> &f,
346 VectorType(
const double,
const double,
const VectorType &)>
347 &id_minus_tau_J_inverse,
350 VectorType &y)
override;
361 const std::function<VectorType(
const double,
const VectorType &)> &f,
390 const std::function<VectorType(
const double,
const VectorType &)> &f,
392 const double delta_t,
394 std::vector<VectorType> &f_stages)
const;
409 template <
typename VectorType>
446 const std::function<VectorType(
const double,
const VectorType &)> &f,
448 VectorType(
const double,
const double,
const VectorType &)>
449 &id_minus_tau_J_inverse,
452 VectorType &y)
override;
465 const std::function<VectorType(
const double,
const VectorType &)> &f,
468 VectorType &solution,
480 std::vector<double> &
b,
481 std::vector<double> &
c)
const;
507 const std::function<VectorType(
const double,
const VectorType &)> &f,
509 const double factor_solution,
510 const double factor_ai,
511 const VectorType &corrent_ri,
513 VectorType &solution,
514 VectorType &next_ri)
const;
528 template <
typename VectorType>
547 const unsigned int max_it = 100,
569 const std::function<VectorType(
const double,
const VectorType &)> &f,
571 VectorType(
const double,
const double,
const VectorType &)>
572 &id_minus_tau_J_inverse,
575 VectorType &y)
override;
614 const std::function<VectorType(
const double,
const VectorType &)> &f,
616 VectorType(
const double,
const double,
const VectorType &)>
617 &id_minus_tau_J_inverse,
621 std::vector<VectorType> &f_stages);
628 const std::function<
void(
const VectorType &, VectorType &)> &get_residual,
629 const std::function<VectorType(
const VectorType &)>
630 &id_minus_tau_J_inverse,
638 const std::function<VectorType(
const double,
const VectorType &)> &f,
641 const VectorType &new_y,
643 VectorType &tendency,
644 VectorType &residual)
const;
668 template <
typename VectorType>
688 const double min_delta = 1e-14,
689 const double max_delta = 1e100,
728 const std::function<VectorType(
const double,
const VectorType &)> &f,
730 VectorType(
const double,
const double,
const VectorType &)>
731 &id_minus_tau_J_inverse,
734 VectorType &y)
override;
745 const std::function<VectorType(
const double,
const VectorType &)> &f,
756 const double min_delta,
757 const double max_delta,
792 const std::function<VectorType(
const double,
const VectorType &)> &f,
794 const double delta_t,
796 std::vector<VectorType> &f_stages);
841 std::vector<double>
b1;
846 std::vector<double>
b2;
double evolve_one_time_step(const std::function< VectorType(const double, const VectorType &)> &f, double t, double delta_t, VectorType &y)
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)
EmbeddedExplicitRungeKutta()=default
~EmbeddedExplicitRungeKutta() override
const Status & get_status() const override
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)
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
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)
const Status & get_status() const override
void initialize(const runge_kutta_method method) override
ExplicitRungeKutta()=default
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) const
double evolve_one_time_step(const std::function< VectorType(const double, const VectorType &)> &f, double t, double delta_t, VectorType &y)
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
ExplicitRungeKutta(const runge_kutta_method method)
void compute_stages(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, std::vector< VectorType > &f_stages)
void newton_solve(const std::function< void(const VectorType &, VectorType &)> &get_residual, const std::function< VectorType(const VectorType &)> &id_minus_tau_J_inverse, VectorType &y)
const Status & get_status() const override
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
void compute_residual(const std::function< VectorType(const double, const VectorType &)> &f, double t, double delta_t, const VectorType &new_y, const VectorType &y, VectorType &tendency, VectorType &residual) const
void set_newton_solver_parameters(const unsigned int max_it, const double tolerance)
ImplicitRungeKutta()=default
ImplicitRungeKutta(const runge_kutta_method method, const unsigned int max_it=100, const double tolerance=1e-6)
void initialize(const runge_kutta_method method) override
LowStorageRungeKutta(const runge_kutta_method method)
void get_coefficients(std::vector< double > &a, std::vector< double > &b, std::vector< double > &c) const
LowStorageRungeKutta()=default
double evolve_one_time_step(const std::function< VectorType(const double, const VectorType &)> &f, double t, double delta_t, VectorType &solution, VectorType &vec_ri, VectorType &vec_ki)
const Status & get_status() const 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
void compute_one_stage(const std::function< VectorType(const double, const VectorType &)> &f, const double t, const double factor_solution, const double factor_ai, const VectorType &corrent_ri, VectorType &vec_ki, VectorType &solution, VectorType &next_ri) const
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
virtual ~RungeKutta() override=default
virtual void initialize(const runge_kutta_method method)=0
std::vector< std::vector< double > > a
virtual 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)=0
virtual 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)=0
virtual ~TimeStepping()=default
virtual const Status & get_status() const =0
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
embedded_runge_kutta_time_step
@ RK_CLASSIC_FOURTH_ORDER
@ LOW_STORAGE_RK_STAGE9_ORDER5
@ LOW_STORAGE_RK_STAGE3_ORDER3
@ LOW_STORAGE_RK_STAGE7_ORDER4
@ LOW_STORAGE_RK_STAGE5_ORDER4
embedded_runge_kutta_time_step exit_delta_t
unsigned int n_iterations
runge_kutta_method method
runge_kutta_method method
unsigned int n_iterations
runge_kutta_method method
runge_kutta_method method