9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
12 #include "Teuchos_DefaultComm.hpp"
14 #include "Tempus_config.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
16 #include "Tempus_StepperTrapezoidal.hpp"
18 #include "../TestModels/SinCosModel.hpp"
19 #include "../TestModels/VanDerPolModel.hpp"
20 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
22 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
23 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
25 #ifdef Tempus_ENABLE_MPI
26 #include "Epetra_MpiComm.h"
28 #include "Epetra_SerialComm.h"
40 using Teuchos::rcp_const_cast;
41 using Teuchos::ParameterList;
42 using Teuchos::sublist;
43 using Teuchos::getParametersFromXmlFile;
55 RCP<ParameterList> pList =
56 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
59 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
62 RCP<ParameterList> tempusPL = sublist(pList,
"Tempus",
true);
66 RCP<Tempus::IntegratorBasic<double> > integrator =
67 Tempus::integratorBasic<double>(tempusPL, model);
69 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
71 RCP<const ParameterList> defaultPL =
72 integrator->getStepper()->getValidParameters();
73 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
75 std::cout << std::endl;
76 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
77 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
84 RCP<Tempus::IntegratorBasic<double> > integrator =
85 Tempus::integratorBasic<double>(model,
"Trapezoidal Method");
87 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
88 RCP<const ParameterList> defaultPL =
89 integrator->getStepper()->getValidParameters();
91 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
93 std::cout << std::endl;
94 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
95 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
109 RCP<ParameterList> pList =
110 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
111 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
114 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
133 stepper->setModel(model);
134 stepper->initialize();
138 ParameterList tscPL = pl->sublist(
"Default Integrator")
139 .sublist(
"Time Step Control");
140 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
141 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
142 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
143 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
144 timeStepControl->setInitTimeStep(dt);
145 timeStepControl->initialize();
148 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
149 stepper->getModel()->getNominalValues();
150 auto icSoln = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
152 rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x_dot());
154 icState->setTime (timeStepControl->getInitTime());
155 icState->setIndex (timeStepControl->getInitIndex());
156 icState->setTimeStep(0.0);
157 icState->setOrder (stepper->getOrder());
162 solutionHistory->setName(
"Forward States");
164 solutionHistory->setStorageLimit(2);
165 solutionHistory->addState(icState);
168 RCP<Tempus::IntegratorBasic<double> > integrator =
169 Tempus::integratorBasic<double>();
170 integrator->setStepperWStepper(stepper);
171 integrator->setTimeStepControl(timeStepControl);
172 integrator->setSolutionHistory(solutionHistory);
174 integrator->initialize();
178 bool integratorStatus = integrator->advanceTime();
179 TEST_ASSERT(integratorStatus)
183 double time = integrator->getTime();
184 double timeFinal =pl->sublist(
"Default Integrator")
185 .sublist(
"Time Step Control").get<
double>(
"Final Time");
186 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
189 RCP<Thyra::VectorBase<double> > x = integrator->getX();
190 RCP<const Thyra::VectorBase<double> > x_exact =
191 model->getExactSolution(time).get_x();
194 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
195 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
198 std::cout <<
" Stepper = Trapezoidal" << std::endl;
199 std::cout <<
" =========================" << std::endl;
200 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
201 << get_ele(*(x_exact), 1) << std::endl;
202 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
203 << get_ele(*(x ), 1) << std::endl;
204 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" "
205 << get_ele(*(xdiff ), 1) << std::endl;
206 std::cout <<
" =========================" << std::endl;
207 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841021, 1.0e-4 );
208 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.541002, 1.0e-4 );
216 RCP<Tempus::IntegratorBasic<double> > integrator;
217 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
218 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
219 std::vector<double> StepSize;
220 std::vector<double> xErrorNorm;
221 std::vector<double> xDotErrorNorm;
222 const int nTimeStepSizes = 7;
225 for (
int n=0; n<nTimeStepSizes; n++) {
228 RCP<ParameterList> pList =
229 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
236 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
243 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
244 pl->sublist(
"Default Integrator")
245 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
246 integrator = Tempus::integratorBasic<double>(pl, model);
253 RCP<Thyra::VectorBase<double> > x0 =
254 model->getNominalValues().get_x()->clone_v();
255 RCP<Thyra::VectorBase<double> > xdot0 =
256 model->getNominalValues().get_x_dot()->clone_v();
257 integrator->initializeSolutionHistory(0.0, x0, xdot0);
261 bool integratorStatus = integrator->advanceTime();
262 TEST_ASSERT(integratorStatus)
265 time = integrator->getTime();
266 double timeFinal =pl->sublist(
"Default Integrator")
267 .sublist(
"Time Step Control").get<
double>(
"Final Time");
268 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
272 RCP<const SolutionHistory<double> > solutionHistory =
273 integrator->getSolutionHistory();
274 writeSolution(
"Tempus_Trapezoidal_SinCos.dat", solutionHistory);
277 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
278 double time_i = (*solutionHistory)[i]->getTime();
280 rcp_const_cast<Thyra::VectorBase<double> > (
281 model->getExactSolution(time_i).get_x()),
282 rcp_const_cast<Thyra::VectorBase<double> > (
283 model->getExactSolution(time_i).get_x_dot()));
284 state->setTime((*solutionHistory)[i]->getTime());
285 solnHistExact->addState(state);
287 writeSolution(
"Tempus_Trapezoidal_SinCos-Ref.dat", solnHistExact);
291 StepSize.push_back(dt);
292 auto solution = Thyra::createMember(model->get_x_space());
293 Thyra::copy(*(integrator->getX()),solution.ptr());
294 solutions.push_back(solution);
295 auto solutionDot = Thyra::createMember(model->get_x_space());
296 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
297 solutionsDot.push_back(solutionDot);
298 if (n == nTimeStepSizes-1) {
299 StepSize.push_back(0.0);
300 auto solutionExact = Thyra::createMember(model->get_x_space());
301 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
302 solutions.push_back(solutionExact);
303 auto solutionDotExact = Thyra::createMember(model->get_x_space());
304 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
305 solutionDotExact.ptr());
306 solutionsDot.push_back(solutionDotExact);
311 double xDotSlope = 0.0;
312 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
313 double order = stepper->getOrder();
316 solutions, xErrorNorm, xSlope,
317 solutionsDot, xDotErrorNorm, xDotSlope);
319 TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
320 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.000832086, 1.0e-4 );
321 TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
322 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.000832086, 1.0e-4 );
324 Teuchos::TimeMonitor::summarize();
332 RCP<Tempus::IntegratorBasic<double> > integrator;
333 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
334 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
335 std::vector<double> StepSize;
336 std::vector<double> xErrorNorm;
337 std::vector<double> xDotErrorNorm;
338 const int nTimeStepSizes = 4;
341 for (
int n=0; n<nTimeStepSizes; n++) {
344 RCP<ParameterList> pList =
345 getParametersFromXmlFile(
"Tempus_Trapezoidal_VanDerPol.xml");
348 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
353 if (n == nTimeStepSizes-1) dt /= 10.0;
356 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
357 pl->sublist(
"Demo Integrator")
358 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
359 integrator = Tempus::integratorBasic<double>(pl, model);
362 bool integratorStatus = integrator->advanceTime();
363 TEST_ASSERT(integratorStatus)
366 time = integrator->getTime();
367 double timeFinal =pl->sublist(
"Demo Integrator")
368 .sublist(
"Time Step Control").get<
double>(
"Final Time");
369 double tol = 100.0 * std::numeric_limits<double>::epsilon();
370 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
373 StepSize.push_back(dt);
374 auto solution = Thyra::createMember(model->get_x_space());
375 Thyra::copy(*(integrator->getX()),solution.ptr());
376 solutions.push_back(solution);
377 auto solutionDot = Thyra::createMember(model->get_x_space());
378 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
379 solutionsDot.push_back(solutionDot);
383 if ((n == 0) or (n == nTimeStepSizes-1)) {
384 std::string fname =
"Tempus_Trapezoidal_VanDerPol-Ref.dat";
385 if (n == 0) fname =
"Tempus_Trapezoidal_VanDerPol.dat";
386 RCP<const SolutionHistory<double> > solutionHistory =
387 integrator->getSolutionHistory();
393 double xDotSlope = 0.0;
394 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
395 double order = stepper->getOrder();
398 solutions, xErrorNorm, xSlope,
399 solutionsDot, xDotErrorNorm, xDotSlope);
401 TEST_FLOATING_EQUALITY( xSlope, order, 0.10 );
402 TEST_FLOATING_EQUALITY( xDotSlope, order, 0.10 );
403 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.00085855, 1.0e-4 );
404 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.00120695, 1.0e-4 );
406 Teuchos::TimeMonitor::summarize();
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
Trapezoidal method time stepper.
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
van der Pol model problem for nonlinear electrical circuit.
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)
@ 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.