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: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
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 "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
20 
21 #include <fstream>
22 #include <limits>
23 #include <sstream>
24 #include <vector>
25 
26 namespace Tempus_Test {
27 
28 using Teuchos::getParametersFromXmlFile;
30 using Teuchos::RCP;
31 using Teuchos::rcp;
32 using Teuchos::rcp_const_cast;
33 using Teuchos::sublist;
34 
38 
39 // ************************************************************
40 // ************************************************************
41 TEUCHOS_UNIT_TEST(BDF2, SinCosAdapt)
42 {
44  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
45  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
46  std::vector<double> StepSize;
47 
48  // Read params from .xml file
49  RCP<ParameterList> pList =
50  getParametersFromXmlFile("Tempus_BDF2_SinCos_AdaptDt.xml");
51  // Set initial time step = 2*dt specified in input file (for convergence
52  // study)
53  double dt = pList->sublist("Tempus")
54  .sublist("Default Integrator")
55  .sublist("Time Step Control")
56  .get<double>("Initial Time Step");
57  dt *= 2.0;
58 
59  // Setup the SinCosModel
60  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
61  const int nTimeStepSizes = scm_pl->get<int>("Number of Time Step Sizes", 7);
62  std::string output_file_string =
63  scm_pl->get<std::string>("Output File Name", "Tempus_BDF2_SinCos");
64  std::string output_file_name = output_file_string + ".dat";
65  std::string err_out_file_name = output_file_string + "-Error.dat";
66  double time = 0.0;
67  for (int n = 0; n < nTimeStepSizes; n++) {
68  auto model = rcp(new SinCosModel<double>(scm_pl));
69 
70  dt /= 2;
71 
72  RCP<ParameterList> tempusPL =
73  getParametersFromXmlFile("Tempus_BDF2_SinCos_AdaptDt.xml");
74  RCP<ParameterList> pl = sublist(tempusPL, "Tempus", true);
75 
76  // Setup the Integrator and reset initial time step
77  pl->sublist("Default Integrator")
78  .sublist("Time Step Control")
79  .set("Initial Time Step", dt / 4.0);
80  // Ensure time step does not get larger than the initial time step size,
81  // as that would mess up the convergence rates.
82  pl->sublist("Default Integrator")
83  .sublist("Time Step Control")
84  .set("Maximum Time Step", dt);
85  // Ensure time step does not get too small and therefore too many steps.
86  pl->sublist("Default Integrator")
87  .sublist("Time Step Control")
88  .set("Minimum Time Step", dt / 4.0);
89  // For the SinCos problem eta is directly related to dt
90  pl->sublist("Default Integrator")
91  .sublist("Time Step Control")
92  .sublist("Time Step Control Strategy")
93  .set("Minimum Value Monitoring Function", dt * 0.99);
94  integrator = Tempus::createIntegratorBasic<double>(pl, model);
95 
96  // Initial Conditions
97  // During the Integrator construction, the initial SolutionState
98  // is set by default to model->getNominalVales().get_x(). However,
99  // the application can set it also by integrator->initializeSolutionHistory.
101  model->getNominalValues().get_x()->clone_v();
102  integrator->initializeSolutionHistory(0.0, x0);
103 
104  // Integrate to timeMax
105  bool integratorStatus = integrator->advanceTime();
106  TEST_ASSERT(integratorStatus)
107 
108  // Test if at 'Final Time'
109  time = integrator->getTime();
110  double timeFinal = pl->sublist("Default Integrator")
111  .sublist("Time Step Control")
112  .get<double>("Final Time");
113  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
114 
115  // Time-integrated solution and the exact solution
116  RCP<Thyra::VectorBase<double>> x = integrator->getX();
118  model->getExactSolution(time).get_x();
119 
120  // Plot sample solution and exact solution
121  if (n == 0) {
122  std::ofstream ftmp(output_file_name);
123  // Warning: the following assumes serial run
124  FILE *gold_file = fopen("Tempus_BDF2_SinCos_AdaptDt_gold.dat", "r");
125  RCP<const SolutionHistory<double>> solutionHistory =
126  integrator->getSolutionHistory();
128  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
129  char time_gold_char[100];
130  fgets(time_gold_char, 100, gold_file);
131  double time_gold;
132  sscanf(time_gold_char, "%lf", &time_gold);
133  RCP<const SolutionState<double>> solutionState = (*solutionHistory)[i];
134  double time_i = solutionState->getTime();
135  // Throw error if time does not match time in gold file to specified
136  // tolerance
137  TEST_FLOATING_EQUALITY(time_i, time_gold, 1.0e-5);
138  RCP<const Thyra::VectorBase<double>> x_plot = solutionState->getX();
139  x_exact_plot = model->getExactSolution(time_i).get_x();
140  ftmp << time_i << " " << get_ele(*(x_plot), 0) << " "
141  << get_ele(*(x_plot), 1) << " " << get_ele(*(x_exact_plot), 0)
142  << " " << get_ele(*(x_exact_plot), 1) << std::endl;
143  }
144  ftmp.close();
145  }
146 
147  // Store off the final solution and step size
148  StepSize.push_back(dt);
149  auto solution = Thyra::createMember(model->get_x_space());
150  Thyra::copy(*(integrator->getX()), solution.ptr());
151  solutions.push_back(solution);
152  auto solutionDot = Thyra::createMember(model->get_x_space());
153  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
154  solutionsDot.push_back(solutionDot);
155  if (n == nTimeStepSizes - 1) { // Add exact solution last in vector.
156  StepSize.push_back(0.0);
157  auto solutionExact = Thyra::createMember(model->get_x_space());
158  Thyra::copy(*(model->getExactSolution(time).get_x()),
159  solutionExact.ptr());
160  solutions.push_back(solutionExact);
161  auto solutionDotExact = Thyra::createMember(model->get_x_space());
162  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
163  solutionDotExact.ptr());
164  solutionsDot.push_back(solutionDotExact);
165  }
166  }
167 
168  // Check the order and intercept
169  if (nTimeStepSizes > 1) {
170  double xSlope = 0.0;
171  double xDotSlope = 0.0;
172  std::vector<double> xErrorNorm;
173  std::vector<double> xDotErrorNorm;
174  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
175  // double order = stepper->getOrder();
176  writeOrderError("Tempus_BDF2_SinCos-Error.dat", stepper, StepSize,
177  solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
178  xDotSlope, out);
179 
180  TEST_FLOATING_EQUALITY(xSlope, 1.932, 0.01);
181  TEST_FLOATING_EQUALITY(xDotSlope, 1.932, 0.01);
182  TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.000192591, 1.0e-4);
183  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.000192591, 1.0e-4);
184  }
185 
187 }
188 
189 } // namespace Tempus_Test
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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...
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.