Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_Test_BDF2_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_StepperBDF2.hpp"
18
19#include "../TestModels/VanDerPolModel.hpp"
20#include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
21
22#include <fstream>
23#include <limits>
24#include <sstream>
25#include <vector>
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(BDF2, VanDerPol)
43{
44 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
45 std::vector<double> StepSize;
46 std::vector<double> ErrorNorm;
47
48 // Read params from .xml file
49 RCP<ParameterList> pList =
50 getParametersFromXmlFile("Tempus_BDF2_VanDerPol.xml");
51 // Set initial time step = 2*dt specified in input file (for convergence
52 // study)
53 //
54 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
55 double dt = pl->sublist("Demo Integrator")
56 .sublist("Time Step Control")
57 .get<double>("Initial Time Step");
58 dt *= 2.0;
59
60 RCP<ParameterList> vdpm_pl = sublist(pList, "VanDerPolModel", true);
61 const int nTimeStepSizes = vdpm_pl->get<int>("Number of Time Step Sizes", 3);
62 // const int nTimeStepSizes = 5;
63 double order = 0.0;
64
65 for (int n = 0; n < nTimeStepSizes; n++) {
66 // Setup the VanDerPolModel
67 auto model = rcp(new VanDerPolModel<double>(vdpm_pl));
68
69 // Set the step size
70 dt /= 2;
71 if (n == nTimeStepSizes - 1) dt /= 10.0;
72
73 // Setup the Integrator and reset initial time step
74 pl->sublist("Demo Integrator")
75 .sublist("Time Step Control")
76 .set("Initial Time Step", dt);
77 RCP<Tempus::IntegratorBasic<double>> integrator =
78 Tempus::createIntegratorBasic<double>(pl, model);
79 order = integrator->getStepper()->getOrder();
80
81 // Integrate to timeMax
82 bool integratorStatus = integrator->advanceTime();
83 TEST_ASSERT(integratorStatus)
84
85 // Test if at 'Final Time'
86 double time = integrator->getTime();
87 double timeFinal = pl->sublist("Demo Integrator")
88 .sublist("Time Step Control")
89 .get<double>("Final Time");
90 double tol = 100.0 * std::numeric_limits<double>::epsilon();
91 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
92
93 // Store off the final solution and step size
94 auto solution = Thyra::createMember(model->get_x_space());
95 Thyra::copy(*(integrator->getX()), solution.ptr());
96 solutions.push_back(solution);
97 StepSize.push_back(dt);
98
99 // Output finest temporal solution for plotting
100 // This only works for ONE MPI process
101 if ((n == 0) || (n == nTimeStepSizes - 1)) {
102 std::string fname = "Tempus_BDF2_VanDerPol-Ref.dat";
103 if (n == 0) fname = "Tempus_BDF2_VanDerPol.dat";
104 std::ofstream ftmp(fname);
105 RCP<const SolutionHistory<double>> solutionHistory =
106 integrator->getSolutionHistory();
107 int nStates = solutionHistory->getNumStates();
108 for (int i = 0; i < nStates; i++) {
109 RCP<const SolutionState<double>> solutionState = (*solutionHistory)[i];
110 RCP<const Thyra::VectorBase<double>> x = solutionState->getX();
111 double ttime = solutionState->getTime();
112 ftmp << ttime << " " << get_ele(*x, 0) << " " << get_ele(*x, 1)
113 << std::endl;
114 }
115 ftmp.close();
116 }
117 }
118
119 // Calculate the error - use the most temporally refined mesh for
120 // the reference solution.
121 auto ref_solution = solutions[solutions.size() - 1];
122 std::vector<double> StepSizeCheck;
123 for (std::size_t i = 0; i < (solutions.size() - 1); ++i) {
124 auto tmp = solutions[i];
125 Thyra::Vp_StV(tmp.ptr(), -1.0, *ref_solution);
126 const double L2norm = Thyra::norm_2(*tmp);
127 StepSizeCheck.push_back(StepSize[i]);
128 ErrorNorm.push_back(L2norm);
129 }
130
131 if (nTimeStepSizes > 2) {
132 // Check the order and intercept
133 double slope =
134 computeLinearRegressionLogLog<double>(StepSizeCheck, ErrorNorm);
135 out << " Stepper = BDF2" << std::endl;
136 out << " =========================" << std::endl;
137 out << " Expected order: " << order << std::endl;
138 out << " Observed order: " << slope << std::endl;
139 out << " =========================" << std::endl;
140 TEST_FLOATING_EQUALITY(slope, order, 0.10);
141 out << "\n\n ** Slope on BDF2 Method = " << slope << "\n"
142 << std::endl;
143 }
144
145 // Write error data
146 {
147 std::ofstream ftmp("Tempus_BDF2_VanDerPol-Error.dat");
148 double error0 = 0.8 * ErrorNorm[0];
149 for (std::size_t n = 0; n < StepSizeCheck.size(); n++) {
150 ftmp << StepSizeCheck[n] << " " << ErrorNorm[n] << " "
151 << error0 * (pow(StepSize[n] / StepSize[0], order)) << std::endl;
152 }
153 ftmp.close();
154 }
155
156 Teuchos::TimeMonitor::summarize();
157}
158
159} // 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.
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)