16#ifndef dealii_time_stepping_h
17#define dealii_time_stepping_h
174 template <
typename VectorType>
195 std::vector<std::function<VectorType(
const double,
const VectorType &)>>
197 std::vector<std::function<
198 VectorType(
const double,
const double,
const VectorType &)>> &J_inverse,
221 template <
typename VectorType>
249 std::vector<std::function<VectorType(
const double,
const VectorType &)>>
251 std::vector<std::function<
252 VectorType(
const double,
const double,
const VectorType &)>> &J_inverse,
255 VectorType &y)
override;
270 const std::function<VectorType(
const double,
const VectorType &)> &f,
272 VectorType(
const double,
const double,
const VectorType &)>
273 & id_minus_tau_J_inverse,
287 std::vector<double>
b;
292 std::vector<double>
c;
297 std::vector<std::vector<double>>
a;
306 template <
typename VectorType>
345 const std::function<VectorType(
const double,
const VectorType &)> &f,
347 VectorType(
const double,
const double,
const VectorType &)>
348 & id_minus_tau_J_inverse,
351 VectorType &y)
override;
362 const std::function<VectorType(
const double,
const VectorType &)> &f,
391 const std::function<VectorType(
const double,
const VectorType &)> &f,
393 const double delta_t,
394 const VectorType & y,
395 std::vector<VectorType> &f_stages)
const;
410 template <
typename VectorType>
447 const std::function<VectorType(
const double,
const VectorType &)> &f,
449 VectorType(
const double,
const double,
const VectorType &)>
450 & id_minus_tau_J_inverse,
453 VectorType &y)
override;
466 const std::function<VectorType(
const double,
const VectorType &)> &f,
469 VectorType &solution,
481 std::vector<double> &
b,
482 std::vector<double> &
c)
const;
508 const std::function<VectorType(
const double,
const VectorType &)> &f,
510 const double factor_solution,
511 const double factor_ai,
512 const VectorType &corrent_ri,
514 VectorType & solution,
515 VectorType & next_ri)
const;
529 template <
typename VectorType>
548 const unsigned int max_it = 100,
570 const std::function<VectorType(
const double,
const VectorType &)> &f,
572 VectorType(
const double,
const double,
const VectorType &)>
573 & id_minus_tau_J_inverse,
576 VectorType &y)
override;
615 const std::function<VectorType(
const double,
const VectorType &)> &f,
617 VectorType(
const double,
const double,
const VectorType &)>
618 & id_minus_tau_J_inverse,
622 std::vector<VectorType> &f_stages);
629 const std::function<
void(
const VectorType &, VectorType &)> &get_residual,
630 const std::function<VectorType(
const VectorType &)>
631 & id_minus_tau_J_inverse,
639 const std::function<VectorType(
const double,
const VectorType &)> &f,
642 const VectorType &new_y,
644 VectorType & tendency,
645 VectorType & residual)
const;
669 template <
typename VectorType>
689 const double min_delta = 1e-14,
690 const double max_delta = 1e100,
729 const std::function<VectorType(
const double,
const VectorType &)> &f,
731 VectorType(
const double,
const double,
const VectorType &)>
732 & id_minus_tau_J_inverse,
735 VectorType &y)
override;
746 const std::function<VectorType(
const double,
const VectorType &)> &f,
757 const double min_delta,
758 const double max_delta,
793 const std::function<VectorType(
const double,
const VectorType &)> &f,
795 const double delta_t,
796 const VectorType & y,
797 std::vector<VectorType> &f_stages);
842 std::vector<double>
b1;
847 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