9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
13 #include "Thyra_VectorStdOps.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
18 #include "../TestModels/SinCosModel.hpp"
19 #include "../TestModels/VanDerPolModel.hpp"
20 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
29 using Teuchos::rcp_const_cast;
30 using Teuchos::rcp_dynamic_cast;
31 using Teuchos::ParameterList;
32 using Teuchos::parameterList;
33 using Teuchos::sublist;
34 using Teuchos::getParametersFromXmlFile;
45 std::vector<std::string> RKMethods;
46 RKMethods.push_back(
"General DIRK");
47 RKMethods.push_back(
"RK Backward Euler");
48 RKMethods.push_back(
"DIRK 1 Stage Theta Method");
49 RKMethods.push_back(
"RK Implicit 1 Stage 1st order Radau IA");
50 RKMethods.push_back(
"RK Implicit Midpoint");
51 RKMethods.push_back(
"SDIRK 2 Stage 2nd order");
52 RKMethods.push_back(
"RK Implicit 2 Stage 2nd order Lobatto IIIB");
53 RKMethods.push_back(
"SDIRK 2 Stage 3rd order");
54 RKMethods.push_back(
"EDIRK 2 Stage 3rd order");
55 RKMethods.push_back(
"EDIRK 2 Stage Theta Method");
56 RKMethods.push_back(
"SDIRK 3 Stage 4th order");
57 RKMethods.push_back(
"SDIRK 5 Stage 4th order");
58 RKMethods.push_back(
"SDIRK 5 Stage 5th order");
59 RKMethods.push_back(
"SDIRK 2(1) Pair");
60 RKMethods.push_back(
"RK Trapezoidal Rule");
61 RKMethods.push_back(
"RK Crank-Nicolson");
63 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
65 std::string RKMethod = RKMethods[m];
68 RCP<ParameterList> pList =
69 getParametersFromXmlFile(
"Tempus_DIRK_SinCos.xml");
72 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
75 RCP<ParameterList> tempusPL = sublist(pList,
"Tempus",
true);
76 tempusPL->sublist(
"Default Stepper").set(
"Stepper Type", RKMethods[m]);
78 if (RKMethods[m] ==
"DIRK 1 Stage Theta Method" ||
79 RKMethods[m] ==
"EDIRK 2 Stage Theta Method") {
81 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
82 RCP<ParameterList> solverPL = parameterList();
83 *solverPL = *(sublist(stepperPL,
"Default Solver",
true));
84 tempusPL->sublist(
"Default Stepper").remove(
"Zero Initial Guess");
85 tempusPL->sublist(
"Default Stepper").remove(
"Default Solver");
86 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Zero Initial Guess", 0);
87 tempusPL->sublist(
"Default Stepper").remove(
"Reset Initial Guess");
88 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Reset Initial Guess", 1);
89 tempusPL->sublist(
"Default Stepper").set(
"Default Solver", *solverPL);
90 tempusPL->sublist(
"Default Stepper").set<
double>(
"theta", 0.5);
91 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 2nd order") {
93 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
94 RCP<ParameterList> solverPL = parameterList();
95 *solverPL = *(sublist(stepperPL,
"Default Solver",
true));
96 tempusPL->sublist(
"Default Stepper").remove(
"Default Solver");
97 tempusPL->sublist(
"Default Stepper").remove(
"Zero Initial Guess");
98 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Zero Initial Guess", 0);
99 tempusPL->sublist(
"Default Stepper").remove(
"Reset Initial Guess");
100 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Reset Initial Guess", 1);
101 tempusPL->sublist(
"Default Stepper").set(
"Default Solver", *solverPL);
102 tempusPL->sublist(
"Default Stepper")
103 .set<
double>(
"gamma", 0.2928932188134524);
104 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 3rd order") {
106 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
107 RCP<ParameterList> solverPL = parameterList();
108 *solverPL = *(sublist(stepperPL,
"Default Solver",
true));
109 tempusPL->sublist(
"Default Stepper").remove(
"Zero Initial Guess");
110 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Zero Initial Guess", 0);
111 tempusPL->sublist(
"Default Stepper").remove(
"Default Solver");
112 tempusPL->sublist(
"Default Stepper").remove(
"Reset Initial Guess");
113 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Reset Initial Guess", 1);
114 tempusPL->sublist(
"Default Stepper").set(
"Default Solver", *solverPL);
115 tempusPL->sublist(
"Default Stepper")
116 .set<std::string>(
"Gamma Type",
"3rd Order A-stable");
117 tempusPL->sublist(
"Default Stepper")
118 .set<
double>(
"gamma", 0.7886751345948128);
119 }
else if (RKMethods[m] ==
"RK Crank-Nicolson") {
121 tempusPL->sublist(
"Default Stepper")
122 .set(
"Stepper Type",
"RK Trapezoidal Rule");
123 }
else if (RKMethods[m] ==
"General DIRK") {
125 Teuchos::RCP<Teuchos::ParameterList> tableauPL = Teuchos::parameterList();
126 tableauPL->set<std::string>(
"A",
"0.2928932188134524 0.0; 0.7071067811865476 0.2928932188134524");
127 tableauPL->set<std::string>(
"b",
"0.7071067811865476 0.2928932188134524");
128 tableauPL->set<std::string>(
"c",
"0.2928932188134524 1.0");
129 tableauPL->set<
int>(
"order", 2);
130 tableauPL->set<std::string>(
"bstar",
"");
131 tempusPL->sublist(
"Default Stepper").set(
"Tableau", *tableauPL);
137 RCP<Tempus::IntegratorBasic<double> > integrator =
138 Tempus::integratorBasic<double>(tempusPL, model);
140 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
141 RCP<ParameterList> defaultPL =
142 Teuchos::rcp_const_cast<Teuchos::ParameterList>(
143 integrator->getStepper()->getValidParameters());
144 defaultPL->remove(
"Description");
146 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
148 std::cout << std::endl;
149 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
150 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
157 RCP<Tempus::IntegratorBasic<double> > integrator =
158 Tempus::integratorBasic<double>(model, RKMethods[m]);
160 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
161 RCP<ParameterList> defaultPL =
162 Teuchos::rcp_const_cast<Teuchos::ParameterList>(
163 integrator->getStepper()->getValidParameters());
164 defaultPL->remove(
"Description");
166 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
168 std::cout << std::endl;
169 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
170 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
185 RCP<ParameterList> pList =
186 getParametersFromXmlFile(
"Tempus_DIRK_SinCos.xml");
187 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
190 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
195 RCP<Tempus::StepperFactory<double> > sf =
197 RCP<Tempus::Stepper<double> > stepper =
198 sf->createStepper(
"SDIRK 2 Stage 2nd order");
199 stepper->setModel(model);
200 stepper->initialize();
204 ParameterList tscPL = pl->sublist(
"Default Integrator")
205 .sublist(
"Time Step Control");
206 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
207 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
208 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
209 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
210 timeStepControl->setInitTimeStep(dt);
211 timeStepControl->initialize();
214 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
215 stepper->getModel()->getNominalValues();
216 auto icSolution = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
218 icState->setTime (timeStepControl->getInitTime());
219 icState->setIndex (timeStepControl->getInitIndex());
220 icState->setTimeStep(0.0);
221 icState->setOrder (stepper->getOrder());
226 solutionHistory->setName(
"Forward States");
228 solutionHistory->setStorageLimit(2);
229 solutionHistory->addState(icState);
232 RCP<Tempus::IntegratorBasic<double> > integrator =
233 Tempus::integratorBasic<double>();
234 integrator->setStepperWStepper(stepper);
235 integrator->setTimeStepControl(timeStepControl);
236 integrator->setSolutionHistory(solutionHistory);
238 integrator->initialize();
242 bool integratorStatus = integrator->advanceTime();
243 TEST_ASSERT(integratorStatus)
247 double time = integrator->getTime();
248 double timeFinal =pl->sublist(
"Default Integrator")
249 .sublist(
"Time Step Control").get<
double>(
"Final Time");
250 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
253 RCP<Thyra::VectorBase<double> > x = integrator->getX();
254 RCP<const Thyra::VectorBase<double> > x_exact =
255 model->getExactSolution(time).get_x();
258 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
259 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
262 std::cout <<
" Stepper = SDIRK 2 Stage 2nd order" << std::endl;
263 std::cout <<
" =========================" << std::endl;
264 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
265 << get_ele(*(x_exact), 1) << std::endl;
266 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
267 << get_ele(*(x ), 1) << std::endl;
268 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" "
269 << get_ele(*(xdiff ), 1) << std::endl;
270 std::cout <<
" =========================" << std::endl;
271 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841470, 1.0e-4 );
272 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.540304, 1.0e-4 );
280 std::vector<std::string> RKMethods;
281 RKMethods.push_back(
"General DIRK");
282 RKMethods.push_back(
"RK Backward Euler");
283 RKMethods.push_back(
"DIRK 1 Stage Theta Method");
284 RKMethods.push_back(
"RK Implicit 1 Stage 1st order Radau IA");
285 RKMethods.push_back(
"RK Implicit Midpoint");
286 RKMethods.push_back(
"SDIRK 2 Stage 2nd order");
287 RKMethods.push_back(
"RK Implicit 2 Stage 2nd order Lobatto IIIB");
288 RKMethods.push_back(
"SDIRK 2 Stage 3rd order");
289 RKMethods.push_back(
"EDIRK 2 Stage 3rd order");
290 RKMethods.push_back(
"EDIRK 2 Stage Theta Method");
291 RKMethods.push_back(
"SDIRK 3 Stage 4th order");
292 RKMethods.push_back(
"SDIRK 5 Stage 4th order");
293 RKMethods.push_back(
"SDIRK 5 Stage 5th order");
294 RKMethods.push_back(
"SDIRK 2(1) Pair");
295 RKMethods.push_back(
"RK Trapezoidal Rule");
296 RKMethods.push_back(
"RK Crank-Nicolson");
297 RKMethods.push_back(
"SSPDIRK22");
298 RKMethods.push_back(
"SSPDIRK32");
299 RKMethods.push_back(
"SSPDIRK23");
300 RKMethods.push_back(
"SSPDIRK33");
301 RKMethods.push_back(
"SDIRK 3 Stage 2nd order");
303 std::vector<double> RKMethodErrors;
304 RKMethodErrors.push_back(2.52738e-05);
305 RKMethodErrors.push_back(0.0124201);
306 RKMethodErrors.push_back(5.20785e-05);
307 RKMethodErrors.push_back(0.0124201);
308 RKMethodErrors.push_back(5.20785e-05);
309 RKMethodErrors.push_back(2.52738e-05);
310 RKMethodErrors.push_back(5.20785e-05);
311 RKMethodErrors.push_back(1.40223e-06);
312 RKMethodErrors.push_back(2.17004e-07);
313 RKMethodErrors.push_back(5.20785e-05);
314 RKMethodErrors.push_back(6.41463e-08);
315 RKMethodErrors.push_back(3.30631e-10);
316 RKMethodErrors.push_back(1.35728e-11);
317 RKMethodErrors.push_back(0.0001041);
318 RKMethodErrors.push_back(5.20785e-05);
319 RKMethodErrors.push_back(5.20785e-05);
320 RKMethodErrors.push_back(1.30205e-05);
321 RKMethodErrors.push_back(5.7869767e-06);
322 RKMethodErrors.push_back(1.00713e-07);
323 RKMethodErrors.push_back(3.94916e-08);
324 RKMethodErrors.push_back(2.52738e-05);
326 TEUCHOS_ASSERT( RKMethods.size() == RKMethodErrors.size() );
328 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
330 std::string RKMethod = RKMethods[m];
331 std::replace(RKMethod.begin(), RKMethod.end(),
' ',
'_');
332 std::replace(RKMethod.begin(), RKMethod.end(),
'/',
'.');
334 RCP<Tempus::IntegratorBasic<double> > integrator;
335 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
336 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
337 std::vector<double> StepSize;
338 std::vector<double> xErrorNorm;
339 std::vector<double> xDotErrorNorm;
341 const int nTimeStepSizes = 2;
344 for (
int n=0; n<nTimeStepSizes; n++) {
347 RCP<ParameterList> pList =
348 getParametersFromXmlFile(
"Tempus_DIRK_SinCos.xml");
351 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
355 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
356 pl->sublist(
"Default Stepper").set(
"Stepper Type", RKMethods[m]);
357 if (RKMethods[m] ==
"DIRK 1 Stage Theta Method" ||
358 RKMethods[m] ==
"EDIRK 2 Stage Theta Method") {
359 pl->sublist(
"Default Stepper").set<
double>(
"theta", 0.5);
360 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 2nd order") {
361 pl->sublist(
"Default Stepper").set(
"gamma", 0.2928932188134524);
362 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 3rd order") {
363 pl->sublist(
"Default Stepper")
364 .set<std::string>(
"Gamma Type",
"3rd Order A-stable");
370 pl->sublist(
"Default Integrator")
371 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
372 integrator = Tempus::integratorBasic<double>(pl, model);
378 RCP<Thyra::VectorBase<double> > x0 =
379 model->getNominalValues().get_x()->clone_v();
380 integrator->initializeSolutionHistory(0.0, x0);
383 bool integratorStatus = integrator->advanceTime();
384 TEST_ASSERT(integratorStatus)
387 time = integrator->getTime();
388 double timeFinal = pl->sublist(
"Default Integrator")
389 .sublist(
"Time Step Control").get<
double>(
"Final Time");
390 double tol = 100.0 * std::numeric_limits<double>::epsilon();
391 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
395 RCP<const SolutionHistory<double> > solutionHistory =
396 integrator->getSolutionHistory();
397 writeSolution(
"Tempus_"+RKMethod+
"_SinCos.dat", solutionHistory);
400 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
401 double time_i = (*solutionHistory)[i]->getTime();
403 rcp_const_cast<Thyra::VectorBase<double> > (
404 model->getExactSolution(time_i).get_x()),
405 rcp_const_cast<Thyra::VectorBase<double> > (
406 model->getExactSolution(time_i).get_x_dot()));
407 state->setTime((*solutionHistory)[i]->getTime());
408 solnHistExact->addState(state);
410 writeSolution(
"Tempus_"+RKMethod+
"_SinCos-Ref.dat", solnHistExact);
414 StepSize.push_back(dt);
415 auto solution = Thyra::createMember(model->get_x_space());
416 Thyra::copy(*(integrator->getX()),solution.ptr());
417 solutions.push_back(solution);
418 auto solutionDot = Thyra::createMember(model->get_x_space());
419 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
420 solutionsDot.push_back(solutionDot);
421 if (n == nTimeStepSizes-1) {
422 StepSize.push_back(0.0);
423 auto solutionExact = Thyra::createMember(model->get_x_space());
424 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
425 solutions.push_back(solutionExact);
426 auto solutionDotExact = Thyra::createMember(model->get_x_space());
427 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
428 solutionDotExact.ptr());
429 solutionsDot.push_back(solutionDotExact);
435 double xDotSlope = 0.0;
436 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
437 double order = stepper->getOrder();
440 solutions, xErrorNorm, xSlope,
441 solutionsDot, xDotErrorNorm, xDotSlope);
443 TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
444 TEST_FLOATING_EQUALITY( xErrorNorm[0], RKMethodErrors[m], 5.0e-4 );
457 std::vector<std::string> RKMethods;
458 RKMethods.push_back(
"SDIRK 2 Stage 2nd order");
460 std::string RKMethod = RKMethods[0];
461 std::replace(RKMethod.begin(), RKMethod.end(),
' ',
'_');
462 std::replace(RKMethod.begin(), RKMethod.end(),
'/',
'.');
464 RCP<Tempus::IntegratorBasic<double> > integrator;
465 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
466 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
467 std::vector<double> StepSize;
468 std::vector<double> xErrorNorm;
469 std::vector<double> xDotErrorNorm;
471 const int nTimeStepSizes = 3;
474 for (
int n=0; n<nTimeStepSizes; n++) {
477 RCP<ParameterList> pList =
478 getParametersFromXmlFile(
"Tempus_DIRK_VanDerPol.xml");
481 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
485 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
486 pl->sublist(
"Default Stepper").set(
"Stepper Type", RKMethods[0]);
487 pl->sublist(
"Default Stepper").set(
"gamma", 0.2928932188134524);
491 if (n == nTimeStepSizes-1) dt /= 10.0;
494 pl->sublist(
"Default Integrator")
495 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
496 integrator = Tempus::integratorBasic<double>(pl, model);
499 bool integratorStatus = integrator->advanceTime();
500 TEST_ASSERT(integratorStatus)
503 time = integrator->getTime();
504 double timeFinal =pl->sublist(
"Default Integrator")
505 .sublist(
"Time Step Control").get<
double>(
"Final Time");
506 double tol = 100.0 * std::numeric_limits<double>::epsilon();
507 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
510 StepSize.push_back(dt);
511 auto solution = Thyra::createMember(model->get_x_space());
512 Thyra::copy(*(integrator->getX()),solution.ptr());
513 solutions.push_back(solution);
514 auto solutionDot = Thyra::createMember(model->get_x_space());
515 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
516 solutionsDot.push_back(solutionDot);
520 if ((n == 0) or (n == nTimeStepSizes-1)) {
521 std::string fname =
"Tempus_"+RKMethod+
"_VanDerPol-Ref.dat";
522 if (n == 0) fname =
"Tempus_"+RKMethod+
"_VanDerPol.dat";
523 RCP<const SolutionHistory<double> > solutionHistory =
524 integrator->getSolutionHistory();
531 double xDotSlope = 0.0;
532 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
533 double order = stepper->getOrder();
536 solutions, xErrorNorm, xSlope,
537 solutionsDot, xDotErrorNorm, xDotSlope);
539 TEST_FLOATING_EQUALITY( xSlope, order, 0.06 );
540 TEST_FLOATING_EQUALITY( xErrorNorm[0], 1.07525e-05, 1.0e-4 );
545 Teuchos::TimeMonitor::summarize();
554 std::vector<std::string> IntegratorList;
555 IntegratorList.push_back(
"Embedded_Integrator_PID");
556 IntegratorList.push_back(
"Embedded_Integrator");
559 const int refIstep = 217;
561 for(
auto integratorChoice : IntegratorList){
563 std::cout <<
"Using Integrator: " << integratorChoice <<
" !!!" << std::endl;
566 RCP<ParameterList> pList =
567 getParametersFromXmlFile(
"Tempus_DIRK_VanDerPol.xml");
571 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
575 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
576 pl->set(
"Integrator Name", integratorChoice);
579 RCP<Tempus::IntegratorBasic<double> > integrator =
580 Tempus::integratorBasic<double>(pl, model);
582 const std::string RKMethod_ =
583 pl->sublist(integratorChoice).get<std::string>(
"Stepper Name");
586 bool integratorStatus = integrator->advanceTime();
587 TEST_ASSERT(integratorStatus);
590 double time = integrator->getTime();
591 double timeFinal = pl->sublist(integratorChoice)
592 .sublist(
"Time Step Control").get<
double>(
"Final Time");
593 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
597 RCP<Thyra::VectorBase<double> > x = integrator->getX();
598 RCP<Thyra::VectorBase<double> > xref = x->clone_v();
599 Thyra::set_ele(0, -1.5484458614405929, xref.ptr());
600 Thyra::set_ele(1, 1.0181127316101317, xref.ptr());
603 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
604 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *xref, -1.0, *(x));
605 const double L2norm = Thyra::norm_2(*xdiff);
608 if (integratorChoice ==
"Embedded_Integrator_PID"){
609 const double absTol = pl->sublist(integratorChoice).
610 sublist(
"Time Step Control").get<
double>(
"Maximum Absolute Error");
611 const double relTol = pl->sublist(integratorChoice).
612 sublist(
"Time Step Control").get<
double>(
"Maximum Relative Error");
618 const int iStep = integrator->getSolutionHistory()->
619 getCurrentState()->getIndex();
624 TEST_FLOATING_EQUALITY(std::log10(L2norm),std::log10(absTol), 0.3 );
625 TEST_FLOATING_EQUALITY(std::log10(L2norm),std::log10(relTol), 0.3 );
627 TEST_COMPARE(iStep, <=, refIstep);
631 std::ofstream ftmp(
"Tempus_"+integratorChoice+RKMethod_+
"_VDP_Example.dat");
632 RCP<const SolutionHistory<double> > solutionHistory =
633 integrator->getSolutionHistory();
634 int nStates = solutionHistory->getNumStates();
636 for (
int i=0; i<nStates; i++) {
637 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
638 double time_i = solutionState->getTime();
639 RCP<const Thyra::VectorBase<double> > x_plot = solutionState->getX();
641 ftmp << time_i <<
" "
642 << Thyra::get_ele(*(x_plot), 0) <<
" "
643 << Thyra::get_ele(*(x_plot), 1) <<
" " << std::endl;
648 Teuchos::TimeMonitor::summarize();
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
van der Pol model problem for nonlinear electrical circuit.
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > solutionHistory)
void writeOrderError(const std::string filename, Teuchos::RCP< Tempus::Stepper< Scalar > > stepper, std::vector< Scalar > &StepSize, std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar >>> &solutions, std::vector< Scalar > &xErrorNorm, Scalar &xSlope, std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar >>> &solutionsDot, std::vector< Scalar > &xDotErrorNorm, Scalar &xDotSlope, std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar >>> &solutionsDotDot, std::vector< Scalar > &xDotDotErrorNorm, Scalar &xDotDotSlope)
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
@ STORAGE_TYPE_STATIC
Keep a fix number of states.
Teuchos::RCP< SolutionState< Scalar > > createSolutionStateX(const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdot=Teuchos::null, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdotdot=Teuchos::null)
Nonmember constructor from non-const solution vectors, x.