Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_BDF2_ASA.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_IntegratorAdjointSensitivity.hpp"
18
19#include "Thyra_VectorStdOps.hpp"
20#include "Thyra_MultiVectorStdOps.hpp"
21
22#include "../TestModels/SinCosModel.hpp"
23#include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
24
25#include "Thyra_DefaultMultiVectorProductVector.hpp"
26
27#include <fstream>
28#include <limits>
29#include <sstream>
30#include <vector>
31
32namespace Tempus_Test {
33
34using Teuchos::getParametersFromXmlFile;
35using Teuchos::ParameterList;
36using Teuchos::RCP;
37using Teuchos::sublist;
38
42
43// ************************************************************
44// ************************************************************
45TEUCHOS_UNIT_TEST(BDF2, SinCos_ASA)
46{
47 std::vector<double> StepSize;
48 std::vector<double> ErrorNorm;
49 const int nTimeStepSizes = 7;
50 double dt = 0.2;
51 double order = 0.0;
52 Teuchos::RCP<const Teuchos::Comm<int> > comm =
53 Teuchos::DefaultComm<int>::getComm();
54 for (int n = 0; n < nTimeStepSizes; n++) {
55 // Read params from .xml file
56 RCP<ParameterList> pList =
57 getParametersFromXmlFile("Tempus_BDF2_SinCos_SA.xml");
58
59 // Setup the SinCosModel
60 RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
61 RCP<SinCosModel<double> > model =
62 Teuchos::rcp(new SinCosModel<double>(scm_pl));
63
64 dt /= 2;
65
66 // Setup sensitivities
67 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
68 // ParameterList& sens_pl = pl->sublist("Sensitivities");
69 ParameterList& interp_pl = pl->sublist("Default Integrator")
70 .sublist("Solution History")
71 .sublist("Interpolator");
72 interp_pl.set("Interpolator Type", "Lagrange");
73 interp_pl.set("Order", 1);
74
75 // Setup the Integrator and reset initial time step
76 pl->sublist("Default Integrator")
77 .sublist("Time Step Control")
78 .set("Initial Time Step", dt);
79 RCP<Tempus::IntegratorAdjointSensitivity<double> > integrator =
80 Tempus::createIntegratorAdjointSensitivity<double>(pl, model);
81 order = integrator->getStepper()->getOrder();
82
83 // Initial Conditions
84 double t0 = pl->sublist("Default Integrator")
85 .sublist("Time Step Control")
86 .get<double>("Initial Time");
87 RCP<const Thyra::VectorBase<double> > x0 =
88 model->getExactSolution(t0).get_x();
89 const int num_param = model->get_p_space(0)->dim();
90 RCP<Thyra::MultiVectorBase<double> > DxDp0 =
91 Thyra::createMembers(model->get_x_space(), num_param);
92 for (int i = 0; i < num_param; ++i)
93 Thyra::assign(DxDp0->col(i).ptr(),
94 *(model->getExactSensSolution(i, t0).get_x()));
95 integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
96 DxDp0, Teuchos::null, Teuchos::null);
97
98 // Integrate to timeMax
99 bool integratorStatus = integrator->advanceTime();
100 TEST_ASSERT(integratorStatus)
101
102 // Test if at 'Final Time'
103 double time = integrator->getTime();
104 double timeFinal = pl->sublist("Default Integrator")
105 .sublist("Time Step Control")
106 .get<double>("Final Time");
107 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
108
109 // Time-integrated solution and the exact solution along with
110 // sensitivities (relying on response g(x) = x). Note we must transpose
111 // dg/dp since the integrator returns it in gradient form.
112 RCP<const Thyra::VectorBase<double> > x = integrator->getX();
113 RCP<const Thyra::MultiVectorBase<double> > DgDp = integrator->getDgDp();
114 RCP<Thyra::MultiVectorBase<double> > DxDp =
115 Thyra::createMembers(model->get_x_space(), num_param);
116 {
117 Thyra::ConstDetachedMultiVectorView<double> dgdp_view(*DgDp);
118 Thyra::DetachedMultiVectorView<double> dxdp_view(*DxDp);
119 const int num_g = DgDp->domain()->dim();
120 for (int i = 0; i < num_g; ++i)
121 for (int j = 0; j < num_param; ++j) dxdp_view(i, j) = dgdp_view(j, i);
122 }
123 RCP<const Thyra::VectorBase<double> > x_exact =
124 model->getExactSolution(time).get_x();
125 RCP<Thyra::MultiVectorBase<double> > DxDp_exact =
126 Thyra::createMembers(model->get_x_space(), num_param);
127 for (int i = 0; i < num_param; ++i)
128 Thyra::assign(DxDp_exact->col(i).ptr(),
129 *(model->getExactSensSolution(i, time).get_x()));
130
131 // Plot sample solution and exact solution
132 if (comm->getRank() == 0 && n == nTimeStepSizes - 1) {
133 typedef Thyra::DefaultProductVector<double> DPV;
134 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
135
136 std::ofstream ftmp("Tempus_BDF2_SinCos_AdjSens.dat");
137 RCP<const SolutionHistory<double> > solutionHistory =
138 integrator->getSolutionHistory();
139 for (int i = 0; i < solutionHistory->getNumStates(); i++) {
140 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
141 const double time_i = solutionState->getTime();
142 RCP<const DPV> x_prod_plot =
143 Teuchos::rcp_dynamic_cast<const DPV>(solutionState->getX());
144 RCP<const Thyra::VectorBase<double> > x_plot =
145 x_prod_plot->getVectorBlock(0);
146 RCP<const DMVPV> adjoint_prod_plot =
147 Teuchos::rcp_dynamic_cast<const DMVPV>(
148 x_prod_plot->getVectorBlock(1));
149 RCP<const Thyra::MultiVectorBase<double> > adjoint_plot =
150 adjoint_prod_plot->getMultiVector();
151 RCP<const Thyra::VectorBase<double> > x_exact_plot =
152 model->getExactSolution(time_i).get_x();
153 ftmp << std::fixed << std::setprecision(7) << time_i << std::setw(11)
154 << get_ele(*(x_plot), 0) << std::setw(11) << get_ele(*(x_plot), 1)
155 << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 0)
156 << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 1)
157 << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 0)
158 << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 1)
159 << std::setw(11) << get_ele(*(x_exact_plot), 0) << std::setw(11)
160 << get_ele(*(x_exact_plot), 1) << std::endl;
161 }
162 ftmp.close();
163 }
164
165 // Calculate the error
166 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
167 RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
168 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
169 Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
170 StepSize.push_back(dt);
171 double L2norm = Thyra::norm_2(*xdiff);
172 L2norm *= L2norm;
173 Teuchos::Array<double> L2norm_DxDp(num_param);
174 Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
175 for (int i = 0; i < num_param; ++i)
176 L2norm += L2norm_DxDp[i] * L2norm_DxDp[i];
177 L2norm = std::sqrt(L2norm);
178 ErrorNorm.push_back(L2norm);
179
180 // out << " n = " << n << " dt = " << dt << " error = " << L2norm
181 // << std::endl;
182 }
183
184 // Check the order and intercept
185 double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
186 out << " Stepper = BDF2" << std::endl;
187 out << " =========================" << std::endl;
188 out << " Expected order: " << order << std::endl;
189 out << " Observed order: " << slope << std::endl;
190 out << " =========================" << std::endl;
191 TEST_FLOATING_EQUALITY(slope, 1.95006, 0.015);
192 TEST_FLOATING_EQUALITY(ErrorNorm[0], 0.0137394, 1.0e-4);
193
194 if (comm->getRank() == 0) {
195 std::ofstream ftmp("Tempus_BDF2_SinCos_AdjSens-Error.dat");
196 double error0 = 0.8 * ErrorNorm[0];
197 for (int n = 0; n < nTimeStepSizes; n++) {
198 ftmp << StepSize[n] << " " << ErrorNorm[n] << " "
199 << error0 * (StepSize[n] / StepSize[0]) << std::endl;
200 }
201 ftmp.close();
202 }
203}
204
205} // 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.
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)