Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_BDF2_FSA.hpp
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_IntegratorForwardSensitivity.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// ************************************************************
45void test_sincos_fsa(const bool use_combined_method,
46 const bool use_dfdp_as_tangent, Teuchos::FancyOStream& out,
47 bool& success)
48{
49 std::vector<double> StepSize;
50 std::vector<double> ErrorNorm;
51 const int nTimeStepSizes = 7;
52 double dt = 0.2;
53 double order = 0.0;
54 Teuchos::RCP<const Teuchos::Comm<int> > comm =
55 Teuchos::DefaultComm<int>::getComm();
56 for (int n = 0; n < nTimeStepSizes; n++) {
57 // Read params from .xml file
58 RCP<ParameterList> pList =
59 getParametersFromXmlFile("Tempus_BDF2_SinCos_SA.xml");
60
61 // Setup the SinCosModel
62 RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
63 scm_pl->set("Use DfDp as Tangent", use_dfdp_as_tangent);
64 RCP<SinCosModel<double> > model =
65 Teuchos::rcp(new SinCosModel<double>(scm_pl));
66
67 dt /= 2;
68
69 // Setup sensitivities
70 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
71 ParameterList& sens_pl = pl->sublist("Sensitivities");
72 if (use_combined_method)
73 sens_pl.set("Sensitivity Method", "Combined");
74 else {
75 sens_pl.set("Sensitivity Method", "Staggered");
76 // sens_pl.set("Reuse State Linear Solver", true);
77 }
78 sens_pl.set("Use DfDp as Tangent", use_dfdp_as_tangent);
79 ParameterList& interp_pl = pl->sublist("Default Integrator")
80 .sublist("Solution History")
81 .sublist("Interpolator");
82 interp_pl.set("Interpolator Type", "Lagrange");
83 interp_pl.set("Order", 1);
84
85 // Setup the Integrator and reset initial time step
86 pl->sublist("Default Integrator")
87 .sublist("Time Step Control")
88 .set("Initial Time Step", dt);
89 RCP<Tempus::IntegratorForwardSensitivity<double> > integrator =
90 Tempus::createIntegratorForwardSensitivity<double>(pl, model);
91 order = integrator->getStepper()->getOrder();
92
93 // Initial Conditions
94 double t0 = pl->sublist("Default Integrator")
95 .sublist("Time Step Control")
96 .get<double>("Initial Time");
97 RCP<const Thyra::VectorBase<double> > x0 =
98 model->getExactSolution(t0).get_x();
99 const int num_param = model->get_p_space(0)->dim();
100 RCP<Thyra::MultiVectorBase<double> > DxDp0 =
101 Thyra::createMembers(model->get_x_space(), num_param);
102 for (int i = 0; i < num_param; ++i)
103 Thyra::assign(DxDp0->col(i).ptr(),
104 *(model->getExactSensSolution(i, t0).get_x()));
105 integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
106 DxDp0, Teuchos::null, Teuchos::null);
107
108 // Integrate to timeMax
109 bool integratorStatus = integrator->advanceTime();
110 TEST_ASSERT(integratorStatus)
111
112 // Test if at 'Final Time'
113 double time = integrator->getTime();
114 double timeFinal = pl->sublist("Default Integrator")
115 .sublist("Time Step Control")
116 .get<double>("Final Time");
117 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
118
119 // Time-integrated solution and the exact solution
120 RCP<const Thyra::VectorBase<double> > x = integrator->getX();
121 RCP<const Thyra::MultiVectorBase<double> > DxDp = integrator->getDxDp();
122 RCP<const Thyra::VectorBase<double> > x_exact =
123 model->getExactSolution(time).get_x();
124 RCP<Thyra::MultiVectorBase<double> > DxDp_exact =
125 Thyra::createMembers(model->get_x_space(), num_param);
126 for (int i = 0; i < num_param; ++i)
127 Thyra::assign(DxDp_exact->col(i).ptr(),
128 *(model->getExactSensSolution(i, time).get_x()));
129
130 // Plot sample solution and exact solution
131 if (comm->getRank() == 0 && n == nTimeStepSizes - 1) {
132 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
133
134 std::ofstream ftmp("Tempus_BDF2_SinCos_Sens.dat");
135 RCP<const SolutionHistory<double> > solutionHistory =
136 integrator->getSolutionHistory();
137 RCP<Thyra::MultiVectorBase<double> > DxDp_exact_plot =
138 Thyra::createMembers(model->get_x_space(), num_param);
139 for (int i = 0; i < solutionHistory->getNumStates(); i++) {
140 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
141 double time_i = solutionState->getTime();
142 RCP<const DMVPV> x_prod_plot =
143 Teuchos::rcp_dynamic_cast<const DMVPV>(solutionState->getX());
144 RCP<const Thyra::VectorBase<double> > x_plot =
145 x_prod_plot->getMultiVector()->col(0);
146 RCP<const Thyra::MultiVectorBase<double> > DxDp_plot =
147 x_prod_plot->getMultiVector()->subView(
148 Teuchos::Range1D(1, num_param));
149 RCP<const Thyra::VectorBase<double> > x_exact_plot =
150 model->getExactSolution(time_i).get_x();
151 for (int j = 0; j < num_param; ++j)
152 Thyra::assign(DxDp_exact_plot->col(j).ptr(),
153 *(model->getExactSensSolution(j, time_i).get_x()));
154 ftmp << std::fixed << std::setprecision(7) << time_i << std::setw(11)
155 << get_ele(*(x_plot), 0) << std::setw(11) << get_ele(*(x_plot), 1);
156 for (int j = 0; j < num_param; ++j)
157 ftmp << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 0)
158 << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 1);
159 ftmp << std::setw(11) << get_ele(*(x_exact_plot), 0) << std::setw(11)
160 << get_ele(*(x_exact_plot), 1);
161 for (int j = 0; j < num_param; ++j)
162 ftmp << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 0)
163 << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 1);
164 ftmp << std::endl;
165 }
166 ftmp.close();
167 }
168
169 // Calculate the error
170 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
171 RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
172 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
173 Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
174 StepSize.push_back(dt);
175 double L2norm = Thyra::norm_2(*xdiff);
176 L2norm *= L2norm;
177 Teuchos::Array<double> L2norm_DxDp(num_param);
178 Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
179 for (int i = 0; i < num_param; ++i)
180 L2norm += L2norm_DxDp[i] * L2norm_DxDp[i];
181 L2norm = std::sqrt(L2norm);
182 ErrorNorm.push_back(L2norm);
183
184 out << " n = " << n << " dt = " << dt << " error = " << L2norm << std::endl;
185 }
186
187 // Check the order and intercept
188 double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
189 out << " Stepper = BDF2" << std::endl;
190 out << " =========================" << std::endl;
191 out << " Expected order: " << order << std::endl;
192 out << " Observed order: " << slope << std::endl;
193 out << " =========================" << std::endl;
194 TEST_FLOATING_EQUALITY(slope, 1.94162, 0.015);
195 TEST_FLOATING_EQUALITY(ErrorNorm[0], 0.0143743, 1.0e-4);
196
197 if (comm->getRank() == 0) {
198 std::ofstream ftmp("Tempus_BDF2_SinCos_Sens-Error.dat");
199 double error0 = 0.8 * ErrorNorm[0];
200 for (int n = 0; n < nTimeStepSizes; n++) {
201 ftmp << StepSize[n] << " " << ErrorNorm[n] << " "
202 << error0 * (StepSize[n] / StepSize[0]) << std::endl;
203 }
204 ftmp.close();
205 }
206}
207
208} // namespace Tempus_Test
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers.
void test_sincos_fsa(const bool use_combined_method, const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out, bool &success)