Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_LeapfrogTest.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 
14 #include "Thyra_VectorStdOps.hpp"
15 
16 #include "Tempus_config.hpp"
17 #include "Tempus_IntegratorBasic.hpp"
18 #include "Tempus_StepperLeapfrog.hpp"
19 
20 #include "../TestModels/HarmonicOscillatorModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
22 
23 #include <fstream>
24 #include <sstream>
25 #include <limits>
26 #include <vector>
27 
28 namespace Tempus_Test {
29 
30 using Teuchos::getParametersFromXmlFile;
32 using Teuchos::RCP;
33 using Teuchos::rcp;
34 using Teuchos::rcp_const_cast;
35 using Teuchos::sublist;
36 
40 
41 // ************************************************************
42 // ************************************************************
43 TEUCHOS_UNIT_TEST(Leapfrog, ConstructingFromDefaults)
44 {
45  double dt = 0.1;
46  std::vector<std::string> options;
47  options.push_back("Default Parameters");
48  options.push_back("ICConsistency and Check");
49 
50  for (const auto& option : options) {
51  // Read params from .xml file
52  RCP<ParameterList> pList =
53  getParametersFromXmlFile("Tempus_Leapfrog_SinCos.xml");
54  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
55 
56  // Setup the HarmonicOscillatorModel
57  RCP<ParameterList> hom_pl = sublist(pList, "HarmonicOscillatorModel", true);
58  auto model = rcp(new HarmonicOscillatorModel<double>(hom_pl));
59 
60  // Setup Stepper for field solve ----------------------------
61  auto stepper = rcp(new Tempus::StepperLeapfrog<double>());
62  stepper->setModel(model);
63  if (option == "ICConsistency and Check") {
64  stepper->setICConsistency("Consistent");
65  stepper->setICConsistencyCheck(true);
66  }
67  stepper->initialize();
68 
69  // Setup TimeStepControl ------------------------------------
70  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
71  ParameterList tscPL =
72  pl->sublist("Default Integrator").sublist("Time Step Control");
73  timeStepControl->setInitIndex(tscPL.get<int>("Initial Time Index"));
74  timeStepControl->setInitTime(tscPL.get<double>("Initial Time"));
75  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
76  timeStepControl->setInitTimeStep(dt);
77  timeStepControl->initialize();
78 
79  // Setup initial condition SolutionState --------------------
80  auto inArgsIC = model->getNominalValues();
81  auto icX = rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x());
82  auto icXDot =
83  rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x_dot());
84  auto icXDotDot =
85  rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x_dot_dot());
86  auto icState = Tempus::createSolutionStateX<double>(icX, icXDot, icXDotDot);
87  icState->setTime(timeStepControl->getInitTime());
88  icState->setIndex(timeStepControl->getInitIndex());
89  icState->setTimeStep(0.0);
90  icState->setOrder(stepper->getOrder());
91  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
92 
93  // Setup SolutionHistory ------------------------------------
94  auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
95  solutionHistory->setName("Forward States");
96  solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
97  solutionHistory->setStorageLimit(2);
98  solutionHistory->addState(icState);
99 
100  // Ensure ICs are consistent and stepper memory is set (e.g., xDot is set).
101  stepper->setInitialConditions(solutionHistory);
102 
103  // Setup Integrator -----------------------------------------
105  Tempus::createIntegratorBasic<double>();
106  integrator->setStepper(stepper);
107  integrator->setTimeStepControl(timeStepControl);
108  integrator->setSolutionHistory(solutionHistory);
109  // integrator->setObserver(...);
110  integrator->initialize();
111 
112  // Integrate to timeMax
113  bool integratorStatus = integrator->advanceTime();
114  TEST_ASSERT(integratorStatus)
115 
116  // Test if at 'Final Time'
117  double time = integrator->getTime();
118  double timeFinal = pl->sublist("Default Integrator")
119  .sublist("Time Step Control")
120  .get<double>("Final Time");
121  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
122 
123  // Time-integrated solution and the exact solution
124  RCP<Thyra::VectorBase<double>> x = integrator->getX();
126  model->getExactSolution(time).get_x();
127 
128  // Calculate the error
129  RCP<Thyra::VectorBase<double>> xdiff = x->clone_v();
130  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
131 
132  // Check the order and intercept
133  out << " Stepper = " << stepper->description() << "\n with "
134  << option << std::endl;
135  out << " =========================" << std::endl;
136  out << " Exact solution : " << get_ele(*(x_exact), 0) << std::endl;
137  out << " Computed solution: " << get_ele(*(x), 0) << std::endl;
138  out << " Difference : " << get_ele(*(xdiff), 0) << std::endl;
139  out << " =========================" << std::endl;
140  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.167158, 1.0e-4);
141  }
142 }
143 
144 // ************************************************************
145 // ************************************************************
146 TEUCHOS_UNIT_TEST(Leapfrog, SinCos)
147 {
149  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
150  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
151  std::vector<double> StepSize;
152  std::vector<double> xErrorNorm;
153  std::vector<double> xDotErrorNorm;
154  const int nTimeStepSizes = 9;
155  double time = 0.0;
156 
157  // Read params from .xml file
158  RCP<ParameterList> pList =
159  getParametersFromXmlFile("Tempus_Leapfrog_SinCos.xml");
160 
161  // Setup the HarmonicOscillatorModel
162  RCP<ParameterList> hom_pl = sublist(pList, "HarmonicOscillatorModel", true);
163  auto model = rcp(new HarmonicOscillatorModel<double>(hom_pl));
164 
165  // Setup the Integrator and reset initial time step
166  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
167 
168  // Set initial time step = 2*dt specified in input file (for convergence
169  // study)
170  double dt = pl->sublist("Default Integrator")
171  .sublist("Time Step Control")
172  .get<double>("Initial Time Step");
173  dt *= 2.0;
174 
175  for (int n = 0; n < nTimeStepSizes; n++) {
176  // Perform time-step refinement
177  dt /= 2;
178  out << "\n \n time step #" << n << " (out of " << nTimeStepSizes - 1
179  << "), dt = " << dt << "\n";
180  pl->sublist("Default Integrator")
181  .sublist("Time Step Control")
182  .set("Initial Time Step", dt);
183  integrator = Tempus::createIntegratorBasic<double>(pl, model);
184 
185  // Integrate to timeMax
186  bool integratorStatus = integrator->advanceTime();
187  TEST_ASSERT(integratorStatus)
188 
189  // Test if at 'Final Time'
190  time = integrator->getTime();
191  double timeFinal = pl->sublist("Default Integrator")
192  .sublist("Time Step Control")
193  .get<double>("Final Time");
194  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
195 
196  // Plot sample solution and exact solution at most-refined resolution
197  if (n == nTimeStepSizes - 1) {
198  RCP<const SolutionHistory<double>> solutionHistory =
199  integrator->getSolutionHistory();
200  writeSolution("Tempus_Leapfrog_SinCos.dat", solutionHistory);
201 
202  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
203  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
204  double time_i = (*solutionHistory)[i]->getTime();
205  auto state = Tempus::createSolutionStateX(
206  rcp_const_cast<Thyra::VectorBase<double>>(
207  model->getExactSolution(time_i).get_x()),
208  rcp_const_cast<Thyra::VectorBase<double>>(
209  model->getExactSolution(time_i).get_x_dot()));
210  state->setTime((*solutionHistory)[i]->getTime());
211  solnHistExact->addState(state);
212  }
213  writeSolution("Tempus_Leapfrog_SinCos-Ref.dat", solnHistExact);
214  }
215 
216  // Store off the final solution and step size
217 
218  StepSize.push_back(dt);
219  auto solution = Thyra::createMember(model->get_x_space());
220  Thyra::copy(*(integrator->getX()), solution.ptr());
221  solutions.push_back(solution);
222  auto solutionDot = Thyra::createMember(model->get_x_space());
223  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
224  solutionsDot.push_back(solutionDot);
225  if (n == nTimeStepSizes - 1) { // Add exact solution last in vector.
226  StepSize.push_back(0.0);
227  auto solutionExact = Thyra::createMember(model->get_x_space());
228  Thyra::copy(*(model->getExactSolution(time).get_x()),
229  solutionExact.ptr());
230  solutions.push_back(solutionExact);
231  auto solutionDotExact = Thyra::createMember(model->get_x_space());
232  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
233  solutionDotExact.ptr());
234  solutionsDot.push_back(solutionDotExact);
235  }
236  }
237 
238  // Check the order and intercept
239  double xSlope = 0.0;
240  double xDotSlope = 0.0;
241  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
242  double order = stepper->getOrder();
243  writeOrderError("Tempus_Leapfrog_SinCos-Error.dat", stepper, StepSize,
244  solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
245  xDotSlope, out);
246 
247  TEST_FLOATING_EQUALITY(xSlope, order, 0.02);
248  TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.0157928, 1.0e-4);
249  TEST_FLOATING_EQUALITY(xDotSlope, 1.09387, 0.01);
250  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.563002, 1.0e-4);
251 
253 }
254 
255 } // namespace Tempus_Test
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.
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
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar >> solutionHistory)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Consider the ODE: where is a constant, is a constant damping parameter, is a constant forcing par...
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)
Ptr< T > ptr() const
TimeStepControl manages the time step size. There several mechanisms that effect the time step size a...
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Keep a fix number of states.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Solution state for integrators and steppers.