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