9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
13 #include "Thyra_VectorStdOps.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperForwardEuler.hpp"
19 #include "../TestModels/SinCosModel.hpp"
20 #include "../TestModels/VanDerPolModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
26 namespace Tempus_Test {
30 using Teuchos::rcp_const_cast;
31 using Teuchos::ParameterList;
32 using Teuchos::sublist;
33 using Teuchos::getParametersFromXmlFile;
41 #define TEST_PARAMETERLIST
42 #define TEST_CONSTRUCTING_FROM_DEFAULTS
44 #define TEST_VANDERPOL
45 #define TEST_NUMBER_TIMESTEPS
48 #ifdef TEST_PARAMETERLIST
54 RCP<ParameterList> pList =
55 getParametersFromXmlFile(
"Tempus_ForwardEuler_SinCos.xml");
62 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
65 RCP<ParameterList> tempusPL = sublist(pList,
"Tempus",
true);
69 RCP<Tempus::IntegratorBasic<double> > integrator =
70 Tempus::integratorBasic<double>(tempusPL, model);
72 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Demo Stepper",
true);
73 RCP<const ParameterList> defaultPL =
74 integrator->getStepper()->getValidParameters();
76 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
78 std::cout << std::endl;
79 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
80 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
87 RCP<Tempus::IntegratorBasic<double> > integrator =
88 Tempus::integratorBasic<double>(model,
"Forward Euler");
90 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Demo Stepper",
true);
91 RCP<const ParameterList> defaultPL =
92 integrator->getStepper()->getValidParameters();
94 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
96 std::cout << std::endl;
97 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
98 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
103 #endif // TEST_PARAMETERLIST
106 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
114 RCP<ParameterList> pList =
115 getParametersFromXmlFile(
"Tempus_ForwardEuler_SinCos.xml");
116 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
119 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
125 stepper->setModel(model);
126 stepper->initialize();
130 ParameterList tscPL = pl->sublist(
"Demo Integrator")
131 .sublist(
"Time Step Control");
132 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
133 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
134 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
135 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
136 timeStepControl->setInitTimeStep(dt);
137 timeStepControl->initialize();
140 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
141 stepper->getModel()->getNominalValues();
142 auto icSolution = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
144 icState->setTime (timeStepControl->getInitTime());
145 icState->setIndex (timeStepControl->getInitIndex());
146 icState->setTimeStep(0.0);
157 RCP<Tempus::IntegratorBasic<double> > integrator =
158 Tempus::integratorBasic<double>();
159 integrator->setStepperWStepper(stepper);
160 integrator->setTimeStepControl(timeStepControl);
163 integrator->initialize();
167 bool integratorStatus = integrator->advanceTime();
168 TEST_ASSERT(integratorStatus)
172 double time = integrator->getTime();
173 double timeFinal =pl->sublist(
"Demo Integrator")
174 .sublist(
"Time Step Control").get<
double>(
"Final Time");
175 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
178 RCP<Thyra::VectorBase<double> > x = integrator->getX();
179 RCP<const Thyra::VectorBase<double> > x_exact =
180 model->getExactSolution(time).get_x();
183 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
184 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
187 std::cout <<
" Stepper = ForwardEuler" << std::endl;
188 std::cout <<
" =========================" << std::endl;
189 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
190 << get_ele(*(x_exact), 1) << std::endl;
191 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
192 << get_ele(*(x ), 1) << std::endl;
193 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" "
194 << get_ele(*(xdiff ), 1) << std::endl;
195 std::cout <<
" =========================" << 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 );
199 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
207 RCP<Tempus::IntegratorBasic<double> > integrator;
208 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
209 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
210 std::vector<double> StepSize;
211 std::vector<double> xErrorNorm;
212 std::vector<double> xDotErrorNorm;
213 const int nTimeStepSizes = 7;
216 for (
int n=0; n<nTimeStepSizes; n++) {
219 RCP<ParameterList> pList =
220 getParametersFromXmlFile(
"Tempus_ForwardEuler_SinCos.xml");
227 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
234 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
235 pl->sublist(
"Demo Integrator")
236 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
237 integrator = Tempus::integratorBasic<double>(pl, model);
243 RCP<Thyra::VectorBase<double> > x0 =
244 model->getNominalValues().get_x()->clone_v();
245 integrator->initializeSolutionHistory(0.0, x0);
248 bool integratorStatus = integrator->advanceTime();
249 TEST_ASSERT(integratorStatus)
252 RCP<Tempus::PhysicsState<double> > physicsState =
253 integrator->getSolutionHistory()->getCurrentState()->getPhysicsState();
254 TEST_EQUALITY(physicsState->getName(),
"Tempus::PhysicsState");
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);
263 RCP<Thyra::VectorBase<double> > x = integrator->getX();
264 RCP<const Thyra::VectorBase<double> > x_exact =
265 model->getExactSolution(time).get_x();
270 integrator->getSolutionHistory();
271 writeSolution(
"Tempus_ForwardEuler_SinCos.dat", solutionHistory);
274 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
275 double time_i = (*solutionHistory)[i]->getTime();
277 model->getExactSolution(time_i).get_x(),
278 model->getExactSolution(time_i).get_x_dot()));
279 state->setTime((*solutionHistory)[i]->getTime());
280 solnHistExact->addState(state);
282 writeSolution(
"Tempus_ForwardEuler_SinCos-Ref.dat", solnHistExact);
286 StepSize.push_back(dt);
287 auto solution = Thyra::createMember(model->get_x_space());
288 Thyra::copy(*(integrator->getX()),solution.ptr());
289 solutions.push_back(solution);
290 auto solutionDot = Thyra::createMember(model->get_x_space());
291 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
292 solutionsDot.push_back(solutionDot);
293 if (n == nTimeStepSizes-1) {
294 StepSize.push_back(0.0);
295 auto solutionExact = Thyra::createMember(model->get_x_space());
296 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
297 solutions.push_back(solutionExact);
298 auto solutionDotExact = Thyra::createMember(model->get_x_space());
299 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
300 solutionDotExact.ptr());
301 solutionsDot.push_back(solutionDotExact);
307 double xDotSlope = 0.0;
308 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
309 double order = stepper->getOrder();
312 solutions, xErrorNorm, xSlope,
313 solutionsDot, xDotErrorNorm, xDotSlope);
315 TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
316 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.051123, 1.0e-4 );
321 Teuchos::TimeMonitor::summarize();
323 #endif // TEST_SINCOS
326 #ifdef TEST_VANDERPOL
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;
339 for (
int n=0; n<nTimeStepSizes; n++) {
342 RCP<ParameterList> pList =
343 getParametersFromXmlFile(
"Tempus_ForwardEuler_VanDerPol.xml");
346 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
351 if (n == nTimeStepSizes-1) dt /= 10.0;
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::integratorBasic<double>(pl, model);
360 bool integratorStatus = integrator->advanceTime();
361 TEST_ASSERT(integratorStatus)
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);
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);
381 if ((n == 0) or (n == nTimeStepSizes-1)) {
382 std::string fname =
"Tempus_ForwardEuler_VanDerPol-Ref.dat";
383 if (n == 0) fname =
"Tempus_ForwardEuler_VanDerPol.dat";
385 integrator->getSolutionHistory();
392 double xDotSlope = 0.0;
393 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
394 double order = stepper->getOrder();
397 solutions, xErrorNorm, xSlope,
398 solutionsDot, xDotErrorNorm, xDotSlope);
400 TEST_FLOATING_EQUALITY( xSlope, order, 0.10 );
401 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.387476, 1.0e-4 );
406 Teuchos::TimeMonitor::summarize();
408 #endif // TEST_VANDERPOL
411 #ifdef TEST_NUMBER_TIMESTEPS
417 std::vector<double> StepSize;
418 std::vector<double> ErrorNorm;
424 RCP<ParameterList> pList =
425 getParametersFromXmlFile(
"Tempus_ForwardEuler_NumberOfTimeSteps.xml");
428 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
432 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
437 const int numTimeSteps = pl->sublist(
"Demo Integrator")
438 .sublist(
"Time Step Control")
439 .get<
int>(
"Number of Time Steps");
440 const std::string integratorStepperType =
441 pl->sublist(
"Demo Integrator")
442 .sublist(
"Time Step Control")
443 .get<std::string>(
"Integrator Step Type");
445 RCP<Tempus::IntegratorBasic<double> > integrator =
446 Tempus::integratorBasic<double>(pl, model);
449 bool integratorStatus = integrator->advanceTime();
450 TEST_ASSERT(integratorStatus)
454 TEST_EQUALITY(numTimeSteps, integrator->getIndex());
456 #endif // TEST_NUMBER_TIMESTEPS
Forward Euler time stepper.
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...