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