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