Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_Test_BackwardEuler_SinCos.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_StepperBackwardEuler.hpp"
17 
18 #include "../TestModels/SinCosModel.hpp"
19 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
20 
21 #include <vector>
22 #include <fstream>
23 #include <sstream>
24 #include <limits>
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(BackwardEuler, SinCos)
42 {
44  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
45  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
46  std::vector<double> StepSize;
47  std::vector<double> xErrorNorm;
48  std::vector<double> xDotErrorNorm;
49  const int nTimeStepSizes = 7;
50  double dt = 0.2;
51  double time = 0.0;
52  for (int n = 0; n < nTimeStepSizes; n++) {
53  // Read params from .xml file
54  RCP<ParameterList> pList =
55  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
56 
57  // std::ofstream ftmp("PL.txt");
58  // pList->print(ftmp);
59  // ftmp.close();
60 
61  // Setup the SinCosModel
62  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
63  // RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
64  auto model = rcp(new SinCosModel<double>(scm_pl));
65 
66  dt /= 2;
67 
68  // Setup the Integrator and reset initial time step
69  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
70  pl->sublist("Default Integrator")
71  .sublist("Time Step Control")
72  .set("Initial Time Step", dt);
73  integrator = Tempus::createIntegratorBasic<double>(pl, model);
74 
75  // Initial Conditions
76  // During the Integrator construction, the initial SolutionState
77  // is set by default to model->getNominalVales().get_x(). However,
78  // the application can set it also by integrator->initializeSolutionHistory.
80  model->getNominalValues().get_x()->clone_v();
81  integrator->initializeSolutionHistory(0.0, x0);
82 
83  // Integrate to timeMax
84  bool integratorStatus = integrator->advanceTime();
85  TEST_ASSERT(integratorStatus)
86 
87  // Test if at 'Final Time'
88  time = integrator->getTime();
89  double timeFinal = pl->sublist("Default Integrator")
90  .sublist("Time Step Control")
91  .get<double>("Final Time");
92  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
93 
94  // Plot sample solution and exact solution
95  if (n == 0) {
96  RCP<const SolutionHistory<double>> solutionHistory =
97  integrator->getSolutionHistory();
98  writeSolution("Tempus_BackwardEuler_SinCos.dat", solutionHistory);
99 
100  auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
101  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
102  double time_i = (*solutionHistory)[i]->getTime();
103  auto state = Tempus::createSolutionStateX(
104  rcp_const_cast<Thyra::VectorBase<double>>(
105  model->getExactSolution(time_i).get_x()),
106  rcp_const_cast<Thyra::VectorBase<double>>(
107  model->getExactSolution(time_i).get_x_dot()));
108  state->setTime((*solutionHistory)[i]->getTime());
109  solnHistExact->addState(state);
110  }
111  writeSolution("Tempus_BackwardEuler_SinCos-Ref.dat", solnHistExact);
112  }
113 
114  // Store off the final solution and step size
115  StepSize.push_back(dt);
116  auto solution = Thyra::createMember(model->get_x_space());
117  Thyra::copy(*(integrator->getX()), solution.ptr());
118  solutions.push_back(solution);
119  auto solutionDot = Thyra::createMember(model->get_x_space());
120  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
121  solutionsDot.push_back(solutionDot);
122  if (n == nTimeStepSizes - 1) { // Add exact solution last in vector.
123  StepSize.push_back(0.0);
124  auto solutionExact = Thyra::createMember(model->get_x_space());
125  Thyra::copy(*(model->getExactSolution(time).get_x()),
126  solutionExact.ptr());
127  solutions.push_back(solutionExact);
128  auto solutionDotExact = Thyra::createMember(model->get_x_space());
129  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
130  solutionDotExact.ptr());
131  solutionsDot.push_back(solutionDotExact);
132  }
133  }
134 
135  // Check the order and intercept
136  double xSlope = 0.0;
137  double xDotSlope = 0.0;
138  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
139  double order = stepper->getOrder();
140  writeOrderError("Tempus_BackwardEuler_SinCos-Error.dat", stepper, StepSize,
141  solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
142  xDotSlope, out);
143 
144  TEST_FLOATING_EQUALITY(xSlope, order, 0.01);
145  TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.0486418, 1.0e-4);
146  TEST_FLOATING_EQUALITY(xDotSlope, order, 0.01);
147  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.0486418, 1.0e-4);
148 
150 }
151 
152 } // 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.
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...
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar >> solutionHistory)
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.