Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_ForwardEulerTest.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 #include "Thyra_DetachedVectorView.hpp"
15 
16 #include "Tempus_IntegratorBasic.hpp"
17 
18 #include "Tempus_StepperForwardEuler.hpp"
19 
20 #include "../TestModels/SinCosModel.hpp"
21 #include "../TestModels/VanDerPolModel.hpp"
22 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
23 
24 #include <fstream>
25 #include <vector>
26 
27 namespace Tempus_Test {
28 
29 using Teuchos::RCP;
30 using Teuchos::rcp;
31 using Teuchos::rcp_const_cast;
32 using Teuchos::ParameterList;
33 using Teuchos::sublist;
34 using Teuchos::getParametersFromXmlFile;
35 
39 
40 
41 // ************************************************************
42 // ************************************************************
43 TEUCHOS_UNIT_TEST(ForwardEuler, ParameterList)
44 {
45  // Read params from .xml file
46  RCP<ParameterList> pList =
47  getParametersFromXmlFile("Tempus_ForwardEuler_SinCos.xml");
48 
49  //std::ofstream ftmp("PL.txt");
50  //pList->print(ftmp);
51  //ftmp.close();
52 
53  // Setup the SinCosModel
54  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
55  auto model = rcp(new SinCosModel<double> (scm_pl));
56 
57  RCP<ParameterList> tempusPL = sublist(pList, "Tempus", true);
58 
59  // Test constructor IntegratorBasic(tempusPL, model)
60  {
61  RCP<Tempus::IntegratorBasic<double> > integrator =
62  Tempus::integratorBasic<double>(tempusPL, model);
63 
64  RCP<ParameterList> stepperPL = sublist(tempusPL, "Demo Stepper", true);
65  RCP<const ParameterList> defaultPL =
66  integrator->getStepper()->getValidParameters();
67 
68  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
69  if (!pass) {
70  std::cout << std::endl;
71  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
72  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
73  }
74  TEST_ASSERT(pass)
75  }
76 
77  // Test constructor IntegratorBasic(model, stepperType)
78  {
79  RCP<Tempus::IntegratorBasic<double> > integrator =
80  Tempus::integratorBasic<double>(model, "Forward Euler");
81 
82  RCP<ParameterList> stepperPL = sublist(tempusPL, "Demo Stepper", true);
83  RCP<const ParameterList> defaultPL =
84  integrator->getStepper()->getValidParameters();
85 
86  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
87  if (!pass) {
88  std::cout << std::endl;
89  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
90  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
91  }
92  TEST_ASSERT(pass)
93  }
94 }
95 
96 
97 // ************************************************************
98 // ************************************************************
99 TEUCHOS_UNIT_TEST(ForwardEuler, ConstructingFromDefaults)
100 {
101  double dt = 0.1;
102 
103  // Read params from .xml file
104  RCP<ParameterList> pList =
105  getParametersFromXmlFile("Tempus_ForwardEuler_SinCos.xml");
106  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
107 
108  // Setup the SinCosModel
109  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
110  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
111  auto model = rcp(new SinCosModel<double>(scm_pl));
112 
113  // Setup Stepper for field solve ----------------------------
114  auto stepper = rcp(new Tempus::StepperForwardEuler<double>());
115  stepper->setModel(model);
116  stepper->initialize();
117 
118  // Setup TimeStepControl ------------------------------------
119  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
120  ParameterList tscPL = pl->sublist("Demo Integrator")
121  .sublist("Time Step Control");
122  timeStepControl->setStepType (tscPL.get<std::string>("Integrator Step Type"));
123  timeStepControl->setInitIndex(tscPL.get<int> ("Initial Time Index"));
124  timeStepControl->setInitTime (tscPL.get<double>("Initial Time"));
125  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
126  timeStepControl->setInitTimeStep(dt);
127  timeStepControl->initialize();
128 
129  // Setup initial condition SolutionState --------------------
130  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
131  stepper->getModel()->getNominalValues();
132  auto icSolution = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
133  auto icState = Tempus::createSolutionStateX(icSolution);
134  icState->setTime (timeStepControl->getInitTime());
135  icState->setIndex (timeStepControl->getInitIndex());
136  icState->setTimeStep(0.0);
137  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
138 
139  // Setup SolutionHistory ------------------------------------
140  auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
141  solutionHistory->setName("Forward States");
142  solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
143  solutionHistory->setStorageLimit(2);
144  solutionHistory->addState(icState);
145 
146  // Setup Integrator -----------------------------------------
147  RCP<Tempus::IntegratorBasic<double> > integrator =
148  Tempus::integratorBasic<double>();
149  integrator->setStepperWStepper(stepper);
150  integrator->setTimeStepControl(timeStepControl);
151  integrator->setSolutionHistory(solutionHistory);
152  //integrator->setObserver(...);
153  integrator->initialize();
154 
155 
156  // Integrate to timeMax
157  bool integratorStatus = integrator->advanceTime();
158  TEST_ASSERT(integratorStatus)
159 
160 
161  // Test if at 'Final Time'
162  double time = integrator->getTime();
163  double timeFinal =pl->sublist("Demo Integrator")
164  .sublist("Time Step Control").get<double>("Final Time");
165  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
166 
167  // Time-integrated solution and the exact solution
168  RCP<Thyra::VectorBase<double> > x = integrator->getX();
169  RCP<const Thyra::VectorBase<double> > x_exact =
170  model->getExactSolution(time).get_x();
171 
172  // Calculate the error
173  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
174  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
175 
176  // Check the order and intercept
177  std::cout << " Stepper = ForwardEuler" << std::endl;
178  std::cout << " =========================" << std::endl;
179  std::cout << " Exact solution : " << get_ele(*(x_exact), 0) << " "
180  << get_ele(*(x_exact), 1) << std::endl;
181  std::cout << " Computed solution: " << get_ele(*(x ), 0) << " "
182  << get_ele(*(x ), 1) << std::endl;
183  std::cout << " Difference : " << get_ele(*(xdiff ), 0) << " "
184  << get_ele(*(xdiff ), 1) << std::endl;
185  std::cout << " =========================" << std::endl;
186  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.882508, 1.0e-4 );
187  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.570790, 1.0e-4 );
188 }
189 
190 
191 // ************************************************************
192 // ************************************************************
193 TEUCHOS_UNIT_TEST(ForwardEuler, SinCos)
194 {
195  RCP<Tempus::IntegratorBasic<double> > integrator;
196  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
197  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
198  std::vector<double> StepSize;
199  std::vector<double> xErrorNorm;
200  std::vector<double> xDotErrorNorm;
201  const int nTimeStepSizes = 7;
202  double dt = 0.2;
203  double time = 0.0;
204  for (int n=0; n<nTimeStepSizes; n++) {
205 
206  // Read params from .xml file
207  RCP<ParameterList> pList =
208  getParametersFromXmlFile("Tempus_ForwardEuler_SinCos.xml");
209 
210  //std::ofstream ftmp("PL.txt");
211  //pList->print(ftmp);
212  //ftmp.close();
213 
214  // Setup the SinCosModel
215  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
216  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
217  auto model = rcp(new SinCosModel<double> (scm_pl));
218 
219  dt /= 2;
220 
221  // Setup the Integrator and reset initial time step
222  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
223  pl->sublist("Demo Integrator")
224  .sublist("Time Step Control").set("Initial Time Step", dt);
225  integrator = Tempus::integratorBasic<double>(pl, model);
226 
227  // Initial Conditions
228  // During the Integrator construction, the initial SolutionState
229  // is set by default to model->getNominalVales().get_x(). However,
230  // the application can set it also by integrator->initializeSolutionHistory.
231  RCP<Thyra::VectorBase<double> > x0 =
232  model->getNominalValues().get_x()->clone_v();
233  integrator->initializeSolutionHistory(0.0, x0);
234 
235  // Integrate to timeMax
236  bool integratorStatus = integrator->advanceTime();
237  TEST_ASSERT(integratorStatus)
238 
239  // Test PhysicsState
240  RCP<Tempus::PhysicsState<double> > physicsState =
241  integrator->getSolutionHistory()->getCurrentState()->getPhysicsState();
242  TEST_EQUALITY(physicsState->getName(), "Tempus::PhysicsState");
243 
244  // Test if at 'Final Time'
245  time = integrator->getTime();
246  double timeFinal = pl->sublist("Demo Integrator")
247  .sublist("Time Step Control").get<double>("Final Time");
248  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
249 
250  // Time-integrated solution and the exact solution
251  RCP<Thyra::VectorBase<double> > x = integrator->getX();
252  RCP<const Thyra::VectorBase<double> > x_exact =
253  model->getExactSolution(time).get_x();
254 
255  // Plot sample solution and exact solution
256  if (n == 0) {
257  RCP<const SolutionHistory<double> > solutionHistory =
258  integrator->getSolutionHistory();
259  writeSolution("Tempus_ForwardEuler_SinCos.dat", solutionHistory);
260 
261  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
262  for (int i=0; i<solutionHistory->getNumStates(); i++) {
263  double time_i = (*solutionHistory)[i]->getTime();
264  auto state = Tempus::createSolutionStateX(
265  rcp_const_cast<Thyra::VectorBase<double> > (
266  model->getExactSolution(time_i).get_x()),
267  rcp_const_cast<Thyra::VectorBase<double> > (
268  model->getExactSolution(time_i).get_x_dot()));
269  state->setTime((*solutionHistory)[i]->getTime());
270  solnHistExact->addState(state);
271  }
272  writeSolution("Tempus_ForwardEuler_SinCos-Ref.dat", solnHistExact);
273  }
274 
275  // Store off the final solution and step size
276  StepSize.push_back(dt);
277  auto solution = Thyra::createMember(model->get_x_space());
278  Thyra::copy(*(integrator->getX()),solution.ptr());
279  solutions.push_back(solution);
280  auto solutionDot = Thyra::createMember(model->get_x_space());
281  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
282  solutionsDot.push_back(solutionDot);
283  if (n == nTimeStepSizes-1) { // Add exact solution last in vector.
284  StepSize.push_back(0.0);
285  auto solutionExact = Thyra::createMember(model->get_x_space());
286  Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
287  solutions.push_back(solutionExact);
288  auto solutionDotExact = Thyra::createMember(model->get_x_space());
289  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
290  solutionDotExact.ptr());
291  solutionsDot.push_back(solutionDotExact);
292  }
293  }
294 
295  // Check the order and intercept
296  double xSlope = 0.0;
297  double xDotSlope = 0.0;
298  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
299  double order = stepper->getOrder();
300  writeOrderError("Tempus_ForwardEuler_SinCos-Error.dat",
301  stepper, StepSize,
302  solutions, xErrorNorm, xSlope,
303  solutionsDot, xDotErrorNorm, xDotSlope);
304 
305  TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
306  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.051123, 1.0e-4 );
307  // xDot not yet available for Forward Euler.
308  //TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
309  //TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.0486418, 1.0e-4 );
310 
311  Teuchos::TimeMonitor::summarize();
312 }
313 
314 
315 // ************************************************************
316 // ************************************************************
317 TEUCHOS_UNIT_TEST(ForwardEuler, VanDerPol)
318 {
319  RCP<Tempus::IntegratorBasic<double> > integrator;
320  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
321  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
322  std::vector<double> StepSize;
323  std::vector<double> xErrorNorm;
324  std::vector<double> xDotErrorNorm;
325  const int nTimeStepSizes = 7; // 8 for Error plot
326  double dt = 0.2;
327  for (int n=0; n<nTimeStepSizes; n++) {
328 
329  // Read params from .xml file
330  RCP<ParameterList> pList =
331  getParametersFromXmlFile("Tempus_ForwardEuler_VanDerPol.xml");
332 
333  // Setup the VanDerPolModel
334  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
335  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
336 
337  // Set the step size
338  dt /= 2;
339  if (n == nTimeStepSizes-1) dt /= 10.0;
340 
341  // Setup the Integrator and reset initial time step
342  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
343  pl->sublist("Demo Integrator")
344  .sublist("Time Step Control").set("Initial Time Step", dt);
345  integrator = Tempus::integratorBasic<double>(pl, model);
346 
347  // Integrate to timeMax
348  bool integratorStatus = integrator->advanceTime();
349  TEST_ASSERT(integratorStatus)
350 
351  // Test if at 'Final Time'
352  double time = integrator->getTime();
353  double timeFinal =pl->sublist("Demo Integrator")
354  .sublist("Time Step Control").get<double>("Final Time");
355  double tol = 100.0 * std::numeric_limits<double>::epsilon();
356  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
357 
358  // Store off the final solution and step size
359  StepSize.push_back(dt);
360  auto solution = Thyra::createMember(model->get_x_space());
361  Thyra::copy(*(integrator->getX()),solution.ptr());
362  solutions.push_back(solution);
363  auto solutionDot = Thyra::createMember(model->get_x_space());
364  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
365  solutionsDot.push_back(solutionDot);
366 
367  // Output finest temporal solution for plotting
368  // This only works for ONE MPI process
369  if ((n == 0) or (n == nTimeStepSizes-1)) {
370  std::string fname = "Tempus_ForwardEuler_VanDerPol-Ref.dat";
371  if (n == 0) fname = "Tempus_ForwardEuler_VanDerPol.dat";
372  RCP<const SolutionHistory<double> > solutionHistory =
373  integrator->getSolutionHistory();
374  writeSolution(fname, solutionHistory);
375  }
376  }
377 
378  // Check the order and intercept
379  double xSlope = 0.0;
380  double xDotSlope = 0.0;
381  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
382  double order = stepper->getOrder();
383  writeOrderError("Tempus_ForwardEuler_VanDerPol-Error.dat",
384  stepper, StepSize,
385  solutions, xErrorNorm, xSlope,
386  solutionsDot, xDotErrorNorm, xDotSlope);
387 
388  TEST_FLOATING_EQUALITY( xSlope, order, 0.10 );
389  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.387476, 1.0e-4 );
390  // xDot not yet available for Forward Euler.
391  //TEST_FLOATING_EQUALITY( xDotSlope, 1.74898, 0.10 );
392  //TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 1.0038, 1.0e-4 );
393 
394  Teuchos::TimeMonitor::summarize();
395 }
396 
397 
398 // ************************************************************
399 // ************************************************************
400 TEUCHOS_UNIT_TEST(ForwardEuler, NumberTimeSteps)
401 {
402 
403  std::vector<double> StepSize;
404  std::vector<double> ErrorNorm;
405  //const int nTimeStepSizes = 7;
406  //double dt = 0.2;
407  //double order = 0.0;
408 
409  // Read params from .xml file
410  RCP<ParameterList> pList =
411  getParametersFromXmlFile("Tempus_ForwardEuler_NumberOfTimeSteps.xml");
412 
413  // Setup the VanDerPolModel
414  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
415  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
416 
417  // Setup the Integrator and reset initial time step
418  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
419 
420  //dt = pl->sublist("Demo Integrator")
421  // .sublist("Time Step Control")
422  // .get<double>("Initial Time Step");
423  const int numTimeSteps = pl->sublist("Demo Integrator")
424  .sublist("Time Step Control")
425  .get<int>("Number of Time Steps");
426  const std::string integratorStepperType =
427  pl->sublist("Demo Integrator")
428  .sublist("Time Step Control")
429  .get<std::string>("Integrator Step Type");
430 
431  RCP<Tempus::IntegratorBasic<double> > integrator =
432  Tempus::integratorBasic<double>(pl, model);
433 
434  // Integrate to timeMax
435  bool integratorStatus = integrator->advanceTime();
436  TEST_ASSERT(integratorStatus)
437 
438  // check that the number of time steps taken is whats is set
439  // in the parameter list
440  TEST_EQUALITY(numTimeSteps, integrator->getIndex());
441 }
442 
443 
444 // ************************************************************
445 // ************************************************************
446 TEUCHOS_UNIT_TEST(ForwardEuler, Variable_TimeSteps)
447 {
448  // Read params from .xml file
449  RCP<ParameterList> pList =
450  getParametersFromXmlFile("Tempus_ForwardEuler_VanDerPol.xml");
451 
452  // Setup the VanDerPolModel
453  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
454  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
455 
456  // Setup the Integrator and reset initial time step
457  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
458 
459  // Set parameters for this test.
460  pl->sublist("Demo Integrator")
461  .sublist("Time Step Control").set("Initial Time Step", 0.01);
462 
463  pl->sublist("Demo Integrator")
464  .sublist("Time Step Control")
465  .sublist("Time Step Control Strategy")
466  .sublist("basic_vs").set("Reduction Factor", 0.9);
467  pl->sublist("Demo Integrator")
468  .sublist("Time Step Control")
469  .sublist("Time Step Control Strategy")
470  .sublist("basic_vs").set("Amplification Factor", 1.15);
471  pl->sublist("Demo Integrator")
472  .sublist("Time Step Control")
473  .sublist("Time Step Control Strategy")
474  .sublist("basic_vs").set("Minimum Value Monitoring Function", 0.05);
475  pl->sublist("Demo Integrator")
476  .sublist("Time Step Control")
477  .sublist("Time Step Control Strategy")
478  .sublist("basic_vs").set("Maximum Value Monitoring Function", 0.1);
479 
480  pl->sublist("Demo Integrator")
481  .sublist("Solution History").set("Storage Type", "Static");
482  pl->sublist("Demo Integrator")
483  .sublist("Solution History").set("Storage Limit", 3);
484 
485  RCP<Tempus::IntegratorBasic<double> > integrator =
486  Tempus::integratorBasic<double>(pl, model);
487 
488  // Integrate to timeMax
489  bool integratorStatus = integrator->advanceTime();
490  TEST_ASSERT(integratorStatus)
491 
492  // Check 'Final Time'
493  double time = integrator->getTime();
494  double timeFinal =pl->sublist("Demo Integrator")
495  .sublist("Time Step Control").get<double>("Final Time");
496  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
497 
498  // Check TimeStep size
499  auto state = integrator->getCurrentState();
500  double dt = state->getTimeStep();
501  TEST_FLOATING_EQUALITY(dt, 0.008310677297208358, 1.0e-12);
502 
503  // Check number of time steps taken
504  const int numTimeSteps = 60;
505  TEST_EQUALITY(numTimeSteps, integrator->getIndex());
506 
507  // Time-integrated solution and the reference solution
508  RCP<Thyra::VectorBase<double> > x = integrator->getX();
509  RCP<Thyra::VectorBase<double> > x_ref = x->clone_v();
510  {
511  Thyra::DetachedVectorView<double> x_ref_view( *x_ref );
512  x_ref_view[0] = -1.931946840284863;
513  x_ref_view[1] = 0.645346748303107;
514  }
515 
516  // Calculate the error
517  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
518  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_ref, -1.0, *(x));
519 
520  // Check the solution
521  std::cout << " Stepper = ForwardEuler" << std::endl;
522  std::cout << " =========================" << std::endl;
523  std::cout << " Reference solution: " << get_ele(*(x_ref), 0) << " "
524  << get_ele(*(x_ref), 1) << std::endl;
525  std::cout << " Computed solution : " << get_ele(*(x ), 0) << " "
526  << get_ele(*(x ), 1) << std::endl;
527  std::cout << " Difference : " << get_ele(*(xdiff), 0) << " "
528  << get_ele(*(xdiff), 1) << std::endl;
529  std::cout << " =========================" << std::endl;
530  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), get_ele(*(x_ref), 0), 1.0e-12);
531  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), get_ele(*(x_ref), 1), 1.0e-12);
532 }
533 
534 
535 } // namespace Tempus_Test
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 mechanicisms that effect the time step size...
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...