Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_Test_NewmarkExplicitAForm_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 
13 #include "Tempus_config.hpp"
14 #include "Tempus_IntegratorBasic.hpp"
15 
16 #include "Tempus_StepperFactory.hpp"
17 #include "Tempus_StepperNewmarkImplicitAForm.hpp"
18 #include "Tempus_StepperNewmarkImplicitDForm.hpp"
19 
20 #include "../TestModels/HarmonicOscillatorModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
22 
23 #include <fstream>
24 #include <limits>
25 #include <sstream>
26 #include <vector>
27 
28 namespace Tempus_Test {
29 
30 using Teuchos::getParametersFromXmlFile;
32 using Teuchos::RCP;
33 using Teuchos::rcp_const_cast;
34 using Teuchos::rcp_dynamic_cast;
35 using Teuchos::sublist;
36 
40 
41 // ************************************************************
42 TEUCHOS_UNIT_TEST(NewmarkExplicitAForm, SinCos)
43 {
45  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
46  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
47  std::vector<double> StepSize;
48  std::vector<double> xErrorNorm;
49  std::vector<double> xDotErrorNorm;
50  const int nTimeStepSizes = 9;
51  double time = 0.0;
52 
53  // Read params from .xml file
54  RCP<ParameterList> pList =
55  getParametersFromXmlFile("Tempus_Test_NewmarkExplicitAForm_SinCos.xml");
56 
57  // Setup the HarmonicOscillatorModel
58  RCP<ParameterList> hom_pl = sublist(pList, "HarmonicOscillatorModel", true);
61 
62  // Setup the Integrator and reset initial time step
63  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
64  RCP<ParameterList> stepperPL = sublist(pl, "Default Stepper", true);
65  stepperPL->remove("Zero Initial Guess");
66 
67  // Set initial time step = 2*dt specified in input file (for convergence
68  // study)
69  //
70  double dt = pl->sublist("Default Integrator")
71  .sublist("Time Step Control")
72  .get<double>("Initial Time Step");
73  dt *= 2.0;
74 
75  for (int n = 0; n < nTimeStepSizes; n++) {
76  // Perform time-step refinement
77  dt /= 2;
78  out << "\n \n time step #" << n << " (out of " << nTimeStepSizes - 1
79  << "), dt = " << dt << "\n";
80  pl->sublist("Default Integrator")
81  .sublist("Time Step Control")
82  .set("Initial Time Step", dt);
83  integrator = Tempus::createIntegratorBasic<double>(pl, model);
84 
85  // Integrate to timeMax
86  bool integratorStatus = integrator->advanceTime();
87  TEST_ASSERT(integratorStatus)
88 
89  // Test if at 'Final Time'
90  time = integrator->getTime();
91  double timeFinal = pl->sublist("Default Integrator")
92  .sublist("Time Step Control")
93  .get<double>("Final Time");
94  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
95 
96  // Plot sample solution and exact solution
97  if (n == 0) {
98  RCP<const SolutionHistory<double>> solutionHistory =
99  integrator->getSolutionHistory();
100  writeSolution("Tempus_Test_NewmarkExplicitAForm_SinCos.dat",
101  solutionHistory);
102 
103  RCP<Tempus::SolutionHistory<double>> solnHistExact =
105  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
106  double time_i = (*solutionHistory)[i]->getTime();
108  rcp_const_cast<Thyra::VectorBase<double>>(
109  model->getExactSolution(time_i).get_x()),
110  rcp_const_cast<Thyra::VectorBase<double>>(
111  model->getExactSolution(time_i).get_x_dot()));
112  state->setTime((*solutionHistory)[i]->getTime());
113  solnHistExact->addState(state);
114  }
115  writeSolution("Tempus_Test_NewmarkExplicitAForm_SinCos-Ref.dat",
116  solnHistExact);
117  }
118 
119  // Store off the final solution and step size
120  StepSize.push_back(dt);
121  auto solution = Thyra::createMember(model->get_x_space());
122  Thyra::copy(*(integrator->getX()), solution.ptr());
123  solutions.push_back(solution);
124  auto solutionDot = Thyra::createMember(model->get_x_space());
125  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
126  solutionsDot.push_back(solutionDot);
127  if (n == nTimeStepSizes - 1) { // Add exact solution last in vector.
128  StepSize.push_back(0.0);
129  auto solutionExact = Thyra::createMember(model->get_x_space());
130  Thyra::copy(*(model->getExactSolution(time).get_x()),
131  solutionExact.ptr());
132  solutions.push_back(solutionExact);
133  auto solutionDotExact = Thyra::createMember(model->get_x_space());
134  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
135  solutionDotExact.ptr());
136  solutionsDot.push_back(solutionDotExact);
137  }
138  }
139 
140  // Check the order and intercept
141  double xSlope = 0.0;
142  double xDotSlope = 0.0;
143  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
144  double order = stepper->getOrder();
145  writeOrderError("Tempus_Test_NewmarkExplicitAForm_SinCos-Error.dat", stepper,
146  StepSize, solutions, xErrorNorm, xSlope, solutionsDot,
147  xDotErrorNorm, xDotSlope, out);
148 
149  TEST_FLOATING_EQUALITY(xSlope, order, 0.02);
150  TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.0157928, 1.0e-4);
151  TEST_FLOATING_EQUALITY(xDotSlope, order, 0.02);
152  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.233045, 1.0e-4);
153 
155 }
156 
157 } // 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)
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
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar >> solutionHistory)
Consider the ODE: where is a constant, is a constant damping parameter, is a constant forcing par...
bool remove(std::string const &name, bool throwIfNotExists=true)
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.