Tempus  Version of the Day
Time Integration
Tempus_BackwardEulerTest.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_StepperBackwardEuler.hpp"
17 
18 #include "../TestModels/SinCosModel.hpp"
19 #include "../TestModels/CDR_Model.hpp"
20 #include "../TestModels/VanDerPolModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
22 
23 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
24 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
25 
26 #ifdef Tempus_ENABLE_MPI
27 #include "Epetra_MpiComm.h"
28 #else
29 #include "Epetra_SerialComm.h"
30 #endif
31 
32 #include <vector>
33 #include <fstream>
34 #include <sstream>
35 #include <limits>
36 
37 namespace Tempus_Test {
38 
39 using Teuchos::RCP;
40 using Teuchos::rcp;
41 using Teuchos::rcp_const_cast;
42 using Teuchos::ParameterList;
43 using Teuchos::sublist;
44 using Teuchos::getParametersFromXmlFile;
45 
49 
50 
51 // ************************************************************
52 // ************************************************************
53 TEUCHOS_UNIT_TEST(BackwardEuler, ParameterList)
54 {
55  // Read params from .xml file
56  RCP<ParameterList> pList =
57  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
58 
59  // Setup the SinCosModel
60  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
61  auto model = rcp(new SinCosModel<double> (scm_pl));
62 
63  RCP<ParameterList> tempusPL = sublist(pList, "Tempus", true);
64 
65  // Test constructor IntegratorBasic(tempusPL, model)
66  {
67  RCP<Tempus::IntegratorBasic<double> > integrator =
68  Tempus::integratorBasic<double>(tempusPL, model);
69 
70  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
71  // Match Predictor for comparison
72  stepperPL->set("Predictor Stepper Type", "None");
73  RCP<const ParameterList> defaultPL =
74  integrator->getStepper()->getValidParameters();
75 
76  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
77  if (!pass) {
78  std::cout << std::endl;
79  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
80  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
81  }
82  TEST_ASSERT(pass)
83  }
84 
85  // Test constructor IntegratorBasic(model, stepperType)
86  {
87  RCP<Tempus::IntegratorBasic<double> > integrator =
88  Tempus::integratorBasic<double>(model, "Backward Euler");
89 
90  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
91  RCP<const ParameterList> defaultPL =
92  integrator->getStepper()->getValidParameters();
93 
94  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
95  if (!pass) {
96  std::cout << std::endl;
97  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
98  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
99  }
100  TEST_ASSERT(pass)
101  }
102 }
103 
104 
105 // ************************************************************
106 // ************************************************************
107 TEUCHOS_UNIT_TEST(BackwardEuler, ConstructingFromDefaults)
108 {
109  double dt = 0.1;
110 
111  // Read params from .xml file
112  RCP<ParameterList> pList =
113  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
114  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
115 
116  // Setup the SinCosModel
117  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
118  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
119  auto model = rcp(new SinCosModel<double>(scm_pl));
120 
121  // Setup Stepper for field solve ----------------------------
122  auto stepper = rcp(new Tempus::StepperBackwardEuler<double>());
123  stepper->setModel(model);
124  stepper->initialize();
125 
126  // Setup TimeStepControl ------------------------------------
127  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
128  ParameterList tscPL = pl->sublist("Default Integrator")
129  .sublist("Time Step Control");
130  timeStepControl->setStepType (tscPL.get<std::string>("Integrator Step Type"));
131  timeStepControl->setInitIndex(tscPL.get<int> ("Initial Time Index"));
132  timeStepControl->setInitTime (tscPL.get<double>("Initial Time"));
133  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
134  timeStepControl->setInitTimeStep(dt);
135  timeStepControl->initialize();
136 
137  // Setup initial condition SolutionState --------------------
138  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
139  stepper->getModel()->getNominalValues();
140  auto icSolution =
141  rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
142  auto icState = Tempus::createSolutionStateX(icSolution);
143  icState->setTime (timeStepControl->getInitTime());
144  icState->setIndex (timeStepControl->getInitIndex());
145  icState->setTimeStep(0.0);
146  icState->setOrder (stepper->getOrder());
147  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
148 
149  // Setup SolutionHistory ------------------------------------
150  auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
151  solutionHistory->setName("Forward States");
152  solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
153  solutionHistory->setStorageLimit(2);
154  solutionHistory->addState(icState);
155 
156  // Setup Integrator -----------------------------------------
157  RCP<Tempus::IntegratorBasic<double> > integrator =
158  Tempus::integratorBasic<double>();
159  integrator->setStepperWStepper(stepper);
160  integrator->setTimeStepControl(timeStepControl);
161  integrator->setSolutionHistory(solutionHistory);
162  //integrator->setObserver(...);
163  integrator->initialize();
164 
165 
166  // Integrate to timeMax
167  bool integratorStatus = integrator->advanceTime();
168  TEST_ASSERT(integratorStatus)
169 
170 
171  // Test if at 'Final Time'
172  double time = integrator->getTime();
173  double timeFinal =pl->sublist("Default Integrator")
174  .sublist("Time Step Control").get<double>("Final Time");
175  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
176 
177  // Time-integrated solution and the exact solution
178  RCP<Thyra::VectorBase<double> > x = integrator->getX();
179  RCP<const Thyra::VectorBase<double> > x_exact =
180  model->getExactSolution(time).get_x();
181 
182  // Calculate the error
183  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
184  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
185 
186  // Check the order and intercept
187  std::cout << " Stepper = BackwardEuler" << std::endl;
188  std::cout << " =========================" << std::endl;
189  std::cout << " Exact solution : " << get_ele(*(x_exact), 0) << " "
190  << get_ele(*(x_exact), 1) << std::endl;
191  std::cout << " Computed solution: " << get_ele(*(x ), 0) << " "
192  << get_ele(*(x ), 1) << std::endl;
193  std::cout << " Difference : " << get_ele(*(xdiff ), 0) << " "
194  << get_ele(*(xdiff ), 1) << std::endl;
195  std::cout << " =========================" << std::endl;
196  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.798923, 1.0e-4 );
197  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.516729, 1.0e-4 );
198 }
199 
200 
201 // ************************************************************
202 // ************************************************************
203 TEUCHOS_UNIT_TEST(BackwardEuler, SinCos)
204 {
205  RCP<Tempus::IntegratorBasic<double> > integrator;
206  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
207  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
208  std::vector<double> StepSize;
209  std::vector<double> xErrorNorm;
210  std::vector<double> xDotErrorNorm;
211  const int nTimeStepSizes = 7;
212  double dt = 0.2;
213  double time = 0.0;
214  for (int n=0; n<nTimeStepSizes; n++) {
215 
216  // Read params from .xml file
217  RCP<ParameterList> pList =
218  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
219 
220  //std::ofstream ftmp("PL.txt");
221  //pList->print(ftmp);
222  //ftmp.close();
223 
224  // Setup the SinCosModel
225  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
226  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
227  auto model = rcp(new SinCosModel<double>(scm_pl));
228 
229  dt /= 2;
230 
231  // Setup the Integrator and reset initial time step
232  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
233  pl->sublist("Default Integrator")
234  .sublist("Time Step Control").set("Initial Time Step", dt);
235  integrator = Tempus::integratorBasic<double>(pl, model);
236 
237  // Initial Conditions
238  // During the Integrator construction, the initial SolutionState
239  // is set by default to model->getNominalVales().get_x(). However,
240  // the application can set it also by integrator->initializeSolutionHistory.
241  RCP<Thyra::VectorBase<double> > x0 =
242  model->getNominalValues().get_x()->clone_v();
243  integrator->initializeSolutionHistory(0.0, x0);
244 
245  // Integrate to timeMax
246  bool integratorStatus = integrator->advanceTime();
247  TEST_ASSERT(integratorStatus)
248 
249  // Test if at 'Final Time'
250  time = integrator->getTime();
251  double timeFinal =pl->sublist("Default Integrator")
252  .sublist("Time Step Control").get<double>("Final Time");
253  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
254 
255  // Plot sample solution and exact solution
256  if (n == 0) {
257  RCP<const SolutionHistory<double> > solutionHistory =
258  integrator->getSolutionHistory();
259  writeSolution("Tempus_BackwardEuler_SinCos.dat", solutionHistory);
260 
261  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
262  for (int i=0; i<solutionHistory->getNumStates(); i++) {
263  double time_i = (*solutionHistory)[i]->getTime();
264  auto state = Tempus::createSolutionStateX(
265  rcp_const_cast<Thyra::VectorBase<double> > (
266  model->getExactSolution(time_i).get_x()),
267  rcp_const_cast<Thyra::VectorBase<double> > (
268  model->getExactSolution(time_i).get_x_dot()));
269  state->setTime((*solutionHistory)[i]->getTime());
270  solnHistExact->addState(state);
271  }
272  writeSolution("Tempus_BackwardEuler_SinCos-Ref.dat", solnHistExact);
273  }
274 
275  // Store off the final solution and step size
276  StepSize.push_back(dt);
277  auto solution = Thyra::createMember(model->get_x_space());
278  Thyra::copy(*(integrator->getX()),solution.ptr());
279  solutions.push_back(solution);
280  auto solutionDot = Thyra::createMember(model->get_x_space());
281  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
282  solutionsDot.push_back(solutionDot);
283  if (n == nTimeStepSizes-1) { // Add exact solution last in vector.
284  StepSize.push_back(0.0);
285  auto solutionExact = Thyra::createMember(model->get_x_space());
286  Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
287  solutions.push_back(solutionExact);
288  auto solutionDotExact = Thyra::createMember(model->get_x_space());
289  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
290  solutionDotExact.ptr());
291  solutionsDot.push_back(solutionDotExact);
292  }
293  }
294 
295  // Check the order and intercept
296  double xSlope = 0.0;
297  double xDotSlope = 0.0;
298  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
299  double order = stepper->getOrder();
300  writeOrderError("Tempus_BackwardEuler_SinCos-Error.dat",
301  stepper, StepSize,
302  solutions, xErrorNorm, xSlope,
303  solutionsDot, xDotErrorNorm, xDotSlope);
304 
305  TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
306  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.0486418, 1.0e-4 );
307  TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
308  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.0486418, 1.0e-4 );
309 
310  Teuchos::TimeMonitor::summarize();
311 }
312 
313 
314 // ************************************************************
315 // ************************************************************
316 TEUCHOS_UNIT_TEST(BackwardEuler, CDR)
317 {
318  // Create a communicator for Epetra objects
319  RCP<Epetra_Comm> comm;
320 #ifdef Tempus_ENABLE_MPI
321  comm = rcp(new Epetra_MpiComm(MPI_COMM_WORLD));
322 #else
323  comm = rcp(new Epetra_SerialComm);
324 #endif
325 
326  RCP<Tempus::IntegratorBasic<double> > integrator;
327  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
328  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
329  std::vector<double> StepSize;
330  std::vector<double> xErrorNorm;
331  std::vector<double> xDotErrorNorm;
332  const int nTimeStepSizes = 5;
333  double dt = 0.2;
334  for (int n=0; n<nTimeStepSizes; n++) {
335 
336  // Read params from .xml file
337  RCP<ParameterList> pList =
338  getParametersFromXmlFile("Tempus_BackwardEuler_CDR.xml");
339 
340  // Create CDR Model
341  RCP<ParameterList> model_pl = sublist(pList, "CDR Model", true);
342  const int num_elements = model_pl->get<int>("num elements");
343  const double left_end = model_pl->get<double>("left end");
344  const double right_end = model_pl->get<double>("right end");
345  const double a_convection = model_pl->get<double>("a (convection)");
346  const double k_source = model_pl->get<double>("k (source)");
347 
348  auto model = rcp(new Tempus_Test::CDR_Model<double>(comm,
349  num_elements,
350  left_end,
351  right_end,
352  a_convection,
353  k_source));
354 
355  // Set the factory
356  ::Stratimikos::DefaultLinearSolverBuilder builder;
357 
358  auto p = rcp(new ParameterList);
359  p->set("Linear Solver Type", "Belos");
360  p->set("Preconditioner Type", "None");
361  builder.setParameterList(p);
362 
363  RCP< ::Thyra::LinearOpWithSolveFactoryBase<double> >
364  lowsFactory = builder.createLinearSolveStrategy("");
365 
366  model->set_W_factory(lowsFactory);
367 
368  // Set the step size
369  dt /= 2;
370 
371  // Setup the Integrator and reset initial time step
372  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
373  pl->sublist("Demo Integrator")
374  .sublist("Time Step Control").set("Initial Time Step", dt);
375  integrator = Tempus::integratorBasic<double>(pl, model);
376 
377  // Integrate to timeMax
378  bool integratorStatus = integrator->advanceTime();
379  TEST_ASSERT(integratorStatus)
380 
381  // Test if at 'Final Time'
382  double time = integrator->getTime();
383  double timeFinal =pl->sublist("Demo Integrator")
384  .sublist("Time Step Control").get<double>("Final Time");
385  double tol = 100.0 * std::numeric_limits<double>::epsilon();
386  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
387 
388  // Store off the final solution and step size
389  StepSize.push_back(dt);
390  auto solution = Thyra::createMember(model->get_x_space());
391  Thyra::copy(*(integrator->getX()),solution.ptr());
392  solutions.push_back(solution);
393  auto solutionDot = Thyra::createMember(model->get_x_space());
394  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
395  solutionsDot.push_back(solutionDot);
396 
397  // Output finest temporal solution for plotting
398  // This only works for ONE MPI process
399  if ((n == nTimeStepSizes-1) && (comm->NumProc() == 1)) {
400  std::ofstream ftmp("Tempus_BackwardEuler_CDR.dat");
401  ftmp << "TITLE=\"Backward Euler Solution to CDR\"\n"
402  << "VARIABLES=\"z\",\"T\"\n";
403  const double dx = std::fabs(left_end-right_end) /
404  static_cast<double>(num_elements);
405  RCP<const SolutionHistory<double> > solutionHistory =
406  integrator->getSolutionHistory();
407  int nStates = solutionHistory->getNumStates();
408  for (int i=0; i<nStates; i++) {
409  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
410  RCP<const Thyra::VectorBase<double> > x = solutionState->getX();
411  double ttime = solutionState->getTime();
412  ftmp << "ZONE T=\"Time="<<ttime<<"\", I="
413  <<num_elements+1<<", F=BLOCK\n";
414  for (int j = 0; j < num_elements+1; j++) {
415  const double x_coord = left_end + static_cast<double>(j) * dx;
416  ftmp << x_coord << " ";
417  }
418  ftmp << std::endl;
419  for (int j=0; j<num_elements+1; j++) ftmp << get_ele(*x, j) << " ";
420  ftmp << std::endl;
421  }
422  ftmp.close();
423  }
424  }
425 
426  // Check the order and intercept
427  double xSlope = 0.0;
428  double xDotSlope = 0.0;
429  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
430  writeOrderError("Tempus_BackwardEuler_CDR-Error.dat",
431  stepper, StepSize,
432  solutions, xErrorNorm, xSlope,
433  solutionsDot, xDotErrorNorm, xDotSlope);
434 
435  TEST_FLOATING_EQUALITY( xSlope, 1.32213, 0.01 );
436  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.116919, 1.0e-4 );
437  TEST_FLOATING_EQUALITY( xDotSlope, 1.32052, 0.01 );
438  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.449888, 1.0e-4 );
439  // At small dt, slopes should be equal to order.
440  //double order = stepper->getOrder();
441  //TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
442  //TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
443 
444  // Write fine mesh solution at final time
445  // This only works for ONE MPI process
446  if (comm->NumProc() == 1) {
447  RCP<ParameterList> pList =
448  getParametersFromXmlFile("Tempus_BackwardEuler_CDR.xml");
449  RCP<ParameterList> model_pl = sublist(pList, "CDR Model", true);
450  const int num_elements = model_pl->get<int>("num elements");
451  const double left_end = model_pl->get<double>("left end");
452  const double right_end = model_pl->get<double>("right end");
453 
454  const Thyra::VectorBase<double>& x = *(solutions[solutions.size()-1]);
455 
456  std::ofstream ftmp("Tempus_BackwardEuler_CDR-Solution.dat");
457  for (int n = 0; n < num_elements+1; n++) {
458  const double dx = std::fabs(left_end-right_end) /
459  static_cast<double>(num_elements);
460  const double x_coord = left_end + static_cast<double>(n) * dx;
461  ftmp << x_coord << " " << Thyra::get_ele(x,n) << std::endl;
462  }
463  ftmp.close();
464  }
465 
466  Teuchos::TimeMonitor::summarize();
467 }
468 
469 
470 // ************************************************************
471 // ************************************************************
472 TEUCHOS_UNIT_TEST(BackwardEuler, VanDerPol)
473 {
474  RCP<Tempus::IntegratorBasic<double> > integrator;
475  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
476  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
477  std::vector<double> StepSize;
478  std::vector<double> xErrorNorm;
479  std::vector<double> xDotErrorNorm;
480  const int nTimeStepSizes = 4;
481  double dt = 0.05;
482  for (int n=0; n<nTimeStepSizes; n++) {
483 
484  // Read params from .xml file
485  RCP<ParameterList> pList =
486  getParametersFromXmlFile("Tempus_BackwardEuler_VanDerPol.xml");
487 
488  // Setup the VanDerPolModel
489  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
490  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
491 
492  // Set the step size
493  dt /= 2;
494  if (n == nTimeStepSizes-1) dt /= 10.0;
495 
496  // Setup the Integrator and reset initial time step
497  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
498  pl->sublist("Demo Integrator")
499  .sublist("Time Step Control").set("Initial Time Step", dt);
500  integrator = Tempus::integratorBasic<double>(pl, model);
501 
502  // Integrate to timeMax
503  bool integratorStatus = integrator->advanceTime();
504  TEST_ASSERT(integratorStatus)
505 
506  // Test if at 'Final Time'
507  double time = integrator->getTime();
508  double timeFinal =pl->sublist("Demo Integrator")
509  .sublist("Time Step Control").get<double>("Final Time");
510  double tol = 100.0 * std::numeric_limits<double>::epsilon();
511  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
512 
513  // Store off the final solution and step size
514  StepSize.push_back(dt);
515  auto solution = Thyra::createMember(model->get_x_space());
516  Thyra::copy(*(integrator->getX()),solution.ptr());
517  solutions.push_back(solution);
518  auto solutionDot = Thyra::createMember(model->get_x_space());
519  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
520  solutionsDot.push_back(solutionDot);
521 
522  // Output finest temporal solution for plotting
523  // This only works for ONE MPI process
524  if ((n == 0) or (n == nTimeStepSizes-1)) {
525  std::string fname = "Tempus_BackwardEuler_VanDerPol-Ref.dat";
526  if (n == 0) fname = "Tempus_BackwardEuler_VanDerPol.dat";
527  RCP<const SolutionHistory<double> > solutionHistory =
528  integrator->getSolutionHistory();
529  writeSolution(fname, solutionHistory);
530  }
531  }
532 
533  // Check the order and intercept
534  double xSlope = 0.0;
535  double xDotSlope = 0.0;
536  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
537  double order = stepper->getOrder();
538  writeOrderError("Tempus_BackwardEuler_VanDerPol-Error.dat",
539  stepper, StepSize,
540  solutions, xErrorNorm, xSlope,
541  solutionsDot, xDotErrorNorm, xDotSlope);
542 
543  TEST_FLOATING_EQUALITY( xSlope, order, 0.10 );
544  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.571031, 1.0e-4 );
545  TEST_FLOATING_EQUALITY( xDotSlope, 1.74898, 0.10 );
546  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 1.0038, 1.0e-4 );
547  // At small dt, slopes should be equal to order.
548  //TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
549 
550  Teuchos::TimeMonitor::summarize();
551 }
552 
553 
554 // ************************************************************
555 // ************************************************************
556 TEUCHOS_UNIT_TEST(BackwardEuler, OptInterface)
557 {
558  // Read params from .xml file
559  RCP<ParameterList> pList =
560  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
561 
562  // Setup the SinCosModel
563  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
564  auto model = rcp(new SinCosModel<double>(scm_pl));
565 
566  // Setup the Integrator
567  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
568  RCP<Tempus::IntegratorBasic<double> >integrator =
569  Tempus::integratorBasic<double>(pl, model);
570 
571  // Integrate to timeMax
572  bool integratorStatus = integrator->advanceTime();
573  TEST_ASSERT(integratorStatus);
574 
575  // Get solution history
576  RCP<const SolutionHistory<double> > solutionHistory =
577  integrator->getSolutionHistory();
578 
579  // Get the stepper and cast to optimization interface
580  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
581  RCP<Tempus::StepperOptimizationInterface<double> > opt_stepper =
582  Teuchos::rcp_dynamic_cast< Tempus::StepperOptimizationInterface<double> >(
583  stepper, true);
584 
585  // Check stencil length
586  TEST_EQUALITY( opt_stepper->stencilLength(), 2);
587 
588  // Create needed vectors/multivectors
589  Teuchos::Array< RCP<const Thyra::VectorBase<double> > > x(2);
590  Teuchos::Array<double> t(2);
591  RCP< const Thyra::VectorBase<double> > p =
592  model->getNominalValues().get_p(0);
593  RCP< Thyra::VectorBase<double> > x_dot =
594  Thyra::createMember(model->get_x_space());
595  RCP< Thyra::VectorBase<double> > f =
596  Thyra::createMember(model->get_f_space());
597  RCP< Thyra::VectorBase<double> > f2 =
598  Thyra::createMember(model->get_f_space());
599  RCP< Thyra::LinearOpBase<double> > dfdx =
600  model->create_W_op();
601  RCP< Thyra::LinearOpBase<double> > dfdx2 =
602  model->create_W_op();
603  RCP< Thyra::MultiVectorBase<double> > dfdx_mv =
604  Teuchos::rcp_dynamic_cast< Thyra::MultiVectorBase<double> >(dfdx,true);
605  RCP< Thyra::MultiVectorBase<double> > dfdx_mv2 =
606  Teuchos::rcp_dynamic_cast< Thyra::MultiVectorBase<double> >(dfdx2,true);
607  const int num_p = p->range()->dim();
608  RCP< Thyra::MultiVectorBase<double> > dfdp =
609  Thyra::createMembers(model->get_f_space(), num_p);
610  RCP< Thyra::MultiVectorBase<double> > dfdp2 =
611  Thyra::createMembers(model->get_f_space(), num_p);
612  RCP< Thyra::LinearOpWithSolveBase<double> > W =
613  model->create_W();
614  RCP< Thyra::LinearOpWithSolveBase<double> > W2 =
615  model->create_W();
616  RCP< Thyra::MultiVectorBase<double> > tmp =
617  Thyra::createMembers(model->get_x_space(), num_p);
618  RCP< Thyra::MultiVectorBase<double> > tmp2 =
619  Thyra::createMembers(model->get_x_space(), num_p);
620  std::vector<double> nrms(num_p);
621  double err;
622 
623  // Loop over states, checking residuals and derivatives
624  const int n = solutionHistory->getNumStates();
625  for (int i=1; i<n; ++i) {
626  RCP<const SolutionState<double> > state = (*solutionHistory)[i];
627  RCP<const SolutionState<double> > prev_state = (*solutionHistory)[i-1];
628 
629  // Fill x, t stencils
630  x[0] = state->getX();
631  x[1] = prev_state->getX();
632  t[0] = state->getTime();
633  t[1] = prev_state->getTime();
634 
635  // Compute x_dot
636  const double dt = t[0]-t[1];
637  Thyra::V_StVpStV(x_dot.ptr(), 1.0/dt, *(x[0]), -1.0/dt, *(x[1]));
638 
639  // Create model inargs
640  typedef Thyra::ModelEvaluatorBase MEB;
641  MEB::InArgs<double> in_args = model->createInArgs();
642  MEB::OutArgs<double> out_args = model->createOutArgs();
643  in_args.set_x(x[0]);
644  in_args.set_x_dot(x_dot);
645  in_args.set_t(t[0]);
646  in_args.set_p(0,p);
647 
648  const double tol = 1.0e-14;
649 
650  // Check residual
651  opt_stepper->computeStepResidual(*f, x, t, *p, 0);
652  out_args.set_f(f2);
653  model->evalModel(in_args, out_args);
654  out_args.set_f(Teuchos::null);
655  Thyra::V_VmV(f.ptr(), *f, *f2);
656  err = Thyra::norm(*f);
657  TEST_FLOATING_EQUALITY(err, 0.0, tol);
658 
659  // Check df/dx_n
660  // df/dx_n = df/dx_dot * dx_dot/dx_n + df/dx_n = 1/dt*df/dx_dot + df/dx_n
661  opt_stepper->computeStepJacobian(*dfdx, x, t, *p, 0, 0);
662  out_args.set_W_op(dfdx2);
663  in_args.set_alpha(1.0/dt);
664  in_args.set_beta(1.0);
665  model->evalModel(in_args, out_args);
666  out_args.set_W_op(Teuchos::null);
667  Thyra::V_VmV(dfdx_mv.ptr(), *dfdx_mv, *dfdx_mv2);
668  Thyra::norms(*dfdx_mv, Teuchos::arrayViewFromVector(nrms));
669  err = 0.0;
670  for (auto nrm : nrms) err += nrm;
671  TEST_FLOATING_EQUALITY(err, 0.0, tol);
672 
673  // Check df/dx_{n-1}
674  // df/dx_{n-1} = df/dx_dot * dx_dot/dx_{n-1} = -1/dt*df/dx_dot
675  opt_stepper->computeStepJacobian(*dfdx, x, t, *p, 0, 1);
676  out_args.set_W_op(dfdx2);
677  in_args.set_alpha(-1.0/dt);
678  in_args.set_beta(0.0);
679  model->evalModel(in_args, out_args);
680  out_args.set_W_op(Teuchos::null);
681  Thyra::V_VmV(dfdx_mv.ptr(), *dfdx_mv, *dfdx_mv2);
682  Thyra::norms(*dfdx_mv, Teuchos::arrayViewFromVector(nrms));
683  err = 0.0;
684  for (auto nrm : nrms) err += nrm;
685  TEST_FLOATING_EQUALITY(err, 0.0, tol);
686 
687  // Check df/dp
688  opt_stepper->computeStepParamDeriv(*dfdp, x, t, *p, 0);
689  out_args.set_DfDp(
690  0, MEB::Derivative<double>(dfdp2, MEB::DERIV_MV_JACOBIAN_FORM));
691  model->evalModel(in_args, out_args);
692  out_args.set_DfDp(0, MEB::Derivative<double>());
693  Thyra::V_VmV(dfdp.ptr(), *dfdp, *dfdp2);
694  Thyra::norms(*dfdp, Teuchos::arrayViewFromVector(nrms));
695  err = 0.0;
696  for (auto nrm : nrms) err += nrm;
697  TEST_FLOATING_EQUALITY(err, 0.0, tol);
698 
699  // Check W
700  opt_stepper->computeStepSolver(*W, x, t, *p, 0);
701  out_args.set_W(W2);
702  in_args.set_alpha(1.0/dt);
703  in_args.set_beta(1.0);
704  model->evalModel(in_args, out_args);
705  out_args.set_W(Teuchos::null);
706  // note: dfdp overwritten above so dfdp != dfdp2
707  Thyra::solve(*W, Thyra::NOTRANS, *dfdp2, tmp.ptr());
708  Thyra::solve(*W2, Thyra::NOTRANS, *dfdp2, tmp2.ptr());
709  Thyra::V_VmV(tmp.ptr(), *tmp, *tmp2);
710  Thyra::norms(*tmp, Teuchos::arrayViewFromVector(nrms));
711  err = 0.0;
712  for (auto nrm : nrms) err += nrm;
713  TEST_FLOATING_EQUALITY(err, 0.0, tol);
714  }
715 
716  Teuchos::TimeMonitor::summarize();
717 }
718 
719 
720 } // 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...
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
1D CGFEM model for convection/diffusion/reaction
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.