Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_DIRKTest.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
12 
13 #include "Thyra_VectorStdOps.hpp"
14 
15 #include "Tempus_IntegratorBasic.hpp"
16 #include "Tempus_StepperDIRK.hpp"
17 
18 #include "../TestModels/SinCosModel.hpp"
19 #include "../TestModels/VanDerPolModel.hpp"
20 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
21 
22 #include <fstream>
23 #include <vector>
24 
25 namespace Tempus_Test {
26 
27 using Teuchos::RCP;
28 using Teuchos::rcp;
29 using Teuchos::rcp_const_cast;
30 using Teuchos::ParameterList;
31 using Teuchos::parameterList;
32 using Teuchos::sublist;
33 using Teuchos::getParametersFromXmlFile;
34 
38 
39 // Comment out any of the following tests to exclude from build/run.
40 #define TEST_PARAMETERLIST
41 #define TEST_CONSTRUCTING_FROM_DEFAULTS
42 #define TEST_SINCOS
43 #define TEST_VANDERPOL
44 #define TEST_EMBEDDED_DIRK
45 
46 
47 #ifdef TEST_PARAMETERLIST
48 // ************************************************************
49 // ************************************************************
50 TEUCHOS_UNIT_TEST(DIRK, ParameterList)
51 {
52  std::vector<std::string> RKMethods;
53  RKMethods.push_back("RK Backward Euler");
54  RKMethods.push_back("IRK 1 Stage Theta Method");
55  RKMethods.push_back("Implicit Midpoint");
56  RKMethods.push_back("SDIRK 1 Stage 1st order");
57  RKMethods.push_back("SDIRK 2 Stage 2nd order");
58  RKMethods.push_back("SDIRK 2 Stage 3rd order");
59  RKMethods.push_back("EDIRK 2 Stage 3rd order");
60  RKMethods.push_back("EDIRK 2 Stage Theta Method");
61  RKMethods.push_back("SDIRK 3 Stage 4th order");
62  RKMethods.push_back("SDIRK 5 Stage 4th order");
63  RKMethods.push_back("SDIRK 5 Stage 5th order");
64  RKMethods.push_back("SDIRK 2(1) Pair");
65 
66  for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
67 
68  std::string RKMethod = RKMethods[m];
69 
70  // Read params from .xml file
71  RCP<ParameterList> pList =
72  getParametersFromXmlFile("Tempus_DIRK_SinCos.xml");
73 
74  // Setup the SinCosModel
75  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
76  auto model = rcp(new SinCosModel<double>(scm_pl));
77 
78  RCP<ParameterList> tempusPL = sublist(pList, "Tempus", true);
79  tempusPL->sublist("Default Stepper").set("Stepper Type", RKMethods[m]);
80 
81  if (RKMethods[m] == "IRK 1 Stage Theta Method" ||
82  RKMethods[m] == "Implicit Midpoint" ||
83  RKMethods[m] == "EDIRK 2 Stage Theta Method") {
84  // Construct in the same order as default.
85  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
86  RCP<ParameterList> solverPL = parameterList();
87  *solverPL = *(sublist(stepperPL, "Default Solver", true));
88  tempusPL->sublist("Default Stepper").remove("Default Solver");
89  tempusPL->sublist("Default Stepper")
90  .set<double>("theta", 0.5);
91  tempusPL->sublist("Default Stepper").remove("Zero Initial Guess");
92  tempusPL->sublist("Default Stepper").set<bool>("Zero Initial Guess", 0);
93  tempusPL->sublist("Default Stepper").set("Default Solver", *solverPL);
94  } else if (RKMethods[m] == "SDIRK 2 Stage 2nd order") {
95  // Construct in the same order as default.
96  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
97  RCP<ParameterList> solverPL = parameterList();
98  *solverPL = *(sublist(stepperPL, "Default Solver", true));
99  tempusPL->sublist("Default Stepper").remove("Default Solver");
100  tempusPL->sublist("Default Stepper")
101  .set<double>("gamma", 0.2928932188134524);
102  tempusPL->sublist("Default Stepper").remove("Zero Initial Guess");
103  tempusPL->sublist("Default Stepper").set<bool>("Zero Initial Guess", 0);
104  tempusPL->sublist("Default Stepper").set("Default Solver", *solverPL);
105  } else if (RKMethods[m] == "SDIRK 2 Stage 3rd order") {
106  // Construct in the same order as default.
107  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
108  RCP<ParameterList> solverPL = parameterList();
109  *solverPL = *(sublist(stepperPL, "Default Solver", true));
110  tempusPL->sublist("Default Stepper").remove("Default Solver");
111  tempusPL->sublist("Default Stepper").set("3rd Order A-stable", true);
112  tempusPL->sublist("Default Stepper").set("2nd Order L-stable", false);
113  tempusPL->sublist("Default Stepper")
114  .set<double>("gamma", 0.7886751345948128);
115  tempusPL->sublist("Default Stepper").remove("Zero Initial Guess");
116  tempusPL->sublist("Default Stepper").set<bool>("Zero Initial Guess", 0);
117  tempusPL->sublist("Default Stepper").set("Default Solver", *solverPL);
118  }
119 
120  // Test constructor IntegratorBasic(tempusPL, model)
121  {
122  RCP<Tempus::IntegratorBasic<double> > integrator =
123  Tempus::integratorBasic<double>(tempusPL, model);
124 
125  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
126  RCP<ParameterList> defaultPL =
127  integrator->getStepper()->getDefaultParameters();
128  defaultPL->remove("Description");
129 
130  // Adjust default parameters for pseudonyms.
131  if ( RKMethods[m] == "Implicit Midpoint" ) {
132  defaultPL->set("Stepper Type", "Implicit Midpoint");
133  }
134 
135  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
136  if (!pass) {
137  std::cout << std::endl;
138  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
139  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
140  }
141  TEST_ASSERT(pass)
142  }
143 
144  // Test constructor IntegratorBasic(model, stepperType)
145  {
146  RCP<Tempus::IntegratorBasic<double> > integrator =
147  Tempus::integratorBasic<double>(model, RKMethods[m]);
148 
149  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
150  RCP<ParameterList> defaultPL =
151  integrator->getStepper()->getDefaultParameters();
152  defaultPL->remove("Description");
153 
154  // Adjust default parameters for pseudonyms.
155  if ( RKMethods[m] == "Implicit Midpoint" ) {
156  defaultPL->set("Stepper Type", "Implicit Midpoint");
157  }
158 
159  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
160  if (!pass) {
161  std::cout << std::endl;
162  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
163  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
164  }
165  TEST_ASSERT(pass)
166  }
167  }
168 }
169 #endif // TEST_PARAMETERLIST
170 
171 
172 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
173 // ************************************************************
174 // ************************************************************
175 TEUCHOS_UNIT_TEST(DIRK, ConstructingFromDefaults)
176 {
177  double dt = 0.025;
178 
179  // Read params from .xml file
180  RCP<ParameterList> pList =
181  getParametersFromXmlFile("Tempus_DIRK_SinCos.xml");
182  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
183 
184  // Setup the SinCosModel
185  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
186  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
187  auto model = rcp(new SinCosModel<double>(scm_pl));
188 
189  // Setup Stepper for field solve ----------------------------
190  auto stepper = rcp(new Tempus::StepperDIRK<double>());
191  stepper->setModel(model);
192  stepper->initialize();
193 
194  // Setup TimeStepControl ------------------------------------
195  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
196  ParameterList tscPL = pl->sublist("Default Integrator")
197  .sublist("Time Step Control");
198  timeStepControl->setStepType (tscPL.get<std::string>("Integrator Step Type"));
199  timeStepControl->setInitIndex(tscPL.get<int> ("Initial Time Index"));
200  timeStepControl->setInitTime (tscPL.get<double>("Initial Time"));
201  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
202  timeStepControl->setInitTimeStep(dt);
203  timeStepControl->initialize();
204 
205  // Setup initial condition SolutionState --------------------
206  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
207  stepper->getModel()->getNominalValues();
208  auto icSolution = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
209  auto icState = rcp(new Tempus::SolutionState<double>(icSolution));
210  icState->setTime (timeStepControl->getInitTime());
211  icState->setIndex (timeStepControl->getInitIndex());
212  icState->setTimeStep(0.0);
213  icState->setOrder (stepper->getOrder());
214  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
215 
216  // Setup SolutionHistory ------------------------------------
218  solutionHistory->setName("Forward States");
220  solutionHistory->setStorageLimit(2);
221  solutionHistory->addState(icState);
222 
223  // Setup Integrator -----------------------------------------
224  RCP<Tempus::IntegratorBasic<double> > integrator =
225  Tempus::integratorBasic<double>();
226  integrator->setStepperWStepper(stepper);
227  integrator->setTimeStepControl(timeStepControl);
228  integrator->setSolutionHistory(solutionHistory);
229  //integrator->setObserver(...);
230  integrator->initialize();
231 
232 
233  // Integrate to timeMax
234  bool integratorStatus = integrator->advanceTime();
235  TEST_ASSERT(integratorStatus)
236 
237 
238  // Test if at 'Final Time'
239  double time = integrator->getTime();
240  double timeFinal =pl->sublist("Default Integrator")
241  .sublist("Time Step Control").get<double>("Final Time");
242  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
243 
244  // Time-integrated solution and the exact solution
245  RCP<Thyra::VectorBase<double> > x = integrator->getX();
246  RCP<const Thyra::VectorBase<double> > x_exact =
247  model->getExactSolution(time).get_x();
248 
249  // Calculate the error
250  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
251  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
252 
253  // Check the order and intercept
254  std::cout << " Stepper = SDIRK 2 Stage 2nd order" << std::endl;
255  std::cout << " =========================" << std::endl;
256  std::cout << " Exact solution : " << get_ele(*(x_exact), 0) << " "
257  << get_ele(*(x_exact), 1) << std::endl;
258  std::cout << " Computed solution: " << get_ele(*(x ), 0) << " "
259  << get_ele(*(x ), 1) << std::endl;
260  std::cout << " Difference : " << get_ele(*(xdiff ), 0) << " "
261  << get_ele(*(xdiff ), 1) << std::endl;
262  std::cout << " =========================" << std::endl;
263  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841470, 1.0e-4 );
264  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.540304, 1.0e-4 );
265 }
266 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
267 
268 
269 #ifdef TEST_SINCOS
270 // ************************************************************
271 // ************************************************************
272 TEUCHOS_UNIT_TEST(DIRK, SinCos)
273 {
274  std::vector<std::string> RKMethods;
275  RKMethods.push_back("RK Backward Euler");
276  RKMethods.push_back("IRK 1 Stage Theta Method");
277  RKMethods.push_back("Implicit Midpoint");
278  RKMethods.push_back("SDIRK 1 Stage 1st order");
279  RKMethods.push_back("SDIRK 2 Stage 2nd order");
280  RKMethods.push_back("SDIRK 2 Stage 3rd order");
281  RKMethods.push_back("EDIRK 2 Stage 3rd order");
282  RKMethods.push_back("EDIRK 2 Stage Theta Method");
283  RKMethods.push_back("SDIRK 3 Stage 4th order");
284  RKMethods.push_back("SDIRK 5 Stage 4th order");
285  RKMethods.push_back("SDIRK 5 Stage 5th order");
286  RKMethods.push_back("SDIRK 2(1) Pair");
287 
288  std::vector<double> RKMethodErrors;
289  RKMethodErrors.push_back(0.0124201);
290  RKMethodErrors.push_back(5.20785e-05);
291  RKMethodErrors.push_back(5.20785e-05);
292  RKMethodErrors.push_back(0.0124201);
293  RKMethodErrors.push_back(2.52738e-05);
294  RKMethodErrors.push_back(1.40223e-06);
295  RKMethodErrors.push_back(2.17004e-07);
296  RKMethodErrors.push_back(5.20785e-05);
297  RKMethodErrors.push_back(6.41463e-08);
298  RKMethodErrors.push_back(3.30631e-10);
299  RKMethodErrors.push_back(1.35728e-11);
300  RKMethodErrors.push_back(0.0001041);
301 
302  for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
303 
304  std::string RKMethod = RKMethods[m];
305  std::replace(RKMethod.begin(), RKMethod.end(), ' ', '_');
306  std::replace(RKMethod.begin(), RKMethod.end(), '/', '.');
307 
308  RCP<Tempus::IntegratorBasic<double> > integrator;
309  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
310  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
311  std::vector<double> StepSize;
312  std::vector<double> xErrorNorm;
313  std::vector<double> xDotErrorNorm;
314 
315  const int nTimeStepSizes = 2; // 7 for error plots
316  double dt = 0.05;
317  double time = 0.0;
318  for (int n=0; n<nTimeStepSizes; n++) {
319 
320  // Read params from .xml file
321  RCP<ParameterList> pList =
322  getParametersFromXmlFile("Tempus_DIRK_SinCos.xml");
323 
324  // Setup the SinCosModel
325  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
326  auto model = rcp(new SinCosModel<double>(scm_pl));
327 
328  // Set the Stepper
329  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
330  pl->sublist("Default Stepper").set("Stepper Type", RKMethods[m]);
331  if (RKMethods[m] == "IRK 1 Stage Theta Method" ||
332  RKMethods[m] == "Implicit Midpoint" ||
333  RKMethods[m] == "EDIRK 2 Stage Theta Method") {
334  pl->sublist("Default Stepper").set<double>("theta", 0.5);
335  } else if (RKMethods[m] == "SDIRK 2 Stage 2nd order") {
336  pl->sublist("Default Stepper").set("gamma", 0.2928932188134524);
337  } else if (RKMethods[m] == "SDIRK 2 Stage 3rd order") {
338  pl->sublist("Default Stepper").set("3rd Order A-stable", true);
339  pl->sublist("Default Stepper").set("2nd Order L-stable", false);
340  pl->sublist("Default Stepper").set("gamma", 0.7886751345948128);
341  }
342 
343  dt /= 2;
344 
345  // Setup the Integrator and reset initial time step
346  pl->sublist("Default Integrator")
347  .sublist("Time Step Control").set("Initial Time Step", dt);
348  integrator = Tempus::integratorBasic<double>(pl, model);
349 
350  // Initial Conditions
351  // During the Integrator construction, the initial SolutionState
352  // is set by default to model->getNominalVales().get_x(). However,
353  // the application can set it also by integrator->initializeSolutionHistory.
354  RCP<Thyra::VectorBase<double> > x0 =
355  model->getNominalValues().get_x()->clone_v();
356  integrator->initializeSolutionHistory(0.0, x0);
357 
358  // Integrate to timeMax
359  bool integratorStatus = integrator->advanceTime();
360  TEST_ASSERT(integratorStatus)
361 
362  // Test if at 'Final Time'
363  time = integrator->getTime();
364  double timeFinal = pl->sublist("Default Integrator")
365  .sublist("Time Step Control").get<double>("Final Time");
366  double tol = 100.0 * std::numeric_limits<double>::epsilon();
367  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
368 
369  // Plot sample solution and exact solution
370  if (n == 0) {
371  RCP<const SolutionHistory<double> > solutionHistory =
372  integrator->getSolutionHistory();
373  writeSolution("Tempus_"+RKMethod+"_SinCos.dat", solutionHistory);
374 
375  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
376  for (int i=0; i<solutionHistory->getNumStates(); i++) {
377  double time_i = (*solutionHistory)[i]->getTime();
378  auto state = rcp(new Tempus::SolutionState<double>(
379  model->getExactSolution(time_i).get_x(),
380  model->getExactSolution(time_i).get_x_dot()));
381  state->setTime((*solutionHistory)[i]->getTime());
382  solnHistExact->addState(state);
383  }
384  writeSolution("Tempus_"+RKMethod+"_SinCos-Ref.dat", solnHistExact);
385  }
386 
387  // Store off the final solution and step size
388  StepSize.push_back(dt);
389  auto solution = Thyra::createMember(model->get_x_space());
390  Thyra::copy(*(integrator->getX()),solution.ptr());
391  solutions.push_back(solution);
392  auto solutionDot = Thyra::createMember(model->get_x_space());
393  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
394  solutionsDot.push_back(solutionDot);
395  if (n == nTimeStepSizes-1) { // Add exact solution last in vector.
396  StepSize.push_back(0.0);
397  auto solutionExact = Thyra::createMember(model->get_x_space());
398  Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
399  solutions.push_back(solutionExact);
400  auto solutionDotExact = Thyra::createMember(model->get_x_space());
401  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
402  solutionDotExact.ptr());
403  solutionsDot.push_back(solutionDotExact);
404  }
405  }
406 
407  // Check the order and intercept
408  double xSlope = 0.0;
409  double xDotSlope = 0.0;
410  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
411  double order = stepper->getOrder();
412  writeOrderError("Tempus_"+RKMethod+"_SinCos-Error.dat",
413  stepper, StepSize,
414  solutions, xErrorNorm, xSlope,
415  solutionsDot, xDotErrorNorm, xDotSlope);
416 
417  TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
418  TEST_FLOATING_EQUALITY( xErrorNorm[0], RKMethodErrors[m], 5.0e-4 );
419  // xDot not yet available for DIRK methods.
420  //TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
421  //TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.0486418, 1.0e-4 );
422 
423  }
424 }
425 #endif // TEST_SINCOS
426 
427 
428 #ifdef TEST_VANDERPOL
429 // ************************************************************
430 // ************************************************************
431 TEUCHOS_UNIT_TEST(DIRK, VanDerPol)
432 {
433  std::vector<std::string> RKMethods;
434  RKMethods.push_back("SDIRK 2 Stage 2nd order");
435 
436  std::string RKMethod = RKMethods[0];
437  std::replace(RKMethod.begin(), RKMethod.end(), ' ', '_');
438  std::replace(RKMethod.begin(), RKMethod.end(), '/', '.');
439 
440  RCP<Tempus::IntegratorBasic<double> > integrator;
441  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
442  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
443  std::vector<double> StepSize;
444  std::vector<double> xErrorNorm;
445  std::vector<double> xDotErrorNorm;
446 
447  const int nTimeStepSizes = 3; // 8 for error plot
448  double dt = 0.05;
449  double time = 0.0;
450  for (int n=0; n<nTimeStepSizes; n++) {
451 
452  // Read params from .xml file
453  RCP<ParameterList> pList =
454  getParametersFromXmlFile("Tempus_DIRK_VanDerPol.xml");
455 
456  // Setup the VanDerPolModel
457  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
458  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
459 
460  // Set the Stepper
461  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
462  pl->sublist("Default Stepper").set("Stepper Type", RKMethods[0]);
463  pl->sublist("Default Stepper").set("gamma", 0.2928932188134524);
464 
465  // Set the step size
466  dt /= 2;
467  if (n == nTimeStepSizes-1) dt /= 10.0;
468 
469  // Setup the Integrator and reset initial time step
470  pl->sublist("Default Integrator")
471  .sublist("Time Step Control").set("Initial Time Step", dt);
472  integrator = Tempus::integratorBasic<double>(pl, model);
473 
474  // Integrate to timeMax
475  bool integratorStatus = integrator->advanceTime();
476  TEST_ASSERT(integratorStatus)
477 
478  // Test if at 'Final Time'
479  time = integrator->getTime();
480  double timeFinal =pl->sublist("Default Integrator")
481  .sublist("Time Step Control").get<double>("Final Time");
482  double tol = 100.0 * std::numeric_limits<double>::epsilon();
483  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
484 
485  // Store off the final solution and step size
486  StepSize.push_back(dt);
487  auto solution = Thyra::createMember(model->get_x_space());
488  Thyra::copy(*(integrator->getX()),solution.ptr());
489  solutions.push_back(solution);
490  auto solutionDot = Thyra::createMember(model->get_x_space());
491  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
492  solutionsDot.push_back(solutionDot);
493 
494  // Output finest temporal solution for plotting
495  // This only works for ONE MPI process
496  if ((n == 0) or (n == nTimeStepSizes-1)) {
497  std::string fname = "Tempus_"+RKMethod+"_VanDerPol-Ref.dat";
498  if (n == 0) fname = "Tempus_"+RKMethod+"_VanDerPol.dat";
499  RCP<const SolutionHistory<double> > solutionHistory =
500  integrator->getSolutionHistory();
501  writeSolution(fname, solutionHistory);
502  }
503  }
504 
505  // Check the order and intercept
506  double xSlope = 0.0;
507  double xDotSlope = 0.0;
508  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
509  double order = stepper->getOrder();
510  writeOrderError("Tempus_"+RKMethod+"_VanDerPol-Error.dat",
511  stepper, StepSize,
512  solutions, xErrorNorm, xSlope,
513  solutionsDot, xDotErrorNorm, xDotSlope);
514 
515  TEST_FLOATING_EQUALITY( xSlope, order, 0.06 );
516  TEST_FLOATING_EQUALITY( xErrorNorm[0], 1.07525e-05, 1.0e-4 );
517  // xDot not yet available for DIRK methods.
518  //TEST_FLOATING_EQUALITY( xDotSlope, 1.74898, 0.10 );
519  //TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 1.0038, 1.0e-4 );
520 
521  Teuchos::TimeMonitor::summarize();
522 }
523 #endif // TEST_VANDERPOL
524 
525 
526 #ifdef TEST_EMBEDDED_DIRK
527 // ************************************************************
528 // ************************************************************
529 TEUCHOS_UNIT_TEST(DIRK, EmbeddedVanDerPol)
530 {
531 
532  std::vector<std::string> IntegratorList;
533  IntegratorList.push_back("Embedded_Integrator_PID");
534  IntegratorList.push_back("Embedded_Integrator");
535 
536  // the embedded solution will test the following:
537  const int refIstep = 217;
538 
539  for(auto integratorChoice : IntegratorList){
540 
541  std::cout << "Using Integrator: " << integratorChoice << " !!!" << std::endl;
542 
543  // Read params from .xml file
544  RCP<ParameterList> pList =
545  getParametersFromXmlFile("Tempus_DIRK_VanDerPol.xml");
546 
547 
548  // Setup the VanDerPolModel
549  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
550  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
551 
552  // Set the Integrator and Stepper
553  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
554  pl->set("Integrator Name", integratorChoice);
555 
556  // Setup the Integrator
557  RCP<Tempus::IntegratorBasic<double> > integrator =
558  Tempus::integratorBasic<double>(pl, model);
559 
560  const std::string RKMethod_ =
561  pl->sublist(integratorChoice).get<std::string>("Stepper Name");
562 
563  // Integrate to timeMax
564  bool integratorStatus = integrator->advanceTime();
565  TEST_ASSERT(integratorStatus);
566 
567  // Test if at 'Final Time'
568  double time = integrator->getTime();
569  double timeFinal = pl->sublist(integratorChoice)
570  .sublist("Time Step Control").get<double>("Final Time");
571  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
572 
573 
574  // Numerical reference solution at timeFinal (for \epsilon = 0.1)
575  RCP<Thyra::VectorBase<double> > x = integrator->getX();
576  RCP<Thyra::VectorBase<double> > xref = x->clone_v();
577  Thyra::set_ele(0, -1.5484458614405929, xref.ptr());
578  Thyra::set_ele(1, 1.0181127316101317, xref.ptr());
579 
580  // Calculate the error
581  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
582  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *xref, -1.0, *(x));
583  const double L2norm = Thyra::norm_2(*xdiff);
584 
585  // Test number of steps, failures, and accuracy
586  if (integratorChoice == "Embedded_Integrator_PID"){
587  const double absTol = pl->sublist(integratorChoice).
588  sublist("Time Step Control").get<double>("Maximum Absolute Error");
589  const double relTol = pl->sublist(integratorChoice).
590  sublist("Time Step Control").get<double>("Maximum Relative Error");
591 
592 
593  // get the number of time steps and number of step failure
594  //const int nFailure_c = integrator->getSolutionHistory()->
595  //getCurrentState()->getMetaData()->getNFailures();
596  const int iStep = integrator->getSolutionHistory()->
597  getCurrentState()->getIndex();
598  //const int nFail = integrator->getSolutionHistory()->
599  // getCurrentState()->getMetaData()->getNRunningFailures();
600 
601  // Should be close to the prescribed tolerance
602  TEST_FLOATING_EQUALITY(std::log10(L2norm),std::log10(absTol), 0.3 );
603  TEST_FLOATING_EQUALITY(std::log10(L2norm),std::log10(relTol), 0.3 );
604  // test for number of steps
605  TEST_COMPARE(iStep, <=, refIstep);
606  }
607 
608  // Plot sample solution and exact solution
609  std::ofstream ftmp("Tempus_"+integratorChoice+RKMethod_+"_VDP_Example.dat");
610  RCP<const SolutionHistory<double> > solutionHistory =
611  integrator->getSolutionHistory();
612  int nStates = solutionHistory->getNumStates();
613  //RCP<const Thyra::VectorBase<double> > x_exact_plot;
614  for (int i=0; i<nStates; i++) {
615  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
616  double time_i = solutionState->getTime();
617  RCP<const Thyra::VectorBase<double> > x_plot = solutionState->getX();
618  //x_exact_plot = model->getExactSolution(time_i).get_x();
619  ftmp << time_i << " "
620  << Thyra::get_ele(*(x_plot), 0) << " "
621  << Thyra::get_ele(*(x_plot), 1) << " " << std::endl;
622  }
623  ftmp.close();
624  }
625 
626  Teuchos::TimeMonitor::summarize();
627 }
628 #endif
629 
630 
631 } // namespace Tempus_Test
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)
Diagonally Implicit Runge-Kutta (DIRK) time stepper.
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...