Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_Test_BDF2_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
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/SinCosModel.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, SinCos)
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
49 // Read params from .xml file
50 RCP<ParameterList> pList = getParametersFromXmlFile("Tempus_BDF2_SinCos.xml");
51 // Set initial time step = 2*dt specified in input file (for convergence
52 // study)
53 double dt = pList->sublist("Tempus")
54 .sublist("Default Integrator")
55 .sublist("Time Step Control")
56 .get<double>("Initial Time Step");
57 dt *= 2.0;
58
59 // Setup the SinCosModel
60 RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
61 const int nTimeStepSizes = scm_pl->get<int>("Number of Time Step Sizes", 7);
62 std::string output_file_string =
63 scm_pl->get<std::string>("Output File Name", "Tempus_BDF2_SinCos");
64 std::string output_file_name = output_file_string + ".dat";
65 std::string ref_out_file_name = output_file_string + "-Ref.dat";
66 std::string err_out_file_name = output_file_string + "-Error.dat";
67 double time = 0.0;
68 for (int n = 0; n < nTimeStepSizes; n++) {
69 auto model = rcp(new SinCosModel<double>(scm_pl));
70
71 dt /= 2;
72
73 // Setup the Integrator and reset initial time step
74 RCP<ParameterList> tempusPL =
75 getParametersFromXmlFile("Tempus_BDF2_SinCos.xml");
76 RCP<ParameterList> pl = sublist(tempusPL, "Tempus", true);
77 pl->sublist("Default Integrator")
78 .sublist("Time Step Control")
79 .set("Initial Time Step", dt);
80 integrator = Tempus::createIntegratorBasic<double>(pl, model);
81
82 // Initial Conditions
83 // During the Integrator construction, the initial SolutionState
84 // is set by default to model->getNominalVales().get_x(). However,
85 // the application can set it also by integrator->initializeSolutionHistory.
86 RCP<Thyra::VectorBase<double>> x0 =
87 model->getNominalValues().get_x()->clone_v();
88 integrator->initializeSolutionHistory(0.0, x0);
89
90 // Integrate to timeMax
91 bool integratorStatus = integrator->advanceTime();
92 TEST_ASSERT(integratorStatus)
93
94 // Test if at 'Final Time'
95 time = integrator->getTime();
96 double timeFinal = pl->sublist("Default Integrator")
97 .sublist("Time Step Control")
98 .get<double>("Final Time");
99 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
100
101 // Plot sample solution and exact solution
102 if (n == 0) {
103 RCP<const SolutionHistory<double>> solutionHistory =
104 integrator->getSolutionHistory();
105 writeSolution(output_file_name, solutionHistory);
106
107 auto solnHistExact = rcp(new Tempus::SolutionHistory<double>());
108 for (int i = 0; i < solutionHistory->getNumStates(); i++) {
109 double time_i = (*solutionHistory)[i]->getTime();
110 auto state = Tempus::createSolutionStateX(
111 rcp_const_cast<Thyra::VectorBase<double>>(
112 model->getExactSolution(time_i).get_x()),
113 rcp_const_cast<Thyra::VectorBase<double>>(
114 model->getExactSolution(time_i).get_x_dot()));
115 state->setTime((*solutionHistory)[i]->getTime());
116 solnHistExact->addState(state);
117 }
118 writeSolution(ref_out_file_name, solnHistExact);
119 }
120
121 // Store off the final solution and step size
122 StepSize.push_back(dt);
123 auto solution = Thyra::createMember(model->get_x_space());
124 Thyra::copy(*(integrator->getX()), solution.ptr());
125 solutions.push_back(solution);
126 auto solutionDot = Thyra::createMember(model->get_x_space());
127 Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
128 solutionsDot.push_back(solutionDot);
129 if (n == nTimeStepSizes - 1) { // Add exact solution last in vector.
130 StepSize.push_back(0.0);
131 auto solutionExact = Thyra::createMember(model->get_x_space());
132 Thyra::copy(*(model->getExactSolution(time).get_x()),
133 solutionExact.ptr());
134 solutions.push_back(solutionExact);
135 auto solutionDotExact = Thyra::createMember(model->get_x_space());
136 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
137 solutionDotExact.ptr());
138 solutionsDot.push_back(solutionDotExact);
139 }
140 }
141
142 // Check the order and intercept
143 if (nTimeStepSizes > 1) {
144 double xSlope = 0.0;
145 double xDotSlope = 0.0;
146 std::vector<double> xErrorNorm;
147 std::vector<double> xDotErrorNorm;
148 RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
149 double order = stepper->getOrder();
150 writeOrderError(err_out_file_name, stepper, StepSize, solutions, xErrorNorm,
151 xSlope, solutionsDot, xDotErrorNorm, xDotSlope, out);
152
153 TEST_FLOATING_EQUALITY(xSlope, order, 0.01);
154 TEST_FLOATING_EQUALITY(xDotSlope, order, 0.01);
155 TEST_FLOATING_EQUALITY(xErrorNorm[0], 5.13425e-05, 1.0e-4);
156 TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 5.13425e-05, 1.0e-4);
157 }
158
159 Teuchos::TimeMonitor::summarize();
160}
161
162} // namespace Tempus_Test
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers.
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
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)
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.