Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_Test_BDF2_SinCosAdapt.cpp
Go to the documentation of this file.
1 //@HEADER
2 // *****************************************************************************
3 // Tempus: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
12 #include "Teuchos_TimeMonitor.hpp"
13 #include "Teuchos_DefaultComm.hpp"
14 
15 #include "Tempus_config.hpp"
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperBDF2.hpp"
18 
19 #include "../TestModels/SinCosModel.hpp"
20 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
21 
22 #include <fstream>
23 #include <limits>
24 #include <sstream>
25 #include <vector>
26 
27 namespace Tempus_Test {
28 
29 using Teuchos::getParametersFromXmlFile;
31 using Teuchos::RCP;
32 using Teuchos::rcp;
33 using Teuchos::rcp_const_cast;
34 using Teuchos::sublist;
35 
39 
40 // ************************************************************
41 // ************************************************************
42 TEUCHOS_UNIT_TEST(BDF2, SinCosAdapt)
43 {
45  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
46  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
47  std::vector<double> StepSize;
48 
49  // Read params from .xml file
50  RCP<ParameterList> pList =
51  getParametersFromXmlFile("Tempus_BDF2_SinCos_AdaptDt.xml");
52  // Set initial time step = 2*dt specified in input file (for convergence
53  // study)
54  double dt = pList->sublist("Tempus")
55  .sublist("Default Integrator")
56  .sublist("Time Step Control")
57  .get<double>("Initial Time Step");
58  dt *= 2.0;
59 
60  // Setup the SinCosModel
61  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
62  const int nTimeStepSizes = scm_pl->get<int>("Number of Time Step Sizes", 7);
63  std::string output_file_string =
64  scm_pl->get<std::string>("Output File Name", "Tempus_BDF2_SinCos");
65  std::string output_file_name = output_file_string + ".dat";
66  std::string err_out_file_name = output_file_string + "-Error.dat";
67  double time = 0.0;
68  for (int n = 0; n < nTimeStepSizes; n++) {
69  auto model = rcp(new SinCosModel<double>(scm_pl));
70 
71  dt /= 2;
72 
73  RCP<ParameterList> tempusPL =
74  getParametersFromXmlFile("Tempus_BDF2_SinCos_AdaptDt.xml");
75  RCP<ParameterList> pl = sublist(tempusPL, "Tempus", true);
76 
77  // Setup the Integrator and reset initial time step
78  pl->sublist("Default Integrator")
79  .sublist("Time Step Control")
80  .set("Initial Time Step", dt / 4.0);
81  // Ensure time step does not get larger than the initial time step size,
82  // as that would mess up the convergence rates.
83  pl->sublist("Default Integrator")
84  .sublist("Time Step Control")
85  .set("Maximum Time Step", dt);
86  // Ensure time step does not get too small and therefore too many steps.
87  pl->sublist("Default Integrator")
88  .sublist("Time Step Control")
89  .set("Minimum Time Step", dt / 4.0);
90  // For the SinCos problem eta is directly related to dt
91  pl->sublist("Default Integrator")
92  .sublist("Time Step Control")
93  .sublist("Time Step Control Strategy")
94  .set("Minimum Value Monitoring Function", dt * 0.99);
95  integrator = Tempus::createIntegratorBasic<double>(pl, model);
96 
97  // Initial Conditions
98  // During the Integrator construction, the initial SolutionState
99  // is set by default to model->getNominalVales().get_x(). However,
100  // the application can set it also by integrator->initializeSolutionHistory.
102  model->getNominalValues().get_x()->clone_v();
103  integrator->initializeSolutionHistory(0.0, x0);
104 
105  // Integrate to timeMax
106  bool integratorStatus = integrator->advanceTime();
107  TEST_ASSERT(integratorStatus)
108 
109  // Test if at 'Final Time'
110  time = integrator->getTime();
111  double timeFinal = pl->sublist("Default Integrator")
112  .sublist("Time Step Control")
113  .get<double>("Final Time");
114  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
115 
116  // Time-integrated solution and the exact solution
117  RCP<Thyra::VectorBase<double>> x = integrator->getX();
119  model->getExactSolution(time).get_x();
120 
121  // Plot sample solution and exact solution
122  if (n == 0) {
123  std::ofstream ftmp(output_file_name);
124  // Warning: the following assumes serial run
125  FILE *gold_file = fopen("Tempus_BDF2_SinCos_AdaptDt_gold.dat", "r");
126  RCP<const SolutionHistory<double>> solutionHistory =
127  integrator->getSolutionHistory();
129  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
130  char time_gold_char[100];
131  fgets(time_gold_char, 100, gold_file);
132  double time_gold;
133  sscanf(time_gold_char, "%lf", &time_gold);
134  RCP<const SolutionState<double>> solutionState = (*solutionHistory)[i];
135  double time_i = solutionState->getTime();
136  // Throw error if time does not match time in gold file to specified
137  // tolerance
138  TEST_FLOATING_EQUALITY(time_i, time_gold, 1.0e-5);
139  RCP<const Thyra::VectorBase<double>> x_plot = solutionState->getX();
140  x_exact_plot = model->getExactSolution(time_i).get_x();
141  ftmp << time_i << " " << get_ele(*(x_plot), 0) << " "
142  << get_ele(*(x_plot), 1) << " " << get_ele(*(x_exact_plot), 0)
143  << " " << get_ele(*(x_exact_plot), 1) << std::endl;
144  }
145  ftmp.close();
146  }
147 
148  // Store off the final solution and step size
149  StepSize.push_back(dt);
150  auto solution = Thyra::createMember(model->get_x_space());
151  Thyra::copy(*(integrator->getX()), solution.ptr());
152  solutions.push_back(solution);
153  auto solutionDot = Thyra::createMember(model->get_x_space());
154  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
155  solutionsDot.push_back(solutionDot);
156  if (n == nTimeStepSizes - 1) { // Add exact solution last in vector.
157  StepSize.push_back(0.0);
158  auto solutionExact = Thyra::createMember(model->get_x_space());
159  Thyra::copy(*(model->getExactSolution(time).get_x()),
160  solutionExact.ptr());
161  solutions.push_back(solutionExact);
162  auto solutionDotExact = Thyra::createMember(model->get_x_space());
163  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
164  solutionDotExact.ptr());
165  solutionsDot.push_back(solutionDotExact);
166  }
167  }
168 
169  // Check the order and intercept
170  if (nTimeStepSizes > 1) {
171  double xSlope = 0.0;
172  double xDotSlope = 0.0;
173  std::vector<double> xErrorNorm;
174  std::vector<double> xDotErrorNorm;
175  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
176  // double order = stepper->getOrder();
177  writeOrderError("Tempus_BDF2_SinCos-Error.dat", stepper, StepSize,
178  solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
179  xDotSlope, out);
180 
181  TEST_FLOATING_EQUALITY(xSlope, 1.932, 0.01);
182  TEST_FLOATING_EQUALITY(xDotSlope, 1.932, 0.01);
183  TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.000192591, 1.0e-4);
184  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.000192591, 1.0e-4);
185  }
186 
188 }
189 
190 } // namespace Tempus_Test
T & get(const std::string &name, T def_value)
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::FancyOStream &out)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
#define TEST_ASSERT(v1)
T * get() const
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation with a...
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Solution state for integrators and steppers.