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