9 #ifndef Tempus_StepperRKButcherTableau_hpp
10 #define Tempus_StepperRKButcherTableau_hpp
14 #pragma clang system_header
17 #include "Tempus_StepperExplicitRK.hpp"
18 #include "Tempus_StepperDIRK.hpp"
39 template<
class Scalar>
52 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
55 std::string ICConsistency,
56 bool ICConsistencyCheck,
61 this->
setup(appModel, obs, useFSAL, ICConsistency,
62 ICConsistencyCheck, useEmbedded);
67 std::ostringstream Description;
72 return Description.str();
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);
111 template<
class Scalar>
124 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
127 std::string ICConsistency,
128 bool ICConsistencyCheck,
133 this->
setup(appModel, obs, useFSAL, ICConsistency,
134 ICConsistencyCheck, useEmbedded);
139 std::ostringstream Description;
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"
151 <<
"b = [ 1/6 1/3 1/3 1/6 ]'";
152 return Description.str();
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);
165 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
166 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
167 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
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;
176 b(0) = onesixth; b(1) = onethird; b(2) = onethird; b(3) = onesixth;
179 c(0) = zero; c(1) = onehalf; c(2) = onehalf; c(3) = one;
211 template<
class Scalar>
224 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
227 std::string ICConsistency,
228 bool ICConsistencyCheck,
233 this->
setup(appModel, obs, useFSAL, ICConsistency,
234 ICConsistencyCheck, useEmbedded);
239 std::ostringstream Description;
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"
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();
258 typedef Teuchos::ScalarTraits<Scalar> ST;
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);
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);
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;
281 b(0) = A(3,0); b(1) = A(3,1); b(2) = A(3,2); b(3) = A(3,3);
284 c(0) = zero; c(1) = onehalf; c(2) = threefourths; c(3) = one;
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));
324 template<
class Scalar>
337 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
340 std::string ICConsistency,
341 bool ICConsistencyCheck,
346 this->
setup(appModel, obs, useFSAL, ICConsistency,
347 ICConsistencyCheck, useEmbedded);
352 std::ostringstream Description;
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"
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();
374 typedef Teuchos::ScalarTraits<Scalar> ST;
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);
382 const Scalar one = ST::one();
383 const Scalar zero = ST::zero();
386 A(1,0) = as<Scalar>(one/(3*one));;
388 A(2,0) = as<Scalar>(one/(6*one));;
389 A(2,1) = as<Scalar>(one/(6*one));;
391 A(3,0) = as<Scalar>(one/(8*one));;
392 A(3,2) = as<Scalar>(3*one/(8*one));;
394 A(4,0) = as<Scalar>(one/(2*one));;
395 A(4,2) = as<Scalar>(-3*one/(2*one));;
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));
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));
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));
444 template<
class Scalar>
458 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
461 std::string ICConsistency,
462 bool ICConsistencyCheck,
467 this->
setup(appModel, obs, useFSAL, ICConsistency,
468 ICConsistencyCheck, useEmbedded);
473 std::ostringstream Description;
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"
484 <<
"b = [ 1/8 3/8 3/8 1/8 ]'";
485 return Description.str();
493 typedef Teuchos::ScalarTraits<Scalar> ST;
496 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
497 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
498 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
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));
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;
514 b(0) =oneeighth; b(1) =threeeighths; b(2) =threeeighths; b(3) =oneeighth;
517 c(0) = zero; c(1) = onethird; c(2) = twothirds; c(3) = one;
548 template<
class Scalar>
561 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
564 std::string ICConsistency,
565 bool ICConsistencyCheck,
570 this->
setup(appModel, obs, useFSAL, ICConsistency,
571 ICConsistencyCheck, useEmbedded);
576 std::ostringstream Description;
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"
587 <<
"b = [ 1/6 2/3 0 1/6 ]'";
588 return Description.str();
594 typedef Teuchos::ScalarTraits<Scalar> ST;
596 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
597 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
598 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
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();
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;
613 b(0) = onesixth; b(1) = twothirds; b(2) = zero; b(3) = onesixth;
616 c(0) = zero; c(1) = onehalf; c(2) = one; c(3) = one;
646 template<
class Scalar>
653 this->
setStepperType(
"RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
659 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
662 std::string ICConsistency,
663 bool ICConsistencyCheck,
666 this->
setStepperType(
"RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
668 this->
setup(appModel, obs, useFSAL, ICConsistency,
669 ICConsistencyCheck, useEmbedded);
674 std::ostringstream Description;
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"
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();
693 typedef Teuchos::ScalarTraits<Scalar> ST;
695 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
696 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
697 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
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();
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;
715 b(0) =onefourth; b(1) =zero; b(2) =zero; b(3) =zero; b(4) =threefourths;
718 c(0) =zero; c(1) =onefifth; c(2) =onefifth; c(3) =onethird; c(4) =twothirds;
744 template<
class Scalar>
757 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
760 std::string ICConsistency,
761 bool ICConsistencyCheck,
766 this->
setup(appModel, obs, useFSAL, ICConsistency,
767 ICConsistencyCheck, useEmbedded);
772 std::ostringstream Description;
774 <<
"c = [ 0 1/2 1 ]'\n"
778 <<
"b = [ 1/6 4/6 1/6 ]'";
779 return Description.str();
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);
795 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
796 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
797 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
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;
805 b(0) = onesixth; b(1) = foursixth; b(2) = onesixth;
808 c(0) = zero; c(1) = onehalf; c(2) = one;
845 template<
class Scalar>
858 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
861 std::string ICConsistency,
862 bool ICConsistencyCheck,
867 this->
setup(appModel, obs, useFSAL, ICConsistency,
868 ICConsistencyCheck, useEmbedded);
873 std::ostringstream Description;
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"
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();
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);
904 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
905 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
906 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
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;
914 b(0) = onesixth; b(1) = onesixth; b(2) = foursixth;
917 c(0) = zero; c(1) = one; c(2) = onehalf;
947 template<
class Scalar>
960 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
963 std::string ICConsistency,
964 bool ICConsistencyCheck,
969 this->
setup(appModel, obs, useFSAL, ICConsistency,
970 ICConsistencyCheck, useEmbedded);
975 std::ostringstream Description;
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"
985 <<
"b = [ 1/4 0 3/4 ]'";
986 return Description.str();
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);
1002 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1003 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1004 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
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;
1012 b(0) = onefourth; b(1) = zero; b(2) = threefourths;
1015 c(0) = zero; c(1) = onethird; c(2) = twothirds;
1044 template<
class Scalar>
1057 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1060 std::string ICConsistency,
1061 bool ICConsistencyCheck,
1066 this->
setup(appModel, obs, useFSAL, ICConsistency,
1067 ICConsistencyCheck, useEmbedded);
1072 std::ostringstream Description;
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"
1082 return Description.str();
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);
1095 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1096 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1097 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1100 A(0,0) = zero; A(0,1) = zero;
1101 A(1,0) = onehalf; A(1,1) = zero;
1104 b(0) = zero; b(1) = one;
1107 c(0) = zero; c(1) = onehalf;
1132 template<
class Scalar>
1145 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1148 std::string ICConsistency,
1149 bool ICConsistencyCheck,
1154 this->
setup(appModel, obs, useFSAL, ICConsistency,
1155 ICConsistencyCheck, useEmbedded);
1160 std::ostringstream Description;
1162 <<
"This Stepper is known as 'RK Explicit Trapezoidal' or 'Heuns Method'.\n"
1166 <<
"b = [ 1/2 1/2 ]'";
1167 return Description.str();
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);
1180 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1181 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1182 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1185 A(0,0) = zero; A(0,1) = zero;
1186 A(1,0) = one; A(1,1) = zero;
1189 b(0) = onehalf; b(1) = onehalf;
1192 c(0) = zero; c(1) = one;
1230 template<
class Scalar>
1243 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1246 std::string ICConsistency,
1247 bool ICConsistencyCheck,
1249 const Teuchos::SerialDenseMatrix<int,Scalar>& A,
1250 const Teuchos::SerialDenseVector<int,Scalar>& b,
1251 const Teuchos::SerialDenseVector<int,Scalar>& c,
1255 const Teuchos::SerialDenseVector<int,Scalar>& bstar)
1258 this->
setTableau(A,b,c,order,orderMin,orderMax,bstar);
1260 TEUCHOS_TEST_FOR_EXCEPTION(
1261 this->
tableau_->isImplicit() ==
true, std::logic_error,
1262 "Error - General ERK received an implicit Butcher Tableau!\n");
1264 this->
setup(appModel, obs, useFSAL, ICConsistency,
1265 ICConsistencyCheck, useEmbedded);
1270 std::stringstream Description;
1272 <<
"The format of the Butcher Tableau parameter list is\n"
1273 <<
" <Parameter name=\"A\" type=\"string\" value=\"# # # ;\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"
1281 <<
"Default tableau is RK4 (order=4):\n"
1282 <<
"c = [ 0 1/2 1/2 1 ]'\n"
1287 <<
"b = [ 1/6 1/3 1/3 1/6 ]'";
1288 return Description.str();
1293 if (this->
tableau_ == Teuchos::null) {
1296 auto t = stepper->getTableau();
1299 t->A(),t->b(),t->c(),
1300 t->order(),t->orderMin(),t->orderMax(),
1306 const Teuchos::SerialDenseVector<int,Scalar>& b,
1307 const Teuchos::SerialDenseVector<int,Scalar>& c,
1311 const Teuchos::SerialDenseVector<int,Scalar>&
1312 bstar = Teuchos::SerialDenseVector<int,Scalar>())
1320 Teuchos::RCP<const Teuchos::ParameterList>
1323 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1325 pl->set<std::string>(
"Initial Condition Consistency",
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);
1358 template<
class Scalar>
1371 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1373 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1375 std::string ICConsistency,
1376 bool ICConsistencyCheck,
1378 bool zeroInitialGuess)
1382 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
1383 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1388 std::ostringstream Description;
1393 return Description.str();
1398 Teuchos::RCP<const Teuchos::ParameterList>
1401 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1403 pl->set<
bool>(
"Initial Condition Consistency Check",
1412 typedef Teuchos::ScalarTraits<Scalar> ST;
1414 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1415 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1416 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1458 template<
class Scalar>
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));
1476 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1478 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1480 std::string ICConsistency,
1481 bool ICConsistencyCheck,
1483 bool zeroInitialGuess,
1484 Scalar gamma = Scalar(0.2928932188134524))
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));
1493 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
1494 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1501 std::ostringstream Description;
1503 <<
"Computer Methods for ODEs and DAEs\n"
1504 <<
"U. M. Ascher and L. R. Petzold\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();
1516 Teuchos::RCP<const Teuchos::ParameterList>
1519 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1521 pl->set<
bool>(
"Initial Condition Consistency Check",
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.");
1536 typedef Teuchos::ScalarTraits<Scalar> ST;
1538 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1539 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1540 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
1542 const Scalar one = ST::one();
1543 const Scalar zero = ST::zero();
1546 A(0,0) =
gamma_; A(0,1) = zero;
1547 A(1,0) = Teuchos::as<Scalar>( one -
gamma_ ); A(1,1) =
gamma_;
1550 b(0) = Teuchos::as<Scalar>( one -
gamma_ ); b(1) =
gamma_;
1553 c(0) =
gamma_; c(1) = one;
1595 template<
class Scalar>
1602 typedef Teuchos::ScalarTraits<Scalar> ST;
1604 const Scalar one = ST::one();
1605 gammaDefault_ = as<Scalar>((3*one+ST::squareroot(3*one))/(6*one));
1616 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1618 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1620 std::string ICConsistency,
1621 bool ICConsistencyCheck,
1623 bool zeroInitialGuess,
1624 std::string gammaType =
"3rd Order A-stable",
1625 Scalar gamma = Scalar(0.7886751345948128))
1627 typedef Teuchos::ScalarTraits<Scalar> ST;
1629 const Scalar one = ST::one();
1630 gammaDefault_ = as<Scalar>((3*one+ST::squareroot(3*one))/(6*one));
1637 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
1638 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
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' "
1663 std::ostringstream Description;
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();
1678 Teuchos::RCP<const Teuchos::ParameterList>
1681 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1684 pl->set<
bool>(
"Initial Condition Consistency Check",
false);
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'.");
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'.");
1703 typedef Teuchos::ScalarTraits<Scalar> ST;
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();
1716 }
else if (
gammaType_ ==
"2nd Order L-stable") {
1718 gamma_ = as<Scalar>( (2*one - ST::squareroot(2*one))/(2*one) );
1724 A(0,0) =
gamma_; A(0,1) = zero;
1728 b(0) = as<Scalar>( one/(2*one) ); b(1) = as<Scalar>( one/(2*one) );
1764 template<
class Scalar>
1777 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1779 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1781 std::string ICConsistency,
1782 bool ICConsistencyCheck,
1784 bool zeroInitialGuess)
1788 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
1789 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
1794 std::ostringstream Description;
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"
1804 <<
"b = [ 1/4 3/4 ]'";
1805 return Description.str();
1810 Teuchos::RCP<const Teuchos::ParameterList>
1813 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1815 pl->set<
bool>(
"Initial Condition Consistency Check",
1824 typedef Teuchos::ScalarTraits<Scalar> ST;
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();
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) );
1838 b(0) = as<Scalar>( one/(4*one) ); b(1) = as<Scalar>( 3*one/(4*one) );
1841 c(0) = zero; c(1) = as<Scalar>( 2*one/(3*one) );
1876 template<
class Scalar>
1883 typedef Teuchos::ScalarTraits<Scalar> ST;
1893 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
1895 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
1897 std::string ICConsistency,
1898 bool ICConsistencyCheck,
1900 bool zeroInitialGuess,
1901 Scalar theta = Scalar(0.5))
1903 typedef Teuchos::ScalarTraits<Scalar> ST;
1909 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
1910 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
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");
1925 std::ostringstream Description;
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"
1934 return Description.str();
1939 Teuchos::RCP<const Teuchos::ParameterList>
1942 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1945 pl->set<
bool>(
"Initial Condition Consistency Check",
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.");
1962 typedef Teuchos::ScalarTraits<Scalar> ST;
1964 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
1965 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
1966 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2009 template<
class Scalar>
2016 typedef Teuchos::ScalarTraits<Scalar> ST;
2026 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2028 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2030 std::string ICConsistency,
2031 bool ICConsistencyCheck,
2033 bool zeroInitialGuess,
2034 Scalar theta = Scalar(0.5))
2036 typedef Teuchos::ScalarTraits<Scalar> ST;
2042 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2043 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
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");
2058 std::ostringstream Description;
2060 <<
"Computer Methods for ODEs and DAEs\n"
2061 <<
"U. M. Ascher and L. R. Petzold\n"
2065 <<
" [ 1-theta theta ]\n"
2066 <<
"b = [ 1-theta theta ]'";
2067 return Description.str();
2072 Teuchos::RCP<const Teuchos::ParameterList>
2075 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2078 pl->set<
bool>(
"Initial Condition Consistency Check",
false);
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.");
2094 typedef Teuchos::ScalarTraits<Scalar> ST;
2095 const Scalar one = ST::one();
2096 const Scalar zero = ST::zero();
2099 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2100 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2101 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2104 A(0,0) = zero; A(0,1) = zero;
2105 A(1,0) = Teuchos::as<Scalar>( one -
theta_ ); A(1,1) =
theta_;
2108 b(0) = Teuchos::as<Scalar>( one -
theta_ );
2147 template<
class Scalar>
2160 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2162 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2164 std::string ICConsistency,
2165 bool ICConsistencyCheck,
2167 bool zeroInitialGuess)
2171 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2172 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2177 std::ostringstream Description;
2179 <<
"Also known as Crank-Nicolson Method.\n"
2183 <<
"b = [ 1/2 1/2 ]'";
2184 return Description.str();
2189 Teuchos::RCP<const Teuchos::ParameterList>
2192 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2194 pl->set<
bool>(
"Initial Condition Consistency Check",
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());
2209 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
2210 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
2211 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
2214 A(0,0) = zero; A(0,1) = zero;
2215 A(1,0) = onehalf; A(1,1) = onehalf;
2258 template<
class Scalar>
2271 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2273 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2275 std::string ICConsistency,
2276 bool ICConsistencyCheck,
2278 bool zeroInitialGuess)
2282 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2283 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2288 std::ostringstream Description;
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"
2303 return Description.str();
2308 Teuchos::RCP<const Teuchos::ParameterList>
2311 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2313 pl->set<
bool>(
"Initial Condition Consistency Check",
2322 typedef Teuchos::ScalarTraits<Scalar> ST;
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();
2366 template<
class Scalar>
2379 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2381 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2383 std::string ICConsistency,
2384 bool ICConsistencyCheck,
2386 bool zeroInitialGuess)
2390 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2391 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2396 std::ostringstream Description;
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"
2407 return Description.str();
2412 Teuchos::RCP<const Teuchos::ParameterList>
2415 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2417 pl->set<
bool>(
"Initial Condition Consistency Check",
2426 typedef Teuchos::ScalarTraits<Scalar> ST;
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();
2438 auto emptyBStar = Teuchos::SerialDenseVector<int,Scalar>();
2440 this->
getStepperType(),A,b,c,order,order,order,emptyBStar,
false));
2466 template<
class Scalar>
2473 this->
setStepperType(
"RK Implicit 2 Stage 2nd order Lobatto IIIB");
2479 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2481 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2483 std::string ICConsistency,
2484 bool ICConsistencyCheck,
2486 bool zeroInitialGuess)
2488 this->
setStepperType(
"RK Implicit 2 Stage 2nd order Lobatto IIIB");
2490 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2491 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2496 std::ostringstream Description;
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"
2505 <<
"A = [ 1/2 0 ]\n"
2507 <<
"b = [ 1/2 1/2 ]'";
2508 return Description.str();
2513 Teuchos::RCP<const Teuchos::ParameterList>
2516 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2518 pl->set<
bool>(
"Initial Condition Consistency Check",
2527 typedef Teuchos::ScalarTraits<Scalar> ST;
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();
2537 A(0,0) = as<Scalar>( one/(2*one) );
2539 A(1,0) = as<Scalar>( one/(2*one) );
2543 b(0) = as<Scalar>( one/(2*one) );
2544 b(1) = as<Scalar>( one/(2*one) );
2551 auto emptyBStar = Teuchos::SerialDenseVector<int,Scalar>();
2553 this->
getStepperType(),A,b,c,order,order,order,emptyBStar,
false));
2583 template<
class Scalar>
2596 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2598 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2600 std::string ICConsistency,
2601 bool ICConsistencyCheck,
2603 bool zeroInitialGuess)
2607 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2608 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2613 std::ostringstream Description;
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"
2621 <<
"c = [ 1/4 3/4 11/20 1/2 1 ]'\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 ]'";
2629 return Description.str();
2634 Teuchos::RCP<const Teuchos::ParameterList>
2637 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2639 pl->set<
bool>(
"Initial Condition Consistency Check",
2648 typedef Teuchos::ScalarTraits<Scalar> ST;
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) );
2659 A(0,0) = onequarter;
2665 A(1,0) = as<Scalar>( one / (2*one) );
2666 A(1,1) = onequarter;
2671 A(2,0) = as<Scalar>( 17*one/(50*one) );
2672 A(2,1) = as<Scalar>( -one/(25*one) );
2673 A(2,2) = onequarter;
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;
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;
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) );
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) );
2743 template<
class Scalar>
2756 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2758 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2760 std::string ICConsistency,
2761 bool ICConsistencyCheck,
2763 bool zeroInitialGuess)
2767 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2768 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2773 std::ostringstream Description;
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"
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();
2793 Teuchos::RCP<const Teuchos::ParameterList>
2796 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2798 pl->set<
bool>(
"Initial Condition Consistency Check",
2807 typedef Teuchos::ScalarTraits<Scalar> ST;
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)) );
2824 A(1,0) = as<Scalar>( one/(2*one) - gamma );
2828 A(2,0) = as<Scalar>( 2*gamma );
2829 A(2,1) = as<Scalar>( one - 4*gamma );
2834 b(1) = as<Scalar>( one-2*delta );
2839 c(1) = as<Scalar>( one/(2*one) );
2840 c(2) = as<Scalar>( one - gamma );
2874 template<
class Scalar>
2887 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
2889 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
2891 std::string ICConsistency,
2892 bool ICConsistencyCheck,
2894 bool zeroInitialGuess)
2898 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
2899 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
2904 std::ostringstream Description;
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"
2911 <<
"c = [ (6-sqrt(6))/10 ]\n"
2912 <<
" [ (6+9*sqrt(6))/35 ]\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"
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"
2929 <<
" [ (6-sqrt(6))/10 ]\n"
2930 <<
" [ (-267+88*sqrt(6))/500 ]\n"
2931 <<
" [ (1329-544*sqrt(6))/2500 ]\n"
2935 <<
" [ (6-sqrt(6))/10 ]\n"
2936 <<
" [ (-96+131*sqrt(6))/625 ]\n"
2941 <<
" [ (6-sqrt(6))/10 ]\n"
2945 <<
" [ (16-sqrt(6))/36 ]\n"
2946 <<
" [ (16+sqrt(6))/36 ]'";
2947 return Description.str();
2952 Teuchos::RCP<const Teuchos::ParameterList>
2955 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
2957 pl->set<
bool>(
"Initial Condition Consistency Check",
2966 typedef Teuchos::ScalarTraits<Scalar> ST;
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) );
2984 A(1,0) = as<Scalar>( (-6*one+5*one*sqrt6)/(14*one) );
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) );
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) );
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) );
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) );
3017 c(1) = as<Scalar>( (6*one+9*one*sqrt6)/(35*one) );
3019 c(3) = as<Scalar>( (4*one-sqrt6)/(10*one) );
3020 c(4) = as<Scalar>( (4*one+sqrt6)/(10*one) );
3047 template<
class Scalar>
3060 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
3062 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
3064 std::string ICConsistency,
3065 bool ICConsistencyCheck,
3067 bool zeroInitialGuess)
3071 this->
setup(appModel, obs, solver, useFSAL, ICConsistency,
3072 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
3077 std::ostringstream Description;
3082 <<
"b = [ 1/2 1/2 ]'\n"
3083 <<
"bstar = [ 1 0 ]'";
3084 return Description.str();
3089 Teuchos::RCP<const Teuchos::ParameterList>
3092 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
3094 pl->set<
bool>(
"Initial Condition Consistency Check",
3103 typedef Teuchos::ScalarTraits<Scalar> ST;
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);
3111 const Scalar one = ST::one();
3112 const Scalar zero = ST::zero();
3115 A(0,0) = one; A(0,1) = zero;
3116 A(1,0) = -one; A(1,1) = one;
3119 b(0) = as<Scalar>(one/(2*one));
3120 b(1) = as<Scalar>(one/(2*one));
3168 template<
class Scalar>
3181 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
3183 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
3185 std::string ICConsistency,
3186 bool ICConsistencyCheck,
3188 bool zeroInitialGuess,
3189 const Teuchos::SerialDenseMatrix<int,Scalar>& A,
3190 const Teuchos::SerialDenseVector<int,Scalar>& b,
3191 const Teuchos::SerialDenseVector<int,Scalar>& c,
3195 const Teuchos::SerialDenseVector<int,Scalar>& bstar)
3198 this->
setTableau(A,b,c,order,orderMin,orderMax,bstar);
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");
3204 this->
setup(appModel, obs, useFSAL, ICConsistency,
3205 ICConsistencyCheck, useEmbedded, zeroInitialGuess);
3210 std::stringstream Description;
3212 <<
"The format of the Butcher Tableau parameter list is\n"
3213 <<
" <Parameter name=\"A\" type=\"string\" value=\"# # # ;\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"
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"
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();
3237 if (this->
tableau_ == Teuchos::null) {
3240 auto t = stepper->getTableau();
3243 t->A(),t->b(),t->c(),
3244 t->order(),t->orderMin(),t->orderMax(),
3250 const Teuchos::SerialDenseVector<int,Scalar>& b,
3251 const Teuchos::SerialDenseVector<int,Scalar>& c,
3255 const Teuchos::SerialDenseVector<int,Scalar>&
3256 bstar = Teuchos::SerialDenseVector<int,Scalar>())
3262 Teuchos::RCP<const Teuchos::ParameterList>
3265 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
3267 pl->set<
bool>(
"Initial Condition Consistency Check",
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);
3289 #endif // Tempus_StepperRKButcherTableau_hpp
General Implicit Runge-Kutta Butcher Tableau.
Teuchos::RCP< RKButcherTableau< Scalar > > tableau_
StepperERK_3Stage3rdOrder()
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperERK_3Stage3rdOrderTVD()
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.
StepperSDIRK_ImplicitMidpoint()
Backward Euler Runge-Kutta Butcher Tableau.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperSDIRK_2Stage2ndOrder()
General Explicit Runge-Kutta Butcher Tableau.
virtual bool getICConsistencyCheckDefault() const
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.
std::string getDescription() const
StepperDIRK_1StageTheta()
std::string getDescription() const
StepperEDIRK_2Stage3rdOrder()
RK Explicit 4 Stage 3rd order by Runge.
std::string getDescription() const
void getValidParametersBasicDIRK(Teuchos::RCP< Teuchos::ParameterList > pl) const
std::string gammaTypeDefault_
std::string getDescription() const
std::string getDescription() const
virtual bool getICConsistencyCheckDefault() 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)
StepperERK_5Stage3rdOrderKandG()
RK Trapezoidal Rule (A.K.A. RK Crank-Nicolson)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperERK_BogackiShampine32()
std::string getDescription() const
std::string getDescription() 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.
void setGamma(Scalar gamma)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual bool getICConsistencyCheckDefault() 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)
virtual bool getICConsistencyCheckDefault() const
Explicit RK 3/8th Rule Butcher Tableau.
std::string getDescription() const
void setTheta(Scalar theta)
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)
void setTheta(Scalar theta)
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)
std::string getDescription() const
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)
std::string getDescription() const
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)
StepperSDIRK_5Stage4thOrder()
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)
std::string getDescription() const
virtual void setupDefault()
Default setup for constructor.
Explicit RK Bogacki-Shampine Butcher Tableau.
EDIRK 2 Stage Theta Method.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
std::string getDescription() const
virtual bool getICConsistencyCheckDefault() const
std::string getDescription() const
virtual bool getICConsistencyCheckDefault() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
std::string getDescription() 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)
std::string getDescription() const
virtual bool getICConsistencyCheckDefault() const
std::string getDescription() const
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
std::string getDescription() const
StepperSDIRK_2Stage3rdOrder()
std::string getDescription() const
virtual bool getICConsistencyCheckDefault() const
std::string getDescription() const
StepperDIRK_1Stage1stOrderRadauIA()
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 >())
RK Explicit 3 Stage 3rd order TVD.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
std::string getDescription() const
StepperEDIRK_2StageTheta()
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)
virtual bool getICConsistencyCheckDefault() const
std::string getDescription() const
std::string getDescription() const
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.
StepperDIRK_BackwardEuler()
virtual std::string getDescription() const
virtual bool getICConsistencyCheckDefault() const
void setGammaType(std::string gammaType)
StepperERK_3Stage3rdOrderHeun()
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.
This observer is a composite observer,.
virtual bool getICConsistencyCheckDefault() const
virtual bool getICConsistencyCheckDefault() const
std::string getDescription() const
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
RK Explicit 3 Stage 3rd order.
StepperERK_4Stage4thOrder()
void setGamma(Scalar gamma)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
std::string getDescription() const
StepperSDIRK_5Stage5thOrder()
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 >())
RK Implicit 1 Stage 1st order Radau IA.
StepperSDIRK_3Stage4thOrder()
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperDIRK_2Stage2ndOrderLobattoIIIB()
virtual bool getICConsistencyCheckDefault() 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))
StepperERK_ForwardEuler()
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
std::string getDescription() 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
std::string getDescription() 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)
virtual bool getICConsistencyCheckDefault() const
StepperEDIRK_TrapezoidalRule()
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))
StepperERK_4Stage3rdOrderRunge()
std::string getDescription() const