Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_BDF2Test.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_StepperBDF2.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 <fstream>
33 #include <limits>
34 #include <sstream>
35 #include <vector>
36 
37 // Comment out any of the following tests to exclude from build/run.
38 #define TEST_PARAMETERLIST
39 #define TEST_CONSTRUCTING_FROM_DEFAULTS
40 #define TEST_SINCOS
41 #define TEST_SINCOS_ADAPT
42 #define TEST_CDR
43 #define TEST_VANDERPOL
44 
45 namespace Tempus_Test {
46 
47 using Teuchos::RCP;
48 using Teuchos::rcp;
49 using Teuchos::rcp_const_cast;
50 using Teuchos::ParameterList;
51 using Teuchos::sublist;
52 using Teuchos::getParametersFromXmlFile;
53 
57 
58 
59 #ifdef TEST_PARAMETERLIST
60 // ************************************************************
61 // ************************************************************
62 TEUCHOS_UNIT_TEST(BDF2, ParameterList)
63 {
64  // Read params from .xml file
65  RCP<ParameterList> pList =
66  getParametersFromXmlFile("Tempus_BDF2_SinCos.xml");
67 
68  // Setup the SinCosModel
69  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
70  auto model = rcp(new SinCosModel<double> (scm_pl));
71 
72  RCP<ParameterList> tempusPL = sublist(pList, "Tempus", true);
73 
74  // Test constructor IntegratorBasic(tempusPL, model)
75  {
76  RCP<Tempus::IntegratorBasic<double> > integrator =
77  Tempus::integratorBasic<double>(tempusPL, model);
78 
79  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
80  // Remove Start Up Stepper for comparison
81  stepperPL->remove("Start Up Stepper Name");
82  stepperPL->remove("Default Start Up Stepper");
83  RCP<ParameterList> defaultPL =
84  integrator->getStepper()->getDefaultParameters();
85  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
86  if (!pass) {
87  std::cout << std::endl;
88  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
89  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
90  }
91  TEST_ASSERT(pass)
92  }
93 
94  // Test constructor IntegratorBasic(model, stepperType)
95  {
96  RCP<Tempus::IntegratorBasic<double> > integrator =
97  Tempus::integratorBasic<double>(model, "BDF2");
98 
99  RCP<ParameterList> stepperPL = sublist(tempusPL, "Default Stepper", true);
100  RCP<ParameterList> defaultPL =
101  integrator->getStepper()->getDefaultParameters();
102 
103  bool pass = haveSameValues(*stepperPL, *defaultPL, true);
104  if (!pass) {
105  std::cout << std::endl;
106  std::cout << "stepperPL -------------- \n" << *stepperPL << std::endl;
107  std::cout << "defaultPL -------------- \n" << *defaultPL << std::endl;
108  }
109  TEST_ASSERT(pass)
110  }
111 }
112 #endif // TEST_PARAMETERLIST
113 
114 
115 #ifdef TEST_CONSTRUCTING_FROM_DEFAULTS
116 // ************************************************************
117 // ************************************************************
118 TEUCHOS_UNIT_TEST(BDF2, ConstructingFromDefaults)
119 {
120  double dt = 0.1;
121 
122  // Read params from .xml file
123  RCP<ParameterList> pList =
124  getParametersFromXmlFile("Tempus_BDF2_SinCos.xml");
125  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
126 
127  // Setup the SinCosModel
128  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
129  //RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
130  auto model = rcp(new SinCosModel<double>(scm_pl));
131 
132  // Setup Stepper for field solve ----------------------------
133  auto stepper = rcp(new Tempus::StepperBDF2<double>());
134  stepper->setModel(model);
135  stepper->initialize();
136 
137  // Setup TimeStepControl ------------------------------------
138  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
139  ParameterList tscPL = pl->sublist("Default Integrator")
140  .sublist("Time Step Control");
141  timeStepControl->setStepType (tscPL.get<std::string>("Integrator Step Type"));
142  timeStepControl->setInitIndex(tscPL.get<int> ("Initial Time Index"));
143  timeStepControl->setInitTime (tscPL.get<double>("Initial Time"));
144  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
145  timeStepControl->setInitTimeStep(dt);
146  timeStepControl->initialize();
147 
148  // Setup initial condition SolutionState --------------------
149  Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
150  stepper->getModel()->getNominalValues();
151  auto icSolution = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
152  auto icState = rcp(new Tempus::SolutionState<double>(icSolution));
153  icState->setTime (timeStepControl->getInitTime());
154  icState->setIndex (timeStepControl->getInitIndex());
155  icState->setTimeStep(0.0);
156  icState->setOrder (stepper->getOrder());
157  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
158 
159  // Setup SolutionHistory ------------------------------------
161  solutionHistory->setName("Forward States");
163  solutionHistory->setStorageLimit(3);
164  solutionHistory->addState(icState);
165 
166  // Setup Integrator -----------------------------------------
167  RCP<Tempus::IntegratorBasic<double> > integrator =
168  Tempus::integratorBasic<double>();
169  integrator->setStepperWStepper(stepper);
170  integrator->setTimeStepControl(timeStepControl);
171  integrator->setSolutionHistory(solutionHistory);
172  //integrator->setObserver(...);
173  integrator->initialize();
174 
175 
176  // Integrate to timeMax
177  bool integratorStatus = integrator->advanceTime();
178  TEST_ASSERT(integratorStatus)
179 
180 
181  // Test if at 'Final Time'
182  double time = integrator->getTime();
183  double timeFinal =pl->sublist("Default Integrator")
184  .sublist("Time Step Control").get<double>("Final Time");
185  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
186 
187  // Time-integrated solution and the exact solution
188  RCP<Thyra::VectorBase<double> > x = integrator->getX();
189  RCP<const Thyra::VectorBase<double> > x_exact =
190  model->getExactSolution(time).get_x();
191 
192  // Calculate the error
193  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
194  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
195 
196  // Check the order and intercept
197  std::cout << " Stepper = BDF2" << std::endl;
198  std::cout << " =========================" << std::endl;
199  std::cout << " Exact solution : " << get_ele(*(x_exact), 0) << " "
200  << get_ele(*(x_exact), 1) << std::endl;
201  std::cout << " Computed solution: " << get_ele(*(x ), 0) << " "
202  << get_ele(*(x ), 1) << std::endl;
203  std::cout << " Difference : " << get_ele(*(xdiff ), 0) << " "
204  << get_ele(*(xdiff ), 1) << std::endl;
205  std::cout << " =========================" << std::endl;
206  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.839732, 1.0e-4 );
207  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.542663, 1.0e-4 );
208 }
209 #endif // TEST_CONSTRUCTING_FROM_DEFAULTS
210 
211 
212 #ifdef TEST_SINCOS
213 // ************************************************************
214 // ************************************************************
215 TEUCHOS_UNIT_TEST(BDF2, SinCos)
216 {
217  RCP<Tempus::IntegratorBasic<double> > integrator;
218  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
219  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
220  std::vector<double> StepSize;
221 
222  // Read params from .xml file
223  RCP<ParameterList> pList = getParametersFromXmlFile("Tempus_BDF2_SinCos.xml");
224  //Set initial time step = 2*dt specified in input file (for convergence study)
225  //
226  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
227  double dt = pl->sublist("Default Integrator")
228  .sublist("Time Step Control").get<double>("Initial Time Step");
229  dt *= 2.0;
230 
231  // Setup the SinCosModel
232  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
233  const int nTimeStepSizes = scm_pl->get<int>("Number of Time Step Sizes", 7);
234  std::string output_file_string =
235  scm_pl->get<std::string>("Output File Name", "Tempus_BDF2_SinCos");
236  std::string output_file_name = output_file_string + ".dat";
237  std::string ref_out_file_name = output_file_string + "-Ref.dat";
238  std::string err_out_file_name = output_file_string + "-Error.dat";
239  double time = 0.0;
240  for (int n=0; n<nTimeStepSizes; n++) {
241 
242  //std::ofstream ftmp("PL.txt");
243  //pList->print(ftmp);
244  //ftmp.close();
245 
246  auto model = rcp(new SinCosModel<double>(scm_pl));
247 
248  dt /= 2;
249 
250  // Setup the Integrator and reset initial time step
251  pl->sublist("Default Integrator")
252  .sublist("Time Step Control").set("Initial Time Step", dt);
253  integrator = Tempus::integratorBasic<double>(pl, model);
254 
255  // Initial Conditions
256  // During the Integrator construction, the initial SolutionState
257  // is set by default to model->getNominalVales().get_x(). However,
258  // the application can set it also by integrator->initializeSolutionHistory.
259  RCP<Thyra::VectorBase<double> > x0 =
260  model->getNominalValues().get_x()->clone_v();
261  integrator->initializeSolutionHistory(0.0, x0);
262 
263  // Integrate to timeMax
264  bool integratorStatus = integrator->advanceTime();
265  TEST_ASSERT(integratorStatus)
266 
267  // Test if at 'Final Time'
268  time = integrator->getTime();
269  double timeFinal =pl->sublist("Default Integrator")
270  .sublist("Time Step Control").get<double>("Final Time");
271  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
272 
273  // Plot sample solution and exact solution
274  if (n == 0) {
275  RCP<const SolutionHistory<double> > solutionHistory =
276  integrator->getSolutionHistory();
277  writeSolution(output_file_name, solutionHistory);
278 
279  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
280  for (int i=0; i<solutionHistory->getNumStates(); i++) {
281  double time_i = (*solutionHistory)[i]->getTime();
282  auto state = rcp(new Tempus::SolutionState<double>(
283  model->getExactSolution(time_i).get_x(),
284  model->getExactSolution(time_i).get_x_dot()));
285  state->setTime((*solutionHistory)[i]->getTime());
286  solnHistExact->addState(state);
287  }
288  writeSolution(ref_out_file_name, solnHistExact);
289  }
290 
291  // Store off the final solution and step size
292  StepSize.push_back(dt);
293  auto solution = Thyra::createMember(model->get_x_space());
294  Thyra::copy(*(integrator->getX()),solution.ptr());
295  solutions.push_back(solution);
296  auto solutionDot = Thyra::createMember(model->get_x_space());
297  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
298  solutionsDot.push_back(solutionDot);
299  if (n == nTimeStepSizes-1) { // Add exact solution last in vector.
300  StepSize.push_back(0.0);
301  auto solutionExact = Thyra::createMember(model->get_x_space());
302  Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
303  solutions.push_back(solutionExact);
304  auto solutionDotExact = Thyra::createMember(model->get_x_space());
305  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
306  solutionDotExact.ptr());
307  solutionsDot.push_back(solutionDotExact);
308  }
309  }
310 
311  // Check the order and intercept
312  if (nTimeStepSizes > 1) {
313  double xSlope = 0.0;
314  double xDotSlope = 0.0;
315  std::vector<double> xErrorNorm;
316  std::vector<double> xDotErrorNorm;
317  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
318  double order = stepper->getOrder();
319  writeOrderError(err_out_file_name,
320  stepper, StepSize,
321  solutions, xErrorNorm, xSlope,
322  solutionsDot, xDotErrorNorm, xDotSlope);
323 
324  TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
325  TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
326  TEST_FLOATING_EQUALITY( xErrorNorm[0], 5.13425e-05, 1.0e-4 );
327  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 5.13425e-05, 1.0e-4 );
328  }
329 
330  Teuchos::TimeMonitor::summarize();
331 }
332 #endif // TEST_SINCOS
333 
334 
335 #ifdef TEST_SINCOS_ADAPT
336 // ************************************************************
337 // ************************************************************
338 TEUCHOS_UNIT_TEST(BDF2, SinCosAdapt)
339 {
340  RCP<Tempus::IntegratorBasic<double> > integrator;
341  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
342  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
343  std::vector<double> StepSize;
344 
345  // Read params from .xml file
346  RCP<ParameterList> pList =
347  getParametersFromXmlFile("Tempus_BDF2_SinCos_AdaptDt.xml");
348  //Set initial time step = 2*dt specified in input file (for convergence study)
349  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
350  double dt = pl->sublist("Default Integrator")
351  .sublist("Time Step Control").get<double>("Initial Time Step");
352  dt *= 2.0;
353 
354  // Setup the SinCosModel
355  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
356  const int nTimeStepSizes = scm_pl->get<int>("Number of Time Step Sizes", 7);
357  std::string output_file_string =
358  scm_pl->get<std::string>("Output File Name", "Tempus_BDF2_SinCos");
359  std::string output_file_name = output_file_string + ".dat";
360  std::string err_out_file_name = output_file_string + "-Error.dat";
361  double time = 0.0;
362  for (int n=0; n<nTimeStepSizes; n++) {
363 
364  //std::ofstream ftmp("PL.txt");
365  //pList->print(ftmp);
366  //ftmp.close();
367 
368  auto model = rcp(new SinCosModel<double>(scm_pl));
369 
370  dt /= 2;
371 
372  // Setup the Integrator and reset initial time step
373  pl->sublist("Default Integrator")
374  .sublist("Time Step Control").set("Initial Time Step", dt/4.0);
375  // Ensure time step does not get larger than the initial time step size,
376  // as that would mess up the convergence rates.
377  pl->sublist("Default Integrator")
378  .sublist("Time Step Control").set("Maximum Time Step", dt);
379  // Ensure time step does not get too small and therefore too many steps.
380  pl->sublist("Default Integrator")
381  .sublist("Time Step Control").set("Minimum Time Step", dt/4.0);
382  // For the SinCos problem eta is directly related to dt
383  pl->sublist("Default Integrator")
384  .sublist("Time Step Control")
385  .sublist("Time Step Control Strategy")
386  .sublist("basic_vs")
387  .set("Minimum Value Monitoring Function", dt*0.99);
388  integrator = Tempus::integratorBasic<double>(pl, model);
389 
390  // Initial Conditions
391  // During the Integrator construction, the initial SolutionState
392  // is set by default to model->getNominalVales().get_x(). However,
393  // the application can set it also by integrator->initializeSolutionHistory.
394  RCP<Thyra::VectorBase<double> > x0 =
395  model->getNominalValues().get_x()->clone_v();
396  integrator->initializeSolutionHistory(0.0, x0);
397 
398  // Integrate to timeMax
399  bool integratorStatus = integrator->advanceTime();
400  TEST_ASSERT(integratorStatus)
401 
402  // Test if at 'Final Time'
403  time = integrator->getTime();
404  double timeFinal =pl->sublist("Default Integrator")
405  .sublist("Time Step Control").get<double>("Final Time");
406  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
407 
408  // Time-integrated solution and the exact solution
409  RCP<Thyra::VectorBase<double> > x = integrator->getX();
410  RCP<const Thyra::VectorBase<double> > x_exact =
411  model->getExactSolution(time).get_x();
412 
413  // Plot sample solution and exact solution
414  if (n == 0) {
415  std::ofstream ftmp(output_file_name);
416  //Warning: the following assumes serial run
417  FILE *gold_file = fopen("Tempus_BDF2_SinCos_AdaptDt_gold.dat", "r");
418  RCP<const SolutionHistory<double> > solutionHistory =
419  integrator->getSolutionHistory();
420  RCP<const Thyra::VectorBase<double> > x_exact_plot;
421  for (int i=0; i<solutionHistory->getNumStates(); i++) {
422  char time_gold_char[100];
423  fgets(time_gold_char, 100, gold_file);
424  double time_gold;
425  sscanf(time_gold_char, "%lf", &time_gold);
426  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
427  double time_i = solutionState->getTime();
428  //Throw error if time does not match time in gold file to specified tolerance
429  TEST_FLOATING_EQUALITY( time_i, time_gold, 1.0e-5 );
430  RCP<const Thyra::VectorBase<double> > x_plot = solutionState->getX();
431  x_exact_plot = model->getExactSolution(time_i).get_x();
432  ftmp << time_i << " "
433  << get_ele(*(x_plot), 0) << " "
434  << get_ele(*(x_plot), 1) << " "
435  << get_ele(*(x_exact_plot), 0) << " "
436  << get_ele(*(x_exact_plot), 1) << std::endl;
437  }
438  ftmp.close();
439  }
440 
441  // Store off the final solution and step size
442  StepSize.push_back(dt);
443  auto solution = Thyra::createMember(model->get_x_space());
444  Thyra::copy(*(integrator->getX()),solution.ptr());
445  solutions.push_back(solution);
446  auto solutionDot = Thyra::createMember(model->get_x_space());
447  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
448  solutionsDot.push_back(solutionDot);
449  if (n == nTimeStepSizes-1) { // Add exact solution last in vector.
450  StepSize.push_back(0.0);
451  auto solutionExact = Thyra::createMember(model->get_x_space());
452  Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
453  solutions.push_back(solutionExact);
454  auto solutionDotExact = Thyra::createMember(model->get_x_space());
455  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
456  solutionDotExact.ptr());
457  solutionsDot.push_back(solutionDotExact);
458  }
459  }
460 
461  // Check the order and intercept
462  if (nTimeStepSizes > 1) {
463  double xSlope = 0.0;
464  double xDotSlope = 0.0;
465  std::vector<double> xErrorNorm;
466  std::vector<double> xDotErrorNorm;
467  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
468  //double order = stepper->getOrder();
469  writeOrderError("Tempus_BDF2_SinCos-Error.dat",
470  stepper, StepSize,
471  solutions, xErrorNorm, xSlope,
472  solutionsDot, xDotErrorNorm, xDotSlope);
473 
474  TEST_FLOATING_EQUALITY( xSlope, 1.95089, 0.01 );
475  TEST_FLOATING_EQUALITY( xDotSlope, 1.95089, 0.01 );
476  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.000197325, 1.0e-4 );
477  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.000197325, 1.0e-4 );
478  }
479 
480  Teuchos::TimeMonitor::summarize();
481 }
482 #endif // TEST_SINCOS_ADAPT
483 
484 
485 #ifdef TEST_CDR
486 // ************************************************************
487 // ************************************************************
489 {
490  // Create a communicator for Epetra objects
491  RCP<Epetra_Comm> comm;
492 #ifdef Tempus_ENABLE_MPI
493  comm = rcp(new Epetra_MpiComm(MPI_COMM_WORLD));
494 #else
495  comm = rcp(new Epetra_SerialComm);
496 #endif
497 
498  RCP<Tempus::IntegratorBasic<double> > integrator;
499  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
500  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
501  std::vector<double> StepSize;
502 
503  // Read params from .xml file
504  RCP<ParameterList> pList =
505  getParametersFromXmlFile("Tempus_BDF2_CDR.xml");
506  //Set initial time step = 2*dt specified in input file (for convergence study)
507  //
508  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
509  double dt = pl->sublist("Demo Integrator")
510  .sublist("Time Step Control").get<double>("Initial Time Step");
511  dt *= 2.0;
512  RCP<ParameterList> model_pl = sublist(pList, "CDR Model", true);
513 
514  const int nTimeStepSizes = model_pl->get<int>("Number of Time Step Sizes", 5);
515 
516  for (int n=0; n<nTimeStepSizes; n++) {
517 
518  // Create CDR Model
519  const int num_elements = model_pl->get<int>("num elements");
520  const double left_end = model_pl->get<double>("left end");
521  const double right_end = model_pl->get<double>("right end");
522  const double a_convection = model_pl->get<double>("a (convection)");
523  const double k_source = model_pl->get<double>("k (source)");
524 
525  auto model = rcp(new Tempus_Test::CDR_Model<double>(comm,
526  num_elements,
527  left_end,
528  right_end,
529  a_convection,
530  k_source));
531 
532  // Set the factory
533  ::Stratimikos::DefaultLinearSolverBuilder builder;
534 
535  auto p = rcp(new ParameterList);
536  p->set("Linear Solver Type", "Belos");
537  p->set("Preconditioner Type", "None");
538  builder.setParameterList(p);
539 
540  RCP< ::Thyra::LinearOpWithSolveFactoryBase<double> >
541  lowsFactory = builder.createLinearSolveStrategy("");
542 
543  model->set_W_factory(lowsFactory);
544 
545  // Set the step size
546  dt /= 2;
547 
548  // Setup the Integrator and reset initial time step
549  pl->sublist("Demo Integrator")
550  .sublist("Time Step Control").set("Initial Time Step", dt);
551  integrator = Tempus::integratorBasic<double>(pl, model);
552 
553  // Integrate to timeMax
554  bool integratorStatus = integrator->advanceTime();
555  TEST_ASSERT(integratorStatus)
556 
557  // Test if at 'Final Time'
558  double time = integrator->getTime();
559  double timeFinal =pl->sublist("Demo Integrator")
560  .sublist("Time Step Control").get<double>("Final Time");
561  double tol = 100.0 * std::numeric_limits<double>::epsilon();
562  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
563 
564  // Store off the final solution and step size
565  StepSize.push_back(dt);
566  auto solution = Thyra::createMember(model->get_x_space());
567  Thyra::copy(*(integrator->getX()),solution.ptr());
568  solutions.push_back(solution);
569  auto solutionDot = Thyra::createMember(model->get_x_space());
570  Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
571  solutionsDot.push_back(solutionDot);
572 
573  // Output finest temporal solution for plotting
574  // This only works for ONE MPI process
575  if ((n == nTimeStepSizes-1) && (comm->NumProc() == 1)) {
576  std::ofstream ftmp("Tempus_BDF2_CDR.dat");
577  ftmp << "TITLE=\"BDF2 Solution to CDR\"\n"
578  << "VARIABLES=\"z\",\"T\"\n";
579  const double dx = std::fabs(left_end-right_end) /
580  static_cast<double>(num_elements);
581  RCP<const SolutionHistory<double> > solutionHistory =
582  integrator->getSolutionHistory();
583  int nStates = solutionHistory->getNumStates();
584  for (int i=0; i<nStates; i++) {
585  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
586  RCP<const Thyra::VectorBase<double> > x = solutionState->getX();
587  double ttime = solutionState->getTime();
588  ftmp << "ZONE T=\"Time="<<ttime<<"\", I="
589  <<num_elements+1<<", F=BLOCK\n";
590  for (int j = 0; j < num_elements+1; j++) {
591  const double x_coord = left_end + static_cast<double>(j) * dx;
592  ftmp << x_coord << " ";
593  }
594  ftmp << std::endl;
595  for (int j=0; j<num_elements+1; j++) ftmp << get_ele(*x, j) << " ";
596  ftmp << std::endl;
597  }
598  ftmp.close();
599  }
600  }
601 
602  // Check the order and intercept
603  if (nTimeStepSizes > 2) {
604  double xSlope = 0.0;
605  double xDotSlope = 0.0;
606  std::vector<double> xErrorNorm;
607  std::vector<double> xDotErrorNorm;
608  RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
609  double order = stepper->getOrder();
610  writeOrderError("Tempus_BDF2_CDR-Error.dat",
611  stepper, StepSize,
612  solutions, xErrorNorm, xSlope,
613  solutionsDot, xDotErrorNorm, xDotSlope);
614  TEST_FLOATING_EQUALITY( xSlope, order, 0.35 );
615  TEST_COMPARE(xSlope, >, 0.95);
616  TEST_FLOATING_EQUALITY( xDotSlope, order, 0.35 );
617  TEST_COMPARE(xDotSlope, >, 0.95);
618 
619  TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.0145747, 1.0e-4 );
620  TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.0563621, 1.0e-4 );
621  }
622 
623  // Write fine mesh solution at final time
624  // This only works for ONE MPI process
625  if (comm->NumProc() == 1) {
626  RCP<ParameterList> pListCDR =
627  getParametersFromXmlFile("Tempus_BDF2_CDR.xml");
628  RCP<ParameterList> model_pl_CDR = sublist(pListCDR, "CDR Model", true);
629  const int num_elements = model_pl_CDR->get<int>("num elements");
630  const double left_end = model_pl_CDR->get<double>("left end");
631  const double right_end = model_pl_CDR->get<double>("right end");
632 
633  const Thyra::VectorBase<double>& x = *(solutions[solutions.size()-1]);
634 
635  std::ofstream ftmp("Tempus_BDF2_CDR-Solution.dat");
636  for (int n = 0; n < num_elements+1; n++) {
637  const double dx = std::fabs(left_end-right_end) /
638  static_cast<double>(num_elements);
639  const double x_coord = left_end + static_cast<double>(n) * dx;
640  ftmp << x_coord << " " << Thyra::get_ele(x,n) << std::endl;
641  }
642  ftmp.close();
643  }
644 
645  Teuchos::TimeMonitor::summarize();
646 }
647 #endif // TEST_CDR
648 
649 
650 #ifdef TEST_VANDERPOL
651 // ************************************************************
652 // ************************************************************
653 TEUCHOS_UNIT_TEST(BDF2, VanDerPol)
654 {
655  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
656  std::vector<double> StepSize;
657  std::vector<double> ErrorNorm;
658 
659  // Read params from .xml file
660  RCP<ParameterList> pList =
661  getParametersFromXmlFile("Tempus_BDF2_VanDerPol.xml");
662  //Set initial time step = 2*dt specified in input file (for convergence study)
663  //
664  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
665  double dt = pl->sublist("Demo Integrator")
666  .sublist("Time Step Control").get<double>("Initial Time Step");
667  dt *= 2.0;
668 
669  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
670  const int nTimeStepSizes = vdpm_pl->get<int>("Number of Time Step Sizes", 3);
671  //const int nTimeStepSizes = 5;
672  double order = 0.0;
673 
674  for (int n=0; n<nTimeStepSizes; n++) {
675 
676  // Setup the VanDerPolModel
677  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
678 
679  // Set the step size
680  dt /= 2;
681  if (n == nTimeStepSizes-1) dt /= 10.0;
682 
683  // Setup the Integrator and reset initial time step
684  pl->sublist("Demo Integrator")
685  .sublist("Time Step Control").set("Initial Time Step", dt);
686  RCP<Tempus::IntegratorBasic<double> > integrator =
687  Tempus::integratorBasic<double>(pl, model);
688  order = integrator->getStepper()->getOrder();
689 
690  // Integrate to timeMax
691  bool integratorStatus = integrator->advanceTime();
692  TEST_ASSERT(integratorStatus)
693 
694  // Test if at 'Final Time'
695  double time = integrator->getTime();
696  double timeFinal =pl->sublist("Demo Integrator")
697  .sublist("Time Step Control").get<double>("Final Time");
698  double tol = 100.0 * std::numeric_limits<double>::epsilon();
699  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
700 
701  // Store off the final solution and step size
702  auto solution = Thyra::createMember(model->get_x_space());
703  Thyra::copy(*(integrator->getX()),solution.ptr());
704  solutions.push_back(solution);
705  StepSize.push_back(dt);
706 
707  // Output finest temporal solution for plotting
708  // This only works for ONE MPI process
709  if ((n == 0) or (n == nTimeStepSizes-1)) {
710  std::string fname = "Tempus_BDF2_VanDerPol-Ref.dat";
711  if (n == 0) fname = "Tempus_BDF2_VanDerPol.dat";
712  std::ofstream ftmp(fname);
713  RCP<const SolutionHistory<double> > solutionHistory =
714  integrator->getSolutionHistory();
715  int nStates = solutionHistory->getNumStates();
716  for (int i=0; i<nStates; i++) {
717  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
718  RCP<const Thyra::VectorBase<double> > x = solutionState->getX();
719  double ttime = solutionState->getTime();
720  ftmp << ttime << " " << get_ele(*x, 0) << " " << get_ele(*x, 1)
721  << std::endl;
722  }
723  ftmp.close();
724  }
725  }
726 
727  // Calculate the error - use the most temporally refined mesh for
728  // the reference solution.
729  auto ref_solution = solutions[solutions.size()-1];
730  std::vector<double> StepSizeCheck;
731  for (std::size_t i=0; i < (solutions.size()-1); ++i) {
732  auto tmp = solutions[i];
733  Thyra::Vp_StV(tmp.ptr(), -1.0, *ref_solution);
734  const double L2norm = Thyra::norm_2(*tmp);
735  StepSizeCheck.push_back(StepSize[i]);
736  ErrorNorm.push_back(L2norm);
737  }
738 
739  if (nTimeStepSizes > 2) {
740  // Check the order and intercept
741  double slope = computeLinearRegressionLogLog<double>(StepSizeCheck,ErrorNorm);
742  std::cout << " Stepper = BDF2" << std::endl;
743  std::cout << " =========================" << std::endl;
744  std::cout << " Expected order: " << order << std::endl;
745  std::cout << " Observed order: " << slope << std::endl;
746  std::cout << " =========================" << std::endl;
747  TEST_FLOATING_EQUALITY( slope, order, 0.10 );
748  out << "\n\n ** Slope on BDF2 Method = " << slope
749  << "\n" << std::endl;
750  }
751 
752  // Write error data
753  {
754  std::ofstream ftmp("Tempus_BDF2_VanDerPol-Error.dat");
755  double error0 = 0.8*ErrorNorm[0];
756  for (std::size_t n = 0; n < StepSizeCheck.size(); n++) {
757  ftmp << StepSizeCheck[n] << " " << ErrorNorm[n] << " "
758  << error0*(pow(StepSize[n]/StepSize[0],order)) << std::endl;
759  }
760  ftmp.close();
761  }
762 
763  Teuchos::TimeMonitor::summarize();
764 }
765 #endif // TEST_VANDERPOL
766 
767 } // namespace Tempus_Test
BDF2 (Backward-Difference-Formula-2) time stepper.
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation with a...
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.
van der Pol model problem for nonlinear electrical circuit.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
1D CGFEM model for convection/diffusion/reaction