Tempus  Version of the Day
Time Integration
Tempus_TrapezoidalTest.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 #include "Teuchos_DefaultComm.hpp"
13 
14 #include "Tempus_config.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
16 #include "Tempus_StepperTrapezoidal.hpp"
17 
18 #include "../TestModels/SinCosModel.hpp"
19 #include "../TestModels/VanDerPolModel.hpp"
20 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
21 
22 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
23 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
24 
25 #ifdef Tempus_ENABLE_MPI
26 #include "Epetra_MpiComm.h"
27 #else
28 #include "Epetra_SerialComm.h"
29 #endif
30 
31 #include <vector>
32 #include <fstream>
33 #include <sstream>
34 #include <limits>
35 
36 namespace Tempus_Test {
37 
38 using Teuchos::RCP;
39 using Teuchos::rcp;
40 using Teuchos::rcp_const_cast;
41 using Teuchos::ParameterList;
42 using Teuchos::sublist;
43 using Teuchos::getParametersFromXmlFile;
44 
48 
49 
50 // ************************************************************
51 // ************************************************************
52 TEUCHOS_UNIT_TEST(Trapezoidal, ParameterList)
53 {
54  // Read params from .xml file
55  RCP<ParameterList> pList =
56  getParametersFromXmlFile("Tempus_Trapezoidal_SinCos.xml");
57 
58  // Setup the SinCosModel
59  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
60  auto model = rcp(new SinCosModel<double> (scm_pl));
61 
62  RCP<ParameterList> tempusPL = sublist(pList, "Tempus", true);
63 
64  // Test constructor IntegratorBasic(tempusPL, model)
65  {
66  RCP<Tempus::IntegratorBasic<double> > integrator =
67  Tempus::integratorBasic<double>(tempusPL, model);
68 
69  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
70 
71  RCP<const ParameterList> defaultPL =
72  integrator->getStepper()->getValidParameters();
73  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
74  if (!pass) {
75  std::cout << std::endl;
76  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
77  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
78  }
79  TEST_ASSERT(pass)
80  }
81 
82  // Test constructor IntegratorBasic(model, stepperType)
83  {
84  RCP<Tempus::IntegratorBasic<double> > integrator =
85  Tempus::integratorBasic<double>(model, "Trapezoidal Method");
86 
87  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
88  RCP<const ParameterList> defaultPL =
89  integrator->getStepper()->getValidParameters();
90 
91  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
92  if (!pass) {
93  std::cout << std::endl;
94  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
95  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
96  }
97  TEST_ASSERT(pass)
98  }
99 }
100 
101 
102 // ************************************************************
103 // ************************************************************
104 TEUCHOS_UNIT_TEST(Trapezoidal, ConstructingFromDefaults)
105 {
106  double dt = 0.1;
107 
108  // Read params from .xml file
109  RCP<ParameterList> pList =
110  getParametersFromXmlFile("Tempus_Trapezoidal_SinCos.xml");
111  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
112 
113  // Setup the SinCosModel
114  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
115  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
116  auto model = rcp(new SinCosModel<double>(scm_pl));
117 
118  // Setup Stepper for field solve ----------------------------
119  auto stepper = rcp(new Tempus::StepperTrapezoidal<double>());
120  //{
121  // // Turn on NOX output.
122  // RCP<ParameterList> sPL = stepper->getNonconstParameterList();
123  // std::string solverName = sPL->get<std::string>("Solver Name");
124  // RCP<ParameterList> solverPL = Teuchos::sublist(sPL, solverName, true);
125  // solverPL->sublist("NOX").sublist("Printing").sublist("Output Information")
126  // .set("Outer Iteration", true);
127  // solverPL->sublist("NOX").sublist("Printing").sublist("Output Information")
128  // .set("Parameters", true);
129  // solverPL->sublist("NOX").sublist("Printing").sublist("Output Information")
130  // .set("Details", true);
131  // stepper->setSolver(solverPL);
132  //}
133  stepper->setModel(model);
134  stepper->initialize();
135 
136  // Setup TimeStepControl ------------------------------------
137  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
138  ParameterList tscPL = pl->sublist("Default Integrator")
139  .sublist("Time Step Control");
140  timeStepControl->setStepType (tscPL.get<std::string>("Integrator Step Type"));
141  timeStepControl->setInitIndex(tscPL.get<int> ("Initial Time Index"));
142  timeStepControl->setInitTime (tscPL.get<double>("Initial Time"));
143  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
144  timeStepControl->setInitTimeStep(dt);
145  timeStepControl->initialize();
146 
147  // Setup initial condition SolutionState --------------------
148  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
149  stepper->getModel()->getNominalValues();
150  auto icSoln = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
151  auto icSolnDot =
152  rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x_dot());
153  auto icState = Tempus::createSolutionStateX(icSoln,icSolnDot);
154  icState->setTime (timeStepControl->getInitTime());
155  icState->setIndex (timeStepControl->getInitIndex());
156  icState->setTimeStep(0.0);
157  icState->setOrder (stepper->getOrder());
158  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
159 
160  // Setup SolutionHistory ------------------------------------
161  auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
162  solutionHistory->setName("Forward States");
163  solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
164  solutionHistory->setStorageLimit(2);
165  solutionHistory->addState(icState);
166 
167  // Setup Integrator -----------------------------------------
168  RCP<Tempus::IntegratorBasic<double> > integrator =
169  Tempus::integratorBasic<double>();
170  integrator->setStepperWStepper(stepper);
171  integrator->setTimeStepControl(timeStepControl);
172  integrator->setSolutionHistory(solutionHistory);
173  //integrator->setObserver(...);
174  integrator->initialize();
175 
176 
177  // Integrate to timeMax
178  bool integratorStatus = integrator->advanceTime();
179  TEST_ASSERT(integratorStatus)
180 
181 
182  // Test if at 'Final Time'
183  double time = integrator->getTime();
184  double timeFinal =pl->sublist("Default Integrator")
185  .sublist("Time Step Control").get<double>("Final Time");
186  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
187 
188  // Time-integrated solution and the exact solution
189  RCP<Thyra::VectorBase<double> > x = integrator->getX();
190  RCP<const Thyra::VectorBase<double> > x_exact =
191  model->getExactSolution(time).get_x();
192 
193  // Calculate the error
194  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
195  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
196 
197  // Check the order and intercept
198  std::cout << " Stepper = Trapezoidal" << std::endl;
199  std::cout << " =========================" << std::endl;
200  std::cout << " Exact solution : " << get_ele(*(x_exact), 0) << " "
201  << get_ele(*(x_exact), 1) << std::endl;
202  std::cout << " Computed solution: " << get_ele(*(x ), 0) << " "
203  << get_ele(*(x ), 1) << std::endl;
204  std::cout << " Difference : " << get_ele(*(xdiff ), 0) << " "
205  << get_ele(*(xdiff ), 1) << std::endl;
206  std::cout << " =========================" << std::endl;
207  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.841021, 1.0e-4 );
208  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.541002, 1.0e-4 );
209 }
210 
211 
212 // ************************************************************
213 // ************************************************************
214 TEUCHOS_UNIT_TEST(Trapezoidal, SinCos)
215 {
216  RCP<Tempus::IntegratorBasic<double> > integrator;
217  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
218  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
219  std::vector<double> StepSize;
220  std::vector<double> xErrorNorm;
221  std::vector<double> xDotErrorNorm;
222  const int nTimeStepSizes = 7;
223  double dt = 0.2;
224  double time = 0.0;
225  for (int n=0; n<nTimeStepSizes; n++) {
226 
227  // Read params from .xml file
228  RCP<ParameterList> pList =
229  getParametersFromXmlFile("Tempus_Trapezoidal_SinCos.xml");
230 
231  //std::ofstream ftmp("PL.txt");
232  //pList->print(ftmp);
233  //ftmp.close();
234 
235  // Setup the SinCosModel
236  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
237  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
238  auto model = rcp(new SinCosModel<double>(scm_pl));
239 
240  dt /= 2;
241 
242  // Setup the Integrator and reset initial time step
243  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
244  pl->sublist("Default Integrator")
245  .sublist("Time Step Control").set("Initial Time Step", dt);
246  integrator = Tempus::integratorBasic<double>(pl, model);
247 
248  // Initial Conditions
249  // During the Integrator construction, the initial SolutionState
250  // is set by default to model->getNominalVales().get_x(). However,
251  // the application can set it also by integrator->initializeSolutionHistory.
252  {
253  RCP<Thyra::VectorBase<double> > x0 =
254  model->getNominalValues().get_x()->clone_v();
255  RCP<Thyra::VectorBase<double> > xdot0 =
256  model->getNominalValues().get_x_dot()->clone_v();
257  integrator->initializeSolutionHistory(0.0, x0, xdot0);
258  }
259 
260  // Integrate to timeMax
261  bool integratorStatus = integrator->advanceTime();
262  TEST_ASSERT(integratorStatus)
263 
264  // Test if at 'Final Time'
265  time = integrator->getTime();
266  double timeFinal =pl->sublist("Default Integrator")
267  .sublist("Time Step Control").get<double>("Final Time");
268  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
269 
270  // Plot sample solution and exact solution
271  if (n == 0) {
272  RCP<const SolutionHistory<double> > solutionHistory =
273  integrator->getSolutionHistory();
274  writeSolution("Tempus_Trapezoidal_SinCos.dat", solutionHistory);
275 
276  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
277  for (int i=0; i<solutionHistory->getNumStates(); i++) {
278  double time_i = (*solutionHistory)[i]->getTime();
279  auto state = Tempus::createSolutionStateX(
280  rcp_const_cast<Thyra::VectorBase<double> > (
281  model->getExactSolution(time_i).get_x()),
282  rcp_const_cast<Thyra::VectorBase<double> > (
283  model->getExactSolution(time_i).get_x_dot()));
284  state->setTime((*solutionHistory)[i]->getTime());
285  solnHistExact->addState(state);
286  }
287  writeSolution("Tempus_Trapezoidal_SinCos-Ref.dat", solnHistExact);
288  }
289 
290  // Store off the final solution and step size
291  StepSize.push_back(dt);
292  auto solution = Thyra::createMember(model->get_x_space());
293  Thyra::copy(*(integrator->getX()),solution.ptr());
294  solutions.push_back(solution);
295  auto solutionDot = Thyra::createMember(model->get_x_space());
296  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
297  solutionsDot.push_back(solutionDot);
298  if (n == nTimeStepSizes-1) { // Add exact solution last in vector.
299  StepSize.push_back(0.0);
300  auto solutionExact = Thyra::createMember(model->get_x_space());
301  Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
302  solutions.push_back(solutionExact);
303  auto solutionDotExact = Thyra::createMember(model->get_x_space());
304  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
305  solutionDotExact.ptr());
306  solutionsDot.push_back(solutionDotExact);
307  }
308  }
309 
310  double xSlope = 0.0;
311  double xDotSlope = 0.0;
312  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
313  double order = stepper->getOrder();
314  writeOrderError("Tempus_Trapezoidal_SinCos-Error.dat",
315  stepper, StepSize,
316  solutions, xErrorNorm, xSlope,
317  solutionsDot, xDotErrorNorm, xDotSlope);
318 
319  TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
320  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.000832086, 1.0e-4 );
321  TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
322  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.000832086, 1.0e-4 );
323 
324  Teuchos::TimeMonitor::summarize();
325 }
326 
327 
328 // ************************************************************
329 // ************************************************************
330 TEUCHOS_UNIT_TEST(Trapezoidal, VanDerPol)
331 {
332  RCP<Tempus::IntegratorBasic<double> > integrator;
333  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
334  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
335  std::vector<double> StepSize;
336  std::vector<double> xErrorNorm;
337  std::vector<double> xDotErrorNorm;
338  const int nTimeStepSizes = 4;
339  double dt = 0.05;
340  double time = 0.0;
341  for (int n=0; n<nTimeStepSizes; n++) {
342 
343  // Read params from .xml file
344  RCP<ParameterList> pList =
345  getParametersFromXmlFile("Tempus_Trapezoidal_VanDerPol.xml");
346 
347  // Setup the VanDerPolModel
348  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
349  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
350 
351  // Set the step size
352  dt /= 2;
353  if (n == nTimeStepSizes-1) dt /= 10.0;
354 
355  // Setup the Integrator and reset initial time step
356  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
357  pl->sublist("Demo Integrator")
358  .sublist("Time Step Control").set("Initial Time Step", dt);
359  integrator = Tempus::integratorBasic<double>(pl, model);
360 
361  // Integrate to timeMax
362  bool integratorStatus = integrator->advanceTime();
363  TEST_ASSERT(integratorStatus)
364 
365  // Test if at 'Final Time'
366  time = integrator->getTime();
367  double timeFinal =pl->sublist("Demo Integrator")
368  .sublist("Time Step Control").get<double>("Final Time");
369  double tol = 100.0 * std::numeric_limits<double>::epsilon();
370  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
371 
372  // Store off the final solution and step size
373  StepSize.push_back(dt);
374  auto solution = Thyra::createMember(model->get_x_space());
375  Thyra::copy(*(integrator->getX()),solution.ptr());
376  solutions.push_back(solution);
377  auto solutionDot = Thyra::createMember(model->get_x_space());
378  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
379  solutionsDot.push_back(solutionDot);
380 
381  // Output finest temporal solution for plotting
382  // This only works for ONE MPI process
383  if ((n == 0) or (n == nTimeStepSizes-1)) {
384  std::string fname = "Tempus_Trapezoidal_VanDerPol-Ref.dat";
385  if (n == 0) fname = "Tempus_Trapezoidal_VanDerPol.dat";
386  RCP<const SolutionHistory<double> > solutionHistory =
387  integrator->getSolutionHistory();
388  writeSolution(fname, solutionHistory);
389  }
390  }
391  // Check the order and intercept
392  double xSlope = 0.0;
393  double xDotSlope = 0.0;
394  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
395  double order = stepper->getOrder();
396  writeOrderError("Tempus_Trapezoidal_VanDerPol-Error.dat",
397  stepper, StepSize,
398  solutions, xErrorNorm, xSlope,
399  solutionsDot, xDotErrorNorm, xDotSlope);
400 
401  TEST_FLOATING_EQUALITY( xSlope, order, 0.10 );
402  TEST_FLOATING_EQUALITY( xDotSlope, order, 0.10 );//=order at samll dt
403  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.00085855, 1.0e-4 );
404  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.00120695, 1.0e-4 );
405 
406  Teuchos::TimeMonitor::summarize();
407 }
408 
409 
410 } // 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...
Trapezoidal method time stepper.
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
van der Pol model problem for nonlinear electrical circuit.
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.