Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_ExplicitRKTest.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 #include "Teuchos_DefaultComm.hpp"
13 
14 #include "Thyra_VectorStdOps.hpp"
15 
16 #include "Tempus_IntegratorBasic.hpp"
18 
19 #include "../TestModels/SinCosModel.hpp"
20 #include "../TestModels/VanDerPolModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
22 
23 #include <fstream>
24 #include <vector>
25 
26 namespace Tempus_Test {
27 
28 using Teuchos::RCP;
29 using Teuchos::rcp;
30 using Teuchos::rcp_const_cast;
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_EMBEDDED_VANDERPOL
44 
45 
46 #ifdef TEST_PARAMETERLIST
47 // ************************************************************
48 // ************************************************************
49 TEUCHOS_UNIT_TEST(ExplicitRK, ParameterList)
50 {
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");
66 
67  for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
68 
69  std::string RKMethod = RKMethods[m];
70  std::replace(RKMethod.begin(), RKMethod.end(), ' ', '_');
71  std::replace(RKMethod.begin(), RKMethod.end(), '/', '.');
72 
73  // Read params from .xml file
74  RCP<ParameterList> pList =
75  getParametersFromXmlFile("Tempus_ExplicitRK_SinCos.xml");
76 
77  // Setup the SinCosModel
78  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
79  auto model = rcp(new SinCosModel<double>(scm_pl));
80 
81  // Set the Stepper
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");
85  } else {
86  tempusPL->sublist("Demo Stepper").set("Stepper Type", RKMethods[m]);
87  }
88 
89  // Set IC consistency to default value.
90  tempusPL->sublist("Demo Stepper")
91  .set("Initial Condition Consistency", "None");
92 
93  // Test constructor IntegratorBasic(tempusPL, model)
94  {
95  RCP<Tempus::IntegratorBasic<double> > integrator =
96  Tempus::integratorBasic<double>(tempusPL, model);
97 
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");
105 
106  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
107  if (!pass) {
108  std::cout << std::endl;
109  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
110  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
111  }
112  TEST_ASSERT(pass)
113  }
114 
115  // Test constructor IntegratorBasic(model, stepperType)
116  {
117  RCP<Tempus::IntegratorBasic<double> > integrator =
118  Tempus::integratorBasic<double>(model, RKMethods[m]);
119 
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");
127 
128  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
129  if (!pass) {
130  std::cout << std::endl;
131  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
132  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
133  }
134  TEST_ASSERT(pass)
135  }
136  }
137 }
138 #endif // TEST_PARAMETERLIST
139 
140 
141 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
142 // ************************************************************
143 // ************************************************************
144 TEUCHOS_UNIT_TEST(ExplicitRK, ConstructingFromDefaults)
145 {
146  double dt = 0.1;
147 
148  // Read params from .xml file
149  RCP<ParameterList> pList =
150  getParametersFromXmlFile("Tempus_ExplicitRK_SinCos.xml");
151  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
152 
153  // Setup the SinCosModel
154  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
155  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
156  auto model = rcp(new SinCosModel<double>(scm_pl));
157 
158  // Setup Stepper for field solve ----------------------------
159  RCP<Tempus::StepperFactory<double> > sf =
160  Teuchos::rcp(new Tempus::StepperFactory<double>());
161  RCP<Tempus::Stepper<double> > stepper =
162  sf->createStepper("RK Explicit 4 Stage");
163  stepper->setModel(model);
164  stepper->initialize();
165 
166  // Setup TimeStepControl ------------------------------------
167  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
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();
176 
177  // Setup initial condition SolutionState --------------------
178  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
179  stepper->getModel()->getNominalValues();
180  auto icSolution = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
181  auto icState = rcp(new Tempus::SolutionState<double>(icSolution));
182  icState->setTime (timeStepControl->getInitTime());
183  icState->setIndex (timeStepControl->getInitIndex());
184  icState->setTimeStep(0.0);
185  icState->setOrder (stepper->getOrder());
186  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
187 
188  // Setup SolutionHistory ------------------------------------
190  solutionHistory->setName("Forward States");
192  solutionHistory->setStorageLimit(2);
193  solutionHistory->addState(icState);
194 
195  // Setup Integrator -----------------------------------------
196  RCP<Tempus::IntegratorBasic<double> > integrator =
197  Tempus::integratorBasic<double>();
198  integrator->setStepperWStepper(stepper);
199  integrator->setTimeStepControl(timeStepControl);
200  integrator->setSolutionHistory(solutionHistory);
201  //integrator->setObserver(...);
202  integrator->initialize();
203 
204 
205  // Integrate to timeMax
206  bool integratorStatus = integrator->advanceTime();
207  TEST_ASSERT(integratorStatus)
208 
209 
210  // Test if at 'Final Time'
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);
215 
216  // Time-integrated solution and the exact solution
217  RCP<Thyra::VectorBase<double> > x = integrator->getX();
218  RCP<const Thyra::VectorBase<double> > x_exact =
219  model->getExactSolution(time).get_x();
220 
221  // Calculate the error
222  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
223  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
224 
225  // Check the order and intercept
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 );
237 }
238 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
239 
240 
241 #ifdef TEST_SINCOS
242 // ************************************************************
243 // ************************************************************
244 TEUCHOS_UNIT_TEST(ExplicitRK, SinCos)
245 {
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"); // slope = 3.87816
261  RKMethods.push_back("General ERK Embedded");
262 
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);
279 
280  for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
281 
282  std::string RKMethod = RKMethods[m];
283  std::replace(RKMethod.begin(), RKMethod.end(), ' ', '_');
284  std::replace(RKMethod.begin(), RKMethod.end(), '/', '.');
285 
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;
292 
293  const int nTimeStepSizes = 7;
294  double dt = 0.2;
295  double time = 0.0;
296  for (int n=0; n<nTimeStepSizes; n++) {
297 
298  // Read params from .xml file
299  RCP<ParameterList> pList =
300  getParametersFromXmlFile("Tempus_ExplicitRK_SinCos.xml");
301 
302  // Setup the SinCosModel
303  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
304  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
305  auto model = rcp(new SinCosModel<double>(scm_pl));
306 
307  // Set the Stepper
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");
313  } else {
314  pl->sublist("Demo Stepper").set("Stepper Type", RKMethods[m]);
315  }
316 
317 
318  dt /= 2;
319 
320  // Setup the Integrator and reset initial time step
321  pl->sublist("Demo Integrator")
322  .sublist("Time Step Control").set("Initial Time Step", dt);
323  integrator = Tempus::integratorBasic<double>(pl, model);
324 
325  // Initial Conditions
326  // During the Integrator construction, the initial SolutionState
327  // is set by default to model->getNominalVales().get_x(). However,
328  // the application can set it also by integrator->initializeSolutionHistory.
329  RCP<Thyra::VectorBase<double> > x0 =
330  model->getNominalValues().get_x()->clone_v();
331  integrator->initializeSolutionHistory(0.0, x0);
332 
333  // Integrate to timeMax
334  bool integratorStatus = integrator->advanceTime();
335  TEST_ASSERT(integratorStatus)
336 
337  // Test if at 'Final Time'
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);
342 
343  // Time-integrated solution and the exact solution
344  RCP<Thyra::VectorBase<double> > x = integrator->getX();
345  RCP<const Thyra::VectorBase<double> > x_exact =
346  model->getExactSolution(time).get_x();
347 
348  // Plot sample solution and exact solution
349  if (n == 0) {
350  RCP<const SolutionHistory<double> > solutionHistory =
351  integrator->getSolutionHistory();
352  writeSolution("Tempus_"+RKMethod+"_SinCos.dat", solutionHistory);
353 
354  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
355  for (int i=0; i<solutionHistory->getNumStates(); i++) {
356  double time_i = (*solutionHistory)[i]->getTime();
357  auto state = rcp(new Tempus::SolutionState<double>(
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);
362  }
363  writeSolution("Tempus_"+RKMethod+"_SinCos-Ref.dat", solnHistExact);
364  }
365 
366  // Store off the final solution and step size
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) { // Add exact solution last in vector.
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);
383  }
384  }
385 
386  // Check the order and intercept
387  double xSlope = 0.0;
388  double xDotSlope = 0.0;
389  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
390  double order = stepper->getOrder();
391  writeOrderError("Tempus_"+RKMethod+"_SinCos-Error.dat",
392  stepper, StepSize,
393  solutions, xErrorNorm, xSlope,
394  solutionsDot, xDotErrorNorm, xDotSlope);
395 
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 );
400  // xDot not yet available for ExplicitRK methods.
401  //TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
402  //TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.0486418, 1.0e-4 );
403 
404  }
405  //Teuchos::TimeMonitor::summarize();
406 }
407 #endif // TEST_SINCOS
408 
409 
410 #ifdef TEST_EMBEDDED_VANDERPOL
411 // ************************************************************
412 // ************************************************************
413 TEUCHOS_UNIT_TEST(ExplicitRK, EmbeddedVanDerPol)
414 {
415 
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");
422 
423  // the embedded solution will test the following:
424  // using the starting stepsize routine, this has now decreased
425  const int refIstep = 45;
426 
427  for(auto integratorChoice : IntegratorList){
428 
429  std::cout << "Using Integrator: " << integratorChoice << " !!!" << std::endl;
430 
431  // Read params from .xml file
432  RCP<ParameterList> pList =
433  getParametersFromXmlFile("Tempus_ExplicitRK_VanDerPol.xml");
434 
435 
436  // Setup the VanDerPolModel
437  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
438  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
439 
440 
441  // Set the Integrator and Stepper
442  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
443  pl->set("Integrator Name", integratorChoice);
444 
445  // Setup the Integrator
446  RCP<Tempus::IntegratorBasic<double> > integrator =
447  Tempus::integratorBasic<double>(pl, model);
448 
449  const std::string RKMethod =
450  pl->sublist(integratorChoice).get<std::string>("Stepper Name");
451 
452  // Integrate to timeMax
453  bool integratorStatus = integrator->advanceTime();
454  TEST_ASSERT(integratorStatus);
455 
456  // Test if at 'Final Time'
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);
461 
462 
463  // Numerical reference solution at timeFinal (for \epsilon = 0.1)
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());
468 
469  // Calculate the error
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);
473 
474  // Test number of steps, failures, and accuracy
475  if ((integratorChoice == "Embedded_Integrator_PID") or
476  (integratorChoice == "Embedded_Integrator_PID_General")) {
477 
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");
482 
483  // Should be close to the prescribed tolerance (magnitude)
484  TEST_COMPARE(std::log10(L2norm), <, std::log10(absTol));
485  TEST_COMPARE(std::log10(L2norm), <, std::log10(relTol));
486 
487  // get the number of time steps and number of step failure
488  //const int nFailure_c = integrator->getSolutionHistory()->
489  //getCurrentState()->getMetaData()->getNFailures();
490  const int iStep = integrator->getSolutionHistory()->
491  getCurrentState()->getIndex();
492  const int nFail = integrator->getSolutionHistory()->
493  getCurrentState()->getMetaData()->getNRunningFailures();
494 
495  // test for number of steps
496  TEST_EQUALITY(iStep, refIstep);
497  std::cout << "Tolerance = " << absTol
498  << " L2norm = " << L2norm
499  << " iStep = " << iStep
500  << " nFail = " << nFail << std::endl;
501  }
502 
503  // Plot sample solution and exact solution
504  std::ofstream ftmp("Tempus_"+integratorChoice+RKMethod+"_VDP_Example.dat");
505  RCP<const SolutionHistory<double> > solutionHistory =
506  integrator->getSolutionHistory();
507  int nStates = solutionHistory->getNumStates();
508  //RCP<const Thyra::VectorBase<double> > x_exact_plot;
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();
513  //x_exact_plot = model->getExactSolution(time_i).get_x();
514  ftmp << time_i << " "
515  << Thyra::get_ele(*(x_plot), 0) << " "
516  << Thyra::get_ele(*(x_plot), 1) << " " << std::endl;
517  }
518  ftmp.close();
519  }
520 
521  Teuchos::TimeMonitor::summarize();
522 }
523 #endif // TEST_EMBEDDED_VANDERPOL
524 
525 } // 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)
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...