MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ode.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_ODE
13 #define MFEM_ODE
14 
15 #include "../config/config.hpp"
16 #include "operator.hpp"
17 
18 namespace mfem
19 {
20 
21 /// Abstract class for solving systems of ODEs: dx/dt = f(x,t)
22 class ODESolver
23 {
24 protected:
25  /// Pointer to the associated TimeDependentOperator.
26  TimeDependentOperator *f; // f(.,t) : R^n --> R^n
28 
29 public:
31 
32  /// Associate a TimeDependentOperator with the ODE solver.
33  /** This method has to be called:
34  - Before the first call to Step().
35  - When the dimensions of the associated TimeDependentOperator change.
36  - When a time stepping sequence has to be restarted.
37  - To change the associated TimeDependentOperator. */
38  virtual void Init(TimeDependentOperator &f);
39 
40  /** @brief Perform a time step from time @a t [in] to time @a t [out] based
41  on the requested step size @a dt [in]. */
42  /** @param[in,out] x Approximate solution.
43  @param[in,out] t Time associated with the approximate solution @a x.
44  @param[in,out] dt Time step size.
45 
46  The following rules describe the common behavior of the method:
47  - The input @a x [in] is the approximate solution for the input time
48  @a t [in].
49  - The input @a dt [in] is the desired time step size, defining the desired
50  target time: t [target] = @a t [in] + @a dt [in].
51  - The output @a x [out] is the approximate solution for the output time
52  @a t [out].
53  - The output @a dt [out] is the last time step taken by the method which
54  may be smaller or larger than the input @a dt [in] value, e.g. because
55  of time step control.
56  - The method may perform more than one time step internally; in this case
57  @a dt [out] is the last internal time step size.
58  - The output value of @a t [out] may be smaller or larger than
59  t [target], however, it is not smaller than @a t [in] + @a dt [out], if
60  at least one internal time step was performed.
61  - The value @a x [out] may be obtained by interpolation using internally
62  stored data.
63  - In some cases, the contents of @a x [in] may not be used, e.g. when
64  @a x [out] from a previous Step() call was obtained by interpolation.
65  - In consecutive calls to this method, the output @a t [out] of one
66  Step() call has to be the same as the input @a t [in] to the next
67  Step() call.
68  - If the previous rule has to be broken, e.g. to restart a time stepping
69  sequence, then the ODE solver must be re-initialized by calling Init()
70  between the two Step() calls. */
71  virtual void Step(Vector &x, double &t, double &dt) = 0;
72 
73  /// Perform time integration from time @a t [in] to time @a tf [in].
74  /** @param[in,out] x Approximate solution.
75  @param[in,out] t Time associated with the approximate solution @a x.
76  @param[in,out] dt Time step size.
77  @param[in] tf Requested final time.
78 
79  The default implementation makes consecutive calls to Step() until
80  reaching @a tf.
81  The following rules describe the common behavior of the method:
82  - The input @a x [in] is the approximate solution for the input time
83  @a t [in].
84  - The input @a dt [in] is the initial time step size.
85  - The output @a dt [out] is the last time step taken by the method which
86  may be smaller or larger than the input @a dt [in] value, e.g. because
87  of time step control.
88  - The output value of @a t [out] is not smaller than @a tf [in]. */
89  virtual void Run(Vector &x, double &t, double &dt, double tf)
90  {
91  while (t < tf) { Step(x, t, dt); }
92  }
93 
94  virtual ~ODESolver() { }
95 };
96 
97 
98 /// The classical forward Euler method
100 {
101 private:
102  Vector dxdt;
103 
104 public:
105  virtual void Init(TimeDependentOperator &_f);
106 
107  virtual void Step(Vector &x, double &t, double &dt);
108 };
109 
110 
111 /** A family of explicit second-order RK2 methods. Some choices for the
112  parameter 'a' are:
113  a = 1/2 - the midpoint method
114  a = 1 - Heun's method
115  a = 2/3 - default, has minimal truncation error. */
116 class RK2Solver : public ODESolver
117 {
118 private:
119  double a;
120  Vector dxdt, x1;
121 
122 public:
123  RK2Solver(const double _a = 2./3.) : a(_a) { }
124 
125  virtual void Init(TimeDependentOperator &_f);
126 
127  virtual void Step(Vector &x, double &t, double &dt);
128 };
129 
130 
131 /// Third-order, strong stability preserving (SSP) Runge-Kutta method
132 class RK3SSPSolver : public ODESolver
133 {
134 private:
135  Vector y, k;
136 
137 public:
138  virtual void Init(TimeDependentOperator &_f);
139 
140  virtual void Step(Vector &x, double &t, double &dt);
141 };
142 
143 
144 /// The classical explicit forth-order Runge-Kutta method, RK4
145 class RK4Solver : public ODESolver
146 {
147 private:
148  Vector y, k, z;
149 
150 public:
151  virtual void Init(TimeDependentOperator &_f);
152 
153  virtual void Step(Vector &x, double &t, double &dt);
154 };
155 
156 
157 /** An explicit Runge-Kutta method corresponding to a general Butcher tableau
158  +--------+----------------------+
159  | c[0] | a[0] |
160  | c[1] | a[1] a[2] |
161  | ... | ... |
162  | c[s-2] | ... a[s(s-1)/2-1] |
163  +--------+----------------------+
164  | | b[0] b[1] ... b[s-1] |
165  +--------+----------------------+ */
167 {
168 private:
169  int s;
170  const double *a, *b, *c;
171  Vector y, *k;
172 
173 public:
174  ExplicitRKSolver(int _s, const double *_a, const double *_b,
175  const double *_c);
176 
177  virtual void Init(TimeDependentOperator &_f);
178 
179  virtual void Step(Vector &x, double &t, double &dt);
180 
181  virtual ~ExplicitRKSolver();
182 };
183 
184 
185 /** An 8-stage, 6th order RK method. From Verner's "efficient" 9-stage 6(5)
186  pair. */
188 {
189 private:
190  static const double a[28], b[8], c[7];
191 
192 public:
193  RK6Solver() : ExplicitRKSolver(8, a, b, c) { }
194 };
195 
196 
197 /** A 12-stage, 8th order RK method. From Verner's "efficient" 13-stage 8(7)
198  pair. */
200 {
201 private:
202  static const double a[66], b[12], c[11];
203 
204 public:
205  RK8Solver() : ExplicitRKSolver(12, a, b, c) { }
206 };
207 
208 
209 /** An explicit Adams-Bashforth method. */
211 {
212 private:
213  int s, smax;
214  const double *a;
215  Vector *k;
216  Array<int> idx;
217  ODESolver *RKsolver;
218 
219 public:
220  AdamsBashforthSolver(int _s, const double *_a);
221 
222  virtual void Init(TimeDependentOperator &_f);
223 
224  virtual void Step(Vector &x, double &t, double &dt);
225 
227  {
228  if (RKsolver) { delete RKsolver; }
229  }
230 };
231 
232 
233 /** A 1-stage, 1st order AB method. */
235 {
236 private:
237  static const double a[1];
238 
239 public:
241 };
242 
243 /** A 2-stage, 2st order AB method. */
245 {
246 private:
247  static const double a[2];
248 
249 public:
251 };
252 
253 /** A 3-stage, 3st order AB method. */
255 {
256 private:
257  static const double a[3];
258 
259 public:
261 };
262 
263 /** A 4-stage, 4st order AB method. */
265 {
266 private:
267  static const double a[4];
268 
269 public:
271 };
272 
273 /** A 5-stage, 5st order AB method. */
275 {
276 private:
277  static const double a[5];
278 
279 public:
281 };
282 
283 
284 /** An implicit Adams-Moulton method. */
286 {
287 private:
288  int s, smax;
289  const double *a;
290  Vector *k;
291  Array<int> idx;
292  ODESolver *RKsolver;
293 
294 public:
295  AdamsMoultonSolver(int _s, const double *_a);
296 
297  virtual void Init(TimeDependentOperator &_f);
298 
299  virtual void Step(Vector &x, double &t, double &dt);
300 
302  {
303  if (RKsolver) { delete RKsolver; }
304  };
305 };
306 
307 /** A 0-stage, 1st order AM method. */
309 {
310 private:
311  static const double a[1];
312 
313 public:
315 };
316 
317 
318 /** A 1-stage, 2st order AM method. */
320 {
321 private:
322  static const double a[2];
323 
324 public:
326 };
327 
328 /** A 2-stage, 3st order AM method. */
330 {
331 private:
332  static const double a[3];
333 
334 public:
336 };
337 
338 /** A 3-stage, 4st order AM method. */
340 {
341 private:
342  static const double a[4];
343 
344 public:
346 };
347 
348 /** A 4-stage, 5st order AM method. */
350 {
351 private:
352  static const double a[5];
353 
354 public:
356 };
357 
358 
359 /// Backward Euler ODE solver. L-stable.
361 {
362 protected:
364 
365 public:
366  virtual void Init(TimeDependentOperator &_f);
367 
368  virtual void Step(Vector &x, double &t, double &dt);
369 };
370 
371 
372 /// Implicit midpoint method. A-stable, not L-stable.
374 {
375 protected:
377 
378 public:
379  virtual void Init(TimeDependentOperator &_f);
380 
381  virtual void Step(Vector &x, double &t, double &dt);
382 };
383 
384 
385 /** Two stage, singly diagonal implicit Runge-Kutta (SDIRK) methods;
386  the choices for gamma_opt are:
387  0 - 3rd order method, not A-stable
388  1 - 3rd order method, A-stable, not L-stable (default)
389  2 - 2nd order method, L-stable
390  3 - 2nd order method, L-stable (has solves outside [t,t+dt]). */
391 class SDIRK23Solver : public ODESolver
392 {
393 protected:
394  double gamma;
396 
397 public:
398  SDIRK23Solver(int gamma_opt = 1);
399 
400  virtual void Init(TimeDependentOperator &_f);
401 
402  virtual void Step(Vector &x, double &t, double &dt);
403 };
404 
405 
406 /** Three stage, singly diagonal implicit Runge-Kutta (SDIRK) method of
407  order 4. A-stable, not L-stable. */
408 class SDIRK34Solver : public ODESolver
409 {
410 protected:
411  Vector k, y, z;
412 
413 public:
414  virtual void Init(TimeDependentOperator &_f);
415 
416  virtual void Step(Vector &x, double &t, double &dt);
417 };
418 
419 
420 /** Three stage, singly diagonal implicit Runge-Kutta (SDIRK) method of
421  order 3. L-stable. */
422 class SDIRK33Solver : public ODESolver
423 {
424 protected:
426 
427 public:
428  virtual void Init(TimeDependentOperator &_f);
429 
430  virtual void Step(Vector &x, double &t, double &dt);
431 };
432 
433 
434 /// Generalized-alpha ODE solver from "A generalized-α method for integrating
435 /// the filtered Navier–Stokes equations with a stabilized finite element
436 /// method" by K.E. Jansen, C.H. Whiting and G.M. Hulbert.
438 {
439 protected:
440  mutable Vector xdot,k,y;
442  bool first;
443 
444  void SetRhoInf(double rho_inf);
445  void PrintProperties(std::ostream &out = mfem::out);
446 public:
447 
448  GeneralizedAlphaSolver(double rho = 1.0) { SetRhoInf(rho); };
449 
450  virtual void Init(TimeDependentOperator &_f);
451 
452  virtual void Step(Vector &x, double &t, double &dt);
453 };
454 
455 
456 /// The SIASolver class is based on the Symplectic Integration Algorithm
457 /// described in "A Symplectic Integration Algorithm for Separable Hamiltonian
458 /// Functions" by J. Candy and W. Rozmus, Journal of Computational Physics,
459 /// Vol. 92, pages 230-256 (1991).
460 
461 /** The Symplectic Integration Algorithm (SIA) is designed for systems of first
462  order ODEs derived from a Hamiltonian.
463  H(q,p,t) = T(p) + V(q,t)
464  Which leads to the equations:
465  dq/dt = dT/dp
466  dp/dt = -dV/dq
467  In the integrator the operators P and F are defined to be:
468  P = dT/dp
469  F = -dV/dq
470  */
472 {
473 public:
474  SIASolver() : F_(NULL), P_(NULL) {}
475 
476  virtual void Init(Operator &P, TimeDependentOperator & F);
477 
478  virtual void Step(Vector &q, Vector &p, double &t, double &dt) = 0;
479 
480  virtual void Run(Vector &q, Vector &p, double &t, double &dt, double tf)
481  {
482  while (t < tf) { Step(q, p, t, dt); }
483  }
484 
485  virtual ~SIASolver() {}
486 
487 protected:
488  TimeDependentOperator * F_; // p_{i+1} = p_{i} + dt F(q_{i})
489  Operator * P_; // q_{i+1} = q_{i} + dt P(p_{i+1})
490 
491  mutable Vector dp_;
492  mutable Vector dq_;
493 };
494 
495 // First Order Symplectic Integration Algorithm
496 class SIA1Solver : public SIASolver
497 {
498 public:
500  void Step(Vector &q, Vector &p, double &t, double &dt);
501 };
502 
503 // Second Order Symplectic Integration Algorithm
504 class SIA2Solver : public SIASolver
505 {
506 public:
508  void Step(Vector &q, Vector &p, double &t, double &dt);
509 };
510 
511 // Variable order Symplectic Integration Algorithm (orders 1-4)
512 class SIAVSolver : public SIASolver
513 {
514 public:
515  SIAVSolver(int order);
516  void Step(Vector &q, Vector &p, double &t, double &dt);
517 
518 private:
519  int order_;
520 
521  Array<double> a_;
522  Array<double> b_;
523 };
524 
525 
526 
527 /// Abstract class for solving systems of ODEs: d2x/dt2 = f(x,dx/dt,t)
529 {
530 protected:
531  /// Pointer to the associated TimeDependentOperator.
532  SecondOrderTimeDependentOperator *f; // f(.,.,t) : R^n x R^n --> R^n
534 
535 public:
537 
538  /// Associate a TimeDependentOperator with the ODE solver.
539  /** This method has to be called:
540  - Before the first call to Step().
541  - When the dimensions of the associated TimeDependentOperator change.
542  - When a time stepping sequence has to be restarted.
543  - To change the associated TimeDependentOperator. */
544  virtual void Init(SecondOrderTimeDependentOperator &f);
545 
546  /** @brief Perform a time step from time @a t [in] to time @a t [out] based
547  on the requested step size @a dt [in]. */
548  /** @param[in,out] x Approximate solution.
549  @param[in,out] dxdt Approximate rate.
550  @param[in,out] t Time associated with the
551  approximate solution @a x and rate @ dxdt
552  @param[in,out] dt Time step size.
553 
554  The following rules describe the common behavior of the method:
555  - The input @a x [in] is the approximate solution for the input time
556  @a t [in].
557  - The input @a dxdt [in] is the approximate rate for the input time
558  @a t [in].
559  - The input @a dt [in] is the desired time step size, defining the desired
560  target time: t [target] = @a t [in] + @a dt [in].
561  - The output @a x [out] is the approximate solution for the output time
562  @a t [out].
563  - The output @a dxdt [out] is the approximate rate for the output time
564  @a t [out].
565  - The output @a dt [out] is the last time step taken by the method which
566  may be smaller or larger than the input @a dt [in] value, e.g. because
567  of time step control.
568  - The method may perform more than one time step internally; in this case
569  @a dt [out] is the last internal time step size.
570  - The output value of @a t [out] may be smaller or larger than
571  t [target], however, it is not smaller than @a t [in] + @a dt [out], if
572  at least one internal time step was performed.
573  - The value @a x [out] may be obtained by interpolation using internally
574  stored data.
575  - In some cases, the contents of @a x [in] may not be used, e.g. when
576  @a x [out] from a previous Step() call was obtained by interpolation.
577  - In consecutive calls to this method, the output @a t [out] of one
578  Step() call has to be the same as the input @a t [in] to the next
579  Step() call.
580  - If the previous rule has to be broken, e.g. to restart a time stepping
581  sequence, then the ODE solver must be re-initialized by calling Init()
582  between the two Step() calls. */
583  virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt) = 0;
584 
585  /// Perform time integration from time @a t [in] to time @a tf [in].
586  /** @param[in,out] x Approximate solution.
587  @param[in,out] dxdt Approximate rate.
588  @param[in,out] t Time associated with the approximate solution @a x.
589  @param[in,out] dt Time step size.
590  @param[in] tf Requested final time.
591 
592  The default implementation makes consecutive calls to Step() until
593  reaching @a tf.
594  The following rules describe the common behavior of the method:
595  - The input @a x [in] is the approximate solution for the input time
596  @a t [in].
597  - The input @a dxdt [in] is the approximate rate for the input time
598  @a t [in].
599  - The input @a dt [in] is the initial time step size.
600  - The output @a dt [out] is the last time step taken by the method which
601  may be smaller or larger than the input @a dt [in] value, e.g. because
602  of time step control.
603  - The output value of @a t [out] is not smaller than @a tf [in]. */
604  virtual void Run(Vector &x, Vector &dxdt, double &t, double &dt, double tf)
605  {
606  while (t < tf) { Step(x, dxdt, t, dt); }
607  }
608 
609  virtual ~SecondOrderODESolver() { }
610 };
611 
612 
613 
614 /// The classical newmark method.
615 /// Newmark, N. M. (1959) A method of computation for structural dynamics.
616 /// Journal of Engineering Mechanics, ASCE, 85 (EM3) 67-94.
618 {
619 private:
620  Vector d2xdt2;
621 
622  double beta, gamma;
623  bool first;
624 
625 public:
626  NewmarkSolver(double beta_ = 0.25, double gamma_ = 0.5) { beta = beta_; gamma = gamma_; };
627 
628  virtual void PrintProperties(std::ostream &out = mfem::out);
629 
630  virtual void Init(SecondOrderTimeDependentOperator &_f);
631 
632  virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt);
633 };
634 
636 {
637 public:
639 };
640 
642 {
643 public:
645 };
646 
648 {
649 public:
650  FoxGoodwinSolver() : NewmarkSolver(1.0/12.0, 0.5) { };
651 };
652 
653 /// Generalized-alpha ODE solver
654 /// A Time Integration Algorithm for Structural Dynamics With Improved
655 /// Numerical Dissipation: The Generalized-α Method
656 /// J.Chung and G.M. Hulbert, J. Appl. Mech 60(2), 371-375, 1993
657 /// https://doi.org/10.1115/1.2900803
658 /// rho_inf in [0,1]
660 {
661 protected:
664  bool first;
665 
666 public:
667  GeneralizedAlpha2Solver(double rho_inf = 1.0)
668  {
669  rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
670  rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
671 
672  alpha_m = (2.0 - rho_inf)/(1.0 + rho_inf);
673  alpha_f = 1.0/(1.0 + rho_inf);
674  beta = 0.25*pow(1.0 + alpha_m - alpha_f,2);
675  gamma = 0.5 + alpha_m - alpha_f;
676  };
677 
678  virtual void PrintProperties(std::ostream &out = mfem::out);
679 
680  virtual void Init(SecondOrderTimeDependentOperator &_f);
681 
682  virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt);
683 };
684 
685 /// The classical midpoint method.
687 {
688 public:
690  {
691  alpha_m = 0.5;
692  alpha_f = 0.5;
693  beta = 0.25;
694  gamma = 0.5;
695  };
696 };
697 
698 /// HHT-alpha ODE solver
699 /// Improved numerical dissipation for time integration algorithms
700 /// in structural dynamics
701 /// H.M. Hilber, T.J.R. Hughes and R.L. Taylor 1977
702 /// https://doi.org/10.1002/eqe.4290050306
703 /// alpha in [2/3,1] --> Defined differently than in paper.
705 {
706 public:
707  HHTAlphaSolver(double alpha = 1.0)
708  {
709  alpha = (alpha > 1.0) ? 1.0 : alpha;
710  alpha = (alpha < 2.0/3.0) ? 2.0/3.0 : alpha;
711 
712  alpha_m = 1.0;
713  alpha_f = alpha;
714  beta = (2-alpha)*(2-alpha)/4;
715  gamma = 0.5 + alpha_m - alpha_f;
716  };
717 
718 };
719 
720 /// WBZ-alpha ODE solver
721 /// An alpha modification of Newmark's method
722 /// W.L. Wood, M. Bossak and O.C. Zienkiewicz 1980
723 /// https://doi.org/10.1002/nme.1620151011
724 /// rho_inf in [0,1]
726 {
727 public:
728  WBZAlphaSolver(double rho_inf = 1.0)
729  {
730  rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
731  rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
732 
733  alpha_f = 1.0;
734  alpha_m = 2.0/(1.0 + rho_inf);
735  beta = 0.25*pow(1.0 + alpha_m - alpha_f,2);
736  gamma = 0.5 + alpha_m - alpha_f;
737  };
738 
739 };
740 
741 }
742 
743 #endif
virtual ~SecondOrderODESolver()
Definition: ode.hpp:609
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:39
void PrintProperties(std::ostream &out=mfem::out)
Definition: ode.cpp:658
virtual void Init(Operator &P, TimeDependentOperator &F)
Definition: ode.cpp:713
void SetRhoInf(double rho_inf)
Definition: ode.cpp:647
void Step(Vector &q, Vector &p, double &t, double &dt)
Definition: ode.cpp:722
virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:927
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:539
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:109
Base abstract class for first order time dependent operators.
Definition: operator.hpp:259
MemoryType mem_type
Definition: ode.hpp:27
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:860
virtual void Init(SecondOrderTimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:887
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:447
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:159
Abstract class for solving systems of ODEs: d2x/dt2 = f(x,dx/dt,t)
Definition: ode.hpp:528
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:360
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:610
virtual ~ODESolver()
Definition: ode.hpp:94
virtual void PrintProperties(std::ostream &out=mfem::out)
Definition: ode.cpp:831
AdamsMoultonSolver(int _s, const double *_a)
Definition: ode.cpp:416
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:574
virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:504
GeneralizedAlpha2Solver(double rho_inf=1.0)
Definition: ode.hpp:667
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:100
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:368
virtual ~ExplicitRKSolver()
Definition: ode.cpp:189
The classical midpoint method.
Definition: ode.hpp:686
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:47
Vector dq_
Definition: ode.hpp:492
virtual void Init(SecondOrderTimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:817
virtual void Init(SecondOrderTimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:823
virtual ~SIASolver()
Definition: ode.hpp:485
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:68
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:510
ExplicitRKSolver(int _s, const double *_a, const double *_b, const double *_c)
Definition: ode.cpp:138
WBZAlphaSolver(double rho_inf=1.0)
Definition: ode.hpp:728
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:495
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:489
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:637
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:566
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:145
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
Definition: device.hpp:240
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:24
TimeDependentOperator * F_
Definition: ode.hpp:488
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:27
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:433
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:132
Base abstract class for second order time dependent operators.
Definition: operator.hpp:451
NewmarkSolver(double beta_=0.25, double gamma_=0.5)
Definition: ode.hpp:626
virtual void PrintProperties(std::ostream &out=mfem::out)
Definition: ode.cpp:898
Vector dp_
Definition: ode.hpp:491
Operator * P_
Definition: ode.hpp:489
RK2Solver(const double _a=2./3.)
Definition: ode.hpp:123
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:373
Host memory; using new[] and delete[].
void Step(Vector &q, Vector &p, double &t, double &dt)
Definition: ode.cpp:735
void Step(Vector &q, Vector &p, double &t, double &dt)
Definition: ode.cpp:792
virtual void Run(Vector &x, Vector &dxdt, double &t, double &dt, double tf)
Perform time integration from time t [in] to time tf [in].
Definition: ode.hpp:604
SIAVSolver(int order)
Definition: ode.cpp:750
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:685
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:76
AdamsBashforthSolver(int _s, const double *_a)
Definition: ode.cpp:347
virtual void Step(Vector &q, Vector &p, double &t, double &dt)=0
GeneralizedAlphaSolver(double rho=1.0)
Definition: ode.hpp:448
virtual void Run(Vector &q, Vector &p, double &t, double &dt, double tf)
Definition: ode.hpp:480
const double alpha
Definition: ex15.cpp:336
Vector data type.
Definition: vector.hpp:48
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:603
virtual void Run(Vector &x, double &t, double &dt, double tf)
Perform time integration from time t [in] to time tf [in].
Definition: ode.hpp:89
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:546
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition: ode.hpp:26
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
The classical forward Euler method.
Definition: ode.hpp:99
Abstract operator.
Definition: operator.hpp:24
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:30
SDIRK23Solver(int gamma_opt=1)
Definition: ode.cpp:519
HHTAlphaSolver(double alpha=1.0)
Definition: ode.hpp:707
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: ode.cpp:381
SecondOrderTimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition: ode.hpp:532
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:148