Tempus  Version of the Day
Time Integration
Tempus_OperatorSplitTest.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"
17 #include "Tempus_StepperOperatorSplit.hpp"
18 #include "Tempus_StepperForwardEuler.hpp"
19 #include "Tempus_StepperBackwardEuler.hpp"
20 
21 #include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
22 #include "../TestModels/VanDerPol_IMEX_ImplicitModel.hpp"
23 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
24 
25 #include <fstream>
26 #include <vector>
27 
28 namespace Tempus_Test {
29 
30 using Teuchos::RCP;
31 using Teuchos::rcp;
32 using Teuchos::rcp_const_cast;
33 using Teuchos::ParameterList;
34 using Teuchos::sublist;
35 using Teuchos::getParametersFromXmlFile;
36 
40 
41 
42 // ************************************************************
43 // ************************************************************
44 TEUCHOS_UNIT_TEST(OperatorSplit, ConstructingFromDefaults)
45 {
46  double dt = 0.05;
47 
48  // Read params from .xml file
49  RCP<ParameterList> pList =
50  getParametersFromXmlFile("Tempus_OperatorSplit_VanDerPol.xml");
51  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
52 
53  // Setup the explicit VanDerPol ModelEvaluator
54  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
55  auto explicitModel = rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL));
56 
57  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
58  auto implicitModel = rcp(new VanDerPol_IMEX_ImplicitModel<double>(vdpmPL));
59 
60 
61  // Setup Stepper for field solve ----------------------------
62  auto stepper = rcp(new Tempus::StepperOperatorSplit<double>());
63 
64  RCP<Tempus::StepperFactory<double> > sf =
65  Teuchos::rcp(new Tempus::StepperFactory<double>());
66 
67  auto subStepper1 =
68  sf->createStepperForwardEuler(explicitModel, Teuchos::null);
69 
70  auto subStepper2 =
71  sf->createStepperBackwardEuler(implicitModel, Teuchos::null);
72 
73  stepper->addStepper(subStepper1);
74  stepper->addStepper(subStepper2);
75  stepper->initialize();
76 
77 
78  // Setup TimeStepControl ------------------------------------
79  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
80  ParameterList tscPL = pl->sublist("Demo Integrator")
81  .sublist("Time Step Control");
82  timeStepControl->setStepType (tscPL.get<std::string>("Integrator Step Type"));
83  timeStepControl->setInitIndex(tscPL.get<int> ("Initial Time Index"));
84  timeStepControl->setInitTime (tscPL.get<double>("Initial Time"));
85  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
86  timeStepControl->setInitTimeStep(dt);
87  timeStepControl->initialize();
88 
89  // Setup initial condition SolutionState --------------------
90  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
91  stepper->getModel()->getNominalValues();
92  auto icX = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
93  auto icXDot = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x_dot());
94  auto icState = Tempus::createSolutionStateX(icX, icXDot);
95  icState->setTime (timeStepControl->getInitTime());
96  icState->setIndex (timeStepControl->getInitIndex());
97  icState->setTimeStep(0.0);
98  icState->setOrder (stepper->getOrder());
99  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
100 
101  // Setup SolutionHistory ------------------------------------
102  auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
103  solutionHistory->setName("Forward States");
104  solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
105  solutionHistory->setStorageLimit(2);
106  solutionHistory->addState(icState);
107 
108  // Setup Integrator -----------------------------------------
109  RCP<Tempus::IntegratorBasic<double> > integrator =
110  Tempus::integratorBasic<double>();
111  integrator->setStepperWStepper(stepper);
112  integrator->setTimeStepControl(timeStepControl);
113  integrator->setSolutionHistory(solutionHistory);
114  //integrator->setObserver(...);
115  integrator->initialize();
116 
117 
118  // Integrate to timeMax
119  bool integratorStatus = integrator->advanceTime();
120  TEST_ASSERT(integratorStatus)
121 
122 
123  // Test if at 'Final Time'
124  double time = integrator->getTime();
125  double timeFinal =pl->sublist("Demo Integrator")
126  .sublist("Time Step Control").get<double>("Final Time");
127  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
128 
129  // Time-integrated solution and the exact solution
130  RCP<Thyra::VectorBase<double> > x = integrator->getX();
131 
132  // Check the order and intercept
133  std::cout << " Stepper = " << stepper->description() << std::endl;
134  std::cout << " =========================" << std::endl;
135  std::cout << " Computed solution: " << get_ele(*(x ), 0) << " "
136  << get_ele(*(x ), 1) << std::endl;
137  std::cout << " =========================" << std::endl;
138  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), -2.223910, 1.0e-4);
139  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.565441, 1.0e-4);
140 }
141 
142 
143 // ************************************************************
144 // ************************************************************
145 TEUCHOS_UNIT_TEST(OperatorSplit, VanDerPol)
146 {
147  RCP<Tempus::IntegratorBasic<double> > integrator;
148  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
149  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
150  std::vector<double> StepSize;
151  std::vector<double> xErrorNorm;
152  std::vector<double> xDotErrorNorm;
153  const int nTimeStepSizes = 4; // 8 for Error plot
154  double dt = 0.1;
155  double time = 0.0;
156  for (int n=0; n<nTimeStepSizes; n++) {
157 
158  // Read params from .xml file
159  RCP<ParameterList> pList =
160  getParametersFromXmlFile("Tempus_OperatorSplit_VanDerPol.xml");
161 
162  // Setup the explicit VanDerPol ModelEvaluator
163  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
164  auto explicitModel = rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL));
165 
166  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
167  auto implicitModel = rcp(new VanDerPol_IMEX_ImplicitModel<double>(vdpmPL));
168 
169  // Setup vector of models
170  std::vector<RCP<const Thyra::ModelEvaluator<double> > > models;
171  models.push_back(explicitModel);
172  models.push_back(implicitModel);
173 
174  // Set the step size
175  dt /= 2;
176  if (n == nTimeStepSizes-1) dt /= 10.0;
177 
178  // Setup the Integrator and reset initial time step
179  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
180  pl->sublist("Demo Integrator")
181  .sublist("Time Step Control").set("Initial Time Step", dt);
182  integrator = Tempus::integratorBasic<double>(pl, models);
183 
184  // Integrate to timeMax
185  bool integratorStatus = integrator->advanceTime();
186  TEST_ASSERT(integratorStatus)
187 
188  // Test if at 'Final Time'
189  time = integrator->getTime();
190  double timeFinal =pl->sublist("Demo Integrator")
191  .sublist("Time Step Control").get<double>("Final Time");
192  double tol = 100.0 * std::numeric_limits<double>::epsilon();
193  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
194 
195  // Store off the final solution and step size
196  StepSize.push_back(dt);
197  auto solution = Thyra::createMember(implicitModel->get_x_space());
198  Thyra::copy(*(integrator->getX()),solution.ptr());
199  solutions.push_back(solution);
200  auto solutionDot = Thyra::createMember(implicitModel->get_x_space());
201  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
202  solutionsDot.push_back(solutionDot);
203 
204  // Output finest temporal solution for plotting
205  // This only works for ONE MPI process
206  if ((n == 0) or (n == nTimeStepSizes-1)) {
207  std::string fname = "Tempus_OperatorSplit_VanDerPol-Ref.dat";
208  if (n == 0) fname = "Tempus_OperatorSplit_VanDerPol.dat";
209  RCP<const SolutionHistory<double> > solutionHistory =
210  integrator->getSolutionHistory();
211  writeSolution(fname, solutionHistory);
212  //solutionHistory->printHistory("medium");
213  }
214  }
215 
216  // Check the order and intercept
217  double xSlope = 0.0;
218  double xDotSlope = 0.0;
219  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
220  double order = stepper->getOrder();
221  writeOrderError("Tempus_OperatorSplit_VanDerPol-Error.dat",
222  stepper, StepSize,
223  solutions, xErrorNorm, xSlope,
224  solutionsDot, xDotErrorNorm, xDotSlope);
225 
226  TEST_FLOATING_EQUALITY( xSlope, order, 0.05 );
227  TEST_FLOATING_EQUALITY( xDotSlope, order, 0.05 );//=order at small dt
228  TEST_FLOATING_EQUALITY( xErrorNorm[0], 1.27294, 1.0e-4 );
229  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 12.7102, 1.0e-4 );
230 
231  Teuchos::TimeMonitor::summarize();
232 }
233 
234 
235 } // 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...
OperatorSplit stepper loops through the Stepper list.
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
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.