Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_Test_BackwardEuler_VanDerPol.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
10#include "Teuchos_UnitTestHarness.hpp"
11#include "Teuchos_XMLParameterListHelpers.hpp"
12#include "Teuchos_TimeMonitor.hpp"
13#include "Teuchos_DefaultComm.hpp"
14
15#include "Tempus_config.hpp"
16#include "Tempus_IntegratorBasic.hpp"
17#include "Tempus_StepperBackwardEuler.hpp"
18
19#include "../TestModels/VanDerPolModel.hpp"
20#include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
21
22#include <vector>
23#include <fstream>
24#include <sstream>
25#include <limits>
26
27namespace Tempus_Test {
28
29using Teuchos::getParametersFromXmlFile;
30using Teuchos::ParameterList;
31using Teuchos::RCP;
32using Teuchos::rcp;
33using Teuchos::rcp_const_cast;
34using Teuchos::sublist;
35
39
40// ************************************************************
41// ************************************************************
42TEUCHOS_UNIT_TEST(BackwardEuler, VanDerPol)
43{
44 RCP<Tempus::IntegratorBasic<double>> integrator;
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 = 4;
51 double dt = 0.05;
52 for (int n = 0; n < nTimeStepSizes; n++) {
53 // Read params from .xml file
54 RCP<ParameterList> pList =
55 getParametersFromXmlFile("Tempus_BackwardEuler_VanDerPol.xml");
56
57 // Setup the VanDerPolModel
58 RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
59 auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
60
61 // Set the step size
62 dt /= 2;
63 if (n == nTimeStepSizes - 1) dt /= 10.0;
64
65 // Setup the Integrator and reset initial time step
66 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
67 pl->sublist("Demo Integrator")
68 .sublist("Time Step Control")
69 .set("Initial Time Step", dt);
70 integrator = Tempus::createIntegratorBasic<double>(pl, model);
71
72 // Integrate to timeMax
73 bool integratorStatus = integrator->advanceTime();
74 TEST_ASSERT(integratorStatus)
75
76 // Test if at 'Final Time'
77 double time = integrator->getTime();
78 double timeFinal = pl->sublist("Demo Integrator")
79 .sublist("Time Step Control")
80 .get<double>("Final Time");
81 double tol = 100.0 * std::numeric_limits<double>::epsilon();
82 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
83
84 // Store off the final solution and step size
85 StepSize.push_back(dt);
86 auto solution = Thyra::createMember(model->get_x_space());
87 Thyra::copy(*(integrator->getX()), solution.ptr());
88 solutions.push_back(solution);
89 auto solutionDot = Thyra::createMember(model->get_x_space());
90 Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
91 solutionsDot.push_back(solutionDot);
92
93 // Output finest temporal solution for plotting
94 // This only works for ONE MPI process
95 if ((n == 0) || (n == nTimeStepSizes - 1)) {
96 std::string fname = "Tempus_BackwardEuler_VanDerPol-Ref.dat";
97 if (n == 0) fname = "Tempus_BackwardEuler_VanDerPol.dat";
98 RCP<const SolutionHistory<double>> solutionHistory =
99 integrator->getSolutionHistory();
100 writeSolution(fname, solutionHistory);
101 }
102 }
103
104 // Check the order and intercept
105 double xSlope = 0.0;
106 double xDotSlope = 0.0;
107 RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
108 double order = stepper->getOrder();
109 writeOrderError("Tempus_BackwardEuler_VanDerPol-Error.dat", stepper, StepSize,
110 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
111 xDotSlope, out);
112
113 TEST_FLOATING_EQUALITY(xSlope, order, 0.10);
114 TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.571031, 1.0e-4);
115 TEST_FLOATING_EQUALITY(xDotSlope, 1.74898, 0.10);
116 TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 1.0038, 1.0e-4);
117 // At small dt, slopes should be equal to order.
118 // TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
119
120 Teuchos::TimeMonitor::summarize();
121}
122
123} // namespace Tempus_Test
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers.
van der Pol model problem for nonlinear electrical circuit.
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > solutionHistory)
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)
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)