Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_Test_BackwardEuler_VanDerPol.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/VanDerPolModel.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, VanDerPol)
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 = 4;
50  double dt = 0.05;
51  for (int n = 0; n < nTimeStepSizes; n++) {
52  // Read params from .xml file
53  RCP<ParameterList> pList =
54  getParametersFromXmlFile("Tempus_BackwardEuler_VanDerPol.xml");
55 
56  // Setup the VanDerPolModel
57  RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
58  auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
59 
60  // Set the step size
61  dt /= 2;
62  if (n == nTimeStepSizes - 1) dt /= 10.0;
63 
64  // Setup the Integrator and reset initial time step
65  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
66  pl->sublist("Demo Integrator")
67  .sublist("Time Step Control")
68  .set("Initial Time Step", dt);
69  integrator = Tempus::createIntegratorBasic<double>(pl, model);
70 
71  // Integrate to timeMax
72  bool integratorStatus = integrator->advanceTime();
73  TEST_ASSERT(integratorStatus)
74 
75  // Test if at 'Final Time'
76  double time = integrator->getTime();
77  double timeFinal = pl->sublist("Demo Integrator")
78  .sublist("Time Step Control")
79  .get<double>("Final Time");
80  double tol = 100.0 * std::numeric_limits<double>::epsilon();
81  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
82 
83  // Store off the final solution and step size
84  StepSize.push_back(dt);
85  auto solution = Thyra::createMember(model->get_x_space());
86  Thyra::copy(*(integrator->getX()), solution.ptr());
87  solutions.push_back(solution);
88  auto solutionDot = Thyra::createMember(model->get_x_space());
89  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
90  solutionsDot.push_back(solutionDot);
91 
92  // Output finest temporal solution for plotting
93  // This only works for ONE MPI process
94  if ((n == 0) || (n == nTimeStepSizes - 1)) {
95  std::string fname = "Tempus_BackwardEuler_VanDerPol-Ref.dat";
96  if (n == 0) fname = "Tempus_BackwardEuler_VanDerPol.dat";
97  RCP<const SolutionHistory<double>> solutionHistory =
98  integrator->getSolutionHistory();
99  writeSolution(fname, solutionHistory);
100  }
101  }
102 
103  // Check the order and intercept
104  double xSlope = 0.0;
105  double xDotSlope = 0.0;
106  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
107  double order = stepper->getOrder();
108  writeOrderError("Tempus_BackwardEuler_VanDerPol-Error.dat", stepper, StepSize,
109  solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
110  xDotSlope, out);
111 
112  TEST_FLOATING_EQUALITY(xSlope, order, 0.10);
113  TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.571031, 1.0e-4);
114  TEST_FLOATING_EQUALITY(xDotSlope, 1.74898, 0.10);
115  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 1.0038, 1.0e-4);
116  // At small dt, slopes should be equal to order.
117  // TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
118 
120 }
121 
122 } // namespace Tempus_Test
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)
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...
van der Pol model problem for nonlinear electrical circuit.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Solution state for integrators and steppers.