Tempus  Version of the Day
Time Integration
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"
17 #include "Tempus_StepperExplicitRK.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::ParameterList;
30 using Teuchos::sublist;
31 using Teuchos::getParametersFromXmlFile;
32 
36 
37 // Comment out any of the following tests to exclude from build/run.
38 #define TEST_PARAMETERLIST
39 #define TEST_CONSTRUCTING_FROM_DEFAULTS
40 #define TEST_SINCOS
41 #define TEST_EMBEDDED_VANDERPOL
42 
43 
44 #ifdef TEST_PARAMETERLIST
45 // ************************************************************
46 // ************************************************************
47 TEUCHOS_UNIT_TEST(ExplicitRK, ParameterList)
48 {
49  std::vector<std::string> RKMethods;
50  RKMethods.push_back("Bogacki-Shampine 3(2) Pair");
51  RKMethods.push_back("Merson 4(5) Pair");
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 2 Stage 2nd order by Runge");
62  RKMethods.push_back("RK Explicit Trapezoidal");
63 
64  for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
65 
66  std::string RKMethod = RKMethods[m];
67  std::replace(RKMethod.begin(), RKMethod.end(), ' ', '_');
68  std::replace(RKMethod.begin(), RKMethod.end(), '/', '.');
69 
70  // Read params from .xml file
71  RCP<ParameterList> pList =
72  getParametersFromXmlFile("Tempus_ExplicitRK_SinCos.xml");
73 
74  // Setup the SinCosModel
75  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
76  RCP<SinCosModel<double> > model =
77  Teuchos::rcp(new SinCosModel<double>(scm_pl));
78 
79  // Set the Stepper
80  RCP<ParameterList> tempusPL = sublist(pList, "Tempus", true);
81  if (RKMethods[m] == "General ERK") {
82  tempusPL->sublist("Demo Integrator").set("Stepper Name", "Demo Stepper 2");
83  } else {
84  tempusPL->sublist("Demo Stepper").set("Stepper Type", RKMethods[m]);
85  }
86 
87  // Test constructor IntegratorBasic(tempusPL, model)
88  {
89  RCP<Tempus::IntegratorBasic<double> > integrator =
90  Tempus::integratorBasic<double>(tempusPL, model);
91 
92  RCP<ParameterList> stepperPL = sublist(tempusPL, "Demo Stepper", true);
93  if (RKMethods[m] == "General ERK")
94  stepperPL = sublist(tempusPL, "Demo Stepper 2", true);
95  RCP<ParameterList> defaultPL =
96  integrator->getStepper()->getDefaultParameters();
97  defaultPL->remove("Description");
98 
99  TEST_ASSERT(haveSameValues(*stepperPL, *defaultPL, true));
100  }
101 
102  // Test constructor IntegratorBasic(model, stepperType)
103  {
104  RCP<Tempus::IntegratorBasic<double> > integrator =
105  Tempus::integratorBasic<double>(model, RKMethods[m]);
106 
107  RCP<ParameterList> stepperPL = sublist(tempusPL, "Demo Stepper", true);
108  if (RKMethods[m] == "General ERK")
109  stepperPL = sublist(tempusPL, "Demo Stepper 2", true);
110  RCP<ParameterList> defaultPL =
111  integrator->getStepper()->getDefaultParameters();
112  defaultPL->remove("Description");
113 
114  TEST_ASSERT(haveSameValues(*stepperPL, *defaultPL, true))
115  }
116  }
117 }
118 #endif // TEST_PARAMETERLIST
119 
120 
121 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
122 // ************************************************************
123 // ************************************************************
124 TEUCHOS_UNIT_TEST(ExplicitRK, ConstructingFromDefaults)
125 {
126  double dt = 0.1;
127 
128  // Read params from .xml file
129  RCP<ParameterList> pList =
130  getParametersFromXmlFile("Tempus_ExplicitRK_SinCos.xml");
131  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
132 
133  // Setup the SinCosModel
134  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
135  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
136  RCP<SinCosModel<double> > model =
137  Teuchos::rcp(new SinCosModel<double>(scm_pl));
138 
139  // Setup Stepper for field solve ----------------------------
140  RCP<Tempus::StepperExplicitRK<double> > stepper =
141  Teuchos::rcp(new Tempus::StepperExplicitRK<double>(model));
142 
143  // Setup TimeStepControl ------------------------------------
144  RCP<Tempus::TimeStepControl<double> > timeStepControl =
145  Teuchos::rcp(new Tempus::TimeStepControl<double>());
146  ParameterList tscPL = pl->sublist("Demo Integrator")
147  .sublist("Time Step Control");
148  timeStepControl->setStepType (tscPL.get<std::string>("Integrator Step Type"));
149  timeStepControl->setInitIndex(tscPL.get<int> ("Initial Time Index"));
150  timeStepControl->setInitTime (tscPL.get<double>("Initial Time"));
151  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
152  timeStepControl->setInitTimeStep(dt);
153  timeStepControl->initialize();
154 
155  // Setup initial condition SolutionState --------------------
156  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
157  stepper->getModel()->getNominalValues();
158  RCP<Thyra::VectorBase<double> > icSolution =
159  Teuchos::rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
160  RCP<Tempus::SolutionState<double> > icState =
161  Teuchos::rcp(new Tempus::SolutionState<double>(icSolution));
162  icState->setTime (timeStepControl->getInitTime());
163  icState->setIndex (timeStepControl->getInitIndex());
164  icState->setTimeStep(0.0);
165  icState->setOrder (stepper->getOrder());
166  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
167 
168  // Setup SolutionHistory ------------------------------------
169  RCP<Tempus::SolutionHistory<double> > solutionHistory =
170  Teuchos::rcp(new Tempus::SolutionHistory<double>());
171  solutionHistory->setName("Forward States");
173  solutionHistory->setStorageLimit(2);
174  solutionHistory->addState(icState);
175 
176  // Setup Integrator -----------------------------------------
177  RCP<Tempus::IntegratorBasic<double> > integrator =
178  Tempus::integratorBasic<double>();
179  integrator->setStepperWStepper(stepper);
180  integrator->setTimeStepControl(timeStepControl);
181  integrator->setSolutionHistory(solutionHistory);
182  //integrator->setObserver(...);
183  integrator->initialize();
184 
185 
186  // Integrate to timeMax
187  bool integratorStatus = integrator->advanceTime();
188  TEST_ASSERT(integratorStatus)
189 
190 
191  // Test if at 'Final Time'
192  double time = integrator->getTime();
193  double timeFinal =pl->sublist("Demo Integrator")
194  .sublist("Time Step Control").get<double>("Final Time");
195  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
196 
197  // Time-integrated solution and the exact solution
198  RCP<Thyra::VectorBase<double> > x = integrator->getX();
199  RCP<const Thyra::VectorBase<double> > x_exact =
200  model->getExactSolution(time).get_x();
201 
202  // Calculate the error
203  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
204  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
205 
206  // Check the order and intercept
207  std::cout << " Stepper = RK Explicit 4 Stage" << std::endl;
208  std::cout << " =========================" << std::endl;
209  std::cout << " Exact solution : " << get_ele(*(x_exact), 0) << " "
210  << get_ele(*(x_exact), 1) << std::endl;
211  std::cout << " Computed solution: " << get_ele(*(x ), 0) << " "
212  << get_ele(*(x ), 1) << std::endl;
213  std::cout << " Difference : " << get_ele(*(xdiff ), 0) << " "
214  << get_ele(*(xdiff ), 1) << std::endl;
215  std::cout << " =========================" << std::endl;
216  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841470, 1.0e-4 );
217  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.540303, 1.0e-4 );
218 }
219 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
220 
221 
222 #ifdef TEST_SINCOS
223 // ************************************************************
224 // ************************************************************
225 TEUCHOS_UNIT_TEST(ExplicitRK, SinCos)
226 {
227  std::vector<std::string> RKMethods;
228  RKMethods.push_back("RK Forward Euler");
229  RKMethods.push_back("RK Explicit 4 Stage");
230  RKMethods.push_back("RK Explicit 3/8 Rule");
231  RKMethods.push_back("RK Explicit 4 Stage 3rd order by Runge");
232  RKMethods.push_back("RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
233  RKMethods.push_back("RK Explicit 3 Stage 3rd order");
234  RKMethods.push_back("RK Explicit 3 Stage 3rd order TVD");
235  RKMethods.push_back("RK Explicit 3 Stage 3rd order by Heun");
236  RKMethods.push_back("RK Explicit 2 Stage 2nd order by Runge");
237  RKMethods.push_back("RK Explicit Trapezoidal");
238  RKMethods.push_back("Bogacki-Shampine 3(2) Pair");
239  RKMethods.push_back("General ERK");
240  RKMethods.push_back("General ERK Embedded");
241  //RKMethods.push_back("Merson 4(5) Pair"); // slope = 3.87816
242  std::vector<double> RKMethodErrors;
243  RKMethodErrors.push_back(0.051123);
244  RKMethodErrors.push_back(8.33251e-07);
245  RKMethodErrors.push_back(8.33251e-07);
246  RKMethodErrors.push_back(4.16897e-05);
247  RKMethodErrors.push_back(8.32108e-06);
248  RKMethodErrors.push_back(4.16603e-05);
249  RKMethodErrors.push_back(4.16603e-05);
250  RKMethodErrors.push_back(4.16603e-05);
251  RKMethodErrors.push_back(0.00166645);
252  RKMethodErrors.push_back(0.00166645);
253  RKMethodErrors.push_back(4.16603e-05);
254  RKMethodErrors.push_back(8.33251e-07);
255  RKMethodErrors.push_back(4.16603e-05);
256  //RKMethodErrors.push_back(1.39383e-07);
257 
258  for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
259 
260  std::string RKMethod = RKMethods[m];
261  std::replace(RKMethod.begin(), RKMethod.end(), ' ', '_');
262  std::replace(RKMethod.begin(), RKMethod.end(), '/', '.');
263 
264  RCP<Tempus::IntegratorBasic<double> > integrator;
265  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
266  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
267  std::vector<double> StepSize;
268  std::vector<double> xErrorNorm;
269  std::vector<double> xDotErrorNorm;
270 
271  const int nTimeStepSizes = 7;
272  double dt = 0.2;
273  double time = 0.0;
274  for (int n=0; n<nTimeStepSizes; n++) {
275 
276  // Read params from .xml file
277  RCP<ParameterList> pList =
278  getParametersFromXmlFile("Tempus_ExplicitRK_SinCos.xml");
279 
280  // Setup the SinCosModel
281  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
282  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
283  RCP<SinCosModel<double> > model =
284  Teuchos::rcp(new SinCosModel<double>(scm_pl));
285 
286  // Set the Stepper
287  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
288  if (RKMethods[m] == "General ERK") {
289  pl->sublist("Demo Integrator").set("Stepper Name", "Demo Stepper 2");
290  } else if (RKMethods[m] == "General ERK Embedded"){
291  pl->sublist("Demo Integrator").set("Stepper Name", "General ERK Embedded Stepper");
292  } else {
293  pl->sublist("Demo Stepper").set("Stepper Type", RKMethods[m]);
294  }
295 
296 
297  dt /= 2;
298 
299  // Setup the Integrator and reset initial time step
300  pl->sublist("Demo Integrator")
301  .sublist("Time Step Control").set("Initial Time Step", dt);
302  integrator = Tempus::integratorBasic<double>(pl, model);
303 
304  // Initial Conditions
305  // During the Integrator construction, the initial SolutionState
306  // is set by default to model->getNominalVales().get_x(). However,
307  // the application can set it also by integrator->setInitialState.
308  RCP<Thyra::VectorBase<double> > x0 =
309  model->getNominalValues().get_x()->clone_v();
310  integrator->setInitialState(0.0, x0);
311 
312  // Integrate to timeMax
313  bool integratorStatus = integrator->advanceTime();
314  TEST_ASSERT(integratorStatus)
315 
316  // Test if at 'Final Time'
317  time = integrator->getTime();
318  double timeFinal = pl->sublist("Demo Integrator")
319  .sublist("Time Step Control").get<double>("Final Time");
320  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
321 
322  // Time-integrated solution and the exact solution
323  RCP<Thyra::VectorBase<double> > x = integrator->getX();
324  RCP<const Thyra::VectorBase<double> > x_exact =
325  model->getExactSolution(time).get_x();
326 
327  // Plot sample solution and exact solution
328  if (n == 0) {
329  RCP<const SolutionHistory<double> > solutionHistory =
330  integrator->getSolutionHistory();
331  writeSolution("Tempus_"+RKMethod+"_SinCos.dat", solutionHistory);
332 
333  RCP<Tempus::SolutionHistory<double> > solnHistExact =
334  Teuchos::rcp(new Tempus::SolutionHistory<double>());
335  for (int i=0; i<solutionHistory->getNumStates(); i++) {
336  double time = (*solutionHistory)[i]->getTime();
337  RCP<Tempus::SolutionState<double> > state =
338  Teuchos::rcp(new Tempus::SolutionState<double>(
339  model->getExactSolution(time).get_x(),
340  model->getExactSolution(time).get_x_dot()));
341  state->setTime((*solutionHistory)[i]->getTime());
342  solnHistExact->addState(state);
343  }
344  writeSolution("Tempus_"+RKMethod+"_SinCos-Ref.dat", solnHistExact);
345  }
346 
347  // Store off the final solution and step size
348  StepSize.push_back(dt);
349  auto solution = Thyra::createMember(model->get_x_space());
350  Thyra::copy(*(integrator->getX()),solution.ptr());
351  solutions.push_back(solution);
352  auto solutionDot = Thyra::createMember(model->get_x_space());
353  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
354  solutionsDot.push_back(solutionDot);
355  if (n == nTimeStepSizes-1) { // Add exact solution last in vector.
356  StepSize.push_back(0.0);
357  auto solution = Thyra::createMember(model->get_x_space());
358  Thyra::copy(*(model->getExactSolution(time).get_x()),solution.ptr());
359  solutions.push_back(solution);
360  auto solutionDot = Thyra::createMember(model->get_x_space());
361  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
362  solutionDot.ptr());
363  solutionsDot.push_back(solutionDot);
364  }
365  }
366 
367  // Check the order and intercept
368  double xSlope = 0.0;
369  double xDotSlope = 0.0;
370  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
371  double order = stepper->getOrder();
372  writeOrderError("Tempus_"+RKMethod+"_SinCos-Error.dat",
373  stepper, StepSize,
374  solutions, xErrorNorm, xSlope,
375  solutionsDot, xDotErrorNorm, xDotSlope);
376 
377  TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
378  TEST_FLOATING_EQUALITY( xErrorNorm[0], RKMethodErrors[m], 1.0e-4 );
379  // xDot not yet available for ExplicitRK methods.
380  //TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
381  //TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.0486418, 1.0e-4 );
382 
383  }
384  //Teuchos::TimeMonitor::summarize();
385 }
386 #endif // TEST_SINCOS
387 
388 
389 #ifdef TEST_EMBEDDED_VANDERPOL
390 // ************************************************************
391 // ************************************************************
392 TEUCHOS_UNIT_TEST(ExplicitRK, EmbeddedVanDerPol)
393 {
394 
395  std::vector<std::string> IntegratorList;
396  IntegratorList.push_back("Embedded_Integrator_PID");
397  IntegratorList.push_back("Demo_Integrator");
398  IntegratorList.push_back("Embedded_Integrator");
399  IntegratorList.push_back("General_Embedded_Integrator");
400  IntegratorList.push_back("Embedded_Integrator_PID_General");
401 
402  // the embedded solution will test the following:
403  // using the starting stepsize routine, this has now decreased
404  const int refIstep = 45;
405 
406  for(auto integratorChoice : IntegratorList){
407 
408  std::cout << "Using Integrator: " << integratorChoice << " !!!" << std::endl;
409 
410  // Read params from .xml file
411  RCP<ParameterList> pList =
412  getParametersFromXmlFile("Tempus_ExplicitRK_VanDerPol.xml");
413 
414 
415  // Setup the VanDerPolModel
416  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
417  RCP<VanDerPolModel<double> > model =
418  Teuchos::rcp(new VanDerPolModel<double>(vdpm_pl));
419 
420 
421  // Set the Integrator and Stepper
422  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
423  pl->set("Integrator Name", integratorChoice);
424 
425  // Setup the Integrator
426  RCP<Tempus::IntegratorBasic<double> > integrator =
427  Tempus::integratorBasic<double>(pl, model);
428 
429  const std::string RKMethod =
430  pl->sublist(integratorChoice).get<std::string>("Stepper Name");
431 
432  // Integrate to timeMax
433  bool integratorStatus = integrator->advanceTime();
434  TEST_ASSERT(integratorStatus);
435 
436  // Test if at 'Final Time'
437  double time = integrator->getTime();
438  double timeFinal = pl->sublist(integratorChoice)
439  .sublist("Time Step Control").get<double>("Final Time");
440  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
441 
442 
443  // Numerical reference solution at timeFinal (for \epsilon = 0.1)
444  RCP<Thyra::VectorBase<double> > x = integrator->getX();
445  RCP<Thyra::VectorBase<double> > xref = x->clone_v();
446  Thyra::set_ele(0, -1.5484458614405929, xref.ptr());
447  Thyra::set_ele(1, 1.0181127316101317, xref.ptr());
448 
449  // Calculate the error
450  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
451  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *xref, -1.0, *(x));
452  const double L2norm = Thyra::norm_2(*xdiff);
453 
454  // Test number of steps, failures, and accuracy
455  if ((integratorChoice == "Embedded_Integrator_PID") or
456  (integratorChoice == "Embedded_Integrator_PID_General")) {
457 
458  const double absTol = pl->sublist(integratorChoice).
459  sublist("Time Step Control").get<double>("Maximum Absolute Error");
460  const double relTol = pl->sublist(integratorChoice).
461  sublist("Time Step Control").get<double>("Maximum Relative Error");
462 
463  // Should be close to the prescribed tolerance (magnitude)
464  TEST_COMPARE(std::log10(L2norm), <, std::log10(absTol));
465  TEST_COMPARE(std::log10(L2norm), <, std::log10(relTol));
466 
467  // get the number of time steps and number of step failure
468  //const int nFailure_c = integrator->getSolutionHistory()->
469  //getCurrentState()->getMetaData()->getNFailures();
470  const int iStep = integrator->getSolutionHistory()->
471  getCurrentState()->getIndex();
472  const int nFail = integrator->getSolutionHistory()->
473  getCurrentState()->getMetaData()->getNRunningFailures();
474 
475  // test for number of steps
476  TEST_EQUALITY(iStep, refIstep);
477  std::cout << "Tolerance = " << absTol
478  << " L2norm = " << L2norm
479  << " iStep = " << iStep
480  << " nFail = " << nFail << std::endl;
481  }
482 
483  // Plot sample solution and exact solution
484  std::ofstream ftmp("Tempus_"+integratorChoice+RKMethod+"_VDP_Example.dat");
485  RCP<const SolutionHistory<double> > solutionHistory =
486  integrator->getSolutionHistory();
487  int nStates = solutionHistory->getNumStates();
488  //RCP<const Thyra::VectorBase<double> > x_exact_plot;
489  for (int i=0; i<nStates; i++) {
490  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
491  double time = solutionState->getTime();
492  RCP<const Thyra::VectorBase<double> > x_plot = solutionState->getX();
493  //x_exact_plot = model->getExactSolution(time).get_x();
494  ftmp << time << " "
495  << Thyra::get_ele(*(x_plot), 0) << " "
496  << Thyra::get_ele(*(x_plot), 1) << " " << std::endl;
497  }
498  ftmp.close();
499  }
500 
501  Teuchos::TimeMonitor::summarize();
502 }
503 #endif // TEST_EMBEDDED_VANDERPOL
504 
505 } // namespace Tempus_Test
Tempus_Test::writeSolution
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > solutionHistory)
Definition: Tempus_ConvergenceTestUtils.hpp:302
Tempus::StepperExplicitRK
Explicit Runge-Kutta time stepper.
Definition: Tempus_StepperExplicitRK_decl.hpp:64
Tempus_Test::writeOrderError
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)
Definition: Tempus_ConvergenceTestUtils.hpp:184
Tempus_Test
Definition: Tempus_BackwardEuler_ASA.cpp:34
Tempus::solutionHistory
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
Definition: Tempus_SolutionHistory_impl.hpp:504
Tempus::STORAGE_TYPE_STATIC
Keep a fix number of states.
Definition: Tempus_SolutionHistory_decl.hpp:30
Tempus::IntegratorBasic
Basic time integrator.
Definition: Tempus_IntegratorBasic_decl.hpp:35
Tempus_Test::SinCosModel
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
Definition: SinCosModel_decl.hpp:91
Tempus::SolutionState
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
Definition: Tempus_SolutionState_decl.hpp:56
Tempus_Test::TEUCHOS_UNIT_TEST
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
Definition: Tempus_BackwardEuler_ASA.cpp:47
Tempus::TimeStepControl
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
Definition: Tempus_Integrator.hpp:26
Tempus::PASSED
Definition: Tempus_Types.hpp:17
Tempus::SolutionHistory
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Definition: Tempus_Integrator.hpp:25
Tempus_Test::VanDerPolModel
van der Pol model problem for nonlinear electrical circuit.
Definition: VanDerPolModel_decl.hpp:111