Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_StepperRKButcherTableau.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_StepperRKButcherTableau_hpp
10 #define Tempus_StepperRKButcherTableau_hpp
11 
12 // disable clang warnings
13 #ifdef __clang__
14 #pragma clang system_header
15 #endif
16 
17 #include "Tempus_StepperExplicitRK.hpp"
18 #include "Tempus_StepperDIRK.hpp"
20 
21 
22 namespace Tempus {
23 
24 
25 // ----------------------------------------------------------------------------
26 /** \brief Forward Euler Runge-Kutta Butcher Tableau
27  *
28  * The tableau for Forward Euler (order=1) is
29  * \f[
30  * \begin{array}{c|c}
31  * c & A \\ \hline
32  * & b^T
33  * \end{array}
34  * \;\;\;\;\mbox{ where }\;\;\;\;
35  * \begin{array}{c|c} 0 & 0 \\ \hline
36  * & 1 \end{array}
37  * \f]
38  */
39 template<class Scalar>
41  virtual public StepperExplicitRK<Scalar>
42 {
43  public:
45  {
46  this->setStepperType("RK Forward Euler");
47  this->setupTableau();
48  this->setupDefault();
49  }
50 
52  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
53  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
54  bool useFSAL,
55  std::string ICConsistency,
56  bool ICConsistencyCheck,
57  bool useEmbedded)
58  {
59  this->setStepperType("RK Forward Euler");
60  this->setupTableau();
61  this->setup(appModel, obs, useFSAL, ICConsistency,
62  ICConsistencyCheck, useEmbedded);
63  }
64 
65  std::string getDescription() const
66  {
67  std::ostringstream Description;
68  Description << this->getStepperType() << "\n"
69  << "c = [ 0 ]'\n"
70  << "A = [ 0 ]\n"
71  << "b = [ 1 ]'";
72  return Description.str();
73  }
74 
75 protected:
76 
77  void setupTableau()
78  {
79  typedef Teuchos::ScalarTraits<Scalar> ST;
80  Teuchos::SerialDenseMatrix<int,Scalar> A(1,1);
81  Teuchos::SerialDenseVector<int,Scalar> b(1);
82  Teuchos::SerialDenseVector<int,Scalar> c(1);
83  A(0,0) = ST::zero();
84  b(0) = ST::one();
85  c(0) = ST::zero();
86  int order = 1;
87 
88  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
89  this->getStepperType(),A,b,c,order,order,order));
90  }
91 };
92 
93 
94 // ----------------------------------------------------------------------------
95 /** \brief Runge-Kutta 4th order Butcher Tableau
96  *
97  * The tableau for RK4 (order=4) is
98  * \f[
99  * \begin{array}{c|c}
100  * c & A \\ \hline
101  * & b^T
102  * \end{array}
103  * \;\;\;\;\mbox{ where }\;\;\;\;
104  * \begin{array}{c|cccc} 0 & 0 & & & \\
105  * 1/2 & 1/2 & 0 & & \\
106  * 1/2 & 0 & 1/2 & 0 & \\
107  * 1 & 0 & 0 & 1 & 0 \\ \hline
108  * & 1/6 & 1/3 & 1/3 & 1/6 \end{array}
109  * \f]
110  */
111 template<class Scalar>
113  virtual public StepperExplicitRK<Scalar>
114 {
115  public:
117  {
118  this->setStepperType("RK Explicit 4 Stage");
119  this->setupTableau();
120  this->setupDefault();
121  }
122 
124  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
125  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
126  bool useFSAL,
127  std::string ICConsistency,
128  bool ICConsistencyCheck,
129  bool useEmbedded)
130  {
131  this->setStepperType("RK Explicit 4 Stage");
132  this->setupTableau();
133  this->setup(appModel, obs, useFSAL, ICConsistency,
134  ICConsistencyCheck, useEmbedded);
135  }
136 
137  std::string getDescription() const
138  {
139  std::ostringstream Description;
140  Description << this->getStepperType() << "\n"
141  << "\"The\" Runge-Kutta Method (explicit):\n"
142  << "Solving Ordinary Differential Equations I:\n"
143  << "Nonstiff Problems, 2nd Revised Edition\n"
144  << "E. Hairer, S.P. Norsett, G. Wanner\n"
145  << "Table 1.2, pg 138\n"
146  << "c = [ 0 1/2 1/2 1 ]'\n"
147  << "A = [ 0 ] \n"
148  << " [ 1/2 0 ]\n"
149  << " [ 0 1/2 0 ]\n"
150  << " [ 0 0 1 0 ]\n"
151  << "b = [ 1/6 1/3 1/3 1/6 ]'";
152  return Description.str();
153  }
154 
156  {
157  typedef Teuchos::ScalarTraits<Scalar> ST;
158  const Scalar one = ST::one();
159  const Scalar zero = ST::zero();
160  const Scalar onehalf = one/(2*one);
161  const Scalar onesixth = one/(6*one);
162  const Scalar onethird = one/(3*one);
163 
164  int NumStages = 4;
165  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
166  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
167  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
168 
169  // Fill A:
170  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero; A(0,3) = zero;
171  A(1,0) = onehalf; A(1,1) = zero; A(1,2) = zero; A(1,3) = zero;
172  A(2,0) = zero; A(2,1) = onehalf; A(2,2) = zero; A(2,3) = zero;
173  A(3,0) = zero; A(3,1) = zero; A(3,2) = one; A(3,3) = zero;
174 
175  // Fill b:
176  b(0) = onesixth; b(1) = onethird; b(2) = onethird; b(3) = onesixth;
177 
178  // fill c:
179  c(0) = zero; c(1) = onehalf; c(2) = onehalf; c(3) = one;
180 
181  int order = 4;
182 
183  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
184  this->getStepperType(),A,b,c,order,order,order));
185  }
186 };
187 
188 
189 // ----------------------------------------------------------------------------
190 /** \brief Explicit RK Bogacki-Shampine Butcher Tableau
191  *
192  * The tableau (order=3(2)) is
193  * \f[
194  * \begin{array}{c|c}
195  * c & A \\ \hline
196  * & b^T \\
197  * & b^{*T}
198  * \end{array}
199  * \;\;\;\;\mbox{ where }\;\;\;\;
200  * \begin{array}{c|cccc} 0 & 0 & & & \\
201  * 1/2 & 1/2 & 0 & & \\
202  * 3/4 & 0 & 3/4 & 0 & \\
203  * 1 & 2/9 & 1/3 & 4/9 & 0 \\ \hline
204  * & 2/9 & 1/3 & 4/9 & 0 \\
205  * & 7/24 & 1/4 & 1/3 & 1/8 \end{array}
206  * \f]
207  * Reference: P. Bogacki and L.F. Shampine.
208  * A 3(2) pair of Runge–Kutta formulas.
209  * Applied Mathematics Letters, 2(4):321 – 325, 1989.
210  */
211 template<class Scalar>
213  virtual public StepperExplicitRK<Scalar>
214 {
215  public:
217  {
218  this->setStepperType("Bogacki-Shampine 3(2) Pair");
219  this->setupTableau();
220  this->setupDefault();
221  }
222 
224  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
225  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
226  bool useFSAL,
227  std::string ICConsistency,
228  bool ICConsistencyCheck,
229  bool useEmbedded)
230  {
231  this->setStepperType("Bogacki-Shampine 3(2) Pair");
232  this->setupTableau();
233  this->setup(appModel, obs, useFSAL, ICConsistency,
234  ICConsistencyCheck, useEmbedded);
235  }
236 
237  std::string getDescription() const
238  {
239  std::ostringstream Description;
240  Description << this->getStepperType() << "\n"
241  << "P. Bogacki and L.F. Shampine.\n"
242  << "A 3(2) pair of Runge–Kutta formulas.\n"
243  << "Applied Mathematics Letters, 2(4):321 – 325, 1989.\n"
244  << "c = [ 0 1/2 3/4 1 ]'\n"
245  << "A = [ 0 ]\n"
246  << " [ 1/2 0 ]\n"
247  << " [ 0 3/4 0 ]\n"
248  << " [ 2/9 1/3 4/9 0 ]\n"
249  << "b = [ 2/9 1/3 4/9 0 ]'\n"
250  << "bstar = [ 7/24 1/4 1/3 1/8 ]'";
251  return Description.str();
252  }
253 
254 protected:
255 
257  {
258  typedef Teuchos::ScalarTraits<Scalar> ST;
259  using Teuchos::as;
260  int NumStages = 4;
261  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
262  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
263  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
264  Teuchos::SerialDenseVector<int,Scalar> bstar(NumStages);
265 
266  const Scalar one = ST::one();
267  const Scalar zero = ST::zero();
268  const Scalar onehalf = one/(2*one);
269  const Scalar onethird = one/(3*one);
270  const Scalar threefourths = (3*one)/(4*one);
271  const Scalar twoninths = (2*one)/(9*one);
272  const Scalar fourninths = (4*one)/(9*one);
273 
274  // Fill A:
275  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero; A(0,3) = zero;
276  A(1,0) = onehalf; A(1,1) = zero; A(1,2) = zero; A(1,3) = zero;
277  A(2,0) = zero; A(2,1) =threefourths; A(2,2) = zero; A(2,3) = zero;
278  A(3,0) =twoninths; A(3,1) = onethird; A(3,2) =fourninths; A(3,3) = zero;
279 
280  // Fill b:
281  b(0) = A(3,0); b(1) = A(3,1); b(2) = A(3,2); b(3) = A(3,3);
282 
283  // Fill c:
284  c(0) = zero; c(1) = onehalf; c(2) = threefourths; c(3) = one;
285 
286  // Fill bstar
287  bstar(0) = as<Scalar>(7*one/(24*one));
288  bstar(1) = as<Scalar>(1*one/(4*one));
289  bstar(2) = as<Scalar>(1*one/(3*one));
290  bstar(3) = as<Scalar>(1*one/(8*one));
291  int order = 3;
292 
293  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
294  this->getStepperType(),A,b,c,order,order,order,bstar));
295  }
296 };
297 
298 
299 // ----------------------------------------------------------------------------
300 /** \brief Explicit RK Merson Butcher Tableau
301  *
302  * The tableau (order=4(5)) is
303  * \f[
304  * \begin{array}{c|c}
305  * c & A \\ \hline
306  * & b^T \\
307  * & b^{*T}
308  * \end{array}
309  * \;\;\;\;\mbox{ where }\;\;\;\;
310  * \begin{array}{c|ccccc} 0 & 0 & & & & \\
311  * 1/3 & 1/3 & 0 & & & \\
312  * 1/3 & 1/6 & 1/6 & 0 & & \\
313  * 1/2 & 1/8 & 0 & 3/8 & & \\
314  * 1 & 1/2 & 0 & -3/2 & 2 & \\ \hline
315  * & 1/6 & 0 & 0 & 2/3 & 1/6 \\
316  * & 1/10 & 0 & 3/10 & 2/5 & 1/5 \end{array}
317  * \f]
318  * Reference: E. Hairer, S.P. Norsett, G. Wanner,
319  * "Solving Ordinary Differential Equations I:
320  * Nonstiff Problems", 2nd Revised Edition,
321  * Table 4.1, pg 167.
322  *
323  */
324 template<class Scalar>
326  virtual public StepperExplicitRK<Scalar>
327 {
328  public:
330  {
331  this->setStepperType("Merson 4(5) Pair");
332  this->setupTableau();
333  this->setupDefault();
334  }
335 
337  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
338  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
339  bool useFSAL,
340  std::string ICConsistency,
341  bool ICConsistencyCheck,
342  bool useEmbedded)
343  {
344  this->setStepperType("Merson 4(5) Pair");
345  this->setupTableau();
346  this->setup(appModel, obs, useFSAL, ICConsistency,
347  ICConsistencyCheck, useEmbedded);
348  }
349 
350  std::string getDescription() const
351  {
352  std::ostringstream Description;
353  Description << this->getStepperType() << "\n"
354  << "Solving Ordinary Differential Equations I:\n"
355  << "Nonstiff Problems, 2nd Revised Edition\n"
356  << "E. Hairer, S.P. Norsett, G. Wanner\n"
357  << "Table 4.1, pg 167\n"
358  << "c = [ 0 1/3 1/3 1/2 1 ]'\n"
359  << "A = [ 0 ]\n"
360  << " [ 1/3 0 ]\n"
361  << " [ 1/6 1/6 0 ]\n"
362  << " [ 1/8 0 3/8 0 ]\n"
363  << " [ 1/2 0 -3/2 2 0 ]\n"
364  << "b = [ 1/6 0 0 2/3 1/6 ]'\n"
365  << "bstar = [ 1/10 0 3/10 2/5 1/5 ]'";
366  return Description.str();
367  }
368 
369 
370 protected:
371 
373  {
374  typedef Teuchos::ScalarTraits<Scalar> ST;
375  using Teuchos::as;
376  int NumStages = 5;
377  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages, true);
378  Teuchos::SerialDenseVector<int,Scalar> b(NumStages, true);
379  Teuchos::SerialDenseVector<int,Scalar> c(NumStages, true);
380  Teuchos::SerialDenseVector<int,Scalar> bstar(NumStages, true);
381 
382  const Scalar one = ST::one();
383  const Scalar zero = ST::zero();
384 
385  // Fill A:
386  A(1,0) = as<Scalar>(one/(3*one));;
387 
388  A(2,0) = as<Scalar>(one/(6*one));;
389  A(2,1) = as<Scalar>(one/(6*one));;
390 
391  A(3,0) = as<Scalar>(one/(8*one));;
392  A(3,2) = as<Scalar>(3*one/(8*one));;
393 
394  A(4,0) = as<Scalar>(one/(2*one));;
395  A(4,2) = as<Scalar>(-3*one/(2*one));;
396  A(4,3) = 2*one;
397 
398  // Fill b:
399  b(0) = as<Scalar>(one/(6*one));
400  b(3) = as<Scalar>(2*one/(3*one));
401  b(4) = as<Scalar>(one/(6*one));
402 
403  // Fill c:
404  c(0) = zero;
405  c(1) = as<Scalar>(1*one/(3*one));
406  c(2) = as<Scalar>(1*one/(3*one));
407  c(3) = as<Scalar>(1*one/(2*one));
408  c(4) = one;
409 
410  // Fill bstar
411  bstar(0) = as<Scalar>(1*one/(10*one));
412  bstar(2) = as<Scalar>(3*one/(10*one));
413  bstar(3) = as<Scalar>(2*one/(5*one));
414  bstar(4) = as<Scalar>(1*one/(5*one));
415  int order = 4;
416 
417  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
418  this->getStepperType(),A,b,c,order,order,order,bstar));
419  }
420 };
421 
422 
423 // ----------------------------------------------------------------------------
424 /** \brief Explicit RK 3/8th Rule Butcher Tableau
425  *
426  * The tableau (order=4) is
427  * \f[
428  * \begin{array}{c|c}
429  * c & A \\ \hline
430  * & b^T
431  * \end{array}
432  * \;\;\;\;\mbox{ where }\;\;\;\;
433  * \begin{array}{c|cccc} 0 & 0 & & & \\
434  * 1/3 & 1/3 & 0 & & \\
435  * 2/3 &-1/3 & 1 & 0 & \\
436  * 1 & 1 & -1 & 1 & 0 \\ \hline
437  * & 1/8 & 3/8 & 3/8 & 1/8 \end{array}
438  * \f]
439  * Reference: E. Hairer, S.P. Norsett, G. Wanner,
440  * "Solving Ordinary Differential Equations I:
441  * Nonstiff Problems", 2nd Revised Edition,
442  * Table 1.2, pg 138.
443  */
444 template<class Scalar>
446  virtual public StepperExplicitRK<Scalar>
447 {
448 public:
449 
451  {
452  this->setStepperType("RK Explicit 3/8 Rule");
453  this->setupTableau();
454  this->setupDefault();
455  }
456 
458  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
459  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
460  bool useFSAL,
461  std::string ICConsistency,
462  bool ICConsistencyCheck,
463  bool useEmbedded)
464  {
465  this->setStepperType("RK Explicit 3/8 Rule");
466  this->setupTableau();
467  this->setup(appModel, obs, useFSAL, ICConsistency,
468  ICConsistencyCheck, useEmbedded);
469  }
470 
471  std::string getDescription() const
472  {
473  std::ostringstream Description;
474  Description << this->getStepperType() << "\n"
475  << "Solving Ordinary Differential Equations I:\n"
476  << "Nonstiff Problems, 2nd Revised Edition\n"
477  << "E. Hairer, S.P. Norsett, G. Wanner\n"
478  << "Table 1.2, pg 138\n"
479  << "c = [ 0 1/3 2/3 1 ]'\n"
480  << "A = [ 0 ]\n"
481  << " [ 1/3 0 ]\n"
482  << " [-1/3 1 0 ]\n"
483  << " [ 1 -1 1 0 ]\n"
484  << "b = [ 1/8 3/8 3/8 1/8 ]'";
485  return Description.str();
486  }
487 
488 
489 protected:
490 
492  {
493  typedef Teuchos::ScalarTraits<Scalar> ST;
494  using Teuchos::as;
495  int NumStages = 4;
496  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
497  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
498  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
499 
500  const Scalar one = ST::one();
501  const Scalar zero = ST::zero();
502  const Scalar onethird = as<Scalar>(one/(3*one));
503  const Scalar twothirds = as<Scalar>(2*one/(3*one));
504  const Scalar oneeighth = as<Scalar>(one/(8*one));
505  const Scalar threeeighths = as<Scalar>(3*one/(8*one));
506 
507  // Fill A:
508  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero; A(0,3) = zero;
509  A(1,0) = onethird; A(1,1) = zero; A(1,2) = zero; A(1,3) = zero;
510  A(2,0) = -onethird; A(2,1) = one; A(2,2) = zero; A(2,3) = zero;
511  A(3,0) = one; A(3,1) = -one; A(3,2) = one; A(3,3) = zero;
512 
513  // Fill b:
514  b(0) =oneeighth; b(1) =threeeighths; b(2) =threeeighths; b(3) =oneeighth;
515 
516  // Fill c:
517  c(0) = zero; c(1) = onethird; c(2) = twothirds; c(3) = one;
518 
519  int order = 4;
520 
521  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
522  this->getStepperType(),A,b,c,order,order,order));
523  }
524 };
525 
526 
527 // ----------------------------------------------------------------------------
528 /** \brief RK Explicit 4 Stage 3rd order by Runge
529  *
530  * The tableau (order=3) is
531  * \f[
532  * \begin{array}{c|c}
533  * c & A \\ \hline
534  * & b^T
535  * \end{array}
536  * \;\;\;\;\mbox{ where }\;\;\;\;
537  * \begin{array}{c|cccc} 0 & 0 & & & \\
538  * 1/2 & 1/2 & 0 & & \\
539  * 1 & 0 & 1 & 0 & \\
540  * 1 & 0 & 0 & 1 & 0 \\ \hline
541  * & 1/6 & 2/3 & 0 & 1/6 \end{array}
542  * \f]
543  * Reference: E. Hairer, S.P. Norsett, G. Wanner,
544  * "Solving Ordinary Differential Equations I:
545  * Nonstiff Problems", 2nd Revised Edition,
546  * Table 1.1, pg 135.
547  */
548 template<class Scalar>
550  virtual public StepperExplicitRK<Scalar>
551 {
552  public:
554  {
555  this->setStepperType("RK Explicit 4 Stage 3rd order by Runge");
556  this->setupTableau();
557  this->setupDefault();
558  }
559 
561  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
562  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
563  bool useFSAL,
564  std::string ICConsistency,
565  bool ICConsistencyCheck,
566  bool useEmbedded)
567  {
568  this->setStepperType("RK Explicit 4 Stage 3rd order by Runge");
569  this->setupTableau();
570  this->setup(appModel, obs, useFSAL, ICConsistency,
571  ICConsistencyCheck, useEmbedded);
572  }
573 
574  std::string getDescription() const
575  {
576  std::ostringstream Description;
577  Description << this->getStepperType() << "\n"
578  << "Solving Ordinary Differential Equations I:\n"
579  << "Nonstiff Problems, 2nd Revised Edition\n"
580  << "E. Hairer, S.P. Norsett, G. Wanner\n"
581  << "Table 1.1, pg 135\n"
582  << "c = [ 0 1/2 1 1 ]'\n"
583  << "A = [ 0 ]\n"
584  << " [ 1/2 0 ]\n"
585  << " [ 0 1 0 ]\n"
586  << " [ 0 0 1 0 ]\n"
587  << "b = [ 1/6 2/3 0 1/6 ]'";
588  return Description.str();
589  }
590 protected:
591 
593  {
594  typedef Teuchos::ScalarTraits<Scalar> ST;
595  int NumStages = 4;
596  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
597  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
598  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
599 
600  const Scalar one = ST::one();
601  const Scalar onehalf = one/(2*one);
602  const Scalar onesixth = one/(6*one);
603  const Scalar twothirds = 2*one/(3*one);
604  const Scalar zero = ST::zero();
605 
606  // Fill A:
607  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero; A(0,3) = zero;
608  A(1,0) = onehalf; A(1,1) = zero; A(1,2) = zero; A(1,3) = zero;
609  A(2,0) = zero; A(2,1) = one; A(2,2) = zero; A(2,3) = zero;
610  A(3,0) = zero; A(3,1) = zero; A(3,2) = one; A(3,3) = zero;
611 
612  // Fill b:
613  b(0) = onesixth; b(1) = twothirds; b(2) = zero; b(3) = onesixth;
614 
615  // Fill c:
616  c(0) = zero; c(1) = onehalf; c(2) = one; c(3) = one;
617 
618  int order = 3;
619 
620  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
621  this->getStepperType(),A,b,c,order,order,order));
622  }
623 };
624 
625 
626 // ----------------------------------------------------------------------------
627 /** \brief RK Explicit 5 Stage 3rd order by Kinnmark and Gray
628  *
629  * The tableau (order=3) is
630  * \f[
631  * \begin{array}{c|c}
632  * c & A \\ \hline
633  * & b^T
634  * \end{array}
635  * \;\;\;\;\mbox{ where }\;\;\;\;
636  * \begin{array}{c|ccccc} 0 & 0 & & & & \\
637  * 1/5 & 1/5 & 0 & & & \\
638  * 1/5 & 0 & 1/5 & 0 & & \\
639  * 1/3 & 0 & 0 & 1/3 & 0 & \\
640  * 2/3 & 0 & 0 & 0 & 2/3 & 0 \\ \hline
641  * & 1/4 & 0 & 0 & 0 & 3/4 \end{array}
642  * \f]
643  * Reference: Modified by P. Ullrich. From the prim_advance_mod.F90
644  * routine in the HOMME atmosphere model code.
645  */
646 template<class Scalar>
648  virtual public StepperExplicitRK<Scalar>
649 {
650  public:
652  {
653  this->setStepperType("RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
654  this->setupTableau();
655  this->setupDefault();
656  }
657 
659  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
660  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
661  bool useFSAL,
662  std::string ICConsistency,
663  bool ICConsistencyCheck,
664  bool useEmbedded)
665  {
666  this->setStepperType("RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
667  this->setupTableau();
668  this->setup(appModel, obs, useFSAL, ICConsistency,
669  ICConsistencyCheck, useEmbedded);
670  }
671 
672  std::string getDescription() const
673  {
674  std::ostringstream Description;
675  Description << this->getStepperType() << "\n"
676  << "Kinnmark & Gray 5 stage, 3rd order scheme \n"
677  << "Modified by P. Ullrich. From the prim_advance_mod.F90 \n"
678  << "routine in the HOMME atmosphere model code.\n"
679  << "c = [ 0 1/5 1/5 1/3 2/3 ]'\n"
680  << "A = [ 0 ]\n"
681  << " [ 1/5 0 ]\n"
682  << " [ 0 1/5 0 ]\n"
683  << " [ 0 0 1/3 0 ]\n"
684  << " [ 0 0 0 2/3 0 ]\n"
685  << "b = [ 1/4 0 0 0 3/4 ]'";
686  return Description.str();
687  }
688 
689 protected:
690 
692  {
693  typedef Teuchos::ScalarTraits<Scalar> ST;
694  int NumStages = 5;
695  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
696  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
697  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
698 
699  const Scalar one = ST::one();
700  const Scalar onefifth = one/(5*one);
701  const Scalar onefourth = one/(4*one);
702  const Scalar onethird = one/(3*one);
703  const Scalar twothirds = 2*one/(3*one);
704  const Scalar threefourths = 3*one/(4*one);
705  const Scalar zero = ST::zero();
706 
707  // Fill A:
708  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero; A(0,3) = zero; A(0,4) = zero;
709  A(1,0) = onefifth; A(1,1) = zero; A(1,2) = zero; A(1,3) = zero; A(1,4) = zero;
710  A(2,0) = zero; A(2,1) = onefifth; A(2,2) = zero; A(2,3) = zero; A(2,4) = zero;
711  A(3,0) = zero; A(3,1) = zero; A(3,2) = onethird; A(3,3) = zero; A(3,4) = zero;
712  A(4,0) = zero; A(4,1) = zero; A(4,2) = zero; A(4,3) = twothirds; A(4,4) = zero;
713 
714  // Fill b:
715  b(0) =onefourth; b(1) =zero; b(2) =zero; b(3) =zero; b(4) =threefourths;
716 
717  // Fill c:
718  c(0) =zero; c(1) =onefifth; c(2) =onefifth; c(3) =onethird; c(4) =twothirds;
719 
720  int order = 3;
721 
722  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
723  this->getStepperType(),A,b,c,order,order,order));
724  }
725 };
726 
727 
728 // ----------------------------------------------------------------------------
729 /** \brief RK Explicit 3 Stage 3rd order
730  *
731  * The tableau (order=3) is
732  * \f[
733  * \begin{array}{c|c}
734  * c & A \\ \hline
735  * & b^T
736  * \end{array}
737  * \;\;\;\;\mbox{ where }\;\;\;\;
738  * \begin{array}{c|ccc} 0 & 0 & & \\
739  * 1/2 & 1/2 & 0 & \\
740  * 1 & -1 & 2 & 0 \\ \hline
741  * & 1/6 & 4/6 & 1/6 \end{array}
742  * \f]
743  */
744 template<class Scalar>
746  virtual public StepperExplicitRK<Scalar>
747 {
748  public:
750  {
751  this->setStepperType("RK Explicit 3 Stage 3rd order");
752  this->setupTableau();
753  this->setupDefault();
754  }
755 
757  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
758  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
759  bool useFSAL,
760  std::string ICConsistency,
761  bool ICConsistencyCheck,
762  bool useEmbedded)
763  {
764  this->setStepperType("RK Explicit 3 Stage 3rd order");
765  this->setupTableau();
766  this->setup(appModel, obs, useFSAL, ICConsistency,
767  ICConsistencyCheck, useEmbedded);
768  }
769 
770  std::string getDescription() const
771  {
772  std::ostringstream Description;
773  Description << this->getStepperType() << "\n"
774  << "c = [ 0 1/2 1 ]'\n"
775  << "A = [ 0 ]\n"
776  << " [ 1/2 0 ]\n"
777  << " [ -1 2 0 ]\n"
778  << "b = [ 1/6 4/6 1/6 ]'";
779  return Description.str();
780  }
781 
782 protected:
783 
785  {
786  typedef Teuchos::ScalarTraits<Scalar> ST;
787  const Scalar one = ST::one();
788  const Scalar two = Teuchos::as<Scalar>(2*one);
789  const Scalar zero = ST::zero();
790  const Scalar onehalf = one/(2*one);
791  const Scalar onesixth = one/(6*one);
792  const Scalar foursixth = 4*one/(6*one);
793 
794  int NumStages = 3;
795  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
796  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
797  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
798 
799  // Fill A:
800  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero;
801  A(1,0) = onehalf; A(1,1) = zero; A(1,2) = zero;
802  A(2,0) = -one; A(2,1) = two; A(2,2) = zero;
803 
804  // Fill b:
805  b(0) = onesixth; b(1) = foursixth; b(2) = onesixth;
806 
807  // fill c:
808  c(0) = zero; c(1) = onehalf; c(2) = one;
809 
810  int order = 3;
811 
812  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
813  this->getStepperType(),A,b,c,order,order,order));
814  }
815 };
816 
817 
818 // ----------------------------------------------------------------------------
819 /** \brief RK Explicit 3 Stage 3rd order TVD
820  *
821  * The tableau (order=3) is
822  * \f[
823  * \begin{array}{c|c}
824  * c & A \\ \hline
825  * & b^T
826  * \end{array}
827  * \;\;\;\;\mbox{ where }\;\;\;\;
828  * \begin{array}{c|ccc} 0 & 0 & & \\
829  * 1 & 1 & 0 & \\
830  * 1/2 & 1/4 & 1/4 & 0 \\ \hline
831  * & 1/6 & 1/6 & 4/6 \end{array}
832  * \f]
833  * Reference: Sigal Gottlieb and Chi-Wang Shu,
834  * 'Total Variation Diminishing Runge-Kutta Schemes',
835  * Mathematics of Computation,
836  * Volume 67, Number 221, January 1998, pp. 73-85.
837  *
838  * This is also written in the following set of updates.
839  \verbatim
840  u1 = u^n + dt L(u^n)
841  u2 = 3 u^n/4 + u1/4 + dt L(u1)/4
842  u^(n+1) = u^n/3 + 2 u2/2 + 2 dt L(u2)/3
843  \endverbatim
844  */
845 template<class Scalar>
847  virtual public StepperExplicitRK<Scalar>
848 {
849  public:
851  {
852  this->setStepperType("RK Explicit 3 Stage 3rd order TVD");
853  this->setupTableau();
854  this->setupDefault();
855  }
856 
858  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
859  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
860  bool useFSAL,
861  std::string ICConsistency,
862  bool ICConsistencyCheck,
863  bool useEmbedded)
864  {
865  this->setStepperType("RK Explicit 3 Stage 3rd order TVD");
866  this->setupTableau();
867  this->setup(appModel, obs, useFSAL, ICConsistency,
868  ICConsistencyCheck, useEmbedded);
869  }
870 
871  std::string getDescription() const
872  {
873  std::ostringstream Description;
874  Description << this->getStepperType() << "\n"
875  << "Sigal Gottlieb and Chi-Wang Shu\n"
876  << "`Total Variation Diminishing Runge-Kutta Schemes'\n"
877  << "Mathematics of Computation\n"
878  << "Volume 67, Number 221, January 1998, pp. 73-85\n"
879  << "c = [ 0 1 1/2 ]'\n"
880  << "A = [ 0 ]\n"
881  << " [ 1 0 ]\n"
882  << " [ 1/4 1/4 0 ]\n"
883  << "b = [ 1/6 1/6 4/6 ]'\n"
884  << "This is also written in the following set of updates.\n"
885  << "u1 = u^n + dt L(u^n)\n"
886  << "u2 = 3 u^n/4 + u1/4 + dt L(u1)/4\n"
887  << "u^(n+1) = u^n/3 + 2 u2/2 + 2 dt L(u2)/3";
888  return Description.str();
889  }
890 
891 protected:
892 
894  {
895  typedef Teuchos::ScalarTraits<Scalar> ST;
896  const Scalar one = ST::one();
897  const Scalar zero = ST::zero();
898  const Scalar onehalf = one/(2*one);
899  const Scalar onefourth = one/(4*one);
900  const Scalar onesixth = one/(6*one);
901  const Scalar foursixth = 4*one/(6*one);
902 
903  int NumStages = 3;
904  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
905  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
906  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
907 
908  // Fill A:
909  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero;
910  A(1,0) = one; A(1,1) = zero; A(1,2) = zero;
911  A(2,0) = onefourth; A(2,1) = onefourth; A(2,2) = zero;
912 
913  // Fill b:
914  b(0) = onesixth; b(1) = onesixth; b(2) = foursixth;
915 
916  // fill c:
917  c(0) = zero; c(1) = one; c(2) = onehalf;
918 
919  int order = 3;
920 
921  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
922  this->getStepperType(),A,b,c,order,order,order));
923  }
924 };
925 
926 
927 // ----------------------------------------------------------------------------
928 /** \brief RK Explicit 3 Stage 3rd order by Heun
929  *
930  * The tableau (order=3) is
931  * \f[
932  * \begin{array}{c|c}
933  * c & A \\ \hline
934  * & b^T
935  * \end{array}
936  * \;\;\;\;\mbox{ where }\;\;\;\;
937  * \begin{array}{c|ccc} 0 & 0 & & \\
938  * 1/3 & 1/3 & 0 & \\
939  * 2/3 & 0 & 2/3 & 0 \\ \hline
940  * & 1/4 & 0 & 3/4 \end{array}
941  * \f]
942  * Reference: E. Hairer, S.P. Norsett, G. Wanner,
943  * "Solving Ordinary Differential Equations I:
944  * Nonstiff Problems", 2nd Revised Edition,
945  * Table 1.1, pg 135.
946  */
947 template<class Scalar>
949  virtual public StepperExplicitRK<Scalar>
950 {
951  public:
953  {
954  this->setStepperType("RK Explicit 3 Stage 3rd order by Heun");
955  this->setupTableau();
956  this->setupDefault();
957  }
958 
960  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
961  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
962  bool useFSAL,
963  std::string ICConsistency,
964  bool ICConsistencyCheck,
965  bool useEmbedded)
966  {
967  this->setStepperType("RK Explicit 3 Stage 3rd order by Heun");
968  this->setupTableau();
969  this->setup(appModel, obs, useFSAL, ICConsistency,
970  ICConsistencyCheck, useEmbedded);
971  }
972 
973  std::string getDescription() const
974  {
975  std::ostringstream Description;
976  Description << this->getStepperType() << "\n"
977  << "Solving Ordinary Differential Equations I:\n"
978  << "Nonstiff Problems, 2nd Revised Edition\n"
979  << "E. Hairer, S.P. Norsett, G. Wanner\n"
980  << "Table 1.1, pg 135\n"
981  << "c = [ 0 1/3 2/3 ]'\n"
982  << "A = [ 0 ] \n"
983  << " [ 1/3 0 ]\n"
984  << " [ 0 2/3 0 ]\n"
985  << "b = [ 1/4 0 3/4 ]'";
986  return Description.str();
987  }
988 
989 protected:
990 
992  {
993  typedef Teuchos::ScalarTraits<Scalar> ST;
994  const Scalar one = ST::one();
995  const Scalar zero = ST::zero();
996  const Scalar onethird = one/(3*one);
997  const Scalar twothirds = 2*one/(3*one);
998  const Scalar onefourth = one/(4*one);
999  const Scalar threefourths = 3*one/(4*one);
1000 
1001  int NumStages = 3;
1002  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1003  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1004  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1005 
1006  // Fill A:
1007  A(0,0) = zero; A(0,1) = zero; A(0,2) = zero;
1008  A(1,0) = onethird; A(1,1) = zero; A(1,2) = zero;
1009  A(2,0) = zero; A(2,1) = twothirds; A(2,2) = zero;
1010 
1011  // Fill b:
1012  b(0) = onefourth; b(1) = zero; b(2) = threefourths;
1013 
1014  // fill c:
1015  c(0) = zero; c(1) = onethird; c(2) = twothirds;
1016 
1017  int order = 3;
1018 
1019  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1020  this->getStepperType(),A,b,c,order,order,order));
1021  }
1022 };
1023 
1024 
1025 // ----------------------------------------------------------------------------
1026 /** \brief RK Explicit Midpoint
1027  *
1028  * The tableau (order=2) is
1029  * \f[
1030  * \begin{array}{c|c}
1031  * c & A \\ \hline
1032  * & b^T
1033  * \end{array}
1034  * \;\;\;\;\mbox{ where }\;\;\;\;
1035  * \begin{array}{c|cc} 0 & 0 & \\
1036  * 1/2 & 1/2 & 0 \\ \hline
1037  * & 0 & 1 \end{array}
1038  * \f]
1039  * Reference: E. Hairer, S.P. Norsett, G. Wanner,
1040  * "Solving Ordinary Differential Equations I:
1041  * Nonstiff Problems", 2nd Revised Edition,
1042  * Table 1.1, pg 135.
1043  */
1044 template<class Scalar>
1046  virtual public StepperExplicitRK<Scalar>
1047 {
1048  public:
1050  {
1051  this->setStepperType("RK Explicit Midpoint");
1052  this->setupTableau();
1053  this->setupDefault();
1054  }
1055 
1057  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1058  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1059  bool useFSAL,
1060  std::string ICConsistency,
1061  bool ICConsistencyCheck,
1062  bool useEmbedded)
1063  {
1064  this->setStepperType("RK Explicit Midpoint");
1065  this->setupTableau();
1066  this->setup(appModel, obs, useFSAL, ICConsistency,
1067  ICConsistencyCheck, useEmbedded);
1068  }
1069 
1070  std::string getDescription() const
1071  {
1072  std::ostringstream Description;
1073  Description << this->getStepperType() << "\n"
1074  << "Solving Ordinary Differential Equations I:\n"
1075  << "Nonstiff Problems, 2nd Revised Edition\n"
1076  << "E. Hairer, S.P. Norsett, G. Wanner\n"
1077  << "Table 1.1, pg 135\n"
1078  << "c = [ 0 1/2 ]'\n"
1079  << "A = [ 0 ]\n"
1080  << " [ 1/2 0 ]\n"
1081  << "b = [ 0 1 ]'";
1082  return Description.str();
1083  }
1084 
1085 protected:
1086 
1088  {
1089  typedef Teuchos::ScalarTraits<Scalar> ST;
1090  const Scalar one = ST::one();
1091  const Scalar zero = ST::zero();
1092  const Scalar onehalf = one/(2*one);
1093 
1094  int NumStages = 2;
1095  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1096  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1097  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1098 
1099  // Fill A:
1100  A(0,0) = zero; A(0,1) = zero;
1101  A(1,0) = onehalf; A(1,1) = zero;
1102 
1103  // Fill b:
1104  b(0) = zero; b(1) = one;
1105 
1106  // fill c:
1107  c(0) = zero; c(1) = onehalf;
1108 
1109  int order = 2;
1110 
1111  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1112  this->getStepperType(),A,b,c,order,order,order));
1113  }
1114 };
1115 
1116 
1117 // ----------------------------------------------------------------------------
1118 /** \brief RK Explicit Trapezoidal
1119  *
1120  * The tableau (order=2) is
1121  * \f[
1122  * \begin{array}{c|c}
1123  * c & A \\ \hline
1124  * & b^T
1125  * \end{array}
1126  * \;\;\;\;\mbox{ where }\;\;\;\;
1127  * \begin{array}{c|cc} 0 & 0 & \\
1128  * 1 & 1 & 0 \\ \hline
1129  * & 1/2 & 1/2 \end{array}
1130  * \f]
1131  */
1132 template<class Scalar>
1134  virtual public StepperExplicitRK<Scalar>
1135 {
1136  public:
1138  {
1139  this->setStepperType("RK Explicit Trapezoidal");
1140  this->setupTableau();
1141  this->setupDefault();
1142  }
1143 
1145  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1146  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1147  bool useFSAL,
1148  std::string ICConsistency,
1149  bool ICConsistencyCheck,
1150  bool useEmbedded)
1151  {
1152  this->setStepperType("RK Explicit Trapezoidal");
1153  this->setupTableau();
1154  this->setup(appModel, obs, useFSAL, ICConsistency,
1155  ICConsistencyCheck, useEmbedded);
1156  }
1157 
1158  std::string getDescription() const
1159  {
1160  std::ostringstream Description;
1161  Description << this->getStepperType() << "\n"
1162  << "This Stepper is known as 'RK Explicit Trapezoidal' or 'Heuns Method'.\n"
1163  << "c = [ 0 1 ]'\n"
1164  << "A = [ 0 ]\n"
1165  << " [ 1 0 ]\n"
1166  << "b = [ 1/2 1/2 ]'";
1167  return Description.str();
1168  }
1169 
1170 protected:
1171 
1173  {
1174  typedef Teuchos::ScalarTraits<Scalar> ST;
1175  const Scalar one = ST::one();
1176  const Scalar zero = ST::zero();
1177  const Scalar onehalf = one/(2*one);
1178 
1179  int NumStages = 2;
1180  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1181  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1182  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1183 
1184  // Fill A:
1185  A(0,0) = zero; A(0,1) = zero;
1186  A(1,0) = one; A(1,1) = zero;
1187 
1188  // Fill b:
1189  b(0) = onehalf; b(1) = onehalf;
1190 
1191  // fill c:
1192  c(0) = zero; c(1) = one;
1193 
1194  int order = 2;
1195 
1196  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1197  this->getStepperType(),A,b,c,order,order,order));
1198  }
1199 };
1200 
1201 
1202 // ----------------------------------------------------------------------------
1203 /** \brief General Explicit Runge-Kutta Butcher Tableau
1204  *
1205  * The format of the Butcher Tableau parameter list is
1206  \verbatim
1207  <Parameter name="A" type="string" value="# # # ;
1208  # # # ;
1209  # # #">
1210  <Parameter name="b" type="string" value="# # #">
1211  <Parameter name="c" type="string" value="# # #">
1212  \endverbatim
1213  * Note the number of stages is implicit in the number of entries.
1214  * The number of stages must be consistent.
1215  *
1216  * Default tableau is RK4 (order=4):
1217  * \f[
1218  * \begin{array}{c|c}
1219  * c & A \\ \hline
1220  * & b^T
1221  * \end{array}
1222  * \;\;\;\;\mbox{ where }\;\;\;\;
1223  * \begin{array}{c|cccc} 0 & 0 & & & \\
1224  * 1/2 & 1/2 & 0 & & \\
1225  * 1/2 & 0 & 1/2 & 0 & \\
1226  * 1 & 0 & 0 & 1 & 0 \\ \hline
1227  * & 1/6 & 1/3 & 1/3 & 1/6 \end{array}
1228  * \f]
1229  */
1230 template<class Scalar>
1232  virtual public StepperExplicitRK<Scalar>
1233 {
1234 public:
1236  {
1237  this->setStepperType("General ERK");
1238  this->setupTableau();
1239  this->setupDefault();
1240  }
1241 
1243  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1244  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1245  bool useFSAL,
1246  std::string ICConsistency,
1247  bool ICConsistencyCheck,
1248  bool useEmbedded,
1249  const Teuchos::SerialDenseMatrix<int,Scalar>& A,
1250  const Teuchos::SerialDenseVector<int,Scalar>& b,
1251  const Teuchos::SerialDenseVector<int,Scalar>& c,
1252  const int order,
1253  const int orderMin,
1254  const int orderMax,
1255  const Teuchos::SerialDenseVector<int,Scalar>& bstar)
1256  {
1257  this->setStepperType("General ERK");
1258  this->setTableau(A,b,c,order,orderMin,orderMax,bstar);
1259 
1260  TEUCHOS_TEST_FOR_EXCEPTION(
1261  this->tableau_->isImplicit() == true, std::logic_error,
1262  "Error - General ERK received an implicit Butcher Tableau!\n");
1263 
1264  this->setup(appModel, obs, useFSAL, ICConsistency,
1265  ICConsistencyCheck, useEmbedded);
1266  }
1267 
1268  virtual std::string getDescription() const
1269  {
1270  std::stringstream Description;
1271  Description << this->getStepperType() << "\n"
1272  << "The format of the Butcher Tableau parameter list is\n"
1273  << " <Parameter name=\"A\" type=\"string\" value=\"# # # ;\n"
1274  << " # # # ;\n"
1275  << " # # #\"/>\n"
1276  << " <Parameter name=\"b\" type=\"string\" value=\"# # #\"/>\n"
1277  << " <Parameter name=\"c\" type=\"string\" value=\"# # #\"/>\n\n"
1278  << "Note the number of stages is implicit in the number of entries.\n"
1279  << "The number of stages must be consistent.\n"
1280  << "\n"
1281  << "Default tableau is RK4 (order=4):\n"
1282  << "c = [ 0 1/2 1/2 1 ]'\n"
1283  << "A = [ 0 ]\n"
1284  << " [ 1/2 0 ]\n"
1285  << " [ 0 1/2 0 ]\n"
1286  << " [ 0 0 1 0 ]\n"
1287  << "b = [ 1/6 1/3 1/3 1/6 ]'";
1288  return Description.str();
1289  }
1290 
1292  {
1293  if (this->tableau_ == Teuchos::null) {
1294  // Set tableau to the default if null, otherwise keep current tableau.
1295  auto stepper = Teuchos::rcp(new StepperERK_4Stage4thOrder<Scalar>());
1296  auto t = stepper->getTableau();
1297  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1298  this->getStepperType(),
1299  t->A(),t->b(),t->c(),
1300  t->order(),t->orderMin(),t->orderMax(),
1301  t->bstar()));
1302  }
1303  }
1304 
1305  void setTableau(const Teuchos::SerialDenseMatrix<int,Scalar>& A,
1306  const Teuchos::SerialDenseVector<int,Scalar>& b,
1307  const Teuchos::SerialDenseVector<int,Scalar>& c,
1308  const int order,
1309  const int orderMin,
1310  const int orderMax,
1311  const Teuchos::SerialDenseVector<int,Scalar>&
1312  bstar = Teuchos::SerialDenseVector<int,Scalar>())
1313  {
1314  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1315  this->getStepperType(),A,b,c,order,orderMin,orderMax,bstar));
1316  }
1317 
1318  virtual std::string getDefaultICConsistency() const { return "Consistent"; }
1319 
1320  Teuchos::RCP<const Teuchos::ParameterList>
1322  {
1323  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1324  this->getValidParametersBasicERK(pl);
1325  pl->set<std::string>("Initial Condition Consistency",
1326  this->getDefaultICConsistency());
1327 
1328  // Tableau ParameterList
1329  Teuchos::RCP<Teuchos::ParameterList> tableauPL = Teuchos::parameterList();
1330  tableauPL->set<std::string>("A",
1331  "0.0 0.0 0.0 0.0; 0.5 0.0 0.0 0.0; 0.0 0.5 0.0 0.0; 0.0 0.0 1.0 0.0");
1332  tableauPL->set<std::string>("b",
1333  "0.166666666666667 0.333333333333333 0.333333333333333 0.166666666666667");
1334  tableauPL->set<std::string>("c", "0.0 0.5 0.5 1.0");
1335  tableauPL->set<int>("order", 4);
1336  tableauPL->set<std::string>("bstar", "");
1337  pl->set("Tableau", *tableauPL);
1338 
1339  return pl;
1340  }
1341 };
1342 
1343 
1344 // ----------------------------------------------------------------------------
1345 /** \brief Backward Euler Runge-Kutta Butcher Tableau
1346  *
1347  * The tableau for Backward Euler (order=1) is
1348  * \f[
1349  * \begin{array}{c|c}
1350  * c & A \\ \hline
1351  * & b^T
1352  * \end{array}
1353  * \;\;\;\;\mbox{ where }\;\;\;\;
1354  * \begin{array}{c|c} 1 & 1 \\ \hline
1355  * & 1 \end{array}
1356  * \f]
1357  */
1358 template<class Scalar>
1360  virtual public StepperDIRK<Scalar>
1361 {
1362  public:
1364  {
1365  this->setStepperType("RK Backward Euler");
1366  this->setupTableau();
1367  this->setupDefault();
1368  }
1369 
1371  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1372  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1373  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1374  bool useFSAL,
1375  std::string ICConsistency,
1376  bool ICConsistencyCheck,
1377  bool useEmbedded,
1378  bool zeroInitialGuess)
1379  {
1380  this->setStepperType("RK Backward Euler");
1381  this->setupTableau();
1382  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
1383  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1384  }
1385 
1386  std::string getDescription() const
1387  {
1388  std::ostringstream Description;
1389  Description << this->getStepperType() << "\n"
1390  << "c = [ 1 ]'\n"
1391  << "A = [ 1 ]\n"
1392  << "b = [ 1 ]'";
1393  return Description.str();
1394  }
1395 
1396  virtual bool getICConsistencyCheckDefault() const { return false; }
1397 
1398  Teuchos::RCP<const Teuchos::ParameterList>
1400  {
1401  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1402  this->getValidParametersBasicDIRK(pl);
1403  pl->set<bool>("Initial Condition Consistency Check",
1405  return pl;
1406  }
1407 
1408 protected:
1409 
1411  {
1412  typedef Teuchos::ScalarTraits<Scalar> ST;
1413  int NumStages = 1;
1414  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1415  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1416  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1417 
1418  // Fill A:
1419  A(0,0) = ST::one();
1420 
1421  // Fill b:
1422  b(0) = ST::one();
1423 
1424  // Fill c:
1425  c(0) = ST::one();
1426 
1427  int order = 1;
1428 
1429  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1430  this->getStepperType(),A,b,c,order,order,order));
1431  }
1432 };
1433 
1434 
1435 // ----------------------------------------------------------------------------
1436 /** \brief SDIRK 2 Stage 2nd order
1437  *
1438  * The tableau (order=1 or 2) is
1439  * \f[
1440  * \begin{array}{c|c}
1441  * c & A \\ \hline
1442  * & b^T
1443  * \end{array}
1444  * \;\;\;\;\mbox{ where }\;\;\;\;
1445  * \begin{array}{c|cc} \gamma & \gamma & \\
1446  * 1 & 1-\gamma & \gamma \\ \hline
1447  * & 1-\gamma & \gamma \end{array}
1448  * \f]
1449  * The default value is \f$\gamma = (2 - \sqrt{2})/2\f$.
1450  * This will produce an L-stable 2nd order method with the stage
1451  * times within the timestep. Other values of gamma will still
1452  * produce an L-stable scheme, but will only be 1st order accurate.
1453  * L-stability is guaranteed because \f$A_{sj} = b_j\f$.
1454  *
1455  * Reference: U. M. Ascher and L. R. Petzold,
1456  * Computer Methods for ODEs and DAEs, p. 106.
1457  */
1458 template<class Scalar>
1460  virtual public StepperDIRK<Scalar>
1461 {
1462  public:
1464  {
1465  typedef Teuchos::ScalarTraits<Scalar> ST;
1466  const Scalar one = ST::one();
1467  gammaDefault_ = Teuchos::as<Scalar>((2*one-ST::squareroot(2*one))/(2*one));
1468 
1469  this->setStepperType("SDIRK 2 Stage 2nd order");
1470  this->setGamma(gammaDefault_);
1471  this->setupTableau();
1472  this->setupDefault();
1473  }
1474 
1476  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1477  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1478  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1479  bool useFSAL,
1480  std::string ICConsistency,
1481  bool ICConsistencyCheck,
1482  bool useEmbedded,
1483  bool zeroInitialGuess,
1484  Scalar gamma = Scalar(0.2928932188134524))
1485  {
1486  typedef Teuchos::ScalarTraits<Scalar> ST;
1487  const Scalar one = ST::one();
1488  gammaDefault_ = Teuchos::as<Scalar>((2*one-ST::squareroot(2*one))/(2*one));
1489 
1490  this->setStepperType("SDIRK 2 Stage 2nd order");
1491  this->setGamma(gamma);
1492  this->setupTableau();
1493  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
1494  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1495  }
1496 
1497  void setGamma(Scalar gamma) { gamma_ = gamma; this->setupTableau(); }
1498 
1499  std::string getDescription() const
1500  {
1501  std::ostringstream Description;
1502  Description << this->getStepperType() << "\n"
1503  << "Computer Methods for ODEs and DAEs\n"
1504  << "U. M. Ascher and L. R. Petzold\n"
1505  << "p. 106\n"
1506  << "gamma = (2+-sqrt(2))/2\n"
1507  << "c = [ gamma 1 ]'\n"
1508  << "A = [ gamma 0 ]\n"
1509  << " [ 1-gamma gamma ]\n"
1510  << "b = [ 1-gamma gamma ]'";
1511  return Description.str();
1512  }
1513 
1514  virtual bool getICConsistencyCheckDefault() const { return false; }
1515 
1516  Teuchos::RCP<const Teuchos::ParameterList>
1518  {
1519  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1520  this->getValidParametersBasicDIRK(pl);
1521  pl->set<bool>("Initial Condition Consistency Check",
1523  pl->set<double>("gamma",gammaDefault_,
1524  "The default value is gamma = (2-sqrt(2))/2. "
1525  "This will produce an L-stable 2nd order method with the stage "
1526  "times within the timestep. Other values of gamma will still "
1527  "produce an L-stable scheme, but will only be 1st order accurate.");
1528 
1529  return pl;
1530  }
1531 
1532 protected:
1533 
1535  {
1536  typedef Teuchos::ScalarTraits<Scalar> ST;
1537  int NumStages = 2;
1538  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1539  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1540  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1541 
1542  const Scalar one = ST::one();
1543  const Scalar zero = ST::zero();
1544 
1545  // Fill A:
1546  A(0,0) = gamma_; A(0,1) = zero;
1547  A(1,0) = Teuchos::as<Scalar>( one - gamma_ ); A(1,1) = gamma_;
1548 
1549  // Fill b:
1550  b(0) = Teuchos::as<Scalar>( one - gamma_ ); b(1) = gamma_;
1551 
1552  // Fill c:
1553  c(0) = gamma_; c(1) = one;
1554 
1555  int order = 1;
1556  if ( std::abs((gamma_-gammaDefault_)/gamma_) < 1.0e-08 ) order = 2;
1557 
1558  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1559  this->getStepperType(),A,b,c,order,order,order));
1560  }
1561 
1562  private:
1564  Scalar gamma_;
1565 };
1566 
1567 
1568 // ----------------------------------------------------------------------------
1569 /** \brief SDIRK 2 Stage 3rd order
1570  *
1571  * The tableau (order=2 or 3) is
1572  * \f[
1573  * \begin{array}{c|c}
1574  * c & A \\ \hline
1575  * & b^T
1576  * \end{array}
1577  * \;\;\;\;\mbox{ where }\;\;\;\;
1578  * \begin{array}{c|cc} \gamma & \gamma & \\
1579  * 1-\gamma & 1-2\gamma & \gamma \\ \hline
1580  * & 1/2 & 1/2 \end{array}
1581  * \f]
1582  * \f[
1583  * \gamma = \left\{ \begin{array}{cc}
1584  * (2\pm \sqrt{2})/2 & \mbox{then 2nd order and L-stable} \\
1585  * (3\pm \sqrt{3})/6 & \mbox{then 3rd order and A-stable}
1586  * \end{array} \right.
1587  * \f]
1588  * The default value is \f$\gamma = (3 + \sqrt{3})/6\f$.
1589  *
1590  * Reference: E. Hairer, S. P. Norsett, and G. Wanner,
1591  * Solving Ordinary Differential Equations I:
1592  * Nonstiff Problems, 2nd Revised Edition,
1593  * Table 7.2, pg 207.
1594  */
1595 template<class Scalar>
1597  virtual public StepperDIRK<Scalar>
1598 {
1599  public:
1601  {
1602  typedef Teuchos::ScalarTraits<Scalar> ST;
1603  using Teuchos::as;
1604  const Scalar one = ST::one();
1605  gammaDefault_ = as<Scalar>((3*one+ST::squareroot(3*one))/(6*one));
1606  gammaTypeDefault_ = "3rd Order A-stable";
1607 
1608  this->setStepperType("SDIRK 2 Stage 3rd order");
1609  this->setGammaType(gammaTypeDefault_);
1610  this->setGamma(gammaDefault_);
1611  this->setupTableau();
1612  this->setupDefault();
1613  }
1614 
1616  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1617  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1618  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1619  bool useFSAL,
1620  std::string ICConsistency,
1621  bool ICConsistencyCheck,
1622  bool useEmbedded,
1623  bool zeroInitialGuess,
1624  std::string gammaType = "3rd Order A-stable",
1625  Scalar gamma = Scalar(0.7886751345948128))
1626  {
1627  typedef Teuchos::ScalarTraits<Scalar> ST;
1628  using Teuchos::as;
1629  const Scalar one = ST::one();
1630  gammaDefault_ = as<Scalar>((3*one+ST::squareroot(3*one))/(6*one));
1631  gammaTypeDefault_ = "3rd Order A-stable";
1632 
1633  this->setStepperType("SDIRK 2 Stage 3rd order");
1634  this->setGammaType(gammaType);
1635  this->setGamma(gamma);
1636  this->setupTableau();
1637  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
1638  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1639  }
1640 
1641  void setGammaType(std::string gammaType)
1642  {
1643  TEUCHOS_TEST_FOR_EXCEPTION(
1644  !(gammaType == "3rd Order A-stable" or
1645  gammaType == "2nd Order L-stable" or
1646  gammaType == "gamma"), std::logic_error,
1647  "gammaType needs to be '3rd Order A-stable', '2nd Order L-stable' "
1648  "or 'gamma'.");
1649 
1650  gammaType_ = gammaType;
1651  this->setupTableau();
1652  }
1653  void setGamma(Scalar gamma)
1654  {
1655  if ( gammaType_ == "gamma" ) {
1656  gamma_ = gamma;
1657  this->setupTableau();
1658  }
1659  }
1660 
1661  std::string getDescription() const
1662  {
1663  std::ostringstream Description;
1664  Description << this->getStepperType() << "\n"
1665  << "Solving Ordinary Differential Equations I:\n"
1666  << "Nonstiff Problems, 2nd Revised Edition\n"
1667  << "E. Hairer, S. P. Norsett, and G. Wanner\n"
1668  << "Table 7.2, pg 207\n"
1669  << "gamma = (3+sqrt(3))/6 -> 3rd order and A-stable\n"
1670  << "gamma = (2-sqrt(2))/2 -> 2nd order and L-stable\n"
1671  << "c = [ gamma 1-gamma ]'\n"
1672  << "A = [ gamma 0 ]\n"
1673  << " [ 1-2*gamma gamma ]\n"
1674  << "b = [ 1/2 1/2 ]'";
1675  return Description.str();
1676  }
1677 
1678  Teuchos::RCP<const Teuchos::ParameterList>
1680  {
1681  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1682  this->getValidParametersBasicDIRK(pl);
1683 
1684  pl->set<bool>("Initial Condition Consistency Check", false);
1685  pl->set<std::string>("Gamma Type", gammaTypeDefault_,
1686  "Valid values are '3rd Order A-stable' ((3+sqrt(3))/6.), "
1687  "'2nd Order L-stable' ((2-sqrt(2))/2) and 'gamma' (user defined). "
1688  "The default value is '3rd Order A-stable'.");
1689  pl->set<double>("gamma", gammaDefault_,
1690  "Equal to (3+sqrt(3))/6 if 'Gamma Type' = '3rd Order A-stable', or "
1691  "(2-sqrt(2))/2 if 'Gamma Type' = '2nd Order L-stable', or "
1692  "user-defined gamma value if 'Gamma Type = 'gamma'. "
1693  "The default value is gamma = (3+sqrt(3))/6, which matches "
1694  "the default 'Gamma Type' = '3rd Order A-stable'.");
1695 
1696  return pl;
1697  }
1698 
1699 protected:
1700 
1702  {
1703  typedef Teuchos::ScalarTraits<Scalar> ST;
1704  using Teuchos::as;
1705  int NumStages = 2;
1706  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1707  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1708  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1709  const Scalar one = ST::one();
1710  const Scalar zero = ST::zero();
1711 
1712  int order = 0;
1713  if (gammaType_ == "3rd Order A-stable") {
1714  order = 3;
1716  } else if (gammaType_ == "2nd Order L-stable") {
1717  order = 2;
1718  gamma_ = as<Scalar>( (2*one - ST::squareroot(2*one))/(2*one) );
1719  } else if (gammaType_ == "gamma") {
1720  order = 2;
1721  }
1722 
1723  // Fill A:
1724  A(0,0) = gamma_; A(0,1) = zero;
1725  A(1,0) = as<Scalar>(one - 2*gamma_); A(1,1) = gamma_;
1726 
1727  // Fill b:
1728  b(0) = as<Scalar>( one/(2*one) ); b(1) = as<Scalar>( one/(2*one) );
1729 
1730  // Fill c:
1731  c(0) = gamma_; c(1) = as<Scalar>( one - gamma_ );
1732 
1733  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1734  this->getStepperType(),A,b,c,order,2,3));
1735  }
1736 
1737  private:
1738  std::string gammaTypeDefault_;
1739  std::string gammaType_;
1741  Scalar gamma_;
1742 };
1743 
1744 
1745 // ----------------------------------------------------------------------------
1746 /** \brief EDIRK 2 Stage 3rd order
1747  *
1748  * The tableau (order=3) is
1749  * \f[
1750  * \begin{array}{c|c}
1751  * c & A \\ \hline
1752  * & b^T
1753  * \end{array}
1754  * \;\;\;\;\mbox{ where }\;\;\;\;
1755  * \begin{array}{c|cc} 0 & 0 & \\
1756  * 2/3 & 1/3 & 1/3 \\ \hline
1757  * & 1/4 & 3/4 \end{array}
1758  * \f]
1759  * Reference: E. Hairer, S. P. Norsett, and G. Wanner,
1760  * Solving Ordinary Differential Equations I:
1761  * Nonstiff Problems, 2nd Revised Edition,
1762  * Table 7.1, pg 205.
1763  */
1764 template<class Scalar>
1766  virtual public StepperDIRK<Scalar>
1767 {
1768  public:
1770  {
1771  this->setStepperType("EDIRK 2 Stage 3rd order");
1772  this->setupTableau();
1773  this->setupDefault();
1774  }
1775 
1777  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1778  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1779  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1780  bool useFSAL,
1781  std::string ICConsistency,
1782  bool ICConsistencyCheck,
1783  bool useEmbedded,
1784  bool zeroInitialGuess)
1785  {
1786  this->setStepperType("EDIRK 2 Stage 3rd order");
1787  this->setupTableau();
1788  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
1789  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1790  }
1791 
1792  std::string getDescription() const
1793  {
1794  std::ostringstream Description;
1795  Description << this->getStepperType() << "\n"
1796  << "Hammer & Hollingsworth method\n"
1797  << "Solving Ordinary Differential Equations I:\n"
1798  << "Nonstiff Problems, 2nd Revised Edition\n"
1799  << "E. Hairer, S. P. Norsett, and G. Wanner\n"
1800  << "Table 7.1, pg 205\n"
1801  << "c = [ 0 2/3 ]'\n"
1802  << "A = [ 0 0 ]\n"
1803  << " [ 1/3 1/3 ]\n"
1804  << "b = [ 1/4 3/4 ]'";
1805  return Description.str();
1806  }
1807 
1808  virtual bool getICConsistencyCheckDefault() const { return false; }
1809 
1810  Teuchos::RCP<const Teuchos::ParameterList>
1812  {
1813  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1814  this->getValidParametersBasicDIRK(pl);
1815  pl->set<bool>("Initial Condition Consistency Check",
1817  return pl;
1818  }
1819 
1820 protected:
1821 
1823  {
1824  typedef Teuchos::ScalarTraits<Scalar> ST;
1825  using Teuchos::as;
1826  int NumStages = 2;
1827  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1828  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1829  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1830  const Scalar one = ST::one();
1831  const Scalar zero = ST::zero();
1832 
1833  // Fill A:
1834  A(0,0) = zero; A(0,1) = zero;
1835  A(1,0) = as<Scalar>( one/(3*one) ); A(1,1) = as<Scalar>( one/(3*one) );
1836 
1837  // Fill b:
1838  b(0) = as<Scalar>( one/(4*one) ); b(1) = as<Scalar>( 3*one/(4*one) );
1839 
1840  // Fill c:
1841  c(0) = zero; c(1) = as<Scalar>( 2*one/(3*one) );
1842  int order = 3;
1843 
1844  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1845  this->getStepperType(),A,b,c,order,order,order));
1846  }
1847 };
1848 
1849 
1850 // ----------------------------------------------------------------------------
1851 /** \brief SDIRK 1 Stage Theta
1852  *
1853  * The tableau (order = 1 or 2) is
1854  * \f[
1855  * \begin{array}{c|c}
1856  * c & A \\ \hline
1857  * & b^T
1858  * \end{array}
1859  * \;\;\;\;\mbox{ where }\;\;\;\;
1860  * \begin{array}{c|c} \theta & \theta \\
1861  * & 1 \end{array}
1862  * \f]
1863  * Valid values are \f$0 < \theta \leq 1\f$, where \f$\theta\f$ = 0
1864  * implies Forward Euler (not avialble with this stepepr as it makes it
1865  * explicit), \f$theta\f$ = 1/2 implies implicit midpoint
1866  * method (default), and \f$theta\f$ = 1 implies Backward Euler.
1867  * For \f$theta\f$ != 1/2, this method is first-order accurate,
1868  * and with \f$theta\f$ = 1/2, it is second-order accurate.
1869  * This method is A-stable, but becomes L-stable with \f$theta\f$ = 1.
1870  * (A.K.A. Generalized Implicit Midpoint Method.)
1871  *
1872  * Reference: Non-standard finite-difference methods
1873  * in dynamical systems, P. Kama,
1874  * Dissertation, University of Pretoria, pg. 49.
1875  */
1876 template<class Scalar>
1878  virtual public StepperDIRK<Scalar>
1879 {
1880  public:
1882  {
1883  typedef Teuchos::ScalarTraits<Scalar> ST;
1884  thetaDefault_ = ST::one()/(2*ST::one());
1885 
1886  this->setStepperType("DIRK 1 Stage Theta Method");
1887  this->setTheta(thetaDefault_);
1888  this->setupTableau();
1889  this->setupDefault();
1890  }
1891 
1893  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
1894  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
1895  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1896  bool useFSAL,
1897  std::string ICConsistency,
1898  bool ICConsistencyCheck,
1899  bool useEmbedded,
1900  bool zeroInitialGuess,
1901  Scalar theta = Scalar(0.5))
1902  {
1903  typedef Teuchos::ScalarTraits<Scalar> ST;
1904  thetaDefault_ = ST::one()/(2*ST::one());
1905 
1906  this->setStepperType("DIRK 1 Stage Theta Method");
1907  this->setTheta(theta);
1908  this->setupTableau();
1909  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
1910  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1911  }
1912 
1913  void setTheta(Scalar theta)
1914  {
1915  TEUCHOS_TEST_FOR_EXCEPTION(
1916  theta == Teuchos::ScalarTraits<Scalar>::zero(), std::logic_error,
1917  "'theta' can not be zero, as it makes this stepper explicit. \n"
1918  "Try using the 'RK Forward Euler' stepper.\n");
1919  theta_ = theta;
1920  this->setupTableau();
1921  }
1922 
1923  std::string getDescription() const
1924  {
1925  std::ostringstream Description;
1926  Description << this->getStepperType() << "\n"
1927  << "Non-standard finite-difference methods\n"
1928  << "in dynamical systems, P. Kama,\n"
1929  << "Dissertation, University of Pretoria, pg. 49.\n"
1930  << "Comment: Generalized Implicit Midpoint Method\n"
1931  << "c = [ theta ]'\n"
1932  << "A = [ theta ]\n"
1933  << "b = [ 1 ]'";
1934  return Description.str();
1935  }
1936 
1937  virtual bool getICConsistencyCheckDefault() const { return false; }
1938 
1939  Teuchos::RCP<const Teuchos::ParameterList>
1941  {
1942  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1943  this->getValidParametersBasicDIRK(pl);
1944 
1945  pl->set<bool>("Initial Condition Consistency Check",
1947  pl->set<double>("theta",thetaDefault_,
1948  "Valid values are 0 <= theta <= 1, where theta = 0 "
1949  "implies Forward Euler, theta = 1/2 implies implicit midpoint "
1950  "method (default), and theta = 1 implies Backward Euler. "
1951  "For theta != 1/2, this method is first-order accurate, "
1952  "and with theta = 1/2, it is second-order accurate. "
1953  "This method is A-stable, but becomes L-stable with theta=1.");
1954 
1955  return pl;
1956  }
1957 
1958 protected:
1959 
1961  {
1962  typedef Teuchos::ScalarTraits<Scalar> ST;
1963  int NumStages = 1;
1964  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1965  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1966  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1967  A(0,0) = theta_;
1968  b(0) = ST::one();
1969  c(0) = theta_;
1970 
1971  int order = 1;
1972  if ( std::abs((theta_-thetaDefault_)/theta_) < 1.0e-08 ) order = 2;
1973 
1974  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
1975  this->getStepperType(),A,b,c,order,1,2));
1976  }
1977 
1978  private:
1980  Scalar theta_;
1981 };
1982 
1983 
1984 // ----------------------------------------------------------------------------
1985 /** \brief EDIRK 2 Stage Theta Method
1986  *
1987  * The tableau (order=3) is
1988  * \f[
1989  * \begin{array}{c|c}
1990  * c & A \\ \hline
1991  * & b^T
1992  * \end{array}
1993  * \;\;\;\;\mbox{ where }\;\;\;\;
1994  * \begin{array}{c|cc} 0 & 0 & \\
1995  * 1 & 1-\theta & \theta \\ \hline
1996  * & 1-\theta & \theta \end{array}
1997  * \f]
1998  * Valid values are \f$0 < \theta <= 1\f$, where \f$\theta\f$ = 0
1999  * implies Forward Euler (not avialble with this stepepr as it makes it
2000  * explicit), \f$\theta\f$ = 1/2 implies trapezoidal
2001  * method (default), and \f$\theta\f$ = 1 implies Backward Euler.
2002  * For \f$\theta\f$ != 1/2, this method is first-order accurate,
2003  * and with \f$\theta\f$ = 1/2, it is second-order accurate.
2004  * This method is A-stable, but becomes L-stable with \f$\theta\f$=1.
2005  *
2006  * Reference: Computer Methods for ODEs and DAEs,
2007  * U. M. Ascher and L. R. Petzold, p. 113.
2008  */
2009 template<class Scalar>
2011  virtual public StepperDIRK<Scalar>
2012 {
2013  public:
2015  {
2016  typedef Teuchos::ScalarTraits<Scalar> ST;
2017  thetaDefault_ = ST::one()/(2*ST::one());
2018 
2019  this->setStepperType("EDIRK 2 Stage Theta Method");
2020  this->setTheta(thetaDefault_);
2021  this->setupTableau();
2022  this->setupDefault();
2023  }
2024 
2026  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2027  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2028  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2029  bool useFSAL,
2030  std::string ICConsistency,
2031  bool ICConsistencyCheck,
2032  bool useEmbedded,
2033  bool zeroInitialGuess,
2034  Scalar theta = Scalar(0.5))
2035  {
2036  typedef Teuchos::ScalarTraits<Scalar> ST;
2037  thetaDefault_ = ST::one()/(2*ST::one());
2038 
2039  this->setStepperType("EDIRK 2 Stage Theta Method");
2040  this->setTheta(theta);
2041  this->setupTableau();
2042  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2043  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2044  }
2045 
2046  void setTheta(Scalar theta)
2047  {
2048  TEUCHOS_TEST_FOR_EXCEPTION(
2049  theta == Teuchos::ScalarTraits<Scalar>::zero(), std::logic_error,
2050  "'theta' can not be zero, as it makes this stepper explicit. \n"
2051  "Try using the 'RK Forward Euler' stepper.\n");
2052  theta_ = theta;
2053  this->setupTableau();
2054  }
2055 
2056  std::string getDescription() const
2057  {
2058  std::ostringstream Description;
2059  Description << this->getStepperType() << "\n"
2060  << "Computer Methods for ODEs and DAEs\n"
2061  << "U. M. Ascher and L. R. Petzold\n"
2062  << "p. 113\n"
2063  << "c = [ 0 1 ]'\n"
2064  << "A = [ 0 0 ]\n"
2065  << " [ 1-theta theta ]\n"
2066  << "b = [ 1-theta theta ]'";
2067  return Description.str();
2068  }
2069 
2070  virtual bool getICConsistencyCheckDefault() const { return false; }
2071 
2072  Teuchos::RCP<const Teuchos::ParameterList>
2074  {
2075  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2076  this->getValidParametersBasicDIRK(pl);
2077 
2078  pl->set<bool>("Initial Condition Consistency Check", false);
2079  pl->set<double>("theta",thetaDefault_,
2080  "Valid values are 0 < theta <= 1, where theta = 0 "
2081  "implies Forward Euler, theta = 1/2 implies trapezoidal "
2082  "method (default), and theta = 1 implies Backward Euler. "
2083  "For theta != 1/2, this method is first-order accurate, "
2084  "and with theta = 1/2, it is second-order accurate. "
2085  "This method is A-stable, but becomes L-stable with theta=1.");
2086 
2087  return pl;
2088  }
2089 
2090 protected:
2091 
2093  {
2094  typedef Teuchos::ScalarTraits<Scalar> ST;
2095  const Scalar one = ST::one();
2096  const Scalar zero = ST::zero();
2097 
2098  int NumStages = 2;
2099  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2100  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2101  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2102 
2103  // Fill A:
2104  A(0,0) = zero; A(0,1) = zero;
2105  A(1,0) = Teuchos::as<Scalar>( one - theta_ ); A(1,1) = theta_;
2106 
2107  // Fill b:
2108  b(0) = Teuchos::as<Scalar>( one - theta_ );
2109  b(1) = theta_;
2110 
2111  // Fill c:
2112  c(0) = zero;
2113  c(1) = one;
2114 
2115  int order = 1;
2116  if ( std::abs((theta_-thetaDefault_)/theta_) < 1.0e-08 ) order = 2;
2117 
2118  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
2119  this->getStepperType(),A,b,c,order,1,2));
2120  }
2121 
2122  private:
2124  Scalar theta_;
2125 };
2126 
2127 
2128 // ----------------------------------------------------------------------------
2129 /** \brief RK Trapezoidal Rule (A.K.A. RK Crank-Nicolson)
2130  *
2131  * The tableau (order=3) is
2132  * \f[
2133  * \begin{array}{c|c}
2134  * c & A \\ \hline
2135  * & b^T
2136  * \end{array}
2137  * \;\;\;\;\mbox{ where }\;\;\;\;
2138  * \begin{array}{c|cc} 0 & 0 & \\
2139  * 1 & 1/2 & 1/2 \\ \hline
2140  * & 1/2 & 1/2 \end{array}
2141  * \f]
2142  * It is second-order accurate and A-stable.
2143  *
2144  * Reference: Computer Methods for ODEs and DAEs,
2145  * U. M. Ascher and L. R. Petzold, p. 113.
2146  */
2147 template<class Scalar>
2149  virtual public StepperDIRK<Scalar>
2150 {
2151  public:
2153  {
2154  this->setStepperType("RK Trapezoidal Rule");
2155  this->setupTableau();
2156  this->setupDefault();
2157  }
2158 
2160  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2161  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2162  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2163  bool useFSAL,
2164  std::string ICConsistency,
2165  bool ICConsistencyCheck,
2166  bool useEmbedded,
2167  bool zeroInitialGuess)
2168  {
2169  this->setStepperType("RK Trapezoidal Rule");
2170  this->setupTableau();
2171  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2172  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2173  }
2174 
2175  std::string getDescription() const
2176  {
2177  std::ostringstream Description;
2178  Description << this->getStepperType() << "\n"
2179  << "Also known as Crank-Nicolson Method.\n"
2180  << "c = [ 0 1 ]'\n"
2181  << "A = [ 0 0 ]\n"
2182  << " [ 1/2 1/2 ]\n"
2183  << "b = [ 1/2 1/2 ]'";
2184  return Description.str();
2185  }
2186 
2187  virtual bool getICConsistencyCheckDefault() const { return false; }
2188 
2189  Teuchos::RCP<const Teuchos::ParameterList>
2191  {
2192  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2193  this->getValidParametersBasicDIRK(pl);
2194  pl->set<bool>("Initial Condition Consistency Check",
2196  return pl;
2197  }
2198 
2199 protected:
2200 
2202  {
2203  typedef Teuchos::ScalarTraits<Scalar> ST;
2204  const Scalar one = ST::one();
2205  const Scalar zero = ST::zero();
2206  const Scalar onehalf = ST::one()/(2*ST::one());
2207 
2208  int NumStages = 2;
2209  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2210  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2211  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2212 
2213  // Fill A:
2214  A(0,0) = zero; A(0,1) = zero;
2215  A(1,0) = onehalf; A(1,1) = onehalf;
2216 
2217  // Fill b:
2218  b(0) = onehalf;
2219  b(1) = onehalf;
2220 
2221  // Fill c:
2222  c(0) = zero;
2223  c(1) = one;
2224 
2225  int order = 2;
2226 
2227  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
2228  this->getStepperType(),A,b,c,order,order,order));
2229  }
2230 };
2231 
2232 
2233 // ----------------------------------------------------------------------------
2234 /** \brief SDIRK Implicit Midpoint
2235  *
2236  * The tableau (order = 1 or 2) is
2237  * \f[
2238  * \begin{array}{c|c}
2239  * c & A \\ \hline
2240  * & b^T
2241  * \end{array}
2242  * \;\;\;\;\mbox{ where }\;\;\;\;
2243  * \begin{array}{c|c} 1/2 & 1/2 \\ \hline
2244  * & 1 \end{array}
2245  * \f]
2246  * Implicit midpoint method is second-order accurate, and is A-stable.
2247  *
2248  * Reference: Solving Ordinary Differential Equations II:
2249  * Stiff and Differential-Algebraic Problems,
2250  * 2nd Revised Edition, E. Hairer and G. Wanner,
2251  * Table 5.2, pg 72.
2252  *
2253  * Solving Ordinary Differential Equations I:
2254  * Nonstiff Problems, 2nd Revised Edition,
2255  * E. Hairer, S. P. Norsett, and G. Wanner,
2256  * Table 7.1, pg 205,
2257  */
2258 template<class Scalar>
2260  virtual public StepperDIRK<Scalar>
2261 {
2262  public:
2264  {
2265  this->setStepperType("RK Implicit Midpoint");
2266  this->setupTableau();
2267  this->setupDefault();
2268  }
2269 
2271  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2272  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2273  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2274  bool useFSAL,
2275  std::string ICConsistency,
2276  bool ICConsistencyCheck,
2277  bool useEmbedded,
2278  bool zeroInitialGuess)
2279  {
2280  this->setStepperType("RK Implicit Midpoint");
2281  this->setupTableau();
2282  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2283  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2284  }
2285 
2286  std::string getDescription() const
2287  {
2288  std::ostringstream Description;
2289  Description << this->getStepperType() << "\n"
2290  << "A-stable\n"
2291  << "Solving Ordinary Differential Equations II:\n"
2292  << "Stiff and Differential-Algebraic Problems,\n"
2293  << "2nd Revised Edition\n"
2294  << "E. Hairer and G. Wanner\n"
2295  << "Table 5.2, pg 72\n"
2296  << "Solving Ordinary Differential Equations I:\n"
2297  << "Nonstiff Problems, 2nd Revised Edition\n"
2298  << "E. Hairer, S. P. Norsett, and G. Wanner\n"
2299  << "Table 7.1, pg 205\n"
2300  << "c = [ 1/2 ]'\n"
2301  << "A = [ 1/2 ]\n"
2302  << "b = [ 1 ]'";
2303  return Description.str();
2304  }
2305 
2306  virtual bool getICConsistencyCheckDefault() const { return false; }
2307 
2308  Teuchos::RCP<const Teuchos::ParameterList>
2310  {
2311  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2312  this->getValidParametersBasicDIRK(pl);
2313  pl->set<bool>("Initial Condition Consistency Check",
2315  return pl;
2316  }
2317 
2318 protected:
2319 
2321  {
2322  typedef Teuchos::ScalarTraits<Scalar> ST;
2323  int NumStages = 1;
2324  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2325  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2326  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2327  const Scalar onehalf = ST::one()/(2*ST::one());
2328  const Scalar one = ST::one();
2329 
2330  // Fill A:
2331  A(0,0) = onehalf;
2332 
2333  // Fill b:
2334  b(0) = one;
2335 
2336  // Fill c:
2337  c(0) = onehalf;
2338 
2339  int order = 2;
2340 
2341  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
2342  this->getStepperType(),A,b,c,order,order,order));
2343  }
2344 };
2345 
2346 
2347 // ----------------------------------------------------------------------------
2348 /** \brief RK Implicit 1 Stage 1st order Radau IA
2349  *
2350  * The tableau (order = 1) is
2351  * \f[
2352  * \begin{array}{c|c}
2353  * c & A \\ \hline
2354  * & b^T
2355  * \end{array}
2356  * \;\;\;\;\mbox{ where }\;\;\;\;
2357  * \begin{array}{c|c} 0 & 1 \\ \hline
2358  * & 1 \end{array}
2359  * \f]
2360  * and is A-stable.
2361  * Reference: Solving Ordinary Differential Equations II:
2362  * Stiff and Differential-Algebraic Problems,
2363  * 2nd Revised Edition, E. Hairer and G. Wanner,
2364  * Table 5.3, pg 73.
2365  */
2366 template<class Scalar>
2368  virtual public StepperDIRK<Scalar>
2369 {
2370  public:
2372  {
2373  this->setStepperType("RK Implicit 1 Stage 1st order Radau IA");
2374  this->setupTableau();
2375  this->setupDefault();
2376  }
2377 
2379  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2380  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2381  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2382  bool useFSAL,
2383  std::string ICConsistency,
2384  bool ICConsistencyCheck,
2385  bool useEmbedded,
2386  bool zeroInitialGuess)
2387  {
2388  this->setStepperType("RK Implicit 1 Stage 1st order Radau IA");
2389  this->setupTableau();
2390  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2391  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2392  }
2393 
2394  std::string getDescription() const
2395  {
2396  std::ostringstream Description;
2397  Description << this->getStepperType() << "\n"
2398  << "A-stable\n"
2399  << "Solving Ordinary Differential Equations II:\n"
2400  << "Stiff and Differential-Algebraic Problems,\n"
2401  << "2nd Revised Edition\n"
2402  << "E. Hairer and G. Wanner\n"
2403  << "Table 5.3, pg 73\n"
2404  << "c = [ 0 ]'\n"
2405  << "A = [ 1 ]\n"
2406  << "b = [ 1 ]'";
2407  return Description.str();
2408  }
2409 
2410  virtual bool getICConsistencyCheckDefault() const { return false; }
2411 
2412  Teuchos::RCP<const Teuchos::ParameterList>
2414  {
2415  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2416  this->getValidParametersBasicDIRK(pl);
2417  pl->set<bool>("Initial Condition Consistency Check",
2419  return pl;
2420  }
2421 
2422 protected:
2423 
2425  {
2426  typedef Teuchos::ScalarTraits<Scalar> ST;
2427  int NumStages = 1;
2428  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2429  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2430  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2431  const Scalar one = ST::one();
2432  const Scalar zero = ST::zero();
2433  A(0,0) = one;
2434  b(0) = one;
2435  c(0) = zero;
2436  int order = 1;
2437 
2438  auto emptyBStar = Teuchos::SerialDenseVector<int,Scalar>();
2439  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
2440  this->getStepperType(),A,b,c,order,order,order,emptyBStar,false));
2441  }
2442 };
2443 
2444 
2445 // ----------------------------------------------------------------------------
2446 /** \brief RK Implicit 2 Stage 2nd order Lobatto IIIB
2447  *
2448  * The tableau (order = 2) is
2449  * \f[
2450  * \begin{array}{c|c}
2451  * c & A \\ \hline
2452  * & b^T
2453  * \end{array}
2454  * \;\;\;\;\mbox{ where }\;\;\;\;
2455  * \begin{array}{c|cc} 0 & 1/2 & 0 \\
2456  * 1 & 1/2 & 0 \\ \hline
2457  * & 1/2 & 1/2 \end{array}
2458  * \f]
2459  * It is second-order accurate and A-stable.
2460  *
2461  * Reference: Solving Ordinary Differential Equations II:
2462  * Stiff and Differential-Algebraic Problems,
2463  * 2nd Revised Edition, E. Hairer and G. Wanner,
2464  * Table 5.9, pg 76.
2465  */
2466 template<class Scalar>
2468  virtual public StepperDIRK<Scalar>
2469 {
2470  public:
2472  {
2473  this->setStepperType("RK Implicit 2 Stage 2nd order Lobatto IIIB");
2474  this->setupTableau();
2475  this->setupDefault();
2476  }
2477 
2479  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2480  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2481  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2482  bool useFSAL,
2483  std::string ICConsistency,
2484  bool ICConsistencyCheck,
2485  bool useEmbedded,
2486  bool zeroInitialGuess)
2487  {
2488  this->setStepperType("RK Implicit 2 Stage 2nd order Lobatto IIIB");
2489  this->setupTableau();
2490  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2491  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2492  }
2493 
2494  std::string getDescription() const
2495  {
2496  std::ostringstream Description;
2497  Description << this->getStepperType() << "\n"
2498  << "A-stable\n"
2499  << "Solving Ordinary Differential Equations II:\n"
2500  << "Stiff and Differential-Algebraic Problems,\n"
2501  << "2nd Revised Edition\n"
2502  << "E. Hairer and G. Wanner\n"
2503  << "Table 5.9, pg 76\n"
2504  << "c = [ 0 1 ]'\n"
2505  << "A = [ 1/2 0 ]\n"
2506  << " [ 1/2 0 ]\n"
2507  << "b = [ 1/2 1/2 ]'";
2508  return Description.str();
2509  }
2510 
2511  virtual bool getICConsistencyCheckDefault() const { return false; }
2512 
2513  Teuchos::RCP<const Teuchos::ParameterList>
2515  {
2516  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2517  this->getValidParametersBasicDIRK(pl);
2518  pl->set<bool>("Initial Condition Consistency Check",
2520  return pl;
2521  }
2522 
2523 protected:
2524 
2526  {
2527  typedef Teuchos::ScalarTraits<Scalar> ST;
2528  using Teuchos::as;
2529  int NumStages = 2;
2530  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2531  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2532  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2533  const Scalar zero = ST::zero();
2534  const Scalar one = ST::one();
2535 
2536  // Fill A:
2537  A(0,0) = as<Scalar>( one/(2*one) );
2538  A(0,1) = zero;
2539  A(1,0) = as<Scalar>( one/(2*one) );
2540  A(1,1) = zero;
2541 
2542  // Fill b:
2543  b(0) = as<Scalar>( one/(2*one) );
2544  b(1) = as<Scalar>( one/(2*one) );
2545 
2546  // Fill c:
2547  c(0) = zero;
2548  c(1) = one;
2549  int order = 2;
2550 
2551  auto emptyBStar = Teuchos::SerialDenseVector<int,Scalar>();
2552  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
2553  this->getStepperType(),A,b,c,order,order,order,emptyBStar,false));
2554  }
2555 
2556 };
2557 
2558 
2559 // ----------------------------------------------------------------------------
2560 /** \brief SDIRK 5 Stage 4th order
2561  *
2562  * The tableau (order = 4) is
2563  * \f[
2564  * \begin{array}{c|c}
2565  * c & A \\ \hline
2566  * & b^T
2567  * \end{array}
2568  * \;\;\;\;\mbox{ where }\;\;\;\;
2569  * \begin{array}{c|ccccc}
2570  * 1/4 & 1/4 & & & & \\
2571  * 3/4 & 1/2 & 1/4 & & & \\
2572  * 11/20 & 17/50 & -1/25 & 1/4 & & \\
2573  * 1/2 & 371/1360 & -137/2720 & 15/544 & 1/4 & \\
2574  * 1 & 25/24 & -49/48 & 125/16 & -85/12 & 1/4 \\ \hline
2575  * & 25/24 & -49/48 & 125/16 & -85/12 & 1/4 \end{array}
2576  * \f]
2577  * and is L-stable.
2578  *
2579  * Reference: Solving Ordinary Differential Equations II:
2580  * Stiff and Differential-Algebraic Problems,
2581  * 2nd Revised Edition, E. Hairer and G. Wanner, pg100.
2582  */
2583 template<class Scalar>
2585  virtual public StepperDIRK<Scalar>
2586 {
2587  public:
2589  {
2590  this->setStepperType("SDIRK 5 Stage 4th order");
2591  this->setupTableau();
2592  this->setupDefault();
2593  }
2594 
2596  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2597  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2598  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2599  bool useFSAL,
2600  std::string ICConsistency,
2601  bool ICConsistencyCheck,
2602  bool useEmbedded,
2603  bool zeroInitialGuess)
2604  {
2605  this->setStepperType("SDIRK 5 Stage 4th order");
2606  this->setupTableau();
2607  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2608  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2609  }
2610 
2611  std::string getDescription() const
2612  {
2613  std::ostringstream Description;
2614  Description << this->getStepperType() << "\n"
2615  << "L-stable\n"
2616  << "Solving Ordinary Differential Equations II:\n"
2617  << "Stiff and Differential-Algebraic Problems,\n"
2618  << "2nd Revised Edition\n"
2619  << "E. Hairer and G. Wanner\n"
2620  << "pg100 \n"
2621  << "c = [ 1/4 3/4 11/20 1/2 1 ]'\n"
2622  << "A = [ 1/4 ]\n"
2623  << " [ 1/2 1/4 ]\n"
2624  << " [ 17/50 -1/25 1/4 ]\n"
2625  << " [ 371/1360 -137/2720 15/544 1/4 ]\n"
2626  << " [ 25/24 -49/48 125/16 -85/12 1/4 ]\n"
2627  << "b = [ 25/24 -49/48 125/16 -85/12 1/4 ]'";
2628  //<< "b = [ 59/48 -17/96 225/32 -85/12 0 ]'";
2629  return Description.str();
2630  }
2631 
2632  virtual bool getICConsistencyCheckDefault() const { return false; }
2633 
2634  Teuchos::RCP<const Teuchos::ParameterList>
2636  {
2637  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2638  this->getValidParametersBasicDIRK(pl);
2639  pl->set<bool>("Initial Condition Consistency Check",
2641  return pl;
2642  }
2643 
2644 protected:
2645 
2647  {
2648  typedef Teuchos::ScalarTraits<Scalar> ST;
2649  using Teuchos::as;
2650  int NumStages = 5;
2651  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2652  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2653  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2654  const Scalar zero = ST::zero();
2655  const Scalar one = ST::one();
2656  const Scalar onequarter = as<Scalar>( one/(4*one) );
2657 
2658  // Fill A:
2659  A(0,0) = onequarter;
2660  A(0,1) = zero;
2661  A(0,2) = zero;
2662  A(0,3) = zero;
2663  A(0,4) = zero;
2664 
2665  A(1,0) = as<Scalar>( one / (2*one) );
2666  A(1,1) = onequarter;
2667  A(1,2) = zero;
2668  A(1,3) = zero;
2669  A(1,4) = zero;
2670 
2671  A(2,0) = as<Scalar>( 17*one/(50*one) );
2672  A(2,1) = as<Scalar>( -one/(25*one) );
2673  A(2,2) = onequarter;
2674  A(2,3) = zero;
2675  A(2,4) = zero;
2676 
2677  A(3,0) = as<Scalar>( 371*one/(1360*one) );
2678  A(3,1) = as<Scalar>( -137*one/(2720*one) );
2679  A(3,2) = as<Scalar>( 15*one/(544*one) );
2680  A(3,3) = onequarter;
2681  A(3,4) = zero;
2682 
2683  A(4,0) = as<Scalar>( 25*one/(24*one) );
2684  A(4,1) = as<Scalar>( -49*one/(48*one) );
2685  A(4,2) = as<Scalar>( 125*one/(16*one) );
2686  A(4,3) = as<Scalar>( -85*one/(12*one) );
2687  A(4,4) = onequarter;
2688 
2689  // Fill b:
2690  b(0) = as<Scalar>( 25*one/(24*one) );
2691  b(1) = as<Scalar>( -49*one/(48*one) );
2692  b(2) = as<Scalar>( 125*one/(16*one) );
2693  b(3) = as<Scalar>( -85*one/(12*one) );
2694  b(4) = onequarter;
2695 
2696  /*
2697  // Alternate version
2698  b(0) = as<Scalar>( 59*one/(48*one) );
2699  b(1) = as<Scalar>( -17*one/(96*one) );
2700  b(2) = as<Scalar>( 225*one/(32*one) );
2701  b(3) = as<Scalar>( -85*one/(12*one) );
2702  b(4) = zero;
2703  */
2704 
2705  // Fill c:
2706  c(0) = onequarter;
2707  c(1) = as<Scalar>( 3*one/(4*one) );
2708  c(2) = as<Scalar>( 11*one/(20*one) );
2709  c(3) = as<Scalar>( one/(2*one) );
2710  c(4) = one;
2711 
2712  int order = 4;
2713 
2714  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
2715  this->getStepperType(),A,b,c,order,order,order));
2716  }
2717 };
2718 
2719 
2720 // ----------------------------------------------------------------------------
2721 /** \brief SDIRK 3 Stage 4th order
2722  *
2723  * The tableau (order = 4) is
2724  * \f[
2725  * \begin{array}{c|c}
2726  * c & A \\ \hline
2727  * & b^T
2728  * \end{array}
2729  * \;\;\;\;\mbox{ where }\;\;\;\;
2730  * \begin{array}{c|ccc}
2731  * \gamma & \gamma & & \\
2732  * 1/2 & 1/2-\gamma & \gamma & \\
2733  * 1-\gamma & 2\gamma & 1-4\gamma & \gamma \\ \hline
2734  * & \delta & 1-2\delta & \delta \end{array}
2735  * \f]
2736  * where \f$\gamma = (1/\sqrt{3})\cos(\pi/18)+1/2\f$ and
2737  * \f$\delta = 1/(6(2\gamma-1)^2)\f$, and is A-stable.
2738  *
2739  * Reference: Solving Ordinary Differential Equations II:
2740  * Stiff and Differential-Algebraic Problems,
2741  * 2nd Revised Edition, E. Hairer and G. Wanner, p. 100.
2742  */
2743 template<class Scalar>
2745  virtual public StepperDIRK<Scalar>
2746 {
2747  public:
2749  {
2750  this->setStepperType("SDIRK 3 Stage 4th order");
2751  this->setupTableau();
2752  this->setupDefault();
2753  }
2754 
2756  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2757  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2758  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2759  bool useFSAL,
2760  std::string ICConsistency,
2761  bool ICConsistencyCheck,
2762  bool useEmbedded,
2763  bool zeroInitialGuess)
2764  {
2765  this->setStepperType("SDIRK 3 Stage 4th order");
2766  this->setupTableau();
2767  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2768  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2769  }
2770 
2771  std::string getDescription() const
2772  {
2773  std::ostringstream Description;
2774  Description << this->getStepperType() << "\n"
2775  << "A-stable\n"
2776  << "Solving Ordinary Differential Equations II:\n"
2777  << "Stiff and Differential-Algebraic Problems,\n"
2778  << "2nd Revised Edition\n"
2779  << "E. Hairer and G. Wanner\n"
2780  << "p. 100 \n"
2781  << "gamma = (1/sqrt(3))*cos(pi/18)+1/2\n"
2782  << "delta = 1/(6*(2*gamma-1)^2)\n"
2783  << "c = [ gamma 1/2 1-gamma ]'\n"
2784  << "A = [ gamma ]\n"
2785  << " [ 1/2-gamma gamma ]\n"
2786  << " [ 2*gamma 1-4*gamma gamma ]\n"
2787  << "b = [ delta 1-2*delta delta ]'";
2788  return Description.str();
2789  }
2790 
2791  virtual bool getICConsistencyCheckDefault() const { return false; }
2792 
2793  Teuchos::RCP<const Teuchos::ParameterList>
2795  {
2796  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2797  this->getValidParametersBasicDIRK(pl);
2798  pl->set<bool>("Initial Condition Consistency Check",
2800  return pl;
2801  }
2802 
2803 protected:
2804 
2806  {
2807  typedef Teuchos::ScalarTraits<Scalar> ST;
2808  using Teuchos::as;
2809  int NumStages = 3;
2810  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2811  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2812  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2813  const Scalar zero = ST::zero();
2814  const Scalar one = ST::one();
2815  const Scalar pi = as<Scalar>(4*one)*std::atan(one);
2816  const Scalar gamma = as<Scalar>( one/ST::squareroot(3*one)*std::cos(pi/(18*one))+one/(2*one) );
2817  const Scalar delta = as<Scalar>( one/(6*one*std::pow(2*gamma-one,2*one)) );
2818 
2819  // Fill A:
2820  A(0,0) = gamma;
2821  A(0,1) = zero;
2822  A(0,2) = zero;
2823 
2824  A(1,0) = as<Scalar>( one/(2*one) - gamma );
2825  A(1,1) = gamma;
2826  A(1,2) = zero;
2827 
2828  A(2,0) = as<Scalar>( 2*gamma );
2829  A(2,1) = as<Scalar>( one - 4*gamma );
2830  A(2,2) = gamma;
2831 
2832  // Fill b:
2833  b(0) = delta;
2834  b(1) = as<Scalar>( one-2*delta );
2835  b(2) = delta;
2836 
2837  // Fill c:
2838  c(0) = gamma;
2839  c(1) = as<Scalar>( one/(2*one) );
2840  c(2) = as<Scalar>( one - gamma );
2841 
2842  int order = 4;
2843 
2844  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
2845  this->getStepperType(),A,b,c,order,order,order));
2846  }
2847 };
2848 
2849 
2850 // ----------------------------------------------------------------------------
2851 /** \brief SDIRK 5 Stage 5th order
2852  *
2853  * The tableau (order = 5) is
2854  * \f[
2855  * \begin{array}{c|c}
2856  * c & A \\ \hline
2857  * & b^T
2858  * \end{array}
2859  * \;\;\;\;\mbox{ where }\;\;\;\;
2860  * \begin{array}{c|ccccc}
2861  * (6-\sqrt{6})/10 & (6-\sqrt{6})/10 & 0 & 0 & 0 & 0 \\
2862  * (6+9\sqrt{6})/35 & (-6+5\sqrt{6})/14 & (6-\sqrt{6})/10 & 0 & 0 & 0 \\
2863  * 1 & (888+607\sqrt{6})/2850 & (126-161\sqrt{6})/1425 & (6-\sqrt{6})/10 & 0 & 0 \\
2864  * (4-\sqrt{6})/10 & (3153-3082\sqrt{6})/14250 & (3213+1148\sqrt{6})/28500 & (-267+88\sqrt{6})/500 & (6-\sqrt{6})/10 & 0 \\
2865  * (4+\sqrt{6})/10 & (-32583+14638\sqrt{6})/71250 & (-17199+364\sqrt{6})/142500 & (1329-544\sqrt{6})/2500 & (-96+131\sqrt{6})/625 & (6-\sqrt{6})/10 \\ \hline
2866  * & 0 & 0 & 1/9 & (16-\sqrt{6})/36 & (16+\sqrt{6})/36
2867  * \end{array}
2868  * \f]
2869  *
2870  * Reference: Solving Ordinary Differential Equations II:
2871  * Stiff and Differential-Algebraic Problems,
2872  * 2nd Revised Edition, E. Hairer and G. Wanner, pg101.
2873  */
2874 template<class Scalar>
2876  virtual public StepperDIRK<Scalar>
2877 {
2878  public:
2880  {
2881  this->setStepperType("SDIRK 5 Stage 5th order");
2882  this->setupTableau();
2883  this->setupDefault();
2884  }
2885 
2887  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
2888  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
2889  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2890  bool useFSAL,
2891  std::string ICConsistency,
2892  bool ICConsistencyCheck,
2893  bool useEmbedded,
2894  bool zeroInitialGuess)
2895  {
2896  this->setStepperType("SDIRK 5 Stage 5th order");
2897  this->setupTableau();
2898  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
2899  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2900  }
2901 
2902  std::string getDescription() const
2903  {
2904  std::ostringstream Description;
2905  Description << this->getStepperType() << "\n"
2906  << "Solving Ordinary Differential Equations II:\n"
2907  << "Stiff and Differential-Algebraic Problems,\n"
2908  << "2nd Revised Edition\n"
2909  << "E. Hairer and G. Wanner\n"
2910  << "pg101 \n"
2911  << "c = [ (6-sqrt(6))/10 ]\n"
2912  << " [ (6+9*sqrt(6))/35 ]\n"
2913  << " [ 1 ]\n"
2914  << " [ (4-sqrt(6))/10 ]\n"
2915  << " [ (4+sqrt(6))/10 ]\n"
2916  << "A = [ A1 A2 A3 A4 A5 ]\n"
2917  << " A1 = [ (6-sqrt(6))/10 ]\n"
2918  << " [ (-6+5*sqrt(6))/14 ]\n"
2919  << " [ (888+607*sqrt(6))/2850 ]\n"
2920  << " [ (3153-3082*sqrt(6))/14250 ]\n"
2921  << " [ (-32583+14638*sqrt(6))/71250 ]\n"
2922  << " A2 = [ 0 ]\n"
2923  << " [ (6-sqrt(6))/10 ]\n"
2924  << " [ (126-161*sqrt(6))/1425 ]\n"
2925  << " [ (3213+1148*sqrt(6))/28500 ]\n"
2926  << " [ (-17199+364*sqrt(6))/142500 ]\n"
2927  << " A3 = [ 0 ]\n"
2928  << " [ 0 ]\n"
2929  << " [ (6-sqrt(6))/10 ]\n"
2930  << " [ (-267+88*sqrt(6))/500 ]\n"
2931  << " [ (1329-544*sqrt(6))/2500 ]\n"
2932  << " A4 = [ 0 ]\n"
2933  << " [ 0 ]\n"
2934  << " [ 0 ]\n"
2935  << " [ (6-sqrt(6))/10 ]\n"
2936  << " [ (-96+131*sqrt(6))/625 ]\n"
2937  << " A5 = [ 0 ]\n"
2938  << " [ 0 ]\n"
2939  << " [ 0 ]\n"
2940  << " [ 0 ]\n"
2941  << " [ (6-sqrt(6))/10 ]\n"
2942  << "b = [ 0 ]\n"
2943  << " [ 0 ]\n"
2944  << " [ 1/9 ]\n"
2945  << " [ (16-sqrt(6))/36 ]\n"
2946  << " [ (16+sqrt(6))/36 ]'";
2947  return Description.str();
2948  }
2949 
2950  virtual bool getICConsistencyCheckDefault() const { return false; }
2951 
2952  Teuchos::RCP<const Teuchos::ParameterList>
2954  {
2955  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2956  this->getValidParametersBasicDIRK(pl);
2957  pl->set<bool>("Initial Condition Consistency Check",
2959  return pl;
2960  }
2961 
2962 protected:
2963 
2965  {
2966  typedef Teuchos::ScalarTraits<Scalar> ST;
2967  using Teuchos::as;
2968  int NumStages = 5;
2969  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2970  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2971  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2972  const Scalar zero = ST::zero();
2973  const Scalar one = ST::one();
2974  const Scalar sqrt6 = ST::squareroot(as<Scalar>(6*one));
2975  const Scalar gamma = as<Scalar>( (6*one - sqrt6) / (10*one) ); // diagonal
2976 
2977  // Fill A:
2978  A(0,0) = gamma;
2979  A(0,1) = zero;
2980  A(0,2) = zero;
2981  A(0,3) = zero;
2982  A(0,4) = zero;
2983 
2984  A(1,0) = as<Scalar>( (-6*one+5*one*sqrt6)/(14*one) );
2985  A(1,1) = gamma;
2986  A(1,2) = zero;
2987  A(1,3) = zero;
2988  A(1,4) = zero;
2989 
2990  A(2,0) = as<Scalar>( (888*one+607*one*sqrt6)/(2850*one) );
2991  A(2,1) = as<Scalar>( (126*one-161*one*sqrt6)/(1425*one) );
2992  A(2,2) = gamma;
2993  A(2,3) = zero;
2994  A(2,4) = zero;
2995 
2996  A(3,0) = as<Scalar>( (3153*one-3082*one*sqrt6)/(14250*one) );
2997  A(3,1) = as<Scalar>( (3213*one+1148*one*sqrt6)/(28500*one) );
2998  A(3,2) = as<Scalar>( (-267*one+88*one*sqrt6)/(500*one) );
2999  A(3,3) = gamma;
3000  A(3,4) = zero;
3001 
3002  A(4,0) = as<Scalar>( (-32583*one+14638*one*sqrt6)/(71250*one) );
3003  A(4,1) = as<Scalar>( (-17199*one+364*one*sqrt6)/(142500*one) );
3004  A(4,2) = as<Scalar>( (1329*one-544*one*sqrt6)/(2500*one) );
3005  A(4,3) = as<Scalar>( (-96*one+131*sqrt6)/(625*one) );
3006  A(4,4) = gamma;
3007 
3008  // Fill b:
3009  b(0) = zero;
3010  b(1) = zero;
3011  b(2) = as<Scalar>( one/(9*one) );
3012  b(3) = as<Scalar>( (16*one-sqrt6)/(36*one) );
3013  b(4) = as<Scalar>( (16*one+sqrt6)/(36*one) );
3014 
3015  // Fill c:
3016  c(0) = gamma;
3017  c(1) = as<Scalar>( (6*one+9*one*sqrt6)/(35*one) );
3018  c(2) = one;
3019  c(3) = as<Scalar>( (4*one-sqrt6)/(10*one) );
3020  c(4) = as<Scalar>( (4*one+sqrt6)/(10*one) );
3021 
3022  int order = 5;
3023 
3024  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
3025  this->getStepperType(),A,b,c,order,order,order));
3026  }
3027 };
3028 
3029 
3030 // ----------------------------------------------------------------------------
3031 /** \brief SDIRK 2(1) pair
3032  *
3033  * The tableau (order=2(1)) is
3034  * \f[
3035  * \begin{array}{c|c}
3036  * c & A \\ \hline
3037  * & b^T \\
3038  * & b^{*T}
3039  * \end{array}
3040  * \;\;\;\;\mbox{ where }\;\;\;\;
3041  * \begin{array}{c|cccc} 0 & 0 & \\
3042  * 1 & -1 & 1 \\ \hline
3043  * & 1/2 & 1/2 \\
3044  * & 1 & 0 \end{array}
3045  * \f]
3046  */
3047 template<class Scalar>
3049  virtual public StepperDIRK<Scalar>
3050 {
3051  public:
3053  {
3054  this->setStepperType("SDIRK 2(1) Pair");
3055  this->setupTableau();
3056  this->setupDefault();
3057  }
3058 
3060  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
3061  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
3062  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
3063  bool useFSAL,
3064  std::string ICConsistency,
3065  bool ICConsistencyCheck,
3066  bool useEmbedded,
3067  bool zeroInitialGuess)
3068  {
3069  this->setStepperType("SDIRK 2(1) Pair");
3070  this->setupTableau();
3071  this->setup(appModel, obs, solver, useFSAL, ICConsistency,
3072  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
3073  }
3074 
3075  std::string getDescription() const
3076  {
3077  std::ostringstream Description;
3078  Description << this->getStepperType() << "\n"
3079  << "c = [ 1 0 ]'\n"
3080  << "A = [ 1 ]\n"
3081  << " [ -1 1 ]\n"
3082  << "b = [ 1/2 1/2 ]'\n"
3083  << "bstar = [ 1 0 ]'";
3084  return Description.str();
3085  }
3086 
3087  virtual bool getICConsistencyCheckDefault() const { return false; }
3088 
3089  Teuchos::RCP<const Teuchos::ParameterList>
3091  {
3092  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
3093  this->getValidParametersBasicDIRK(pl);
3094  pl->set<bool>("Initial Condition Consistency Check",
3096  return pl;
3097  }
3098 
3099 protected:
3100 
3102  {
3103  typedef Teuchos::ScalarTraits<Scalar> ST;
3104  using Teuchos::as;
3105  int NumStages = 2;
3106  Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
3107  Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
3108  Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
3109  Teuchos::SerialDenseVector<int,Scalar> bstar(NumStages);
3110 
3111  const Scalar one = ST::one();
3112  const Scalar zero = ST::zero();
3113 
3114  // Fill A:
3115  A(0,0) = one; A(0,1) = zero;
3116  A(1,0) = -one; A(1,1) = one;
3117 
3118  // Fill b:
3119  b(0) = as<Scalar>(one/(2*one));
3120  b(1) = as<Scalar>(one/(2*one));
3121 
3122  // Fill c:
3123  c(0) = one;
3124  c(1) = zero;
3125 
3126  // Fill bstar
3127  bstar(0) = one;
3128  bstar(1) = zero;
3129  int order = 2;
3130 
3131  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
3132  this->getStepperType(),A,b,c,order,order,order,bstar));
3133  }
3134 };
3135 
3136 
3137 // ----------------------------------------------------------------------------
3138 /** \brief General Implicit Runge-Kutta Butcher Tableau
3139  *
3140  * The format of the Butcher Tableau parameter list is
3141  \verbatim
3142  <Parameter name="A" type="string" value="# # # ;
3143  # # # ;
3144  # # #">
3145  <Parameter name="b" type="string" value="# # #">
3146  <Parameter name="c" type="string" value="# # #">
3147  \endverbatim
3148  * Note the number of stages is implicit in the number of entries.
3149  * The number of stages must be consistent.
3150  *
3151  * Default tableau is "SDIRK 2 Stage 2nd order":
3152  * \f[
3153  * \begin{array}{c|c}
3154  * c & A \\ \hline
3155  * & b^T
3156  * \end{array}
3157  * \;\;\;\;\mbox{ where }\;\;\;\;
3158  * \begin{array}{c|cc} \gamma & \gamma & \\
3159  * 1 & 1-\gamma & \gamma \\ \hline
3160  * & 1-\gamma & \gamma \end{array}
3161  * \f]
3162  * where \f$\gamma = (2\pm \sqrt{2})/2\f$. This will produce an
3163  * L-stable 2nd order method.
3164  *
3165  * Reference: U. M. Ascher and L. R. Petzold,
3166  * Computer Methods for ODEs and DAEs, p. 106.
3167  */
3168 template<class Scalar>
3170  virtual public StepperDIRK<Scalar>
3171 {
3172  public:
3174  {
3175  this->setStepperType("General DIRK");
3176  this->setupTableau();
3177  this->setupDefault();
3178  }
3179 
3181  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
3182  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
3183  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
3184  bool useFSAL,
3185  std::string ICConsistency,
3186  bool ICConsistencyCheck,
3187  bool useEmbedded,
3188  bool zeroInitialGuess,
3189  const Teuchos::SerialDenseMatrix<int,Scalar>& A,
3190  const Teuchos::SerialDenseVector<int,Scalar>& b,
3191  const Teuchos::SerialDenseVector<int,Scalar>& c,
3192  const int order,
3193  const int orderMin,
3194  const int orderMax,
3195  const Teuchos::SerialDenseVector<int,Scalar>& bstar)
3196  {
3197  this->setStepperType("General DIRK");
3198  this->setTableau(A,b,c,order,orderMin,orderMax,bstar);
3199 
3200  TEUCHOS_TEST_FOR_EXCEPTION(
3201  this->tableau_->isImplicit() != true, std::logic_error,
3202  "Error - General DIRK did not receive a DIRK Butcher Tableau!\n");
3203 
3204  this->setup(appModel, obs, useFSAL, ICConsistency,
3205  ICConsistencyCheck, useEmbedded, zeroInitialGuess);
3206  }
3207 
3208  std::string getDescription() const
3209  {
3210  std::stringstream Description;
3211  Description << this->getStepperType() << "\n"
3212  << "The format of the Butcher Tableau parameter list is\n"
3213  << " <Parameter name=\"A\" type=\"string\" value=\"# # # ;\n"
3214  << " # # # ;\n"
3215  << " # # #\"/>\n"
3216  << " <Parameter name=\"b\" type=\"string\" value=\"# # #\"/>\n"
3217  << " <Parameter name=\"c\" type=\"string\" value=\"# # #\"/>\n\n"
3218  << "Note the number of stages is implicit in the number of entries.\n"
3219  << "The number of stages must be consistent.\n"
3220  << "\n"
3221  << "Default tableau is 'SDIRK 2 Stage 2nd order':\n"
3222  << " Computer Methods for ODEs and DAEs\n"
3223  << " U. M. Ascher and L. R. Petzold\n"
3224  << " p. 106\n"
3225  << " gamma = (2-sqrt(2))/2\n"
3226  << " c = [ gamma 1 ]'\n"
3227  << " A = [ gamma 0 ]\n"
3228  << " [ 1-gamma gamma ]\n"
3229  << " b = [ 1-gamma gamma ]'";
3230  return Description.str();
3231  }
3232 
3233  virtual bool getICConsistencyCheckDefault() const { return false; }
3234 
3236  {
3237  if (this->tableau_ == Teuchos::null) {
3238  // Set tableau to the default if null, otherwise keep current tableau.
3239  auto stepper = Teuchos::rcp(new StepperSDIRK_2Stage2ndOrder<Scalar>());
3240  auto t = stepper->getTableau();
3241  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
3242  this->getStepperType(),
3243  t->A(),t->b(),t->c(),
3244  t->order(),t->orderMin(),t->orderMax(),
3245  t->bstar()));
3246  }
3247  }
3248 
3249  void setTableau(const Teuchos::SerialDenseMatrix<int,Scalar>& A,
3250  const Teuchos::SerialDenseVector<int,Scalar>& b,
3251  const Teuchos::SerialDenseVector<int,Scalar>& c,
3252  const int order,
3253  const int orderMin,
3254  const int orderMax,
3255  const Teuchos::SerialDenseVector<int,Scalar>&
3256  bstar = Teuchos::SerialDenseVector<int,Scalar>())
3257  {
3258  this->tableau_ = Teuchos::rcp(new RKButcherTableau<Scalar>(
3259  this->getStepperType(),A,b,c,order,orderMin,orderMax,bstar));
3260  }
3261 
3262  Teuchos::RCP<const Teuchos::ParameterList>
3264  {
3265  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
3266  this->getValidParametersBasicDIRK(pl);
3267  pl->set<bool>("Initial Condition Consistency Check",
3269 
3270  // Tableau ParameterList
3271  Teuchos::RCP<Teuchos::ParameterList> tableauPL = Teuchos::parameterList();
3272  tableauPL->set<std::string>("A",
3273  "0.2928932188134524 0.0; 0.7071067811865476 0.2928932188134524");
3274  tableauPL->set<std::string>("b",
3275  "0.7071067811865476 0.2928932188134524");
3276  tableauPL->set<std::string>("c", "0.2928932188134524 1.0");
3277  tableauPL->set<int>("order", 2);
3278  tableauPL->set<std::string>("bstar", "");
3279  pl->set("Tableau", *tableauPL);
3280 
3281  return pl;
3282  }
3283 };
3284 
3285 
3286 } // namespace Tempus
3287 
3288 
3289 #endif // Tempus_StepperRKButcherTableau_hpp
General Implicit Runge-Kutta Butcher Tableau.
Teuchos::RCP< RKButcherTableau< Scalar > > tableau_
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
RK Implicit 2 Stage 2nd order Lobatto IIIB.
virtual std::string getDefaultICConsistency() const
std::string getStepperType() const
Explicit Runge-Kutta time stepper.
virtual void setupDefault()
Default setup for constructor.
Backward Euler Runge-Kutta Butcher Tableau.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
General Explicit Runge-Kutta Butcher Tableau.
RK Explicit 5 Stage 3rd order by Kinnmark and Gray.
virtual void setup(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &wrapperModel, const Teuchos::RCP< StepperRKObserver< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
Setup for constructor.
RK Explicit 4 Stage 3rd order by Runge.
void getValidParametersBasicDIRK(Teuchos::RCP< Teuchos::ParameterList > pl) const
StepperERK_Trapezoidal(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
StepperDIRK_BackwardEuler(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
RK Trapezoidal Rule (A.K.A. RK Crank-Nicolson)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperSDIRK_ImplicitMidpoint(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
virtual void setup(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
Setup for constructor.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperERK_ForwardEuler(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
Explicit RK 3/8th Rule Butcher Tableau.
Diagonally Implicit Runge-Kutta (DIRK) time stepper.
StepperERK_4Stage4thOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
StepperSDIRK_21Pair(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
StepperEDIRK_2Stage3rdOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
StepperERK_BogackiShampine32(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
StepperDIRK_1Stage1stOrderRadauIA(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
StepperEDIRK_2StageTheta(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess, Scalar theta=Scalar(0.5))
StepperERK_3_8Rule(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
virtual void setupDefault()
Default setup for constructor.
Explicit RK Bogacki-Shampine Butcher Tableau.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperERK_General(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, const Teuchos::SerialDenseMatrix< int, Scalar > &A, const Teuchos::SerialDenseVector< int, Scalar > &b, const Teuchos::SerialDenseVector< int, Scalar > &c, const int order, const int orderMin, const int orderMax, const Teuchos::SerialDenseVector< int, Scalar > &bstar)
StepperERK_Midpoint(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
StepperERK_4Stage3rdOrderRunge(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
void getValidParametersBasicERK(Teuchos::RCP< Teuchos::ParameterList > pl) const
StepperERK_3Stage3rdOrderHeun(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
void setTableau(const Teuchos::SerialDenseMatrix< int, Scalar > &A, const Teuchos::SerialDenseVector< int, Scalar > &b, const Teuchos::SerialDenseVector< int, Scalar > &c, const int order, const int orderMin, const int orderMax, const Teuchos::SerialDenseVector< int, Scalar > &bstar=Teuchos::SerialDenseVector< int, Scalar >())
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< RKButcherTableau< Scalar > > tableau_
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperERK_3Stage3rdOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
StepperSDIRK_5Stage4thOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
Forward Euler Runge-Kutta Butcher Tableau.
virtual std::string getDescription() const
Explicit RK Merson Butcher Tableau.
StepperERK_Merson45(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
Runge-Kutta 4th order Butcher Tableau.
StepperSDIRK_3Stage4thOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
StepperEDIRK_TrapezoidalRule(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
void setTableau(const Teuchos::SerialDenseMatrix< int, Scalar > &A, const Teuchos::SerialDenseVector< int, Scalar > &b, const Teuchos::SerialDenseVector< int, Scalar > &c, const int order, const int orderMin, const int orderMax, const Teuchos::SerialDenseVector< int, Scalar > &bstar=Teuchos::SerialDenseVector< int, Scalar >())
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperERK_5Stage3rdOrderKandG(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
StepperSDIRK_2Stage3rdOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess, std::string gammaType="3rd Order A-stable", Scalar gamma=Scalar(0.7886751345948128))
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperDIRK_2Stage2ndOrderLobattoIIIB(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperSDIRK_5Stage5thOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
StepperERK_3Stage3rdOrderTVD(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
void setStepperType(std::string s)
RK Explicit 3 Stage 3rd order by Heun.
StepperDIRK_General(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess, const Teuchos::SerialDenseMatrix< int, Scalar > &A, const Teuchos::SerialDenseVector< int, Scalar > &b, const Teuchos::SerialDenseVector< int, Scalar > &c, const int order, const int orderMin, const int orderMax, const Teuchos::SerialDenseVector< int, Scalar > &bstar)
StepperSDIRK_2Stage2ndOrder(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess, Scalar gamma=Scalar(0.2928932188134524))
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperDIRK_1StageTheta(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess, Scalar theta=Scalar(0.5))