13 #include "Teuchos_DefaultComm.hpp"
15 #include "Tempus_config.hpp"
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperTrapezoidal.hpp"
19 #include "../TestModels/SinCosModel.hpp"
20 #include "../TestModels/VanDerPolModel.hpp"
21 #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_Trapezoidal_SinCos.xml");
58 Tempus::createIntegratorBasic<double>(tempusPL, model);
63 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>(
79 model, std::string(
"Trapezoidal Method"));
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(
"Default Parameters");
104 options.push_back(
"ICConsistency and Check");
106 for (
const auto& option : options) {
109 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
132 stepper->setModel(model);
133 if (option ==
"ICConsistency and Check") {
134 stepper->setICConsistency(
"Consistent");
135 stepper->setICConsistencyCheck(
true);
137 stepper->initialize();
143 timeStepControl->setInitIndex(tscPL.
get<
int>(
"Initial Time Index"));
144 timeStepControl->setInitTime(tscPL.
get<
double>(
"Initial Time"));
145 timeStepControl->setFinalTime(tscPL.
get<
double>(
"Final Time"));
146 timeStepControl->setInitTimeStep(dt);
147 timeStepControl->initialize();
150 auto inArgsIC = model->getNominalValues();
155 icState->setTime(timeStepControl->getInitTime());
156 icState->setIndex(timeStepControl->getInitIndex());
157 icState->setTimeStep(0.0);
158 icState->setOrder(stepper->getOrder());
163 solutionHistory->setName(
"Forward States");
165 solutionHistory->setStorageLimit(2);
166 solutionHistory->addState(icState);
170 Tempus::createIntegratorBasic<double>();
171 integrator->setStepper(stepper);
172 integrator->setTimeStepControl(timeStepControl);
173 integrator->setSolutionHistory(solutionHistory);
175 integrator->initialize();
178 bool integratorStatus = integrator->advanceTime();
182 double time = integrator->getTime();
183 double timeFinal = pl->
sublist(
"Default Integrator")
185 .
get<
double>(
"Final Time");
191 model->getExactSolution(time).get_x();
195 Thyra::V_StVpStV(xdiff.
ptr(), 1.0, *x_exact, -1.0, *(x));
198 out <<
" Stepper = " << stepper->description() <<
"\n with "
199 << option << std::endl;
200 out <<
" =========================" << std::endl;
201 out <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
202 << get_ele(*(x_exact), 1) << std::endl;
203 out <<
" Computed solution: " << get_ele(*(x), 0) <<
" "
204 << get_ele(*(x), 1) << std::endl;
205 out <<
" Difference : " << get_ele(*(xdiff), 0) <<
" "
206 << get_ele(*(xdiff), 1) << std::endl;
207 out <<
" =========================" << std::endl;
208 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841021, 1.0e-4);
209 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.541002, 1.0e-4);
218 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
219 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
220 std::vector<double> StepSize;
221 std::vector<double> xErrorNorm;
222 std::vector<double> xDotErrorNorm;
223 const int nTimeStepSizes = 7;
226 for (
int n = 0; n < nTimeStepSizes; n++) {
229 getParametersFromXmlFile(
"Tempus_Trapezoidal_SinCos.xml");
244 pl->
sublist(
"Default Integrator")
246 .
set(
"Initial Time Step", dt);
247 integrator = Tempus::createIntegratorBasic<double>(pl, model);
255 model->getNominalValues().get_x()->clone_v();
257 model->getNominalValues().get_x_dot()->clone_v();
258 integrator->initializeSolutionHistory(0.0, x0, xdot0);
262 bool integratorStatus = integrator->advanceTime();
266 time = integrator->getTime();
267 double timeFinal = pl->sublist(
"Default Integrator")
268 .sublist(
"Time Step Control")
269 .
get<
double>(
"Final Time");
275 integrator->getSolutionHistory();
276 writeSolution(
"Tempus_Trapezoidal_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_Trapezoidal_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);
314 double xDotSlope = 0.0;
316 double order = stepper->getOrder();
317 writeOrderError(
"Tempus_Trapezoidal_SinCos-Error.dat", stepper, StepSize,
318 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
334 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
335 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
336 std::vector<double> StepSize;
337 std::vector<double> xErrorNorm;
338 std::vector<double> xDotErrorNorm;
339 const int nTimeStepSizes = 4;
342 for (
int n = 0; n < nTimeStepSizes; n++) {
345 getParametersFromXmlFile(
"Tempus_Trapezoidal_VanDerPol.xml");
353 if (n == nTimeStepSizes - 1) dt /= 10.0;
359 .
set(
"Initial Time Step", dt);
360 integrator = Tempus::createIntegratorBasic<double>(pl, model);
363 bool integratorStatus = integrator->advanceTime();
367 time = integrator->getTime();
368 double timeFinal = pl->sublist(
"Demo Integrator")
369 .sublist(
"Time Step Control")
370 .
get<
double>(
"Final Time");
371 double tol = 100.0 * std::numeric_limits<double>::epsilon();
375 StepSize.push_back(dt);
376 auto solution = Thyra::createMember(model->get_x_space());
377 Thyra::copy(*(integrator->getX()), solution.ptr());
378 solutions.push_back(solution);
379 auto solutionDot = Thyra::createMember(model->get_x_space());
380 Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
381 solutionsDot.push_back(solutionDot);
385 if ((n == 0) || (n == nTimeStepSizes - 1)) {
386 std::string fname =
"Tempus_Trapezoidal_VanDerPol-Ref.dat";
387 if (n == 0) fname =
"Tempus_Trapezoidal_VanDerPol.dat";
389 integrator->getSolutionHistory();
395 double xDotSlope = 0.0;
397 double order = stepper->getOrder();
398 writeOrderError(
"Tempus_Trapezoidal_VanDerPol-Error.dat", stepper, StepSize,
399 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
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.
Trapezoidal method time stepper.
T & get(const std::string &name, T def_value)
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="")
Solution state for integrators and steppers.