14 #include "Thyra_VectorStdOps.hpp"
15 #include "Thyra_DetachedVectorView.hpp"
17 #include "Tempus_IntegratorBasic.hpp"
19 #include "Tempus_StepperForwardEuler.hpp"
21 #include "../TestModels/SinCosModel.hpp"
22 #include "../TestModels/VanDerPolModel.hpp"
23 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
28 namespace Tempus_Test {
30 using Teuchos::getParametersFromXmlFile;
34 using Teuchos::rcp_const_cast;
35 using Teuchos::sublist;
47 getParametersFromXmlFile(
"Tempus_ForwardEuler_SinCos.xml");
58 Tempus::createIntegratorBasic<double>(tempusPL, model);
62 integrator->getStepper()->getValidParameters();
64 bool pass = haveSameValuesSorted(*stepperPL, *defaultPL,
true);
67 out <<
"stepperPL -------------- \n"
68 << *stepperPL << std::endl;
69 out <<
"defaultPL -------------- \n"
70 << *defaultPL << std::endl;
78 Tempus::createIntegratorBasic<double>(model,
79 std::string(
"Forward Euler"));
83 integrator->getStepper()->getValidParameters();
85 bool pass = haveSameValuesSorted(*stepperPL, *defaultPL,
true);
88 out <<
"stepperPL -------------- \n"
89 << *stepperPL << std::endl;
90 out <<
"defaultPL -------------- \n"
91 << *defaultPL << std::endl;
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");
107 for (
const auto& option : options) {
110 getParametersFromXmlFile(
"Tempus_ForwardEuler_SinCos.xml");
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);
129 stepper->initialize();
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();
142 auto inArgsIC = model()->getNominalValues();
146 icState->setTime(timeStepControl->getInitTime());
147 icState->setIndex(timeStepControl->getInitIndex());
148 icState->setTimeStep(0.0);
153 solutionHistory->setName(
"Forward States");
155 solutionHistory->setStorageLimit(2);
156 solutionHistory->addState(icState);
159 stepper->setInitialConditions(solutionHistory);
163 Tempus::createIntegratorBasic<double>();
164 integrator->setStepper(stepper);
165 integrator->setTimeStepControl(timeStepControl);
166 integrator->setSolutionHistory(solutionHistory);
168 integrator->initialize();
171 bool integratorStatus = integrator->advanceTime();
175 double time = integrator->getTime();
176 double timeFinal = pl->
sublist(
"Demo Integrator")
178 .
get<
double>(
"Final Time");
184 model->getExactSolution(time).get_x();
188 Thyra::V_StVpStV(xdiff.
ptr(), 1.0, *x_exact, -1.0, *(x));
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);
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;
219 for (
int n = 0; n < nTimeStepSizes; n++) {
222 getParametersFromXmlFile(
"Tempus_ForwardEuler_SinCos.xml");
239 .
set(
"Initial Time Step", dt);
240 integrator = Tempus::createIntegratorBasic<double>(pl, model);
247 model->getNominalValues().get_x()->clone_v();
248 integrator->initializeSolutionHistory(0.0, x0);
249 integrator->initialize();
252 bool integratorStatus = integrator->advanceTime();
257 integrator->getSolutionHistory()->getCurrentState()->getPhysicsState();
258 TEST_EQUALITY(physicsState->getName(),
"Tempus::PhysicsState");
261 time = integrator->getTime();
262 double timeFinal = pl->sublist(
"Demo Integrator")
263 .sublist(
"Time Step Control")
264 .
get<
double>(
"Final Time");
270 model->getExactSolution(time).get_x();
275 integrator->getSolutionHistory();
276 writeSolution(
"Tempus_ForwardEuler_SinCos.dat", solutionHistory);
279 for (
int i = 0; i < solutionHistory->getNumStates(); i++) {
280 double time_i = (*solutionHistory)[i]->getTime();
283 model->getExactSolution(time_i).get_x()),
285 model->getExactSolution(time_i).get_x_dot()));
286 state->setTime((*solutionHistory)[i]->getTime());
287 solnHistExact->addState(state);
289 writeSolution(
"Tempus_ForwardEuler_SinCos-Ref.dat", solnHistExact);
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) {
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);
315 double xDotSlope = 0.0;
317 double order = stepper->getOrder();
318 writeOrderError(
"Tempus_ForwardEuler_SinCos-Error.dat", stepper, StepSize,
319 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
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;
343 for (
int n = 0; n < nTimeStepSizes; n++) {
346 getParametersFromXmlFile(
"Tempus_ForwardEuler_VanDerPol.xml");
354 if (n == nTimeStepSizes - 1) dt /= 10.0;
360 .
set(
"Initial Time Step", dt);
361 integrator = Tempus::createIntegratorBasic<double>(pl, model);
364 bool integratorStatus = integrator->advanceTime();
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();
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);
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";
390 integrator->getSolutionHistory();
397 double xDotSlope = 0.0;
399 double order = stepper->getOrder();
400 writeOrderError(
"Tempus_ForwardEuler_VanDerPol-Error.dat", stepper, StepSize,
401 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
417 std::vector<double> StepSize;
418 std::vector<double> ErrorNorm;
425 getParametersFromXmlFile(
"Tempus_ForwardEuler_NumberOfTimeSteps.xml");
437 const int numTimeSteps = pl->
sublist(
"Demo Integrator")
439 .
get<
int>(
"Number of Time Steps");
442 Tempus::createIntegratorBasic<double>(pl, model);
445 bool integratorStatus = integrator->advanceTime();
459 getParametersFromXmlFile(
"Tempus_ForwardEuler_VanDerPol.xml");
471 .
set(
"Initial Time Step", 0.01);
475 .
sublist(
"Time Step Control Strategy")
476 .
set(
"Reduction Factor", 0.9);
479 .
sublist(
"Time Step Control Strategy")
480 .
set(
"Amplification Factor", 1.15);
483 .
sublist(
"Time Step Control Strategy")
484 .
set(
"Minimum Value Monitoring Function", 0.05);
487 .
sublist(
"Time Step Control Strategy")
488 .
set(
"Maximum Value Monitoring Function", 0.1);
492 .
set(
"Storage Type",
"Static");
495 .
set(
"Storage Limit", 3);
498 Tempus::createIntegratorBasic<double>(pl, model);
501 bool integratorStatus = integrator->advanceTime();
505 double time = integrator->getTime();
506 double timeFinal = pl->sublist(
"Demo Integrator")
507 .sublist(
"Time Step Control")
508 .
get<
double>(
"Final Time");
512 auto state = integrator->getCurrentState();
513 double dt = state->getTimeStep();
514 TEST_FLOATING_EQUALITY(dt, 0.008310677297208358, 1.0e-12);
517 const int numTimeSteps = 60;
525 x_ref_view[0] = -1.931946840284863;
526 x_ref_view[1] = 0.645346748303107;
531 Thyra::V_StVpStV(xdiff.
ptr(), 1.0, *x_ref, -1.0, *(x));
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);
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)
Forward Euler time stepper.
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)
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)
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.