Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_IMEX_RKTest.cpp
Go to the documentation of this file.
1//@HEADER
2// *****************************************************************************
3// Tempus: Time Integration and Sensitivity Analysis Package
4//
5// Copyright 2017 NTESS and the Tempus contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8//@HEADER
9
10#include "Teuchos_UnitTestHarness.hpp"
11#include "Teuchos_XMLParameterListHelpers.hpp"
12#include "Teuchos_TimeMonitor.hpp"
13
14#include "Thyra_VectorStdOps.hpp"
15
16#include "Tempus_IntegratorBasic.hpp"
17#include "Tempus_WrapperModelEvaluatorPairIMEX_Basic.hpp"
18#include "Tempus_StepperIMEX_RK.hpp"
19
20#include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
21#include "../TestModels/VanDerPol_IMEX_ImplicitModel.hpp"
22#include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
23
24#include <fstream>
25#include <vector>
26
27namespace Tempus_Test {
28
29using Teuchos::getParametersFromXmlFile;
30using Teuchos::ParameterList;
31using Teuchos::RCP;
32using Teuchos::rcp;
33using Teuchos::rcp_const_cast;
34using Teuchos::sublist;
35
39
40// ************************************************************
41// ************************************************************
42TEUCHOS_UNIT_TEST(IMEX_RK, ConstructingFromDefaults)
43{
44 double dt = 0.025;
45
46 // Read params from .xml file
47 RCP<ParameterList> pList =
48 getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
49 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
50
51 // Setup the explicit VanDerPol ModelEvaluator
52 RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
53 auto explicitModel = rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL));
54
55 // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
56 auto implicitModel = rcp(new VanDerPol_IMEX_ImplicitModel<double>(vdpmPL));
57
58 // Setup the IMEX Pair ModelEvaluator
60 explicitModel, implicitModel));
61
62 // Setup Stepper for field solve ----------------------------
63 auto stepper = rcp(new Tempus::StepperIMEX_RK<double>());
64 stepper->setModel(model);
65 stepper->initialize();
66
67 // Setup TimeStepControl ------------------------------------
68 auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
69 ParameterList tscPL =
70 pl->sublist("Default Integrator").sublist("Time Step Control");
71 timeStepControl->setInitIndex(tscPL.get<int>("Initial Time Index"));
72 timeStepControl->setInitTime(tscPL.get<double>("Initial Time"));
73 timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
74 timeStepControl->setInitTimeStep(dt);
75 timeStepControl->initialize();
76
77 // Setup initial condition SolutionState --------------------
78 auto inArgsIC = model->getNominalValues();
79 auto icSolution = rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x());
80 auto icState = Tempus::createSolutionStateX(icSolution);
81 icState->setTime(timeStepControl->getInitTime());
82 icState->setIndex(timeStepControl->getInitIndex());
83 icState->setTimeStep(0.0);
84 icState->setOrder(stepper->getOrder());
85 icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
86
87 // Setup SolutionHistory ------------------------------------
88 auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
89 solutionHistory->setName("Forward States");
90 solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
91 solutionHistory->setStorageLimit(2);
92 solutionHistory->addState(icState);
93
94 // Setup Integrator -----------------------------------------
95 RCP<Tempus::IntegratorBasic<double>> integrator =
96 Tempus::createIntegratorBasic<double>();
97 integrator->setStepper(stepper);
98 integrator->setTimeStepControl(timeStepControl);
99 integrator->setSolutionHistory(solutionHistory);
100 integrator->initialize();
101
102 // Integrate to timeMax
103 bool integratorStatus = integrator->advanceTime();
104 TEST_ASSERT(integratorStatus)
105
106 // Test if at 'Final Time'
107 double time = integrator->getTime();
108 double timeFinal = pl->sublist("Default Integrator")
109 .sublist("Time Step Control")
110 .get<double>("Final Time");
111 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
112
113 // Time-integrated solution and the exact solution
114 RCP<Thyra::VectorBase<double>> x = integrator->getX();
115
116 // Check the order and intercept
117 out << " Stepper = " << stepper->description() << std::endl;
118 out << " =========================" << std::endl;
119 out << " Computed solution: " << get_ele(*(x), 0) << " "
120 << get_ele(*(x), 1) << std::endl;
121 out << " =========================" << std::endl;
122 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 1.810210, 1.0e-4);
123 TEST_FLOATING_EQUALITY(get_ele(*(x), 1), -0.754602, 1.0e-4);
124}
125
126// ************************************************************
127// ************************************************************
128TEUCHOS_UNIT_TEST(IMEX_RK, VanDerPol)
129{
130 std::vector<std::string> stepperTypes;
131 stepperTypes.push_back("IMEX RK 1st order");
132 stepperTypes.push_back("SSP1_111");
133 stepperTypes.push_back("IMEX RK SSP2");
134 stepperTypes.push_back("SSP2_222");
135 stepperTypes.push_back("IMEX RK ARS 233");
136 stepperTypes.push_back("General IMEX RK");
137 stepperTypes.push_back("IMEX RK SSP3");
138
139 std::vector<double> stepperOrders;
140 stepperOrders.push_back(1.07964);
141 stepperOrders.push_back(1.07964); // SSP1_111
142 stepperOrders.push_back(2.00408);
143 stepperOrders.push_back(2.76941); // SSP2_222
144 stepperOrders.push_back(2.70655);
145 stepperOrders.push_back(2.00211);
146 stepperOrders.push_back(2.00211);
147
148 std::vector<double> stepperErrors;
149 stepperErrors.push_back(0.0046423);
150 stepperErrors.push_back(0.103569); // SSP1_111
151 stepperErrors.push_back(0.0154534);
152 stepperErrors.push_back(0.000533759); // SSP2_222
153 stepperErrors.push_back(0.000298908);
154 stepperErrors.push_back(0.0071546);
155 stepperErrors.push_back(0.0151202);
156
157 std::vector<double> stepperInitDt;
158 stepperInitDt.push_back(0.0125);
159 stepperInitDt.push_back(0.0125);
160 stepperInitDt.push_back(0.05);
161 stepperInitDt.push_back(0.05);
162 stepperInitDt.push_back(0.05);
163 stepperInitDt.push_back(0.05);
164 stepperInitDt.push_back(0.05);
165
166 TEUCHOS_ASSERT(stepperTypes.size() == stepperOrders.size());
167 TEUCHOS_ASSERT(stepperTypes.size() == stepperErrors.size());
168 TEUCHOS_ASSERT(stepperTypes.size() == stepperInitDt.size());
169
170 std::vector<std::string>::size_type m;
171 for (m = 0; m != stepperTypes.size(); m++) {
172 std::string stepperType = stepperTypes[m];
173 std::string stepperName = stepperTypes[m];
174 std::replace(stepperName.begin(), stepperName.end(), ' ', '_');
175 std::replace(stepperName.begin(), stepperName.end(), '/', '.');
176
177 RCP<Tempus::IntegratorBasic<double>> integrator;
178 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
179 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
180 std::vector<double> StepSize;
181 std::vector<double> xErrorNorm;
182 std::vector<double> xDotErrorNorm;
183
184 const int nTimeStepSizes = 3; // 6 for error plot
185 double dt = stepperInitDt[m];
186 double time = 0.0;
187 for (int n = 0; n < nTimeStepSizes; n++) {
188 // Read params from .xml file
189 RCP<ParameterList> pList =
190 getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
191
192 // Setup the explicit VanDerPol ModelEvaluator
193 RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
194 auto explicitModel =
196
197 // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
198 auto implicitModel =
200
201 // Setup the IMEX Pair ModelEvaluator
203 explicitModel, implicitModel));
204
205 // Set the Stepper
206 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
207 if (stepperType == "General IMEX RK") {
208 // use the appropriate stepper sublist
209 pl->sublist("Default Integrator")
210 .set("Stepper Name", "General IMEX RK");
211 }
212 else {
213 pl->sublist("Default Stepper").set("Stepper Type", stepperType);
214 }
215
216 // Set the step size
217 if (n == nTimeStepSizes - 1)
218 dt /= 10.0;
219 else
220 dt /= 2;
221
222 // Setup the Integrator and reset initial time step
223 pl->sublist("Default Integrator")
224 .sublist("Time Step Control")
225 .set("Initial Time Step", dt);
226 integrator = Tempus::createIntegratorBasic<double>(pl, model);
227
228 // Integrate to timeMax
229 bool integratorStatus = integrator->advanceTime();
230 TEST_ASSERT(integratorStatus)
231
232 // Test if at 'Final Time'
233 time = integrator->getTime();
234 double timeFinal = pl->sublist("Default Integrator")
235 .sublist("Time Step Control")
236 .get<double>("Final Time");
237 double tol = 100.0 * std::numeric_limits<double>::epsilon();
238 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
239
240 // Store off the final solution and step size
241 StepSize.push_back(dt);
242 auto solution = Thyra::createMember(model->get_x_space());
243 Thyra::copy(*(integrator->getX()), solution.ptr());
244 solutions.push_back(solution);
245 auto solutionDot = Thyra::createMember(model->get_x_space());
246 Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
247 solutionsDot.push_back(solutionDot);
248
249 // Output finest temporal solution for plotting
250 // This only works for ONE MPI process
251 if ((n == 0) || (n == nTimeStepSizes - 1)) {
252 std::string fname = "Tempus_" + stepperName + "_VanDerPol-Ref.dat";
253 if (n == 0) fname = "Tempus_" + stepperName + "_VanDerPol.dat";
254 RCP<const SolutionHistory<double>> solutionHistory =
255 integrator->getSolutionHistory();
256 writeSolution(fname, solutionHistory);
257 }
258 }
259
260 // Check the order and intercept
261 double xSlope = 0.0;
262 double xDotSlope = 0.0;
263 RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
264 // double order = stepper->getOrder();
265
266 // xDot not yet available for IMEX-RK methods, e.g., are not calc. and zero.
267 solutionsDot.clear();
268
269 writeOrderError("Tempus_" + stepperName + "_VanDerPol-Error.dat", stepper,
270 StepSize, solutions, xErrorNorm, xSlope, solutionsDot,
271 xDotErrorNorm, xDotSlope, out);
272
273 TEST_FLOATING_EQUALITY(xSlope, stepperOrders[m], 0.02);
274 TEST_FLOATING_EQUALITY(xErrorNorm[0], stepperErrors[m], 1.0e-4);
275 // TEST_FLOATING_EQUALITY( xDotSlope, 1.74898, 0.10 );
276 // TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 1.0038, 1.0e-4 );
277 }
278 // Teuchos::TimeMonitor::summarize();
279}
280
281} // namespace Tempus_Test
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers.
Implicit-Explicit Runge-Kutta (IMEX-RK) time stepper.
TimeStepControl manages the time step size. There several mechanisms that effect the time step size a...
ModelEvaluator pair for implicit and explicit (IMEX) evaulations.
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::FancyOStream &out)
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.