9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
12 #include "Teuchos_DefaultComm.hpp"
14 #include "Thyra_VectorStdOps.hpp"
16 #include "Tempus_IntegratorBasic.hpp"
19 #include "../TestModels/SinCosModel.hpp"
20 #include "../TestModels/VanDerPolModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
26 namespace Tempus_Test {
30 using Teuchos::rcp_const_cast;
31 using Teuchos::ParameterList;
32 using Teuchos::sublist;
33 using Teuchos::getParametersFromXmlFile;
40 #define TEST_PARAMETERLIST
41 #define TEST_CONSTRUCTING_FROM_DEFAULTS
43 #define TEST_EMBEDDED_VANDERPOL
46 #ifdef TEST_PARAMETERLIST
51 std::vector<std::string> RKMethods;
52 RKMethods.push_back(
"General ERK");
53 RKMethods.push_back(
"RK Forward Euler");
54 RKMethods.push_back(
"RK Explicit 4 Stage");
55 RKMethods.push_back(
"RK Explicit 3/8 Rule");
56 RKMethods.push_back(
"RK Explicit 4 Stage 3rd order by Runge");
57 RKMethods.push_back(
"RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
58 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order");
59 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order TVD");
60 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order by Heun");
61 RKMethods.push_back(
"RK Explicit Midpoint");
62 RKMethods.push_back(
"RK Explicit Trapezoidal");
63 RKMethods.push_back(
"Heuns Method");
64 RKMethods.push_back(
"Bogacki-Shampine 3(2) Pair");
65 RKMethods.push_back(
"Merson 4(5) Pair");
67 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
69 std::string RKMethod = RKMethods[m];
70 std::replace(RKMethod.begin(), RKMethod.end(),
' ',
'_');
71 std::replace(RKMethod.begin(), RKMethod.end(),
'/',
'.');
74 RCP<ParameterList> pList =
75 getParametersFromXmlFile(
"Tempus_ExplicitRK_SinCos.xml");
78 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
82 RCP<ParameterList> tempusPL = sublist(pList,
"Tempus",
true);
83 if (RKMethods[m] ==
"General ERK") {
84 tempusPL->sublist(
"Demo Integrator").set(
"Stepper Name",
"Demo Stepper 2");
86 tempusPL->sublist(
"Demo Stepper").set(
"Stepper Type", RKMethods[m]);
90 tempusPL->sublist(
"Demo Stepper")
91 .set(
"Initial Condition Consistency",
"None");
95 RCP<Tempus::IntegratorBasic<double> > integrator =
96 Tempus::integratorBasic<double>(tempusPL, model);
98 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Demo Stepper",
true);
99 if (RKMethods[m] ==
"General ERK")
100 stepperPL = sublist(tempusPL,
"Demo Stepper 2",
true);
101 RCP<ParameterList> defaultPL =
102 Teuchos::rcp_const_cast<Teuchos::ParameterList>(
103 integrator->getStepper()->getValidParameters());
104 defaultPL->remove(
"Description");
106 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
108 std::cout << std::endl;
109 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
110 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
117 RCP<Tempus::IntegratorBasic<double> > integrator =
118 Tempus::integratorBasic<double>(model, RKMethods[m]);
120 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Demo Stepper",
true);
121 if (RKMethods[m] ==
"General ERK")
122 stepperPL = sublist(tempusPL,
"Demo Stepper 2",
true);
123 RCP<ParameterList> defaultPL =
124 Teuchos::rcp_const_cast<Teuchos::ParameterList>(
125 integrator->getStepper()->getValidParameters());
126 defaultPL->remove(
"Description");
128 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
130 std::cout << std::endl;
131 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
132 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
138 #endif // TEST_PARAMETERLIST
141 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
149 RCP<ParameterList> pList =
150 getParametersFromXmlFile(
"Tempus_ExplicitRK_SinCos.xml");
151 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
154 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
159 RCP<Tempus::StepperFactory<double> > sf =
161 RCP<Tempus::Stepper<double> > stepper =
162 sf->createStepper(
"RK Explicit 4 Stage");
163 stepper->setModel(model);
164 stepper->initialize();
168 ParameterList tscPL = pl->sublist(
"Demo Integrator")
169 .sublist(
"Time Step Control");
170 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
171 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
172 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
173 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
174 timeStepControl->setInitTimeStep(dt);
175 timeStepControl->initialize();
178 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
179 stepper->getModel()->getNominalValues();
180 auto icSolution = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
182 icState->setTime (timeStepControl->getInitTime());
183 icState->setIndex (timeStepControl->getInitIndex());
184 icState->setTimeStep(0.0);
185 icState->setOrder (stepper->getOrder());
196 RCP<Tempus::IntegratorBasic<double> > integrator =
197 Tempus::integratorBasic<double>();
198 integrator->setStepperWStepper(stepper);
199 integrator->setTimeStepControl(timeStepControl);
202 integrator->initialize();
206 bool integratorStatus = integrator->advanceTime();
207 TEST_ASSERT(integratorStatus)
211 double time = integrator->getTime();
212 double timeFinal =pl->sublist(
"Demo Integrator")
213 .sublist(
"Time Step Control").get<
double>(
"Final Time");
214 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
217 RCP<Thyra::VectorBase<double> > x = integrator->getX();
218 RCP<const Thyra::VectorBase<double> > x_exact =
219 model->getExactSolution(time).get_x();
222 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
223 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
226 std::cout <<
" Stepper = RK Explicit 4 Stage" << std::endl;
227 std::cout <<
" =========================" << std::endl;
228 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
229 << get_ele(*(x_exact), 1) << std::endl;
230 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
231 << get_ele(*(x ), 1) << std::endl;
232 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" "
233 << get_ele(*(xdiff ), 1) << std::endl;
234 std::cout <<
" =========================" << std::endl;
235 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841470, 1.0e-4 );
236 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.540303, 1.0e-4 );
238 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
246 std::vector<std::string> RKMethods;
247 RKMethods.push_back(
"General ERK");
248 RKMethods.push_back(
"RK Forward Euler");
249 RKMethods.push_back(
"RK Explicit 4 Stage");
250 RKMethods.push_back(
"RK Explicit 3/8 Rule");
251 RKMethods.push_back(
"RK Explicit 4 Stage 3rd order by Runge");
252 RKMethods.push_back(
"RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
253 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order");
254 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order TVD");
255 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order by Heun");
256 RKMethods.push_back(
"RK Explicit Midpoint");
257 RKMethods.push_back(
"RK Explicit Trapezoidal");
258 RKMethods.push_back(
"Heuns Method");
259 RKMethods.push_back(
"Bogacki-Shampine 3(2) Pair");
260 RKMethods.push_back(
"Merson 4(5) Pair");
261 RKMethods.push_back(
"General ERK Embedded");
263 std::vector<double> RKMethodErrors;
264 RKMethodErrors.push_back(8.33251e-07);
265 RKMethodErrors.push_back(0.051123);
266 RKMethodErrors.push_back(8.33251e-07);
267 RKMethodErrors.push_back(8.33251e-07);
268 RKMethodErrors.push_back(4.16897e-05);
269 RKMethodErrors.push_back(8.32108e-06);
270 RKMethodErrors.push_back(4.16603e-05);
271 RKMethodErrors.push_back(4.16603e-05);
272 RKMethodErrors.push_back(4.16603e-05);
273 RKMethodErrors.push_back(0.00166645);
274 RKMethodErrors.push_back(0.00166645);
275 RKMethodErrors.push_back(0.00166645);
276 RKMethodErrors.push_back(4.16603e-05);
277 RKMethodErrors.push_back(1.39383e-07);
278 RKMethodErrors.push_back(4.16603e-05);
280 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
282 std::string RKMethod = RKMethods[m];
283 std::replace(RKMethod.begin(), RKMethod.end(),
' ',
'_');
284 std::replace(RKMethod.begin(), RKMethod.end(),
'/',
'.');
286 RCP<Tempus::IntegratorBasic<double> > integrator;
287 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
288 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
289 std::vector<double> StepSize;
290 std::vector<double> xErrorNorm;
291 std::vector<double> xDotErrorNorm;
293 const int nTimeStepSizes = 7;
296 for (
int n=0; n<nTimeStepSizes; n++) {
299 RCP<ParameterList> pList =
300 getParametersFromXmlFile(
"Tempus_ExplicitRK_SinCos.xml");
303 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
308 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
309 if (RKMethods[m] ==
"General ERK") {
310 pl->sublist(
"Demo Integrator").set(
"Stepper Name",
"Demo Stepper 2");
311 }
else if (RKMethods[m] ==
"General ERK Embedded"){
312 pl->sublist(
"Demo Integrator").set(
"Stepper Name",
"General ERK Embedded Stepper");
314 pl->sublist(
"Demo Stepper").set(
"Stepper Type", RKMethods[m]);
321 pl->sublist(
"Demo Integrator")
322 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
323 integrator = Tempus::integratorBasic<double>(pl, model);
329 RCP<Thyra::VectorBase<double> > x0 =
330 model->getNominalValues().get_x()->clone_v();
331 integrator->initializeSolutionHistory(0.0, x0);
334 bool integratorStatus = integrator->advanceTime();
335 TEST_ASSERT(integratorStatus)
338 time = integrator->getTime();
339 double timeFinal = pl->sublist(
"Demo Integrator")
340 .sublist(
"Time Step Control").get<
double>(
"Final Time");
341 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
344 RCP<Thyra::VectorBase<double> > x = integrator->getX();
345 RCP<const Thyra::VectorBase<double> > x_exact =
346 model->getExactSolution(time).get_x();
351 integrator->getSolutionHistory();
352 writeSolution(
"Tempus_"+RKMethod+
"_SinCos.dat", solutionHistory);
355 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
356 double time_i = (*solutionHistory)[i]->getTime();
358 model->getExactSolution(time_i).get_x(),
359 model->getExactSolution(time_i).get_x_dot()));
360 state->setTime((*solutionHistory)[i]->getTime());
361 solnHistExact->addState(state);
363 writeSolution(
"Tempus_"+RKMethod+
"_SinCos-Ref.dat", solnHistExact);
367 StepSize.push_back(dt);
368 auto solution = Thyra::createMember(model->get_x_space());
369 Thyra::copy(*(integrator->getX()),solution.ptr());
370 solutions.push_back(solution);
371 auto solutionDot = Thyra::createMember(model->get_x_space());
372 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
373 solutionsDot.push_back(solutionDot);
374 if (n == nTimeStepSizes-1) {
375 StepSize.push_back(0.0);
376 auto solutionExact = Thyra::createMember(model->get_x_space());
377 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
378 solutions.push_back(solutionExact);
379 auto solutionDotExact = Thyra::createMember(model->get_x_space());
380 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
381 solutionDotExact.ptr());
382 solutionsDot.push_back(solutionDotExact);
388 double xDotSlope = 0.0;
389 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
390 double order = stepper->getOrder();
393 solutions, xErrorNorm, xSlope,
394 solutionsDot, xDotErrorNorm, xDotSlope);
396 double order_tol = 0.01;
397 if (RKMethods[m] ==
"Merson 4(5) Pair") order_tol = 0.04;
398 TEST_FLOATING_EQUALITY( xSlope, order, order_tol );
399 TEST_FLOATING_EQUALITY( xErrorNorm[0], RKMethodErrors[m], 1.0e-4 );
407 #endif // TEST_SINCOS
410 #ifdef TEST_EMBEDDED_VANDERPOL
416 std::vector<std::string> IntegratorList;
417 IntegratorList.push_back(
"Embedded_Integrator_PID");
418 IntegratorList.push_back(
"Demo_Integrator");
419 IntegratorList.push_back(
"Embedded_Integrator");
420 IntegratorList.push_back(
"General_Embedded_Integrator");
421 IntegratorList.push_back(
"Embedded_Integrator_PID_General");
425 const int refIstep = 45;
427 for(
auto integratorChoice : IntegratorList){
429 std::cout <<
"Using Integrator: " << integratorChoice <<
" !!!" << std::endl;
432 RCP<ParameterList> pList =
433 getParametersFromXmlFile(
"Tempus_ExplicitRK_VanDerPol.xml");
437 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
442 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
443 pl->set(
"Integrator Name", integratorChoice);
446 RCP<Tempus::IntegratorBasic<double> > integrator =
447 Tempus::integratorBasic<double>(pl, model);
449 const std::string RKMethod =
450 pl->sublist(integratorChoice).get<std::string>(
"Stepper Name");
453 bool integratorStatus = integrator->advanceTime();
454 TEST_ASSERT(integratorStatus);
457 double time = integrator->getTime();
458 double timeFinal = pl->sublist(integratorChoice)
459 .sublist(
"Time Step Control").get<
double>(
"Final Time");
460 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
464 RCP<Thyra::VectorBase<double> > x = integrator->getX();
465 RCP<Thyra::VectorBase<double> > xref = x->clone_v();
466 Thyra::set_ele(0, -1.5484458614405929, xref.ptr());
467 Thyra::set_ele(1, 1.0181127316101317, xref.ptr());
470 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
471 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *xref, -1.0, *(x));
472 const double L2norm = Thyra::norm_2(*xdiff);
475 if ((integratorChoice ==
"Embedded_Integrator_PID") or
476 (integratorChoice ==
"Embedded_Integrator_PID_General")) {
478 const double absTol = pl->sublist(integratorChoice).
479 sublist(
"Time Step Control").get<
double>(
"Maximum Absolute Error");
480 const double relTol = pl->sublist(integratorChoice).
481 sublist(
"Time Step Control").get<
double>(
"Maximum Relative Error");
484 TEST_COMPARE(std::log10(L2norm), <, std::log10(absTol));
485 TEST_COMPARE(std::log10(L2norm), <, std::log10(relTol));
490 const int iStep = integrator->getSolutionHistory()->
491 getCurrentState()->getIndex();
492 const int nFail = integrator->getSolutionHistory()->
493 getCurrentState()->getMetaData()->getNRunningFailures();
496 TEST_EQUALITY(iStep, refIstep);
497 std::cout <<
"Tolerance = " << absTol
498 <<
" L2norm = " << L2norm
499 <<
" iStep = " << iStep
500 <<
" nFail = " << nFail << std::endl;
504 std::ofstream ftmp(
"Tempus_"+integratorChoice+RKMethod+
"_VDP_Example.dat");
506 integrator->getSolutionHistory();
507 int nStates = solutionHistory->getNumStates();
509 for (
int i=0; i<nStates; i++) {
510 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
511 double time_i = solutionState->getTime();
512 RCP<const Thyra::VectorBase<double> > x_plot = solutionState->getX();
514 ftmp << time_i <<
" "
515 << Thyra::get_ele(*(x_plot), 0) <<
" "
516 << Thyra::get_ele(*(x_plot), 1) <<
" " << std::endl;
521 Teuchos::TimeMonitor::summarize();
523 #endif // TEST_EMBEDDED_VANDERPOL
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 mechanicisms that effect the time step size...
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
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...