12 #include "Teuchos_DefaultComm.hpp"
14 #include "Tempus_config.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
16 #include "Tempus_StepperBackwardEuler.hpp"
18 #include "../TestModels/SinCosModel.hpp"
19 #include "../TestModels/CDR_Model.hpp"
20 #include "../TestModels/VanDerPolModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
23 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
24 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
26 #ifdef Tempus_ENABLE_MPI
27 #include "Epetra_MpiComm.h"
29 #include "Epetra_SerialComm.h"
37 namespace Tempus_Test {
41 using Teuchos::rcp_const_cast;
43 using Teuchos::sublist;
44 using Teuchos::getParametersFromXmlFile;
57 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
68 Tempus::createIntegratorBasic<double>(tempusPL, model);
72 integrator->getStepper()->getValidParameters();
74 bool pass = haveSameValuesSorted(*stepperPL, *defaultPL,
true);
77 out <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
78 out <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
86 Tempus::createIntegratorBasic<double>(model, std::string(
"Backward Euler"));
90 stepperPL->set(
"Predictor Stepper Type",
"None");
92 integrator->getStepper()->getValidParameters();
94 bool pass = haveSameValuesSorted(*stepperPL, *defaultPL,
true);
97 out <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
98 out <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
110 std::vector<std::string> options;
111 options.push_back(
"Default Parameters");
112 options.push_back(
"ICConsistency and Check");
114 for(
const auto& option: options) {
118 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
128 stepper->setModel(model);
129 if ( option ==
"ICConsistency and Check") {
130 stepper->setICConsistency(
"Consistent");
131 stepper->setICConsistencyCheck(
true);
133 stepper->initialize();
139 timeStepControl->setInitIndex(tscPL.
get<
int> (
"Initial Time Index"));
140 timeStepControl->setInitTime (tscPL.
get<
double>(
"Initial Time"));
141 timeStepControl->setFinalTime(tscPL.
get<
double>(
"Final Time"));
142 timeStepControl->setInitTimeStep(dt);
143 timeStepControl->initialize();
146 auto inArgsIC = model->getNominalValues();
150 icState->setTime (timeStepControl->getInitTime());
151 icState->setIndex (timeStepControl->getInitIndex());
152 icState->setTimeStep(0.0);
153 icState->setOrder (stepper->getOrder());
158 solutionHistory->setName(
"Forward States");
160 solutionHistory->setStorageLimit(2);
161 solutionHistory->addState(icState);
164 stepper->setInitialConditions(solutionHistory);
168 Tempus::createIntegratorBasic<double>();
169 integrator->setStepper(stepper);
170 integrator->setTimeStepControl(timeStepControl);
171 integrator->setSolutionHistory(solutionHistory);
173 integrator->initialize();
177 bool integratorStatus = integrator->advanceTime();
182 double time = integrator->getTime();
183 double timeFinal =pl->
sublist(
"Default Integrator")
184 .
sublist(
"Time Step Control").
get<
double>(
"Final Time");
190 model->getExactSolution(time).get_x();
194 Thyra::V_StVpStV(xdiff.
ptr(), 1.0, *x_exact, -1.0, *(x));
197 out <<
" Stepper = " << stepper->description()
198 <<
" with " << option << std::endl;
199 out <<
" =========================" << std::endl;
200 out <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
201 << get_ele(*(x_exact), 1) << std::endl;
202 out <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
203 << get_ele(*(x ), 1) << std::endl;
204 out <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" "
205 << get_ele(*(xdiff ), 1) << std::endl;
206 out <<
" =========================" << std::endl;
207 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.798923, 1.0e-4 );
208 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.516729, 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++) {
230 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
245 pl->
sublist(
"Default Integrator")
246 .
sublist(
"Time Step Control").
set(
"Initial Time Step", dt);
247 integrator = Tempus::createIntegratorBasic<double>(pl, model);
254 model->getNominalValues().get_x()->clone_v();
255 integrator->initializeSolutionHistory(0.0, x0);
258 bool integratorStatus = integrator->advanceTime();
262 time = integrator->getTime();
263 double timeFinal =pl->sublist(
"Default Integrator")
264 .sublist(
"Time Step Control").
get<
double>(
"Final Time");
270 integrator->getSolutionHistory();
271 writeSolution(
"Tempus_BackwardEuler_SinCos.dat", solutionHistory);
274 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
275 double time_i = (*solutionHistory)[i]->getTime();
278 model->getExactSolution(time_i).get_x()),
280 model->getExactSolution(time_i).get_x_dot()));
281 state->setTime((*solutionHistory)[i]->getTime());
282 solnHistExact->addState(state);
284 writeSolution(
"Tempus_BackwardEuler_SinCos-Ref.dat", solnHistExact);
288 StepSize.push_back(dt);
289 auto solution = Thyra::createMember(model->get_x_space());
290 Thyra::copy(*(integrator->getX()),solution.ptr());
291 solutions.push_back(solution);
292 auto solutionDot = Thyra::createMember(model->get_x_space());
293 Thyra::copy(*(integrator->getXDot()),solutionDot.ptr());
294 solutionsDot.push_back(solutionDot);
295 if (n == nTimeStepSizes-1) {
296 StepSize.push_back(0.0);
297 auto solutionExact = Thyra::createMember(model->get_x_space());
298 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
299 solutions.push_back(solutionExact);
300 auto solutionDotExact = Thyra::createMember(model->get_x_space());
301 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
302 solutionDotExact.ptr());
303 solutionsDot.push_back(solutionDotExact);
309 double xDotSlope = 0.0;
311 double order = stepper->getOrder();
314 solutions, xErrorNorm, xSlope,
315 solutionsDot, xDotErrorNorm, xDotSlope);
332 #ifdef Tempus_ENABLE_MPI
333 comm =
rcp(
new Epetra_MpiComm(MPI_COMM_WORLD));
335 comm =
rcp(
new Epetra_SerialComm);
339 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
340 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
341 std::vector<double> StepSize;
342 std::vector<double> xErrorNorm;
343 std::vector<double> xDotErrorNorm;
344 const int nTimeStepSizes = 5;
346 for (
int n=0; n<nTimeStepSizes; n++) {
350 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
354 const int num_elements = model_pl->
get<
int>(
"num elements");
355 const double left_end = model_pl->
get<
double>(
"left end");
356 const double right_end = model_pl->
get<
double>(
"right end");
357 const double a_convection = model_pl->
get<
double>(
"a (convection)");
358 const double k_source = model_pl->
get<
double>(
"k (source)");
368 ::Stratimikos::DefaultLinearSolverBuilder builder;
371 p->set(
"Linear Solver Type",
"Belos");
372 p->set(
"Preconditioner Type",
"None");
373 builder.setParameterList(p);
376 lowsFactory = builder.createLinearSolveStrategy(
"");
378 model->set_W_factory(lowsFactory);
386 .
sublist(
"Time Step Control").
set(
"Initial Time Step", dt);
387 integrator = Tempus::createIntegratorBasic<double>(pl, model);
390 bool integratorStatus = integrator->advanceTime();
394 double time = integrator->getTime();
395 double timeFinal =pl->sublist(
"Demo Integrator")
396 .sublist(
"Time Step Control").
get<
double>(
"Final Time");
397 double tol = 100.0 * std::numeric_limits<double>::epsilon();
401 StepSize.push_back(dt);
402 auto solution = Thyra::createMember(model->get_x_space());
403 Thyra::copy(*(integrator->getX()),solution.ptr());
404 solutions.push_back(solution);
405 auto solutionDot = Thyra::createMember(model->get_x_space());
406 Thyra::copy(*(integrator->getXDot()),solutionDot.ptr());
407 solutionsDot.push_back(solutionDot);
411 if ((n == nTimeStepSizes-1) && (comm->NumProc() == 1)) {
412 std::ofstream ftmp(
"Tempus_BackwardEuler_CDR.dat");
413 ftmp <<
"TITLE=\"Backward Euler Solution to CDR\"\n"
414 <<
"VARIABLES=\"z\",\"T\"\n";
415 const double dx = std::fabs(left_end-right_end) /
416 static_cast<double>(num_elements);
418 integrator->getSolutionHistory();
419 int nStates = solutionHistory->getNumStates();
420 for (
int i=0; i<nStates; i++) {
423 double ttime = solutionState->getTime();
424 ftmp <<
"ZONE T=\"Time="<<ttime<<
"\", I="
425 <<num_elements+1<<
", F=BLOCK\n";
426 for (
int j = 0; j < num_elements+1; j++) {
427 const double x_coord = left_end +
static_cast<double>(j) * dx;
428 ftmp << x_coord <<
" ";
431 for (
int j=0; j<num_elements+1; j++) ftmp << get_ele(*x, j) <<
" ";
440 double xDotSlope = 0.0;
444 solutions, xErrorNorm, xSlope,
445 solutionsDot, xDotErrorNorm, xDotSlope);
458 if (comm->NumProc() == 1) {
460 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
462 const int num_elements = model_pl->
get<
int>(
"num elements");
463 const double left_end = model_pl->
get<
double>(
"left end");
464 const double right_end = model_pl->
get<
double>(
"right end");
468 std::ofstream ftmp(
"Tempus_BackwardEuler_CDR-Solution.dat");
469 for (
int n = 0; n < num_elements+1; n++) {
470 const double dx = std::fabs(left_end-right_end) /
471 static_cast<double>(num_elements);
472 const double x_coord = left_end +
static_cast<double>(n) * dx;
473 ftmp << x_coord <<
" " << Thyra::get_ele(x,n) << std::endl;
487 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
488 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
489 std::vector<double> StepSize;
490 std::vector<double> xErrorNorm;
491 std::vector<double> xDotErrorNorm;
492 const int nTimeStepSizes = 4;
494 for (
int n=0; n<nTimeStepSizes; n++) {
498 getParametersFromXmlFile(
"Tempus_BackwardEuler_VanDerPol.xml");
506 if (n == nTimeStepSizes-1) dt /= 10.0;
511 .
sublist(
"Time Step Control").
set(
"Initial Time Step", dt);
512 integrator = Tempus::createIntegratorBasic<double>(pl, model);
515 bool integratorStatus = integrator->advanceTime();
519 double time = integrator->getTime();
520 double timeFinal =pl->sublist(
"Demo Integrator")
521 .sublist(
"Time Step Control").
get<
double>(
"Final Time");
522 double tol = 100.0 * std::numeric_limits<double>::epsilon();
526 StepSize.push_back(dt);
527 auto solution = Thyra::createMember(model->get_x_space());
528 Thyra::copy(*(integrator->getX()),solution.ptr());
529 solutions.push_back(solution);
530 auto solutionDot = Thyra::createMember(model->get_x_space());
531 Thyra::copy(*(integrator->getXDot()),solutionDot.ptr());
532 solutionsDot.push_back(solutionDot);
536 if ((n == 0) || (n == nTimeStepSizes-1)) {
537 std::string fname =
"Tempus_BackwardEuler_VanDerPol-Ref.dat";
538 if (n == 0) fname =
"Tempus_BackwardEuler_VanDerPol.dat";
540 integrator->getSolutionHistory();
547 double xDotSlope = 0.0;
549 double order = stepper->getOrder();
552 solutions, xErrorNorm, xSlope,
553 solutionsDot, xDotErrorNorm, xDotSlope);
572 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
581 Tempus::createIntegratorBasic<double>(pl, model);
584 bool integratorStatus = integrator->advanceTime();
589 integrator->getSolutionHistory();
604 model->getNominalValues().get_p(0);
606 Thyra::createMember(model->get_x_space());
608 Thyra::createMember(model->get_f_space());
610 Thyra::createMember(model->get_f_space());
612 model->create_W_op();
614 model->create_W_op();
619 const int num_p = p->range()->dim();
621 Thyra::createMembers(model->get_f_space(), num_p);
623 Thyra::createMembers(model->get_f_space(), num_p);
629 Thyra::createMembers(model->get_x_space(), num_p);
631 Thyra::createMembers(model->get_x_space(), num_p);
632 std::vector<double> nrms(num_p);
636 const int n = solutionHistory->getNumStates();
637 for (
int i=1; i<n; ++i) {
642 x[0] = state->getX();
643 x[1] = prev_state->getX();
644 t[0] = state->getTime();
645 t[1] = prev_state->getTime();
648 const double dt = t[0]-t[1];
649 Thyra::V_StVpStV(x_dot.
ptr(), 1.0/dt, *(x[0]), -1.0/dt, *(x[1]));
653 MEB::InArgs<double> in_args = model->createInArgs();
654 MEB::OutArgs<double> out_args = model->createOutArgs();
656 in_args.set_x_dot(x_dot);
660 const double tol = 1.0e-14;
663 opt_stepper->computeStepResidual(*f, x, t, *p, 0);
665 model->evalModel(in_args, out_args);
666 out_args.set_f(Teuchos::null);
667 Thyra::V_VmV(f.ptr(), *f, *f2);
668 err = Thyra::norm(*f);
673 opt_stepper->computeStepJacobian(*dfdx, x, t, *p, 0, 0);
674 out_args.set_W_op(dfdx2);
675 in_args.set_alpha(1.0/dt);
676 in_args.set_beta(1.0);
677 model->evalModel(in_args, out_args);
678 out_args.set_W_op(Teuchos::null);
679 Thyra::V_VmV(dfdx_mv.
ptr(), *dfdx_mv, *dfdx_mv2);
680 Thyra::norms(*dfdx_mv, Teuchos::arrayViewFromVector(nrms));
682 for (
auto nrm : nrms) err += nrm;
687 opt_stepper->computeStepJacobian(*dfdx, x, t, *p, 0, 1);
688 out_args.set_W_op(dfdx2);
689 in_args.set_alpha(-1.0/dt);
690 in_args.set_beta(0.0);
691 model->evalModel(in_args, out_args);
692 out_args.set_W_op(Teuchos::null);
693 Thyra::V_VmV(dfdx_mv.
ptr(), *dfdx_mv, *dfdx_mv2);
694 Thyra::norms(*dfdx_mv, Teuchos::arrayViewFromVector(nrms));
696 for (
auto nrm : nrms) err += nrm;
700 opt_stepper->computeStepParamDeriv(*dfdp, x, t, *p, 0);
702 0, MEB::Derivative<double>(dfdp2, MEB::DERIV_MV_JACOBIAN_FORM));
703 model->evalModel(in_args, out_args);
704 out_args.set_DfDp(0, MEB::Derivative<double>());
705 Thyra::V_VmV(dfdp.
ptr(), *dfdp, *dfdp2);
706 Thyra::norms(*dfdp, Teuchos::arrayViewFromVector(nrms));
708 for (
auto nrm : nrms) err += nrm;
712 opt_stepper->computeStepSolver(*W, x, t, *p, 0);
714 in_args.set_alpha(1.0/dt);
715 in_args.set_beta(1.0);
716 model->evalModel(in_args, out_args);
717 out_args.set_W(Teuchos::null);
721 Thyra::V_VmV(tmp.
ptr(), *tmp, *tmp2);
722 Thyra::norms(*tmp, Teuchos::arrayViewFromVector(nrms));
724 for (
auto nrm : nrms) err += nrm;
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)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#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)
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)
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...
Backward Euler time stepper.
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.
Stepper interface to support full-space optimization.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
#define TEST_EQUALITY(v1, v2)
Solution state for integrators and steppers.
1D CGFEM model for convection/diffusion/reaction