9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
13 #include "Thyra_VectorStdOps.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperOperatorSplit.hpp"
18 #include "Tempus_StepperForwardEuler.hpp"
19 #include "Tempus_StepperBackwardEuler.hpp"
21 #include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
22 #include "../TestModels/VanDerPol_IMEX_ImplicitModel.hpp"
23 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
28 namespace Tempus_Test {
32 using Teuchos::rcp_const_cast;
33 using Teuchos::ParameterList;
34 using Teuchos::sublist;
35 using Teuchos::getParametersFromXmlFile;
49 RCP<ParameterList> pList =
50 getParametersFromXmlFile(
"Tempus_OperatorSplit_VanDerPol.xml");
51 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
54 RCP<ParameterList> vdpmPL = sublist(pList,
"VanDerPolModel",
true);
64 RCP<Tempus::StepperFactory<double> > sf =
68 sf->createStepperForwardEuler(explicitModel, Teuchos::null);
71 sf->createStepperBackwardEuler(implicitModel, Teuchos::null);
73 stepper->addStepper(subStepper1);
74 stepper->addStepper(subStepper2);
75 stepper->initialize();
80 ParameterList tscPL = pl->sublist(
"Demo Integrator")
81 .sublist(
"Time Step Control");
82 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
83 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
84 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
85 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
86 timeStepControl->setInitTimeStep(dt);
87 timeStepControl->initialize();
90 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
91 stepper->getModel()->getNominalValues();
92 auto icX = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
93 auto icXDot = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x_dot());
95 icState->setTime (timeStepControl->getInitTime());
96 icState->setIndex (timeStepControl->getInitIndex());
97 icState->setTimeStep(0.0);
98 icState->setOrder (stepper->getOrder());
103 solutionHistory->setName(
"Forward States");
105 solutionHistory->setStorageLimit(2);
106 solutionHistory->addState(icState);
109 RCP<Tempus::IntegratorBasic<double> > integrator =
110 Tempus::integratorBasic<double>();
111 integrator->setStepperWStepper(stepper);
112 integrator->setTimeStepControl(timeStepControl);
113 integrator->setSolutionHistory(solutionHistory);
115 integrator->initialize();
119 bool integratorStatus = integrator->advanceTime();
120 TEST_ASSERT(integratorStatus)
124 double time = integrator->getTime();
125 double timeFinal =pl->sublist(
"Demo Integrator")
126 .sublist(
"Time Step Control").get<
double>(
"Final Time");
127 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
130 RCP<Thyra::VectorBase<double> > x = integrator->getX();
133 std::cout <<
" Stepper = " << stepper->description() << std::endl;
134 std::cout <<
" =========================" << std::endl;
135 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) <<
" "
136 << get_ele(*(x ), 1) << std::endl;
137 std::cout <<
" =========================" << std::endl;
138 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), -2.223910, 1.0e-4);
139 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.565441, 1.0e-4);
147 RCP<Tempus::IntegratorBasic<double> > integrator;
148 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
149 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
150 std::vector<double> StepSize;
151 std::vector<double> xErrorNorm;
152 std::vector<double> xDotErrorNorm;
153 const int nTimeStepSizes = 4;
156 for (
int n=0; n<nTimeStepSizes; n++) {
159 RCP<ParameterList> pList =
160 getParametersFromXmlFile(
"Tempus_OperatorSplit_VanDerPol.xml");
163 RCP<ParameterList> vdpmPL = sublist(pList,
"VanDerPolModel",
true);
170 std::vector<RCP<const Thyra::ModelEvaluator<double> > > models;
171 models.push_back(explicitModel);
172 models.push_back(implicitModel);
176 if (n == nTimeStepSizes-1) dt /= 10.0;
179 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
180 pl->sublist(
"Demo Integrator")
181 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
182 integrator = Tempus::integratorBasic<double>(pl, models);
185 bool integratorStatus = integrator->advanceTime();
186 TEST_ASSERT(integratorStatus)
189 time = integrator->getTime();
190 double timeFinal =pl->sublist(
"Demo Integrator")
191 .sublist(
"Time Step Control").get<
double>(
"Final Time");
192 double tol = 100.0 * std::numeric_limits<double>::epsilon();
193 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
196 StepSize.push_back(dt);
197 auto solution = Thyra::createMember(implicitModel->get_x_space());
198 Thyra::copy(*(integrator->getX()),solution.ptr());
199 solutions.push_back(solution);
200 auto solutionDot = Thyra::createMember(implicitModel->get_x_space());
201 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
202 solutionsDot.push_back(solutionDot);
206 if ((n == 0) or (n == nTimeStepSizes-1)) {
207 std::string fname =
"Tempus_OperatorSplit_VanDerPol-Ref.dat";
208 if (n == 0) fname =
"Tempus_OperatorSplit_VanDerPol.dat";
209 RCP<const SolutionHistory<double> > solutionHistory =
210 integrator->getSolutionHistory();
218 double xDotSlope = 0.0;
219 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
220 double order = stepper->getOrder();
223 solutions, xErrorNorm, xSlope,
224 solutionsDot, xDotErrorNorm, xDotSlope);
226 TEST_FLOATING_EQUALITY( xSlope, order, 0.05 );
227 TEST_FLOATING_EQUALITY( xDotSlope, order, 0.05 );
228 TEST_FLOATING_EQUALITY( xErrorNorm[0], 1.27294, 1.0e-4 );
229 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 12.7102, 1.0e-4 );
231 Teuchos::TimeMonitor::summarize();
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.
van der Pol model formulated for IMEX-RK.
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > solutionHistory)
OperatorSplit stepper loops through the Stepper list.
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...
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Keep a fix number of states.
van der Pol model formulated for IMEX.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...