Reference documentation for deal.II version 9.3.2
\(\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\}}\)
time_stepping.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2014 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_time_stepping_h
17 #define dealii_time_stepping_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 
24 #include <functional>
25 #include <vector>
26 
28 
33 namespace TimeStepping
34 {
61  {
143  invalid
144  };
145 
146 
147 
153  {
166  };
167 
168 
169 
174  template <typename VectorType>
176  {
177  public:
181  virtual ~TimeStepping() = default;
182 
193  virtual double
195  std::vector<std::function<VectorType(const double, const VectorType &)>>
196  & F,
197  std::vector<std::function<
198  VectorType(const double, const double, const VectorType &)>> &J_inverse,
199  double t,
200  double delta_t,
201  VectorType & y) = 0;
202 
206  struct Status
207  {};
208 
212  virtual const Status &
213  get_status() const = 0;
214  };
215 
216 
217 
221  template <typename VectorType>
222  class RungeKutta : public TimeStepping<VectorType>
223  {
224  public:
228  virtual ~RungeKutta() override = default;
229 
233  virtual void
234  initialize(const runge_kutta_method method) = 0;
235 
247  double
249  std::vector<std::function<VectorType(const double, const VectorType &)>>
250  & F,
251  std::vector<std::function<
252  VectorType(const double, const double, const VectorType &)>> &J_inverse,
253  double t,
254  double delta_t,
255  VectorType &y) override;
256 
268  virtual double
270  const std::function<VectorType(const double, const VectorType &)> &f,
271  const std::function<
272  VectorType(const double, const double, const VectorType &)>
273  & id_minus_tau_J_inverse,
274  double t,
275  double delta_t,
276  VectorType &y) = 0;
277 
278  protected:
282  unsigned int n_stages;
283 
287  std::vector<double> b;
288 
292  std::vector<double> c;
293 
297  std::vector<std::vector<double>> a;
298  };
299 
300 
301 
306  template <typename VectorType>
307  class ExplicitRungeKutta : public RungeKutta<VectorType>
308  {
309  public:
311 
317  ExplicitRungeKutta() = default;
318 
323 
327  void
328  initialize(const runge_kutta_method method) override;
329 
341  double
343  const std::function<VectorType(const double, const VectorType &)> &f,
344  const std::function<
345  VectorType(const double, const double, const VectorType &)>
346  & id_minus_tau_J_inverse,
347  double t,
348  double delta_t,
349  VectorType &y) override;
350 
358  double
360  const std::function<VectorType(const double, const VectorType &)> &f,
361  double t,
362  double delta_t,
363  VectorType &y);
364 
368  struct Status : public TimeStepping<VectorType>::Status
369  {
371  : method(invalid)
372  {}
373 
375  };
376 
380  const Status &
381  get_status() const override;
382 
383  private:
387  void
389  const std::function<VectorType(const double, const VectorType &)> &f,
390  const double t,
391  const double delta_t,
392  const VectorType & y,
393  std::vector<VectorType> &f_stages) const;
394 
399  };
400 
401 
402 
408  template <typename VectorType>
409  class LowStorageRungeKutta : public RungeKutta<VectorType>
410  {
411  public:
413 
419  LowStorageRungeKutta() = default;
420 
425 
429  void
430  initialize(const runge_kutta_method method) override;
431 
443  double
445  const std::function<VectorType(const double, const VectorType &)> &f,
446  const std::function<
447  VectorType(const double, const double, const VectorType &)>
448  & id_minus_tau_J_inverse,
449  double t,
450  double delta_t,
451  VectorType &y) override;
452 
462  double
464  const std::function<VectorType(const double, const VectorType &)> &f,
465  double t,
466  double delta_t,
467  VectorType &solution,
468  VectorType &vec_ri,
469  VectorType &vec_ki);
470 
477  void
478  get_coefficients(std::vector<double> &a,
479  std::vector<double> &b,
480  std::vector<double> &c) const;
481 
485  struct Status : public TimeStepping<VectorType>::Status
486  {
488  : method(invalid)
489  {}
490 
492  };
493 
497  const Status &
498  get_status() const override;
499 
500  private:
504  void
506  const std::function<VectorType(const double, const VectorType &)> &f,
507  const double t,
508  const double factor_solution,
509  const double factor_ai,
510  const VectorType &corrent_ri,
511  VectorType & vec_ki,
512  VectorType & solution,
513  VectorType & next_ri) const;
514 
519  };
520 
521 
522 
527  template <typename VectorType>
528  class ImplicitRungeKutta : public RungeKutta<VectorType>
529  {
530  public:
532 
538  ImplicitRungeKutta() = default;
539 
546  const unsigned int max_it = 100,
547  const double tolerance = 1e-6);
548 
552  void
553  initialize(const runge_kutta_method method) override;
554 
566  double
568  const std::function<VectorType(const double, const VectorType &)> &f,
569  const std::function<
570  VectorType(const double, const double, const VectorType &)>
571  & id_minus_tau_J_inverse,
572  double t,
573  double delta_t,
574  VectorType &y) override;
575 
580  void
582  const double tolerance);
583 
588  struct Status : public TimeStepping<VectorType>::Status
589  {
591  : method(invalid)
594  {}
595 
597  unsigned int n_iterations;
599  };
600 
604  const Status &
605  get_status() const override;
606 
607  private:
611  void
613  const std::function<VectorType(const double, const VectorType &)> &f,
614  const std::function<
615  VectorType(const double, const double, const VectorType &)>
616  & id_minus_tau_J_inverse,
617  double t,
618  double delta_t,
619  VectorType & y,
620  std::vector<VectorType> &f_stages);
621 
625  void
627  const std::function<void(const VectorType &, VectorType &)> &get_residual,
628  const std::function<VectorType(const VectorType &)>
629  & id_minus_tau_J_inverse,
630  VectorType &y);
631 
635  void
637  const std::function<VectorType(const double, const VectorType &)> &f,
638  double t,
639  double delta_t,
640  const VectorType &new_y,
641  const VectorType &y,
642  VectorType & tendency,
643  VectorType & residual) const;
644 
651 
655  unsigned int max_it;
656 
660  double tolerance;
661 
666  };
667 
668 
669 
674  template <typename VectorType>
675  class EmbeddedExplicitRungeKutta : public RungeKutta<VectorType>
676  {
677  public:
679 
686 
692  const double coarsen_param = 1.2,
693  const double refine_param = 0.8,
694  const double min_delta = 1e-14,
695  const double max_delta = 1e100,
696  const double refine_tol = 1e-8,
697  const double coarsen_tol = 1e-12);
698 
703  {
704  free_memory();
705  }
706 
710  void
712 
716  void
717  initialize(const runge_kutta_method method) override;
718 
730  double
732  const std::function<VectorType(const double, const VectorType &)> &f,
733  const std::function<
734  VectorType(const double, const double, const VectorType &)>
735  & id_minus_tau_J_inverse,
736  double t,
737  double delta_t,
738  VectorType &y) override;
739 
747  double
749  const std::function<VectorType(const double, const VectorType &)> &f,
750  double t,
751  double delta_t,
752  VectorType &y);
753 
757  void
759  const double refine_param,
760  const double min_delta,
761  const double max_delta,
762  const double refine_tol,
763  const double coarsen_tol);
764 
771  struct Status : public TimeStepping<VectorType>::Status
772  {
775  unsigned int n_iterations;
777  double error_norm;
778  };
779 
783  const Status &
784  get_status() const override;
785 
786  private:
790  void
792  const std::function<VectorType(const double, const VectorType &)> &f,
793  const double t,
794  const double delta_t,
795  const VectorType & y,
796  std::vector<VectorType> &f_stages);
797 
803 
808  double refine_param;
809 
813  double min_delta_t;
814 
818  double max_delta_t;
819 
824  double refine_tol;
825 
830  double coarsen_tol;
831 
837 
841  std::vector<double> b1;
842 
846  std::vector<double> b2;
847 
852  VectorType *last_stage;
853 
858  };
859 } // namespace TimeStepping
860 
862 
863 #endif
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)
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
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)
void initialize(const runge_kutta_method method) override
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_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(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
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)
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
const Status & get_status() const override
std::vector< double > b
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
std::vector< double > c
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
Definition: config.h:396
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
embedded_runge_kutta_time_step
@ RK_CLASSIC_FOURTH_ORDER
Definition: time_stepping.h:79
@ LOW_STORAGE_RK_STAGE9_ORDER5
@ LOW_STORAGE_RK_STAGE3_ORDER3
Definition: time_stepping.h:86
@ LOW_STORAGE_RK_STAGE7_ORDER4
Definition: time_stepping.h:96
@ LOW_STORAGE_RK_STAGE5_ORDER4
Definition: time_stepping.h:91
static const unsigned int invalid_unsigned_int
Definition: types.h:196
T signaling_nan()
embedded_runge_kutta_time_step exit_delta_t