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