Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_IMEX_RK_PartitionedTest.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 
12 #include "Teuchos_TimeMonitor.hpp"
13 
14 #include "Thyra_VectorStdOps.hpp"
15 
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_WrapperModelEvaluatorPairPartIMEX_Basic.hpp"
18 #include "Tempus_StepperIMEX_RK_Partition.hpp"
19 
20 #include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
21 #include "../TestModels/VanDerPol_IMEXPart_ImplicitModel.hpp"
22 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
23 
24 #include <fstream>
25 #include <vector>
26 
27 namespace Tempus_Test {
28 
29 using Teuchos::getParametersFromXmlFile;
31 using Teuchos::RCP;
32 using Teuchos::rcp;
33 using Teuchos::rcp_const_cast;
34 using Teuchos::sublist;
35 
39 
40 // ************************************************************
41 // ************************************************************
42 TEUCHOS_UNIT_TEST(IMEX_RK_Partitioned, 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  const bool useProductVector = true;
54  auto explicitModel =
55  rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL, useProductVector));
56 
57  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
58  auto implicitModel =
60 
61  // Setup the IMEX Pair ModelEvaluator
62  const int numExplicitBlocks = 1;
63  const int parameterIndex = 4;
65  explicitModel, implicitModel, numExplicitBlocks, parameterIndex));
66 
67  // Setup Stepper for field solve ----------------------------
68  auto stepper = rcp(new Tempus::StepperIMEX_RK_Partition<double>());
69  stepper->setModel(model);
70  stepper->initialize();
71 
72  // Setup TimeStepControl ------------------------------------
73  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
74  ParameterList tscPL =
75  pl->sublist("Default Integrator").sublist("Time Step Control");
76  timeStepControl->setInitIndex(tscPL.get<int>("Initial Time Index"));
77  timeStepControl->setInitTime(tscPL.get<double>("Initial Time"));
78  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
79  timeStepControl->setInitTimeStep(dt);
80  timeStepControl->initialize();
81 
82  // Setup initial condition SolutionState --------------------
83  auto inArgsIC = model->getNominalValues();
84  auto icSolution = rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x());
85  auto icState = Tempus::createSolutionStateX(icSolution);
86  icState->setTime(timeStepControl->getInitTime());
87  icState->setIndex(timeStepControl->getInitIndex());
88  icState->setTimeStep(0.0);
89  icState->setOrder(stepper->getOrder());
90  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
91 
92  // Setup SolutionHistory ------------------------------------
93  auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
94  solutionHistory->setName("Forward States");
95  solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
96  solutionHistory->setStorageLimit(2);
97  solutionHistory->addState(icState);
98 
99  // Setup Integrator -----------------------------------------
101  Tempus::createIntegratorBasic<double>();
102  integrator->setStepper(stepper);
103  integrator->setTimeStepControl(timeStepControl);
104  integrator->setSolutionHistory(solutionHistory);
105  // integrator->setObserver(...);
106  integrator->initialize();
107 
108  // Integrate to timeMax
109  bool integratorStatus = integrator->advanceTime();
110  TEST_ASSERT(integratorStatus)
111 
112  // Test if at 'Final Time'
113  double time = integrator->getTime();
114  double timeFinal = pl->sublist("Default Integrator")
115  .sublist("Time Step Control")
116  .get<double>("Final Time");
117  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
118 
119  // Time-integrated solution and the exact solution
120  RCP<Thyra::VectorBase<double>> x = integrator->getX();
121 
122  // Check the order and intercept
123  out << " Stepper = " << stepper->description() << std::endl;
124  out << " =========================" << std::endl;
125  out << " Computed solution: " << get_ele(*(x), 0) << " "
126  << get_ele(*(x), 1) << std::endl;
127  out << " =========================" << std::endl;
128  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 1.810210, 1.0e-4);
129  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), -0.754602, 1.0e-4);
130 }
131 
132 // ************************************************************
133 // ************************************************************
134 TEUCHOS_UNIT_TEST(IMEX_RK_Partitioned, VanDerPol)
135 {
136  std::vector<std::string> stepperTypes;
137  stepperTypes.push_back("Partitioned IMEX RK 1st order");
138  stepperTypes.push_back("Partitioned IMEX RK SSP2");
139  stepperTypes.push_back("Partitioned IMEX RK ARS 233");
140  stepperTypes.push_back("General Partitioned IMEX RK");
141 
142  std::vector<double> stepperOrders;
143  stepperOrders.push_back(1.07964);
144  stepperOrders.push_back(2.00408);
145  stepperOrders.push_back(2.70655);
146  stepperOrders.push_back(2.00211);
147 
148  std::vector<double> stepperErrors;
149  stepperErrors.push_back(0.0046423);
150  stepperErrors.push_back(0.0154534);
151  stepperErrors.push_back(0.000298908);
152  stepperErrors.push_back(0.0071546);
153 
154  std::vector<double> stepperInitDt;
155  stepperInitDt.push_back(0.0125);
156  stepperInitDt.push_back(0.05);
157  stepperInitDt.push_back(0.05);
158  stepperInitDt.push_back(0.05);
159 
160  std::vector<std::string>::size_type m;
161  for (m = 0; m != stepperTypes.size(); m++) {
162  std::string stepperType = stepperTypes[m];
163  std::string stepperName = stepperTypes[m];
164  std::replace(stepperName.begin(), stepperName.end(), ' ', '_');
165  std::replace(stepperName.begin(), stepperName.end(), '/', '.');
166 
168  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
169  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
170  std::vector<double> StepSize;
171  std::vector<double> xErrorNorm;
172  std::vector<double> xDotErrorNorm;
173 
174  const int nTimeStepSizes = 3; // 6 for error plot
175  double dt = stepperInitDt[m];
176  double time = 0.0;
177  for (int n = 0; n < nTimeStepSizes; n++) {
178  // Read params from .xml file
179  RCP<ParameterList> pList =
180  getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
181 
182  // Setup the explicit VanDerPol ModelEvaluator
183  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
184  const bool useProductVector = true;
185  auto explicitModel = rcp(
186  new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL, useProductVector));
187 
188  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
189  auto implicitModel =
191 
192  // Setup the IMEX Pair ModelEvaluator
193  const int numExplicitBlocks = 1;
194  const int parameterIndex = 4;
195  auto model =
197  explicitModel, implicitModel, numExplicitBlocks, parameterIndex));
198 
199  // Set the Stepper
200  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
201 
202  if (stepperType == "General Partitioned IMEX RK") {
203  // use the appropriate stepper sublist
204  pl->sublist("Default Integrator")
205  .set("Stepper Name", "General IMEX RK");
206  }
207  else {
208  pl->sublist("Default Stepper").set("Stepper Type", stepperType);
209  }
210 
211  // Set the step size
212  if (n == nTimeStepSizes - 1)
213  dt /= 10.0;
214  else
215  dt /= 2;
216 
217  // Setup the Integrator and reset initial time step
218  pl->sublist("Default Integrator")
219  .sublist("Time Step Control")
220  .set("Initial Time Step", dt);
221  integrator = Tempus::createIntegratorBasic<double>(pl, model);
222 
223  // Integrate to timeMax
224  bool integratorStatus = integrator->advanceTime();
225  TEST_ASSERT(integratorStatus)
226 
227  // Test if at 'Final Time'
228  time = integrator->getTime();
229  double timeFinal = pl->sublist("Default Integrator")
230  .sublist("Time Step Control")
231  .get<double>("Final Time");
232  double tol = 100.0 * std::numeric_limits<double>::epsilon();
233  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
234 
235  // Store off the final solution and step size
236  StepSize.push_back(dt);
237  auto solution = Thyra::createMember(model->get_x_space());
238  Thyra::copy(*(integrator->getX()), solution.ptr());
239  solutions.push_back(solution);
240  auto solutionDot = Thyra::createMember(model->get_x_space());
241  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
242  solutionsDot.push_back(solutionDot);
243 
244  // Output finest temporal solution for plotting
245  // This only works for ONE MPI process
246  if ((n == 0) || (n == nTimeStepSizes - 1)) {
247  std::string fname = "Tempus_" + stepperName + "_VanDerPol-Ref.dat";
248  if (n == 0) fname = "Tempus_" + stepperName + "_VanDerPol.dat";
249  RCP<const SolutionHistory<double>> solutionHistory =
250  integrator->getSolutionHistory();
251  writeSolution(fname, solutionHistory);
252  }
253  }
254 
255  // Check the order and intercept
256  double xSlope = 0.0;
257  double xDotSlope = 0.0;
258  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
259  // double order = stepper->getOrder();
260 
261  // xDot not yet available for DIRK methods, e.g., are not calc. and zero.
262  solutionsDot.clear();
263 
264  writeOrderError("Tempus_" + stepperName + "_VanDerPol-Error.dat", stepper,
265  StepSize, solutions, xErrorNorm, xSlope, solutionsDot,
266  xDotErrorNorm, xDotSlope, out);
267 
268  TEST_FLOATING_EQUALITY(xSlope, stepperOrders[m], 0.02);
269  TEST_FLOATING_EQUALITY(xErrorNorm[0], stepperErrors[m], 1.0e-4);
270  // TEST_FLOATING_EQUALITY( xDotSlope, 1.74898, 0.02 );
271  // TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 1.0038, 1.0e-4 );
272  }
274 }
275 
276 } // namespace Tempus_Test
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 the partitioned IMEX-RK.
T & get(const std::string &name, T def_value)
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)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
#define TEST_ASSERT(v1)
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar >> solutionHistory)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Partitioned Implicit-Explicit Runge-Kutta (IMEX-RK) time stepper.
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
ModelEvaluator pair for implicit and explicit (IMEX) evaulations.
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...
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Keep a fix number of states.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Solution state for integrators and steppers.