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_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"
41 using Teuchos::rcp_const_cast;
42 using Teuchos::ParameterList;
43 using Teuchos::sublist;
44 using Teuchos::getParametersFromXmlFile;
56 RCP<ParameterList> pList =
57 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
60 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
63 RCP<ParameterList> tempusPL = sublist(pList,
"Tempus",
true);
67 RCP<Tempus::IntegratorBasic<double> > integrator =
68 Tempus::integratorBasic<double>(tempusPL, model);
70 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
72 stepperPL->set(
"Predictor Stepper Type",
"None");
73 RCP<const ParameterList> defaultPL =
74 integrator->getStepper()->getValidParameters();
76 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
78 std::cout << std::endl;
79 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
80 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
87 RCP<Tempus::IntegratorBasic<double> > integrator =
88 Tempus::integratorBasic<double>(model,
"Backward Euler");
90 RCP<ParameterList> stepperPL = sublist(tempusPL,
"Default Stepper",
true);
91 RCP<const ParameterList> defaultPL =
92 integrator->getStepper()->getValidParameters();
94 bool pass = haveSameValues(*stepperPL, *defaultPL,
true);
96 std::cout << std::endl;
97 std::cout <<
"stepperPL -------------- \n" << *stepperPL << std::endl;
98 std::cout <<
"defaultPL -------------- \n" << *defaultPL << std::endl;
112 RCP<ParameterList> pList =
113 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
114 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
117 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
123 stepper->setModel(model);
124 stepper->initialize();
128 ParameterList tscPL = pl->sublist(
"Default Integrator")
129 .sublist(
"Time Step Control");
130 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
131 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
132 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
133 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
134 timeStepControl->setInitTimeStep(dt);
135 timeStepControl->initialize();
138 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
139 stepper->getModel()->getNominalValues();
141 rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
143 icState->setTime (timeStepControl->getInitTime());
144 icState->setIndex (timeStepControl->getInitIndex());
145 icState->setTimeStep(0.0);
146 icState->setOrder (stepper->getOrder());
151 solutionHistory->setName(
"Forward States");
153 solutionHistory->setStorageLimit(2);
154 solutionHistory->addState(icState);
157 RCP<Tempus::IntegratorBasic<double> > integrator =
158 Tempus::integratorBasic<double>();
159 integrator->setStepperWStepper(stepper);
160 integrator->setTimeStepControl(timeStepControl);
161 integrator->setSolutionHistory(solutionHistory);
163 integrator->initialize();
167 bool integratorStatus = integrator->advanceTime();
168 TEST_ASSERT(integratorStatus)
172 double time = integrator->getTime();
173 double timeFinal =pl->sublist(
"Default Integrator")
174 .sublist(
"Time Step Control").get<
double>(
"Final Time");
175 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
178 RCP<Thyra::VectorBase<double> > x = integrator->getX();
179 RCP<const Thyra::VectorBase<double> > x_exact =
180 model->getExactSolution(time).get_x();
183 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
184 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
187 std::cout <<
" Stepper = BackwardEuler" << std::endl;
188 std::cout <<
" =========================" << std::endl;
189 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) <<
" "
190 << get_ele(*(x_exact), 1) << std::endl;
191 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
192 << get_ele(*(x ), 1) << std::endl;
193 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) <<
" "
194 << get_ele(*(xdiff ), 1) << std::endl;
195 std::cout <<
" =========================" << std::endl;
196 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.798923, 1.0e-4 );
197 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.516729, 1.0e-4 );
205 RCP<Tempus::IntegratorBasic<double> > integrator;
206 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
207 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
208 std::vector<double> StepSize;
209 std::vector<double> xErrorNorm;
210 std::vector<double> xDotErrorNorm;
211 const int nTimeStepSizes = 7;
214 for (
int n=0; n<nTimeStepSizes; n++) {
217 RCP<ParameterList> pList =
218 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
225 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
232 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
233 pl->sublist(
"Default Integrator")
234 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
235 integrator = Tempus::integratorBasic<double>(pl, model);
241 RCP<Thyra::VectorBase<double> > x0 =
242 model->getNominalValues().get_x()->clone_v();
243 integrator->initializeSolutionHistory(0.0, x0);
246 bool integratorStatus = integrator->advanceTime();
247 TEST_ASSERT(integratorStatus)
250 time = integrator->getTime();
251 double timeFinal =pl->sublist(
"Default Integrator")
252 .sublist(
"Time Step Control").get<
double>(
"Final Time");
253 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
257 RCP<const SolutionHistory<double> > solutionHistory =
258 integrator->getSolutionHistory();
259 writeSolution(
"Tempus_BackwardEuler_SinCos.dat", solutionHistory);
262 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
263 double time_i = (*solutionHistory)[i]->getTime();
265 rcp_const_cast<Thyra::VectorBase<double> > (
266 model->getExactSolution(time_i).get_x()),
267 rcp_const_cast<Thyra::VectorBase<double> > (
268 model->getExactSolution(time_i).get_x_dot()));
269 state->setTime((*solutionHistory)[i]->getTime());
270 solnHistExact->addState(state);
272 writeSolution(
"Tempus_BackwardEuler_SinCos-Ref.dat", solnHistExact);
276 StepSize.push_back(dt);
277 auto solution = Thyra::createMember(model->get_x_space());
278 Thyra::copy(*(integrator->getX()),solution.ptr());
279 solutions.push_back(solution);
280 auto solutionDot = Thyra::createMember(model->get_x_space());
281 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
282 solutionsDot.push_back(solutionDot);
283 if (n == nTimeStepSizes-1) {
284 StepSize.push_back(0.0);
285 auto solutionExact = Thyra::createMember(model->get_x_space());
286 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
287 solutions.push_back(solutionExact);
288 auto solutionDotExact = Thyra::createMember(model->get_x_space());
289 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
290 solutionDotExact.ptr());
291 solutionsDot.push_back(solutionDotExact);
297 double xDotSlope = 0.0;
298 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
299 double order = stepper->getOrder();
302 solutions, xErrorNorm, xSlope,
303 solutionsDot, xDotErrorNorm, xDotSlope);
305 TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
306 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.0486418, 1.0e-4 );
307 TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
308 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.0486418, 1.0e-4 );
310 Teuchos::TimeMonitor::summarize();
319 RCP<Epetra_Comm> comm;
320 #ifdef Tempus_ENABLE_MPI
321 comm = rcp(
new Epetra_MpiComm(MPI_COMM_WORLD));
323 comm = rcp(
new Epetra_SerialComm);
326 RCP<Tempus::IntegratorBasic<double> > integrator;
327 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
328 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
329 std::vector<double> StepSize;
330 std::vector<double> xErrorNorm;
331 std::vector<double> xDotErrorNorm;
332 const int nTimeStepSizes = 5;
334 for (
int n=0; n<nTimeStepSizes; n++) {
337 RCP<ParameterList> pList =
338 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
341 RCP<ParameterList> model_pl = sublist(pList,
"CDR Model",
true);
342 const int num_elements = model_pl->get<
int>(
"num elements");
343 const double left_end = model_pl->get<
double>(
"left end");
344 const double right_end = model_pl->get<
double>(
"right end");
345 const double a_convection = model_pl->get<
double>(
"a (convection)");
346 const double k_source = model_pl->get<
double>(
"k (source)");
356 ::Stratimikos::DefaultLinearSolverBuilder builder;
358 auto p = rcp(
new ParameterList);
359 p->set(
"Linear Solver Type",
"Belos");
360 p->set(
"Preconditioner Type",
"None");
361 builder.setParameterList(p);
363 RCP< ::Thyra::LinearOpWithSolveFactoryBase<double> >
364 lowsFactory = builder.createLinearSolveStrategy(
"");
366 model->set_W_factory(lowsFactory);
372 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
373 pl->sublist(
"Demo Integrator")
374 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
375 integrator = Tempus::integratorBasic<double>(pl, model);
378 bool integratorStatus = integrator->advanceTime();
379 TEST_ASSERT(integratorStatus)
382 double time = integrator->getTime();
383 double timeFinal =pl->sublist(
"Demo Integrator")
384 .sublist(
"Time Step Control").get<
double>(
"Final Time");
385 double tol = 100.0 * std::numeric_limits<double>::epsilon();
386 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
389 StepSize.push_back(dt);
390 auto solution = Thyra::createMember(model->get_x_space());
391 Thyra::copy(*(integrator->getX()),solution.ptr());
392 solutions.push_back(solution);
393 auto solutionDot = Thyra::createMember(model->get_x_space());
394 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
395 solutionsDot.push_back(solutionDot);
399 if ((n == nTimeStepSizes-1) && (comm->NumProc() == 1)) {
400 std::ofstream ftmp(
"Tempus_BackwardEuler_CDR.dat");
401 ftmp <<
"TITLE=\"Backward Euler Solution to CDR\"\n"
402 <<
"VARIABLES=\"z\",\"T\"\n";
403 const double dx = std::fabs(left_end-right_end) /
404 static_cast<double>(num_elements);
405 RCP<const SolutionHistory<double> > solutionHistory =
406 integrator->getSolutionHistory();
407 int nStates = solutionHistory->getNumStates();
408 for (
int i=0; i<nStates; i++) {
409 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
410 RCP<const Thyra::VectorBase<double> > x = solutionState->getX();
411 double ttime = solutionState->getTime();
412 ftmp <<
"ZONE T=\"Time="<<ttime<<
"\", I="
413 <<num_elements+1<<
", F=BLOCK\n";
414 for (
int j = 0; j < num_elements+1; j++) {
415 const double x_coord = left_end +
static_cast<double>(j) * dx;
416 ftmp << x_coord <<
" ";
419 for (
int j=0; j<num_elements+1; j++) ftmp << get_ele(*x, j) <<
" ";
428 double xDotSlope = 0.0;
429 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
432 solutions, xErrorNorm, xSlope,
433 solutionsDot, xDotErrorNorm, xDotSlope);
435 TEST_FLOATING_EQUALITY( xSlope, 1.32213, 0.01 );
436 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.116919, 1.0e-4 );
437 TEST_FLOATING_EQUALITY( xDotSlope, 1.32052, 0.01 );
438 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.449888, 1.0e-4 );
446 if (comm->NumProc() == 1) {
447 RCP<ParameterList> pList =
448 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
449 RCP<ParameterList> model_pl = sublist(pList,
"CDR Model",
true);
450 const int num_elements = model_pl->get<
int>(
"num elements");
451 const double left_end = model_pl->get<
double>(
"left end");
452 const double right_end = model_pl->get<
double>(
"right end");
454 const Thyra::VectorBase<double>& x = *(solutions[solutions.size()-1]);
456 std::ofstream ftmp(
"Tempus_BackwardEuler_CDR-Solution.dat");
457 for (
int n = 0; n < num_elements+1; n++) {
458 const double dx = std::fabs(left_end-right_end) /
459 static_cast<double>(num_elements);
460 const double x_coord = left_end +
static_cast<double>(n) * dx;
461 ftmp << x_coord <<
" " << Thyra::get_ele(x,n) << std::endl;
466 Teuchos::TimeMonitor::summarize();
474 RCP<Tempus::IntegratorBasic<double> > integrator;
475 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
476 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
477 std::vector<double> StepSize;
478 std::vector<double> xErrorNorm;
479 std::vector<double> xDotErrorNorm;
480 const int nTimeStepSizes = 4;
482 for (
int n=0; n<nTimeStepSizes; n++) {
485 RCP<ParameterList> pList =
486 getParametersFromXmlFile(
"Tempus_BackwardEuler_VanDerPol.xml");
489 RCP<ParameterList> vdpm_pl = sublist(pList,
"VanDerPolModel",
true);
494 if (n == nTimeStepSizes-1) dt /= 10.0;
497 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
498 pl->sublist(
"Demo Integrator")
499 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
500 integrator = Tempus::integratorBasic<double>(pl, model);
503 bool integratorStatus = integrator->advanceTime();
504 TEST_ASSERT(integratorStatus)
507 double time = integrator->getTime();
508 double timeFinal =pl->sublist(
"Demo Integrator")
509 .sublist(
"Time Step Control").get<
double>(
"Final Time");
510 double tol = 100.0 * std::numeric_limits<double>::epsilon();
511 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
514 StepSize.push_back(dt);
515 auto solution = Thyra::createMember(model->get_x_space());
516 Thyra::copy(*(integrator->getX()),solution.ptr());
517 solutions.push_back(solution);
518 auto solutionDot = Thyra::createMember(model->get_x_space());
519 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
520 solutionsDot.push_back(solutionDot);
524 if ((n == 0) or (n == nTimeStepSizes-1)) {
525 std::string fname =
"Tempus_BackwardEuler_VanDerPol-Ref.dat";
526 if (n == 0) fname =
"Tempus_BackwardEuler_VanDerPol.dat";
527 RCP<const SolutionHistory<double> > solutionHistory =
528 integrator->getSolutionHistory();
535 double xDotSlope = 0.0;
536 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
537 double order = stepper->getOrder();
540 solutions, xErrorNorm, xSlope,
541 solutionsDot, xDotErrorNorm, xDotSlope);
543 TEST_FLOATING_EQUALITY( xSlope, order, 0.10 );
544 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.571031, 1.0e-4 );
545 TEST_FLOATING_EQUALITY( xDotSlope, 1.74898, 0.10 );
546 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 1.0038, 1.0e-4 );
550 Teuchos::TimeMonitor::summarize();
559 RCP<ParameterList> pList =
560 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
563 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
567 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
568 RCP<Tempus::IntegratorBasic<double> >integrator =
569 Tempus::integratorBasic<double>(pl, model);
572 bool integratorStatus = integrator->advanceTime();
573 TEST_ASSERT(integratorStatus);
576 RCP<const SolutionHistory<double> > solutionHistory =
577 integrator->getSolutionHistory();
580 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
581 RCP<Tempus::StepperOptimizationInterface<double> > opt_stepper =
582 Teuchos::rcp_dynamic_cast< Tempus::StepperOptimizationInterface<double> >(
586 TEST_EQUALITY( opt_stepper->stencilLength(), 2);
589 Teuchos::Array< RCP<const Thyra::VectorBase<double> > > x(2);
590 Teuchos::Array<double> t(2);
591 RCP< const Thyra::VectorBase<double> > p =
592 model->getNominalValues().get_p(0);
593 RCP< Thyra::VectorBase<double> > x_dot =
594 Thyra::createMember(model->get_x_space());
595 RCP< Thyra::VectorBase<double> > f =
596 Thyra::createMember(model->get_f_space());
597 RCP< Thyra::VectorBase<double> > f2 =
598 Thyra::createMember(model->get_f_space());
599 RCP< Thyra::LinearOpBase<double> > dfdx =
600 model->create_W_op();
601 RCP< Thyra::LinearOpBase<double> > dfdx2 =
602 model->create_W_op();
603 RCP< Thyra::MultiVectorBase<double> > dfdx_mv =
604 Teuchos::rcp_dynamic_cast< Thyra::MultiVectorBase<double> >(dfdx,
true);
605 RCP< Thyra::MultiVectorBase<double> > dfdx_mv2 =
606 Teuchos::rcp_dynamic_cast< Thyra::MultiVectorBase<double> >(dfdx2,
true);
607 const int num_p = p->range()->dim();
608 RCP< Thyra::MultiVectorBase<double> > dfdp =
609 Thyra::createMembers(model->get_f_space(), num_p);
610 RCP< Thyra::MultiVectorBase<double> > dfdp2 =
611 Thyra::createMembers(model->get_f_space(), num_p);
612 RCP< Thyra::LinearOpWithSolveBase<double> > W =
614 RCP< Thyra::LinearOpWithSolveBase<double> > W2 =
616 RCP< Thyra::MultiVectorBase<double> > tmp =
617 Thyra::createMembers(model->get_x_space(), num_p);
618 RCP< Thyra::MultiVectorBase<double> > tmp2 =
619 Thyra::createMembers(model->get_x_space(), num_p);
620 std::vector<double> nrms(num_p);
624 const int n = solutionHistory->getNumStates();
625 for (
int i=1; i<n; ++i) {
626 RCP<const SolutionState<double> > state = (*solutionHistory)[i];
627 RCP<const SolutionState<double> > prev_state = (*solutionHistory)[i-1];
630 x[0] = state->getX();
631 x[1] = prev_state->getX();
632 t[0] = state->getTime();
633 t[1] = prev_state->getTime();
636 const double dt = t[0]-t[1];
637 Thyra::V_StVpStV(x_dot.ptr(), 1.0/dt, *(x[0]), -1.0/dt, *(x[1]));
640 typedef Thyra::ModelEvaluatorBase MEB;
641 MEB::InArgs<double> in_args = model->createInArgs();
642 MEB::OutArgs<double> out_args = model->createOutArgs();
644 in_args.set_x_dot(x_dot);
648 const double tol = 1.0e-14;
651 opt_stepper->computeStepResidual(*f, x, t, *p, 0);
653 model->evalModel(in_args, out_args);
654 out_args.set_f(Teuchos::null);
655 Thyra::V_VmV(f.ptr(), *f, *f2);
656 err = Thyra::norm(*f);
657 TEST_FLOATING_EQUALITY(err, 0.0, tol);
661 opt_stepper->computeStepJacobian(*dfdx, x, t, *p, 0, 0);
662 out_args.set_W_op(dfdx2);
663 in_args.set_alpha(1.0/dt);
664 in_args.set_beta(1.0);
665 model->evalModel(in_args, out_args);
666 out_args.set_W_op(Teuchos::null);
667 Thyra::V_VmV(dfdx_mv.ptr(), *dfdx_mv, *dfdx_mv2);
668 Thyra::norms(*dfdx_mv, Teuchos::arrayViewFromVector(nrms));
670 for (
auto nrm : nrms) err += nrm;
671 TEST_FLOATING_EQUALITY(err, 0.0, tol);
675 opt_stepper->computeStepJacobian(*dfdx, x, t, *p, 0, 1);
676 out_args.set_W_op(dfdx2);
677 in_args.set_alpha(-1.0/dt);
678 in_args.set_beta(0.0);
679 model->evalModel(in_args, out_args);
680 out_args.set_W_op(Teuchos::null);
681 Thyra::V_VmV(dfdx_mv.ptr(), *dfdx_mv, *dfdx_mv2);
682 Thyra::norms(*dfdx_mv, Teuchos::arrayViewFromVector(nrms));
684 for (
auto nrm : nrms) err += nrm;
685 TEST_FLOATING_EQUALITY(err, 0.0, tol);
688 opt_stepper->computeStepParamDeriv(*dfdp, x, t, *p, 0);
690 0, MEB::Derivative<double>(dfdp2, MEB::DERIV_MV_JACOBIAN_FORM));
691 model->evalModel(in_args, out_args);
692 out_args.set_DfDp(0, MEB::Derivative<double>());
693 Thyra::V_VmV(dfdp.ptr(), *dfdp, *dfdp2);
694 Thyra::norms(*dfdp, Teuchos::arrayViewFromVector(nrms));
696 for (
auto nrm : nrms) err += nrm;
697 TEST_FLOATING_EQUALITY(err, 0.0, tol);
700 opt_stepper->computeStepSolver(*W, x, t, *p, 0);
702 in_args.set_alpha(1.0/dt);
703 in_args.set_beta(1.0);
704 model->evalModel(in_args, out_args);
705 out_args.set_W(Teuchos::null);
707 Thyra::solve(*W, Thyra::NOTRANS, *dfdp2, tmp.ptr());
708 Thyra::solve(*W2, Thyra::NOTRANS, *dfdp2, tmp2.ptr());
709 Thyra::V_VmV(tmp.ptr(), *tmp, *tmp2);
710 Thyra::norms(*tmp, Teuchos::arrayViewFromVector(nrms));
712 for (
auto nrm : nrms) err += nrm;
713 TEST_FLOATING_EQUALITY(err, 0.0, tol);
716 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...
Backward Euler time stepper.
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
1D CGFEM model for convection/diffusion/reaction
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.