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
time_stepping.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2014 - 2023 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
33namespace 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
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:
310 using RungeKutta<VectorType>::evolve_one_time_step;
311
318
323
327 void
328 initialize(const runge_kutta_method method) override;
329
343 double
345 const std::function<VectorType(const double, const VectorType &)> &f,
346 const std::function<
347 VectorType(const double, const double, const VectorType &)>
348 & id_minus_tau_J_inverse,
349 double t,
350 double delta_t,
351 VectorType &y) override;
352
360 double
362 const std::function<VectorType(const double, const VectorType &)> &f,
363 double t,
364 double delta_t,
365 VectorType &y);
366
370 struct Status : public TimeStepping<VectorType>::Status
371 {
373 : method(invalid)
374 {}
375
377 };
378
382 const Status &
383 get_status() const override;
384
385 private:
389 void
391 const std::function<VectorType(const double, const VectorType &)> &f,
392 const double t,
393 const double delta_t,
394 const VectorType & y,
395 std::vector<VectorType> &f_stages) const;
396
401 };
402
403
404
410 template <typename VectorType>
411 class LowStorageRungeKutta : public RungeKutta<VectorType>
412 {
413 public:
414 using RungeKutta<VectorType>::evolve_one_time_step;
415
422
427
431 void
432 initialize(const runge_kutta_method method) override;
433
445 double
447 const std::function<VectorType(const double, const VectorType &)> &f,
448 const std::function<
449 VectorType(const double, const double, const VectorType &)>
450 & id_minus_tau_J_inverse,
451 double t,
452 double delta_t,
453 VectorType &y) override;
454
464 double
466 const std::function<VectorType(const double, const VectorType &)> &f,
467 double t,
468 double delta_t,
469 VectorType &solution,
470 VectorType &vec_ri,
471 VectorType &vec_ki);
472
479 void
480 get_coefficients(std::vector<double> &a,
481 std::vector<double> &b,
482 std::vector<double> &c) const;
483
487 struct Status : public TimeStepping<VectorType>::Status
488 {
490 : method(invalid)
491 {}
492
494 };
495
499 const Status &
500 get_status() const override;
501
502 private:
506 void
508 const std::function<VectorType(const double, const VectorType &)> &f,
509 const double t,
510 const double factor_solution,
511 const double factor_ai,
512 const VectorType &corrent_ri,
513 VectorType & vec_ki,
514 VectorType & solution,
515 VectorType & next_ri) const;
516
521 };
522
523
524
529 template <typename VectorType>
530 class ImplicitRungeKutta : public RungeKutta<VectorType>
531 {
532 public:
533 using RungeKutta<VectorType>::evolve_one_time_step;
534
541
548 const unsigned int max_it = 100,
549 const double tolerance = 1e-6);
550
554 void
555 initialize(const runge_kutta_method method) override;
556
568 double
570 const std::function<VectorType(const double, const VectorType &)> &f,
571 const std::function<
572 VectorType(const double, const double, const VectorType &)>
573 & id_minus_tau_J_inverse,
574 double t,
575 double delta_t,
576 VectorType &y) override;
577
582 void
584 const double tolerance);
585
590 struct Status : public TimeStepping<VectorType>::Status
591 {
593 : method(invalid)
594 , n_iterations(numbers::invalid_unsigned_int)
595 , norm_residual(numbers::signaling_nan<double>())
596 {}
597
599 unsigned int n_iterations;
601 };
602
606 const Status &
607 get_status() const override;
608
609 private:
613 void
615 const std::function<VectorType(const double, const VectorType &)> &f,
616 const std::function<
617 VectorType(const double, const double, const VectorType &)>
618 & id_minus_tau_J_inverse,
619 double t,
620 double delta_t,
621 VectorType & y,
622 std::vector<VectorType> &f_stages);
623
627 void
629 const std::function<void(const VectorType &, VectorType &)> &get_residual,
630 const std::function<VectorType(const VectorType &)>
631 & id_minus_tau_J_inverse,
632 VectorType &y);
633
637 void
639 const std::function<VectorType(const double, const VectorType &)> &f,
640 double t,
641 double delta_t,
642 const VectorType &new_y,
643 const VectorType &y,
644 VectorType & tendency,
645 VectorType & residual) const;
646
650 unsigned int max_it;
651
655 double tolerance;
656
661 };
662
663
664
669 template <typename VectorType>
670 class EmbeddedExplicitRungeKutta : public RungeKutta<VectorType>
671 {
672 public:
673 using RungeKutta<VectorType>::evolve_one_time_step;
674
681
687 const double coarsen_param = 1.2,
688 const double refine_param = 0.8,
689 const double min_delta = 1e-14,
690 const double max_delta = 1e100,
691 const double refine_tol = 1e-8,
692 const double coarsen_tol = 1e-12);
693
698 {
699 free_memory();
700 }
701
705 void
707
711 void
712 initialize(const runge_kutta_method method) override;
713
727 double
729 const std::function<VectorType(const double, const VectorType &)> &f,
730 const std::function<
731 VectorType(const double, const double, const VectorType &)>
732 & id_minus_tau_J_inverse,
733 double t,
734 double delta_t,
735 VectorType &y) override;
736
744 double
746 const std::function<VectorType(const double, const VectorType &)> &f,
747 double t,
748 double delta_t,
749 VectorType &y);
750
754 void
756 const double refine_param,
757 const double min_delta,
758 const double max_delta,
759 const double refine_tol,
760 const double coarsen_tol);
761
768 struct Status : public TimeStepping<VectorType>::Status
769 {
771 : method(invalid)
772 {}
773
776 unsigned int n_iterations;
779 };
780
784 const Status &
785 get_status() const override;
786
787 private:
791 void
793 const std::function<VectorType(const double, const VectorType &)> &f,
794 const double t,
795 const double delta_t,
796 const VectorType & y,
797 std::vector<VectorType> &f_stages);
798
804
810
815
820
826
832
837 bool last_same_as_first = false;
838
842 std::vector<double> b1;
843
847 std::vector<double> b2;
848
853 VectorType *last_stage = nullptr;
854
859 };
860} // namespace TimeStepping
861
863
864#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)
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(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)
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
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:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
@ 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