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