Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_Test_NewmarkImplicitAForm_HarmonicOscillator_Damped_FirstOrder.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(NewmarkImplicitAForm, HarmonicOscillatorDamped_FirstOrder)
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 = 10;
52  double time = 0.0;
53 
54  // Read params from .xml file
55  RCP<ParameterList> pList = getParametersFromXmlFile(
56  "Tempus_Test_NewmarkImplicitAForm_HarmonicOscillator_Damped_FirstOrder."
57  "xml");
58 
59  // Setup the HarmonicOscillatorModel
60  RCP<ParameterList> hom_pl = sublist(pList, "HarmonicOscillatorModel", true);
63 
64  // Setup the Integrator and reset initial time step
65  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
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();
101  "Tempus_Test_NewmarkImplicitAForm_HarmonicOscillator_Damped_"
102  "FirstOrder.dat",
103  solutionHistory);
104 
105  RCP<Tempus::SolutionHistory<double>> solnHistExact =
107  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
108  double time_i = (*solutionHistory)[i]->getTime();
110  rcp_const_cast<Thyra::VectorBase<double>>(
111  model->getExactSolution(time_i).get_x()),
112  rcp_const_cast<Thyra::VectorBase<double>>(
113  model->getExactSolution(time_i).get_x_dot()));
114  state->setTime((*solutionHistory)[i]->getTime());
115  solnHistExact->addState(state);
116  }
118  "Tempus_Test_NewmarkImplicitAForm_HarmonicOscillator_Damped_"
119  "FirstOrder-Ref.dat",
120  solnHistExact);
121  }
122 
123  // Store off the final solution and step size
124  StepSize.push_back(dt);
125  auto solution = Thyra::createMember(model->get_x_space());
126  Thyra::copy(*(integrator->getX()), solution.ptr());
127  solutions.push_back(solution);
128  auto solutionDot = Thyra::createMember(model->get_x_space());
129  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
130  solutionsDot.push_back(solutionDot);
131  if (n == nTimeStepSizes - 1) { // Add exact solution last in vector.
132  StepSize.push_back(0.0);
133  auto solutionExact = Thyra::createMember(model->get_x_space());
134  Thyra::copy(*(model->getExactSolution(time).get_x()),
135  solutionExact.ptr());
136  solutions.push_back(solutionExact);
137  auto solutionDotExact = Thyra::createMember(model->get_x_space());
138  Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
139  solutionDotExact.ptr());
140  solutionsDot.push_back(solutionDotExact);
141  }
142  }
143 
144  // Check the order and intercept
145  double xSlope = 0.0;
146  double xDotSlope = 0.0;
147  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
148  double order = stepper->getOrder();
150  "Tempus_Test_NewmarkImplicitAForm_HarmonicOscillator_Damped_FirstOrder-"
151  "Error.dat",
152  stepper, StepSize, solutions, xErrorNorm, xSlope, solutionsDot,
153  xDotErrorNorm, xDotSlope, out);
154 
155  TEST_FLOATING_EQUALITY(xSlope, order, 0.02);
156  TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.0224726, 1.0e-4);
157  TEST_FLOATING_EQUALITY(xDotSlope, order, 0.02);
158  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.0122223, 1.0e-4);
159 
161 }
162 
163 } // 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)
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.