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" 16 #include "Tempus_StepperFactory.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 if (RKMethods[m] ==
"EDIRK 2 Stage Theta Method")
85 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Use FSAL", 1);
86 tempusPL->sublist(
"Default Stepper").remove(
"Zero Initial Guess");
87 tempusPL->sublist(
"Default Stepper").remove(
"Default Solver");
88 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Zero Initial Guess", 0);
89 tempusPL->sublist(
"Default Stepper").remove(
"Reset Initial Guess");
90 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Reset Initial Guess", 1);
91 tempusPL->sublist(
"Default Stepper").set(
"Default Solver", *solverPL);
92 tempusPL->sublist(
"Default Stepper").set<
double>(
"theta", 0.5);
93 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 2nd order") {
95 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
96 RCP<ParameterList> solverPL = parameterList();
97 *solverPL = *(sublist(stepperPL,
"Default Solver",
true));
98 tempusPL->sublist(
"Default Stepper").remove(
"Default Solver");
99 tempusPL->sublist(
"Default Stepper").remove(
"Zero Initial Guess");
100 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Zero Initial Guess", 0);
101 tempusPL->sublist(
"Default Stepper").remove(
"Reset Initial Guess");
102 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Reset Initial Guess", 1);
103 tempusPL->sublist(
"Default Stepper").set(
"Default Solver", *solverPL);
104 tempusPL->sublist(
"Default Stepper")
105 .set<
double>(
"gamma", 0.2928932188134524);
106 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 3rd order") {
108 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
109 RCP<ParameterList> solverPL = parameterList();
110 *solverPL = *(sublist(stepperPL,
"Default Solver",
true));
111 tempusPL->sublist(
"Default Stepper").remove(
"Zero Initial Guess");
112 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Zero Initial Guess", 0);
113 tempusPL->sublist(
"Default Stepper").remove(
"Default Solver");
114 tempusPL->sublist(
"Default Stepper").remove(
"Reset Initial Guess");
115 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Reset Initial Guess", 1);
116 tempusPL->sublist(
"Default Stepper").set(
"Default Solver", *solverPL);
117 tempusPL->sublist(
"Default Stepper")
118 .set<std::string>(
"Gamma Type",
"3rd Order A-stable");
119 tempusPL->sublist(
"Default Stepper")
120 .set<
double>(
"gamma", 0.7886751345948128);
121 }
else if (RKMethods[m] ==
"RK Trapezoidal Rule") {
122 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Use FSAL", 1);
123 }
else if (RKMethods[m] ==
"RK Crank-Nicolson") {
124 tempusPL->sublist(
"Default Stepper").set<
bool>(
"Use FSAL", 1);
126 tempusPL->sublist(
"Default Stepper")
127 .set(
"Stepper Type",
"RK Trapezoidal Rule");
128 }
else if (RKMethods[m] ==
"General DIRK") {
130 Teuchos::RCP<Teuchos::ParameterList> tableauPL = Teuchos::parameterList();
131 tableauPL->set<std::string>(
"A",
"0.292893218813452 0; 0.707106781186548 0.292893218813452");
132 tableauPL->set<std::string>(
"b",
"0.707106781186548 0.292893218813452");
133 tableauPL->set<std::string>(
"c",
"0.292893218813452 1");
134 tableauPL->set<
int>(
"order", 2);
135 tableauPL->set<std::string>(
"bstar",
"");
136 tempusPL->sublist(
"Default Stepper").set(
"Tableau", *tableauPL);
142 RCP<Tempus::IntegratorBasic<double> > integrator =
143 Tempus::createIntegratorBasic<double>(tempusPL, model);
145 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
146 RCP<ParameterList> defaultPL =
147 Teuchos::rcp_const_cast<Teuchos::ParameterList>(
148 integrator->getStepper()->getValidParameters());
150 bool pass = haveSameValuesSorted(*stepperPL, *defaultPL,
true);
152 std::cout << std::endl;
153 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
154 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
161 RCP<Tempus::IntegratorBasic<double> > integrator =
162 Tempus::createIntegratorBasic<double>(model, RKMethods[m]);
164 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
165 RCP<ParameterList> defaultPL =
166 Teuchos::rcp_const_cast<Teuchos::ParameterList>(
167 integrator->getStepper()->getValidParameters());
173 if (RKMethods[m] ==
"EDIRK 2 Stage Theta Method" ||
174 RKMethods[m] ==
"RK Trapezoidal Rule" ||
175 RKMethods[m] ==
"RK Crank-Nicolson") {
176 stepperPL->set<std::string>(
"Initial Condition Consistency",
"Consistent");
177 stepperPL->remove(
"Default Solver");
178 defaultPL->remove(
"Default Solver");
181 bool pass = haveSameValuesSorted(*stepperPL, *defaultPL,
true);
183 std::cout << std::endl;
184 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
185 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
198 std::vector<std::string> options;
199 options.push_back(
"Default Parameters");
200 options.push_back(
"ICConsistency and Check");
202 for(
const auto& option: options) {
205 RCP<ParameterList> pList =
206 getParametersFromXmlFile(
"Tempus_DIRK_SinCos.xml");
207 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
210 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
215 RCP<Tempus::StepperFactory<double> > sf =
217 RCP<Tempus::Stepper<double> > stepper =
218 sf->createStepper(
"SDIRK 2 Stage 2nd order");
219 stepper->setModel(model);
220 if ( option ==
"ICConsistency and Check") {
221 stepper->setICConsistency(
"Consistent");
222 stepper->setICConsistencyCheck(
true);
224 stepper->initialize();
228 ParameterList tscPL = pl->sublist(
"Default Integrator")
229 .sublist(
"Time Step Control");
230 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
231 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
232 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
233 timeStepControl->setInitTimeStep(dt);
234 timeStepControl->initialize();
237 auto inArgsIC = model->getNominalValues();
240 icState->setTime (timeStepControl->getInitTime());
241 icState->setIndex (timeStepControl->getInitIndex());
242 icState->setTimeStep(0.0);
243 icState->setOrder (stepper->getOrder());
248 solutionHistory->setName(
"Forward States");
250 solutionHistory->setStorageLimit(2);
251 solutionHistory->addState(icState);
254 RCP<Tempus::IntegratorBasic<double> > integrator =
255 Tempus::createIntegratorBasic<double>();
256 integrator->setStepper(stepper);
257 integrator->setTimeStepControl(timeStepControl);
258 integrator->setSolutionHistory(solutionHistory);
260 integrator->initialize();
264 bool integratorStatus = integrator->advanceTime();
265 TEST_ASSERT(integratorStatus)
269 double time = integrator->getTime();
270 double timeFinal =pl->sublist(
"Default Integrator")
271 .sublist(
"Time Step Control").get<
double>(
"Final Time");
272 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
275 RCP<Thyra::VectorBase<double> > x = integrator->getX();
276 RCP<const Thyra::VectorBase<double> > x_exact =
277 model->getExactSolution(time).get_x();
280 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
281 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
284 std::cout <<
" Stepper = " << stepper->description()
285 <<
" with " << option << std::endl;
286 std::cout <<
" =========================" << std::endl;
287 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" " 288 << get_ele(*(x_exact), 1) << std::endl;
289 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" " 290 << get_ele(*(x ), 1) << std::endl;
291 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" " 292 << get_ele(*(xdiff ), 1) << std::endl;
293 std::cout <<
" =========================" << std::endl;
294 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841470, 1.0e-4 );
295 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.540304, 1.0e-4 );
305 std::vector<std::string> RKMethods;
306 RKMethods.push_back(
"EDIRK 2 Stage Theta Method");
307 RKMethods.push_back(
"RK Trapezoidal Rule");
309 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
316 auto stepper = sf->createStepper(RKMethods[m]);
317 stepper->setModel(model);
318 stepper->setUseFSAL(
false);
319 stepper->initialize();
323 timeStepControl->setInitTime (0.0);
324 timeStepControl->setFinalTime(1.0);
325 timeStepControl->setInitTimeStep(dt);
326 timeStepControl->initialize();
329 auto inArgsIC = model->getNominalValues();
332 icState->setTime (timeStepControl->getInitTime());
333 icState->setIndex (timeStepControl->getInitIndex());
334 icState->setTimeStep(0.0);
335 icState->setOrder (stepper->getOrder());
340 solutionHistory->setName(
"Forward States");
342 solutionHistory->setStorageLimit(2);
343 solutionHistory->addState(icState);
346 auto integrator = Tempus::createIntegratorBasic<double>();
347 integrator->setStepper(stepper);
348 integrator->setTimeStepControl(timeStepControl);
349 integrator->setSolutionHistory(solutionHistory);
350 integrator->initialize();
354 bool integratorStatus = integrator->advanceTime();
355 TEST_ASSERT(integratorStatus)
359 double time = integrator->getTime();
360 TEST_FLOATING_EQUALITY(time, 1.0, 1.0e-14);
363 auto x = integrator->getX();
364 auto x_exact = model->getExactSolution(time).get_x();
367 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
368 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
371 std::cout <<
" Stepper = " << stepper->description()
372 <<
"\n with " <<
"useFSAL=false" << std::endl;
373 std::cout <<
" =========================" << std::endl;
374 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" " 375 << get_ele(*(x_exact), 1) << std::endl;
376 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" " 377 << get_ele(*(x ), 1) << std::endl;
378 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" " 379 << get_ele(*(xdiff ), 1) << std::endl;
380 std::cout <<
" =========================" << std::endl;
381 if (RKMethods[m] ==
"EDIRK 2 Stage Theta Method") {
382 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841021, 1.0e-4 );
383 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.541002, 1.0e-4 );
384 }
else if (RKMethods[m] ==
"RK Trapezoidal Rule") {
385 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841021, 1.0e-4 );
386 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.541002, 1.0e-4 );
396 std::vector<std::string> RKMethods;
397 RKMethods.push_back(
"General DIRK");
398 RKMethods.push_back(
"RK Backward Euler");
399 RKMethods.push_back(
"DIRK 1 Stage Theta Method");
400 RKMethods.push_back(
"RK Implicit 1 Stage 1st order Radau IA");
401 RKMethods.push_back(
"RK Implicit Midpoint");
402 RKMethods.push_back(
"SDIRK 2 Stage 2nd order");
403 RKMethods.push_back(
"RK Implicit 2 Stage 2nd order Lobatto IIIB");
404 RKMethods.push_back(
"SDIRK 2 Stage 3rd order");
405 RKMethods.push_back(
"EDIRK 2 Stage 3rd order");
406 RKMethods.push_back(
"EDIRK 2 Stage Theta Method");
407 RKMethods.push_back(
"SDIRK 3 Stage 4th order");
408 RKMethods.push_back(
"SDIRK 5 Stage 4th order");
409 RKMethods.push_back(
"SDIRK 5 Stage 5th order");
410 RKMethods.push_back(
"SDIRK 2(1) Pair");
411 RKMethods.push_back(
"RK Trapezoidal Rule");
412 RKMethods.push_back(
"RK Crank-Nicolson");
413 RKMethods.push_back(
"SSPDIRK22");
414 RKMethods.push_back(
"SSPDIRK32");
415 RKMethods.push_back(
"SSPDIRK23");
416 RKMethods.push_back(
"SSPDIRK33");
417 RKMethods.push_back(
"SDIRK 3 Stage 2nd order");
419 std::vector<double> RKMethodErrors;
420 RKMethodErrors.push_back(2.52738e-05);
421 RKMethodErrors.push_back(0.0124201);
422 RKMethodErrors.push_back(5.20785e-05);
423 RKMethodErrors.push_back(0.0124201);
424 RKMethodErrors.push_back(5.20785e-05);
425 RKMethodErrors.push_back(2.52738e-05);
426 RKMethodErrors.push_back(5.20785e-05);
427 RKMethodErrors.push_back(1.40223e-06);
428 RKMethodErrors.push_back(2.17004e-07);
429 RKMethodErrors.push_back(5.20785e-05);
430 RKMethodErrors.push_back(6.41463e-08);
431 RKMethodErrors.push_back(3.30631e-10);
432 RKMethodErrors.push_back(1.35728e-11);
433 RKMethodErrors.push_back(0.0001041);
434 RKMethodErrors.push_back(5.20785e-05);
435 RKMethodErrors.push_back(5.20785e-05);
436 RKMethodErrors.push_back(1.30205e-05);
437 RKMethodErrors.push_back(5.7869767e-06);
438 RKMethodErrors.push_back(1.00713e-07);
439 RKMethodErrors.push_back(3.94916e-08);
440 RKMethodErrors.push_back(2.52738e-05);
442 TEUCHOS_ASSERT( RKMethods.size() == RKMethodErrors.size() );
444 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
446 std::string RKMethod = RKMethods[m];
447 std::replace(RKMethod.begin(), RKMethod.end(),
' ',
'_');
448 std::replace(RKMethod.begin(), RKMethod.end(),
'/',
'.');
450 RCP<Tempus::IntegratorBasic<double> > integrator;
451 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
452 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
453 std::vector<double> StepSize;
454 std::vector<double> xErrorNorm;
455 std::vector<double> xDotErrorNorm;
457 const int nTimeStepSizes = 2;
460 for (
int n=0; n<nTimeStepSizes; n++) {
463 RCP<ParameterList> pList =
464 getParametersFromXmlFile(
"Tempus_DIRK_SinCos.xml");
467 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
471 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
472 pl->sublist(
"Default Stepper").set(
"Stepper Type", RKMethods[m]);
473 if (RKMethods[m] ==
"DIRK 1 Stage Theta Method" ||
474 RKMethods[m] ==
"EDIRK 2 Stage Theta Method") {
475 pl->sublist(
"Default Stepper").set<
double>(
"theta", 0.5);
476 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 2nd order") {
477 pl->sublist(
"Default Stepper").set(
"gamma", 0.2928932188134524);
478 }
else if (RKMethods[m] ==
"SDIRK 2 Stage 3rd order") {
479 pl->sublist(
"Default Stepper")
480 .set<std::string>(
"Gamma Type",
"3rd Order A-stable");
486 pl->sublist(
"Default Integrator")
487 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
488 integrator = Tempus::createIntegratorBasic<double>(pl, model);
494 RCP<Thyra::VectorBase<double> > x0 =
495 model->getNominalValues().get_x()->clone_v();
496 integrator->initializeSolutionHistory(0.0, x0);
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);
511 RCP<const SolutionHistory<double> > solutionHistory =
512 integrator->getSolutionHistory();
513 writeSolution(
"Tempus_"+RKMethod+
"_SinCos.dat", solutionHistory);
516 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
517 double time_i = (*solutionHistory)[i]->getTime();
520 model->getExactSolution(time_i).get_x()),
522 model->getExactSolution(time_i).get_x_dot()));
523 state->setTime((*solutionHistory)[i]->getTime());
524 solnHistExact->addState(state);
526 writeSolution(
"Tempus_"+RKMethod+
"_SinCos-Ref.dat", solnHistExact);
530 StepSize.push_back(dt);
531 auto solution = Thyra::createMember(model->get_x_space());
532 Thyra::copy(*(integrator->getX()),solution.ptr());
533 solutions.push_back(solution);
534 auto solutionDot = Thyra::createMember(model->get_x_space());
535 Thyra::copy(*(integrator->getXDot()),solutionDot.ptr());
536 solutionsDot.push_back(solutionDot);
537 if (n == nTimeStepSizes-1) {
538 StepSize.push_back(0.0);
539 auto solutionExact = Thyra::createMember(model->get_x_space());
540 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
541 solutions.push_back(solutionExact);
542 auto solutionDotExact = Thyra::createMember(model->get_x_space());
543 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
544 solutionDotExact.ptr());
545 solutionsDot.push_back(solutionDotExact);
551 double xDotSlope = 0.0;
552 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
553 double order = stepper->getOrder();
556 solutions, xErrorNorm, xSlope,
557 solutionsDot, xDotErrorNorm, xDotSlope);
559 TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
560 TEST_FLOATING_EQUALITY( xErrorNorm[0], RKMethodErrors[m], 5.0e-4 );
573 std::vector<std::string> RKMethods;
574 RKMethods.push_back(
"SDIRK 2 Stage 2nd order");
576 std::string RKMethod = RKMethods[0];
577 std::replace(RKMethod.begin(), RKMethod.end(),
' ',
'_');
578 std::replace(RKMethod.begin(), RKMethod.end(),
'/',
'.');
580 RCP<Tempus::IntegratorBasic<double> > integrator;
581 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
582 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
583 std::vector<double> StepSize;
584 std::vector<double> xErrorNorm;
585 std::vector<double> xDotErrorNorm;
587 const int nTimeStepSizes = 3;
590 for (
int n=0; n<nTimeStepSizes; n++) {
593 RCP<ParameterList> pList =
594 getParametersFromXmlFile(
"Tempus_DIRK_VanDerPol.xml");
597 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
601 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
602 pl->sublist(
"Default Stepper").set(
"Stepper Type", RKMethods[0]);
603 pl->sublist(
"Default Stepper").set(
"gamma", 0.2928932188134524);
607 if (n == nTimeStepSizes-1) dt /= 10.0;
610 pl->sublist(
"Default Integrator")
611 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
612 integrator = Tempus::createIntegratorBasic<double>(pl, model);
615 bool integratorStatus = integrator->advanceTime();
616 TEST_ASSERT(integratorStatus)
619 time = integrator->getTime();
620 double timeFinal =pl->sublist(
"Default Integrator")
621 .sublist(
"Time Step Control").get<
double>(
"Final Time");
622 double tol = 100.0 * std::numeric_limits<double>::epsilon();
623 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
626 StepSize.push_back(dt);
627 auto solution = Thyra::createMember(model->get_x_space());
628 Thyra::copy(*(integrator->getX()),solution.ptr());
629 solutions.push_back(solution);
630 auto solutionDot = Thyra::createMember(model->get_x_space());
631 Thyra::copy(*(integrator->getXDot()),solutionDot.ptr());
632 solutionsDot.push_back(solutionDot);
636 if ((n == 0) || (n == nTimeStepSizes-1)) {
637 std::string fname =
"Tempus_"+RKMethod+
"_VanDerPol-Ref.dat";
638 if (n == 0) fname =
"Tempus_"+RKMethod+
"_VanDerPol.dat";
639 RCP<const SolutionHistory<double> > solutionHistory =
640 integrator->getSolutionHistory();
647 double xDotSlope = 0.0;
648 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
649 double order = stepper->getOrder();
652 solutions, xErrorNorm, xSlope,
653 solutionsDot, xDotErrorNorm, xDotSlope);
655 TEST_FLOATING_EQUALITY( xSlope, order, 0.06 );
656 TEST_FLOATING_EQUALITY( xErrorNorm[0], 1.07525e-05, 1.0e-4 );
661 Teuchos::TimeMonitor::summarize();
670 std::vector<std::string> IntegratorList;
671 IntegratorList.push_back(
"Embedded_Integrator_PID");
672 IntegratorList.push_back(
"Embedded_Integrator");
675 const int refIstep = 217;
677 for(
auto integratorChoice : IntegratorList){
679 std::cout <<
"Using Integrator: " << integratorChoice <<
" !!!" << std::endl;
682 RCP<ParameterList> pList =
683 getParametersFromXmlFile(
"Tempus_DIRK_VanDerPol.xml");
687 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
691 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
692 pl->set(
"Integrator Name", integratorChoice);
695 RCP<Tempus::IntegratorBasic<double> > integrator =
696 Tempus::createIntegratorBasic<double>(pl, model);
698 const std::string RKMethod_ =
699 pl->sublist(integratorChoice).get<std::string>(
"Stepper Name");
702 bool integratorStatus = integrator->advanceTime();
703 TEST_ASSERT(integratorStatus);
706 double time = integrator->getTime();
707 double timeFinal = pl->sublist(integratorChoice)
708 .sublist(
"Time Step Control").get<
double>(
"Final Time");
709 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
713 RCP<Thyra::VectorBase<double> > x = integrator->getX();
714 RCP<Thyra::VectorBase<double> > xref = x->clone_v();
715 Thyra::set_ele(0, -1.5484458614405929, xref.ptr());
716 Thyra::set_ele(1, 1.0181127316101317, xref.ptr());
719 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
720 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *xref, -1.0, *(x));
721 const double L2norm = Thyra::norm_2(*xdiff);
724 if (integratorChoice ==
"Embedded_Integrator_PID"){
725 const double absTol = pl->sublist(integratorChoice).
726 sublist(
"Time Step Control").get<
double>(
"Maximum Absolute Error");
727 const double relTol = pl->sublist(integratorChoice).
728 sublist(
"Time Step Control").get<
double>(
"Maximum Relative Error");
734 const int iStep = integrator->getSolutionHistory()->
735 getCurrentState()->getIndex();
740 TEST_FLOATING_EQUALITY(std::log10(L2norm),std::log10(absTol), 0.3 );
741 TEST_FLOATING_EQUALITY(std::log10(L2norm),std::log10(relTol), 0.3 );
743 TEST_COMPARE(iStep, <=, refIstep);
747 std::ofstream ftmp(
"Tempus_"+integratorChoice+RKMethod_+
"_VDP_Example.dat");
748 RCP<const SolutionHistory<double> > solutionHistory =
749 integrator->getSolutionHistory();
750 int nStates = solutionHistory->getNumStates();
752 for (
int i=0; i<nStates; i++) {
753 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
754 double time_i = solutionState->getTime();
755 RCP<const Thyra::VectorBase<double> > x_plot = solutionState->getX();
757 ftmp << time_i <<
" " 758 << Thyra::get_ele(*(x_plot), 0) <<
" " 759 << Thyra::get_ele(*(x_plot), 1) <<
" " << std::endl;
764 Teuchos::TimeMonitor::summarize();
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.
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation with a...
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)
TimeStepControl manages the time step size. There several mechanisms that effect the time step size a...
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Keep a fix number of states.
van der Pol model problem for nonlinear electrical circuit.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...