/usr/include/deal.II/base/time_stepping.h is in libdeal.ii-dev 8.4.2-2+b1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 | // ---------------------------------------------------------------------
//
// Copyright (C) 2014 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------
#ifndef dealii__time_stepping_h
#define dealii__time_stepping_h
#include <deal.II/base/config.h>
#include <deal.II/base/std_cxx11/function.h>
#include <vector>
DEAL_II_NAMESPACE_OPEN
/**
* Namespace containing the time stepping methods.
*
* @author Bruno Turcksin
* @date 2014
*/
namespace TimeStepping
{
/**
* The following Runge-Kutta methods are available:
* - Explicit methods (see ExplicitRungeKutta::initialize):
* - FORWARD_EULER (first order)
* - RK_THIRD_ORDER (third order Runge-Kutta)
* - RK_CLASSIC_FOURTH_ORDER (classical fourth order Runge-Kutta)
* - Implicit methods (see ImplicitRungeKutta::initialize):
* - BACKWARD_EULER (first order)
* - IMPLICIT_MIDPOINT (second order)
* - CRANK_NICOLSON (second order)
* - SDIRK_TWO_STAGES (second order)
* - Embedded explicit methods (see EmbeddedExplicitRungeKutta::initialize):
* - HEUN_EULER (second order)
* - BOGACKI_SHAMPINE (third order)
* - DOPRI: Dormand-Prince (fifth order, method used by ode45 in
* MATLAB)
* - FEHLBERG (fifth order)
* - CASH_KARP (firth order)
*/
enum runge_kutta_method { FORWARD_EULER, RK_THIRD_ORDER, RK_CLASSIC_FOURTH_ORDER,
BACKWARD_EULER, IMPLICIT_MIDPOINT, CRANK_NICOLSON,
SDIRK_TWO_STAGES, HEUN_EULER, BOGACKI_SHAMPINE, DOPRI,
FEHLBERG, CASH_KARP
};
/**
* Reason for exiting evolve_one_time_step when using an embedded method:
* DELTA_T (the time step is in the valid range), MIN_DELTA_T (the time step
* was increased to the minimum acceptable time step), MAX_DELTA_T (the time
* step was reduced to the maximum acceptable time step).
*/
enum embedded_runge_kutta_time_step { DELTA_T, MIN_DELTA_T, MAX_DELTA_T };
/**
* Abstract class for time stepping methods. These methods assume that the
* equation has the form: $ \frac{\partial y}{\partial t} = f(t,y) $.
*/
template <typename VectorType>
class TimeStepping
{
public:
/**
* Virtual destructor.
*/
virtual ~TimeStepping() {}
/**
* Purely virtual function. This function is used to advance from time @p
* t to t+ @p delta_t. @p 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. @p 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.
*/
virtual double evolve_one_time_step
(std::vector<std_cxx11::function<VectorType (const double, const VectorType &)> > &F,
std::vector<std_cxx11::function<VectorType (const double, const double, const VectorType &)> > &J_inverse,
double t,
double delta_t,
VectorType &y) = 0;
/**
* Empty structure used to store information.
*/
struct Status {};
/**
* Purely virtual function that return Status.
*/
virtual const Status &get_status() const = 0;
};
/**
* Base class for the Runge-Kutta method
*
* @author Damien Lebrun-Grandie, Bruno Turcksin
* @date 2014
*/
template <typename VectorType>
class RungeKutta : public TimeStepping<VectorType>
{
public:
/**
* Virtual destructor.
*/
virtual ~RungeKutta() {}
/**
* Purely virtual method used to initialize the Runge-Kutta method.
*/
virtual void initialize(runge_kutta_method method) = 0;
/**
* This function is used to advance from time @p t to t+ @p delta_t. @p 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. @p 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, @p F and @ J_inverse can only contain one element.
*/
double evolve_one_time_step
(std::vector<std_cxx11::function<VectorType (const double, const VectorType &)> > &F,
std::vector<std_cxx11::function<VectorType (const double, const double, const VectorType &)> > &J_inverse,
double t,
double delta_t,
VectorType &y);
/**
* Purely virtual function. This function is used to advance from time @p
* t to t+ @p delta_t. @p 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. @p 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
* J}{\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.
*/
virtual double evolve_one_time_step
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
std_cxx11::function<VectorType (const double, const double, const VectorType &)> id_minus_tau_J_inverse,
double t,
double delta_t,
VectorType &y) = 0;
protected:
/**
* Number of stages of the Runge-Kutta method.
*/
unsigned int n_stages;
/**
* Butcher tableau coefficients.
*/
std::vector<double> b;
/**
* Butcher tableau coefficients.
*/
std::vector<double> c;
/**
* Butcher tableau coefficients.
*/
std::vector<std::vector<double> > a;
};
/**
* ExplicitRungeKutta is derived from RungeKutta and implement the explicit
* methods.
*/
template <typename VectorType>
class ExplicitRungeKutta : public RungeKutta<VectorType>
{
public:
using RungeKutta<VectorType>::evolve_one_time_step;
/**
* Default constructor. initialize(runge_kutta_method) needs to be called
* before the object can be used.
*/
ExplicitRungeKutta() {}
/**
* Constructor. This function calls initialize(runge_kutta_method).
*/
ExplicitRungeKutta(runge_kutta_method method);
/**
* Initialize the explicit Runge-Kutta method.
*/
void initialize(runge_kutta_method method);
/**
* This function is used to advance from time @p t to t+ @p delta_t. @p 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. @p 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 J}{\partial y} $. The input
* parameter 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.
*/
double evolve_one_time_step
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
std_cxx11::function<VectorType (const double, const double, const VectorType &)> id_minus_tau_J_inverse,
double t,
double delta_t,
VectorType &y);
/**
* This function is used to advance from time @p t to t+ @p delta_t. This
* function is similar to the one derived from RungeKutta, 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.
*/
double evolve_one_time_step
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
double t,
double delta_t,
VectorType &y);
/**
* This structure stores the name of the method used.
*/
struct Status : public TimeStepping<VectorType>::Status
{
runge_kutta_method method;
};
/**
* Return the status of the current object.
*/
const Status &get_status() const;
private:
/**
* Compute the different stages needed.
*/
void compute_stages
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
const double t,
const double delta_t,
const VectorType &y,
std::vector<VectorType> &f_stages) const;
/**
* Status structure of the object.
*/
Status status;
};
/**
* This class is derived from RungeKutta and implement the implicit methods.
* This class works only for Diagonal Implicit Runge-Kutta (DIRK) methods.
*/
template <typename VectorType>
class ImplicitRungeKutta : public RungeKutta<VectorType>
{
public:
using RungeKutta<VectorType>::evolve_one_time_step;
/**
* Default constructor. initialize(runge_kutta_method) and
* set_newton_solver_parameters(unsigned int,double) need to be called
* before the object can be used.
*/
ImplicitRungeKutta() {}
/**
* Constructor. This function calls initialize(runge_kutta_method) and
* initialize the maximum number of iterations and the tolerance of the
* Newton solver.
*/
ImplicitRungeKutta(runge_kutta_method method, unsigned int max_it=100, double tolerance=1e-6);
/**
* Initialize the implicit Runge-Kutta method.
*/
void initialize(runge_kutta_method method);
/**
* This function is used to advance from time @p t to t+ @p delta_t. @p 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. @p id_minus_tau_J_inverse is a function that computes $
* (I-\tau J)^{-1}$ where $ I $ is the identity matrix, $ \tau $ is given,
* and $ J $ is the Jacobian $ \frac{\partial J}{\partial y} $. The input
* parameters this function receives 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.
*/
double evolve_one_time_step
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
std_cxx11::function<VectorType (const double, const double, const VectorType &)> id_minus_tau_J_inverse,
double t,
double delta_t,
VectorType &y);
/**
* Set the maximum number of iterations and the tolerance used by the
* Newton solver.
*/
void set_newton_solver_parameters(unsigned int max_it, double tolerance);
/**
* Structure that stores the name of the method, the number of Newton
* iterations and the norm of the residual when exiting the Newton solver.
*/
struct Status : public TimeStepping<VectorType>::Status
{
runge_kutta_method method;
unsigned int n_iterations;
double norm_residual;
};
/**
* Return the status of the current object.
*/
const Status &get_status() const;
private:
/**
* Compute the different stages needed.
*/
void compute_stages
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
std_cxx11::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);
/**
* Newton solver used for the implicit stages.
*/
void newton_solve(std_cxx11::function<void (const VectorType &,VectorType &)> get_residual,
std_cxx11::function<VectorType (const VectorType &)> id_minus_tau_J_inverse,
VectorType &y);
/**
* Compute the residual needed by the Newton solver.
*/
void compute_residual(std_cxx11::function<VectorType (const double, const VectorType &)> f,
double t,
double delta_t,
const VectorType &old_y,
const VectorType &y,
VectorType &tendency,
VectorType &residual) const;
/**
* When using SDIRK, there is no need to compute the linear combination of
* the stages. Thus, when this flag is true, the linear combination is
* skipped.
*/
bool skip_linear_combi;
/**
* Maximum number of iterations of the Newton solver.
*/
unsigned int max_it;
/**
* Tolerance of the Newton solver.
*/
double tolerance;
/**
* Status structure of the object.
*/
Status status;
};
/**
* This is class is derived from RungeKutta and implement embedded explicit
* methods.
*/
template <typename VectorType>
class EmbeddedExplicitRungeKutta : public RungeKutta<VectorType>
{
public:
using RungeKutta<VectorType>::evolve_one_time_step;
/**
* 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() {}
/**
* Constructor. This function calls initialize(runge_kutta_method) and
* initialize the parameters needed for time adaptation.
*/
EmbeddedExplicitRungeKutta(runge_kutta_method method,
double coarsen_param = 1.2,
double refine_param = 0.8,
double min_delta = 1e-14,
double max_delta = 1e100,
double refine_tol = 1e-8,
double coarsen_tol = 1e-12);
/**
* Destructor.
*/
~EmbeddedExplicitRungeKutta()
{
free_memory();
}
/**
* If necessary, deallocate memory allocated by the object.
*/
void free_memory();
/**
* Initialize the embedded explicit Runge-Kutta method.
*/
void initialize(runge_kutta_method method);
/**
* This function is used to advance from time @p t to t+ @p delta_t. @p 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. @p 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 J}{\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.
*/
double evolve_one_time_step
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
std_cxx11::function<VectorType (const double, const double, const VectorType &)> id_minus_tau_J_inverse,
double t,
double delta_t,
VectorType &y);
/**
* This function is used to advance from time @p t to t+ @p 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.
*/
double evolve_one_time_step
(std_cxx11::function<VectorType (const double, const VectorType &)> f,
double t,
double delta_t,
VectorType &y);
/**
* Set the parameters necessary for the time adaptation.
*/
void set_time_adaptation_parameters(double coarsen_param,
double refine_param,
double min_delta,
double max_delta,
double refine_tol,
double coarsen_tol);
/**
* Structure that stores the name of the method, the reason to exit
* evolve_one_time_step, the number of iteration inside n_iterations, a
* guess of what the next time step should be, and an estimate of the norm
* of the error.
*/
struct Status : public TimeStepping<VectorType>::Status
{
runge_kutta_method method;
embedded_runge_kutta_time_step exit_delta_t;
unsigned int n_iterations;
double delta_t_guess;
double error_norm;
};
/**
* Return the status of the current object.
*/
const Status &get_status() const;
private:
/**
* Compute the different stages needed.
*/
void compute_stages(std_cxx11::function<VectorType (const double, const VectorType &)> f,
const double t,
const double delta_t,
const VectorType &y,
std::vector<VectorType> &f_stages);
/**
* This parameter is the factor (>1) by which the time step is multiplied
* when the time stepping can be coarsen.
*/
double coarsen_param;
/**
* This parameter is the factor (<1) by which the time step is multiplied
* when the time stepping must be refined.
*/
double refine_param;
/**
* Smallest time step allowed.
*/
double min_delta_t;
/**
* Largest time step allowed.
*/
double max_delta_t;
/**
* Refinement tolerance: if the error estimate is larger than refine_tol,
* the time step is refined.
*/
double refine_tol;
/**
* Coarsening tolerance: if the error estimate is smaller than coarse_tol,
* the time step is coarsen.
*/
double coarsen_tol;
/**
* If the flag is true, the last stage is the same as the first stage and
* one evaluation of f can be saved.
*/
bool last_same_as_first;
/**
* Butcher tableau coefficients.
*/
std::vector<double> b1;
/**
* Butcher tableau coefficients.
*/
std::vector<double> b2;
/**
* 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.
*/
VectorType *last_stage;
/**
* Status structure of the object.
*/
Status status;
};
}
DEAL_II_NAMESPACE_CLOSE
#endif
|