46 std::vector<std::string> options;
47 options.push_back(
"Default Parameters");
48 options.push_back(
"ICConsistency and Check");
50 for (
const auto& option : options) {
52 RCP<ParameterList> pList =
53 getParametersFromXmlFile(
"Tempus_Leapfrog_SinCos.xml");
54 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
57 RCP<ParameterList> hom_pl = sublist(pList,
"HarmonicOscillatorModel",
true);
62 stepper->setModel(model);
63 if (option ==
"ICConsistency and Check") {
64 stepper->setICConsistency(
"Consistent");
65 stepper->setICConsistencyCheck(
true);
67 stepper->initialize();
72 pl->sublist(
"Default Integrator").sublist(
"Time Step Control");
73 timeStepControl->setInitIndex(tscPL.get<
int>(
"Initial Time Index"));
74 timeStepControl->setInitTime(tscPL.get<
double>(
"Initial Time"));
75 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
76 timeStepControl->setInitTimeStep(dt);
77 timeStepControl->initialize();
80 auto inArgsIC = model->getNominalValues();
81 auto icX = rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x());
83 rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x_dot());
85 rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x_dot_dot());
86 auto icState = Tempus::createSolutionStateX<double>(icX, icXDot, icXDotDot);
87 icState->setTime(timeStepControl->getInitTime());
88 icState->setIndex(timeStepControl->getInitIndex());
89 icState->setTimeStep(0.0);
90 icState->setOrder(stepper->getOrder());
95 solutionHistory->setName(
"Forward States");
97 solutionHistory->setStorageLimit(2);
98 solutionHistory->addState(icState);
101 stepper->setInitialConditions(solutionHistory);
104 RCP<Tempus::IntegratorBasic<double>> integrator =
105 Tempus::createIntegratorBasic<double>();
106 integrator->setStepper(stepper);
107 integrator->setTimeStepControl(timeStepControl);
108 integrator->setSolutionHistory(solutionHistory);
110 integrator->initialize();
113 bool integratorStatus = integrator->advanceTime();
114 TEST_ASSERT(integratorStatus)
117 double time = integrator->getTime();
118 double timeFinal = pl->sublist(
"Default Integrator")
119 .sublist(
"Time Step Control")
120 .get<
double>(
"Final Time");
121 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
124 RCP<Thyra::VectorBase<double>> x = integrator->getX();
125 RCP<const Thyra::VectorBase<double>> x_exact =
126 model->getExactSolution(time).get_x();
129 RCP<Thyra::VectorBase<double>> xdiff = x->clone_v();
130 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
133 out <<
" Stepper = " << stepper->description() <<
"\n with "
134 << option << std::endl;
135 out <<
" =========================" << std::endl;
136 out <<
" Exact solution : " << get_ele(*(x_exact), 0) << std::endl;
137 out <<
" Computed solution: " << get_ele(*(x), 0) << std::endl;
138 out <<
" Difference : " << get_ele(*(xdiff), 0) << std::endl;
139 out <<
" =========================" << std::endl;
140 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.167158, 1.0e-4);
148 RCP<Tempus::IntegratorBasic<double>> integrator;
149 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
150 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
151 std::vector<double> StepSize;
152 std::vector<double> xErrorNorm;
153 std::vector<double> xDotErrorNorm;
154 const int nTimeStepSizes = 9;
158 RCP<ParameterList> pList =
159 getParametersFromXmlFile(
"Tempus_Leapfrog_SinCos.xml");
162 RCP<ParameterList> hom_pl = sublist(pList,
"HarmonicOscillatorModel",
true);
166 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
170 double dt = pl->sublist(
"Default Integrator")
171 .sublist(
"Time Step Control")
172 .get<
double>(
"Initial Time Step");
175 for (
int n = 0; n < nTimeStepSizes; n++) {
178 out <<
"\n \n time step #" << n <<
" (out of " << nTimeStepSizes - 1
179 <<
"), dt = " << dt <<
"\n";
180 pl->sublist(
"Default Integrator")
181 .sublist(
"Time Step Control")
182 .set(
"Initial Time Step", dt);
183 integrator = Tempus::createIntegratorBasic<double>(pl, model);
186 bool integratorStatus = integrator->advanceTime();
187 TEST_ASSERT(integratorStatus)
190 time = integrator->getTime();
191 double timeFinal = pl->sublist(
"Default Integrator")
192 .sublist(
"Time Step Control")
193 .get<
double>(
"Final Time");
194 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
197 if (n == nTimeStepSizes - 1) {
198 RCP<const SolutionHistory<double>> solutionHistory =
199 integrator->getSolutionHistory();
200 writeSolution(
"Tempus_Leapfrog_SinCos.dat", solutionHistory);
203 for (
int i = 0; i < solutionHistory->getNumStates(); i++) {
204 double time_i = (*solutionHistory)[i]->getTime();
207 model->getExactSolution(time_i).get_x()),
209 model->getExactSolution(time_i).get_x_dot()));
210 state->setTime((*solutionHistory)[i]->getTime());
211 solnHistExact->addState(state);
213 writeSolution(
"Tempus_Leapfrog_SinCos-Ref.dat", solnHistExact);
218 StepSize.push_back(dt);
219 auto solution = Thyra::createMember(model->get_x_space());
220 Thyra::copy(*(integrator->getX()), solution.ptr());
221 solutions.push_back(solution);
222 auto solutionDot = Thyra::createMember(model->get_x_space());
223 Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
224 solutionsDot.push_back(solutionDot);
225 if (n == nTimeStepSizes - 1) {
226 StepSize.push_back(0.0);
227 auto solutionExact = Thyra::createMember(model->get_x_space());
228 Thyra::copy(*(model->getExactSolution(time).get_x()),
229 solutionExact.ptr());
230 solutions.push_back(solutionExact);
231 auto solutionDotExact = Thyra::createMember(model->get_x_space());
232 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
233 solutionDotExact.ptr());
234 solutionsDot.push_back(solutionDotExact);
240 double xDotSlope = 0.0;
241 RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
242 double order = stepper->getOrder();
244 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
247 TEST_FLOATING_EQUALITY(xSlope, order, 0.02);
248 TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.0157928, 1.0e-4);
249 TEST_FLOATING_EQUALITY(xDotSlope, 1.09387, 0.01);
250 TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.563002, 1.0e-4);
252 Teuchos::TimeMonitor::summarize();
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
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)