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