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"
36 namespace Tempus_Test {
40 using Teuchos::rcp_const_cast;
41 using Teuchos::ParameterList;
42 using Teuchos::sublist;
43 using Teuchos::getParametersFromXmlFile;
50 #define TEST_PARAMETERLIST
51 #define TEST_CONSTRUCTING_FROM_DEFAULTS
53 #define TEST_VANDERPOL
56 #ifdef TEST_PARAMETERLIST
62 RCP<ParameterList> pList =
63 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
66 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
69 RCP<ParameterList> tempusPL = sublist(pList,
"Tempus",
true);
73 RCP<Tempus::IntegratorBasic<double> > integrator =
74 Tempus::integratorBasic<double>(tempusPL, model);
76 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
78 RCP<const ParameterList> defaultPL =
79 integrator->getStepper()->getValidParameters();
80 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
82 std::cout << std::endl;
83 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
84 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
91 RCP<Tempus::IntegratorBasic<double> > integrator =
92 Tempus::integratorBasic<double>(model,
"Trapezoidal Method");
94 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
95 RCP<const ParameterList> defaultPL =
96 integrator->getStepper()->getValidParameters();
98 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
100 std::cout << std::endl;
101 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
102 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
107 #endif // TEST_PARAMETERLIST
110 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
118 RCP<ParameterList> pList =
119 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
120 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
123 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
142 stepper->setModel(model);
143 stepper->setSolver();
144 stepper->initialize();
148 ParameterList tscPL = pl->sublist(
"Default Integrator")
149 .sublist(
"Time Step Control");
150 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
151 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
152 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
153 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
154 timeStepControl->setInitTimeStep(dt);
155 timeStepControl->initialize();
158 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
159 stepper->getModel()->getNominalValues();
160 auto icSoln = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
162 rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x_dot());
164 icState->setTime (timeStepControl->getInitTime());
165 icState->setIndex (timeStepControl->getInitIndex());
166 icState->setTimeStep(0.0);
167 icState->setOrder (stepper->getOrder());
178 RCP<Tempus::IntegratorBasic<double> > integrator =
179 Tempus::integratorBasic<double>();
180 integrator->setStepperWStepper(stepper);
181 integrator->setTimeStepControl(timeStepControl);
184 integrator->initialize();
188 bool integratorStatus = integrator->advanceTime();
189 TEST_ASSERT(integratorStatus)
193 double time = integrator->getTime();
194 double timeFinal =pl->sublist(
"Default Integrator")
195 .sublist(
"Time Step Control").get<
double>(
"Final Time");
196 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
199 RCP<Thyra::VectorBase<double> > x = integrator->getX();
200 RCP<const Thyra::VectorBase<double> > x_exact =
201 model->getExactSolution(time).get_x();
204 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
205 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
208 std::cout <<
" Stepper = Trapezoidal" << std::endl;
209 std::cout <<
" =========================" << std::endl;
210 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
211 << get_ele(*(x_exact), 1) << std::endl;
212 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
213 << get_ele(*(x ), 1) << std::endl;
214 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" "
215 << get_ele(*(xdiff ), 1) << std::endl;
216 std::cout <<
" =========================" << std::endl;
217 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841021, 1.0e-4 );
218 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.541002, 1.0e-4 );
220 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
228 RCP<Tempus::IntegratorBasic<double> > integrator;
229 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
230 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
231 std::vector<double> StepSize;
232 std::vector<double> xErrorNorm;
233 std::vector<double> xDotErrorNorm;
234 const int nTimeStepSizes = 7;
237 for (
int n=0; n<nTimeStepSizes; n++) {
240 RCP<ParameterList> pList =
241 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
248 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
255 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
256 pl->sublist(
"Default Integrator")
257 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
258 integrator = Tempus::integratorBasic<double>(pl, model);
265 RCP<Thyra::VectorBase<double> > x0 =
266 model->getNominalValues().get_x()->clone_v();
267 RCP<Thyra::VectorBase<double> > xdot0 =
268 model->getNominalValues().get_x_dot()->clone_v();
269 integrator->initializeSolutionHistory(0.0, x0, xdot0);
273 bool integratorStatus = integrator->advanceTime();
274 TEST_ASSERT(integratorStatus)
277 time = integrator->getTime();
278 double timeFinal =pl->sublist(
"Default Integrator")
279 .sublist(
"Time Step Control").get<
double>(
"Final Time");
280 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
285 integrator->getSolutionHistory();
286 writeSolution(
"Tempus_Trapezoidal_SinCos.dat", solutionHistory);
289 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
290 double time_i = (*solutionHistory)[i]->getTime();
292 model->getExactSolution(time_i).get_x(),
293 model->getExactSolution(time_i).get_x_dot()));
294 state->setTime((*solutionHistory)[i]->getTime());
295 solnHistExact->addState(state);
297 writeSolution(
"Tempus_Trapezoidal_SinCos-Ref.dat", solnHistExact);
301 StepSize.push_back(dt);
302 auto solution = Thyra::createMember(model->get_x_space());
303 Thyra::copy(*(integrator->getX()),solution.ptr());
304 solutions.push_back(solution);
305 auto solutionDot = Thyra::createMember(model->get_x_space());
306 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
307 solutionsDot.push_back(solutionDot);
308 if (n == nTimeStepSizes-1) {
309 StepSize.push_back(0.0);
310 auto solutionExact = Thyra::createMember(model->get_x_space());
311 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
312 solutions.push_back(solutionExact);
313 auto solutionDotExact = Thyra::createMember(model->get_x_space());
314 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
315 solutionDotExact.ptr());
316 solutionsDot.push_back(solutionDotExact);
321 double xDotSlope = 0.0;
322 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
323 double order = stepper->getOrder();
326 solutions, xErrorNorm, xSlope,
327 solutionsDot, xDotErrorNorm, xDotSlope);
329 TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
330 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.000832086, 1.0e-4 );
331 TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
332 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.000832086, 1.0e-4 );
334 Teuchos::TimeMonitor::summarize();
336 #endif // TEST_SINCOS
339 #ifdef TEST_VANDERPOL
344 RCP<Tempus::IntegratorBasic<double> > integrator;
345 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
346 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
347 std::vector<double> StepSize;
348 std::vector<double> xErrorNorm;
349 std::vector<double> xDotErrorNorm;
350 const int nTimeStepSizes = 4;
353 for (
int n=0; n<nTimeStepSizes; n++) {
356 RCP<ParameterList> pList =
357 getParametersFromXmlFile(
"Tempus_Trapezoidal_VanDerPol.xml");
360 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
365 if (n == nTimeStepSizes-1) dt /= 10.0;
368 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
369 pl->sublist(
"Demo Integrator")
370 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
371 integrator = Tempus::integratorBasic<double>(pl, model);
374 bool integratorStatus = integrator->advanceTime();
375 TEST_ASSERT(integratorStatus)
378 time = integrator->getTime();
379 double timeFinal =pl->sublist(
"Demo Integrator")
380 .sublist(
"Time Step Control").get<
double>(
"Final Time");
381 double tol = 100.0 * std::numeric_limits<double>::epsilon();
382 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
385 StepSize.push_back(dt);
386 auto solution = Thyra::createMember(model->get_x_space());
387 Thyra::copy(*(integrator->getX()),solution.ptr());
388 solutions.push_back(solution);
389 auto solutionDot = Thyra::createMember(model->get_x_space());
390 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
391 solutionsDot.push_back(solutionDot);
395 if ((n == 0) or (n == nTimeStepSizes-1)) {
396 std::string fname =
"Tempus_Trapezoidal_VanDerPol-Ref.dat";
397 if (n == 0) fname =
"Tempus_Trapezoidal_VanDerPol.dat";
399 integrator->getSolutionHistory();
405 double xDotSlope = 0.0;
406 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
407 double order = stepper->getOrder();
410 solutions, xErrorNorm, xSlope,
411 solutionsDot, xDotErrorNorm, xDotSlope);
413 TEST_FLOATING_EQUALITY( xSlope, order, 0.10 );
414 TEST_FLOATING_EQUALITY( xDotSlope, order, 0.10 );
415 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.00085855, 1.0e-4 );
416 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.00120695, 1.0e-4 );
418 Teuchos::TimeMonitor::summarize();
420 #endif // TEST_VANDERPOL
Trapezoidal method 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...