Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_IMEX_RKTest.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_WrapperModelEvaluatorPairIMEX_Basic.hpp"
17 #include "Tempus_StepperIMEX_RK.hpp"
18 
19 #include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
20 #include "../TestModels/VanDerPol_IMEX_ImplicitModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
22 
23 #include <fstream>
24 #include <vector>
25 
26 namespace Tempus_Test {
27 
28 using Teuchos::RCP;
29 using Teuchos::rcp;
30 using Teuchos::rcp_const_cast;
31 using Teuchos::ParameterList;
32 using Teuchos::sublist;
33 using Teuchos::getParametersFromXmlFile;
34 
38 
39 // Comment out any of the following tests to exclude from build/run.
40 #define TEST_CONSTRUCTING_FROM_DEFAULTS
41 #define TEST_VANDERPOL
42 
43 
44 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
45 // ************************************************************
46 // ************************************************************
47 TEUCHOS_UNIT_TEST(IMEX_RK, ConstructingFromDefaults)
48 {
49  double dt = 0.025;
50 
51  // Read params from .xml file
52  RCP<ParameterList> pList =
53  getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
54  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
55 
56  // Setup the explicit VanDerPol ModelEvaluator
57  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
58  auto explicitModel = rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL));
59 
60  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
61  auto implicitModel = rcp(new VanDerPol_IMEX_ImplicitModel<double>(vdpmPL));
62 
63  // Setup the IMEX Pair ModelEvaluator
65  explicitModel, implicitModel));
66 
67 
68  // Setup Stepper for field solve ----------------------------
69  auto stepper = rcp(new Tempus::StepperIMEX_RK<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 = rcp(new Tempus::SolutionState<double>(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 ------------------------------------
97  solutionHistory->setName("Forward States");
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 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
136 
137 
138 #ifdef TEST_VANDERPOL
139 // ************************************************************
140 // ************************************************************
141 TEUCHOS_UNIT_TEST(IMEX_RK, VanDerPol)
142 {
143  std::vector<std::string> stepperTypes;
144  stepperTypes.push_back("IMEX RK 1st order");
145  stepperTypes.push_back("IMEX RK SSP2" );
146  stepperTypes.push_back("IMEX RK ARS 233" );
147  stepperTypes.push_back("General IMEX RK" );
148 
149  std::vector<double> stepperOrders;
150  stepperOrders.push_back(1.07964);
151  stepperOrders.push_back(2.00408);
152  stepperOrders.push_back(2.70655);
153  stepperOrders.push_back(2.00211);
154 
155  std::vector<double> stepperErrors;
156  stepperErrors.push_back(0.0046423);
157  stepperErrors.push_back(0.0154534);
158  stepperErrors.push_back(0.000298908);
159  stepperErrors.push_back(0.0071546);
160 
161  std::vector<double> stepperInitDt;
162  stepperInitDt.push_back(0.0125);
163  stepperInitDt.push_back(0.05);
164  stepperInitDt.push_back(0.05);
165  stepperInitDt.push_back(0.05);
166 
167  std::vector<std::string>::size_type m;
168  for(m = 0; m != stepperTypes.size(); m++) {
169 
170  std::string stepperType = stepperTypes[m];
171  std::string stepperName = stepperTypes[m];
172  std::replace(stepperName.begin(), stepperName.end(), ' ', '_');
173  std::replace(stepperName.begin(), stepperName.end(), '/', '.');
174 
175  RCP<Tempus::IntegratorBasic<double> > integrator;
176  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
177  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
178  std::vector<double> StepSize;
179  std::vector<double> xErrorNorm;
180  std::vector<double> xDotErrorNorm;
181 
182  const int nTimeStepSizes = 3; // 6 for error plot
183  double dt = stepperInitDt[m];
184  double time = 0.0;
185  for (int n=0; n<nTimeStepSizes; n++) {
186 
187  // Read params from .xml file
188  RCP<ParameterList> pList =
189  getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
190 
191  // Setup the explicit VanDerPol ModelEvaluator
192  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
193  auto explicitModel = rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL));
194 
195  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
196  auto implicitModel = rcp(new VanDerPol_IMEX_ImplicitModel<double>(vdpmPL));
197 
198  // Setup the IMEX Pair ModelEvaluator
200  explicitModel, implicitModel));
201 
202  // Set the Stepper
203  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
204  if (stepperType == "General IMEX RK"){
205  // use the appropriate stepper sublist
206  pl->sublist("Default Integrator").set("Stepper Name", "General IMEX RK");
207  } else {
208  pl->sublist("Default Stepper").set("Stepper Type", stepperType);
209  }
210 
211  // Set the step size
212  if (n == nTimeStepSizes-1) dt /= 10.0;
213  else dt /= 2;
214 
215  // Setup the Integrator and reset initial time step
216  pl->sublist("Default Integrator")
217  .sublist("Time Step Control").set("Initial Time Step", dt);
218  integrator = Tempus::integratorBasic<double>(pl, model);
219 
220  // Integrate to timeMax
221  bool integratorStatus = integrator->advanceTime();
222  TEST_ASSERT(integratorStatus)
223 
224  // Test if at 'Final Time'
225  time = integrator->getTime();
226  double timeFinal =pl->sublist("Default Integrator")
227  .sublist("Time Step Control").get<double>("Final Time");
228  double tol = 100.0 * std::numeric_limits<double>::epsilon();
229  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
230 
231  // Store off the final solution and step size
232  StepSize.push_back(dt);
233  auto solution = Thyra::createMember(model->get_x_space());
234  Thyra::copy(*(integrator->getX()),solution.ptr());
235  solutions.push_back(solution);
236  auto solutionDot = Thyra::createMember(model->get_x_space());
237  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
238  solutionsDot.push_back(solutionDot);
239 
240  // Output finest temporal solution for plotting
241  // This only works for ONE MPI process
242  if ((n == 0) or (n == nTimeStepSizes-1)) {
243  std::string fname = "Tempus_"+stepperName+"_VanDerPol-Ref.dat";
244  if (n == 0) fname = "Tempus_"+stepperName+"_VanDerPol.dat";
245  RCP<const SolutionHistory<double> > solutionHistory =
246  integrator->getSolutionHistory();
247  writeSolution(fname, solutionHistory);
248  }
249  }
250 
251  // Check the order and intercept
252  double xSlope = 0.0;
253  double xDotSlope = 0.0;
254  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
255  //double order = stepper->getOrder();
256  writeOrderError("Tempus_"+stepperName+"_VanDerPol-Error.dat",
257  stepper, StepSize,
258  solutions, xErrorNorm, xSlope,
259  solutionsDot, xDotErrorNorm, xDotSlope);
260 
261  TEST_FLOATING_EQUALITY( xSlope, stepperOrders[m], 0.02 );
262  TEST_FLOATING_EQUALITY( xErrorNorm[0], stepperErrors[m], 1.0e-4 );
263  // xDot not yet available for IMEX_RK.
264  //TEST_FLOATING_EQUALITY( xDotSlope, 1.74898, 0.10 );
265  //TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 1.0038, 1.0e-4 );
266 
267  }
268  //Teuchos::TimeMonitor::summarize();
269 }
270 #endif // TEST_VANDERPOL
271 
272 
273 } // namespace Tempus_Test
ModelEvaluator pair for implicit and explicit (IMEX) evaulations.
Implicit-Explicit Runge-Kutta (IMEX-RK) time stepper.
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)
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Keep a fix number of states.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...