9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
12 #include "Teuchos_DefaultComm.hpp"
14 #include "Tempus_config.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
16 #include "Tempus_IntegratorAdjointSensitivity.hpp"
18 #include "Thyra_VectorStdOps.hpp"
19 #include "Thyra_MultiVectorStdOps.hpp"
21 #include "../TestModels/SinCosModel.hpp"
22 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
24 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
25 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
26 #include "Thyra_DefaultMultiVectorProductVector.hpp"
27 #include "Thyra_DefaultProductVector.hpp"
37 using Teuchos::ParameterList;
38 using Teuchos::sublist;
39 using Teuchos::getParametersFromXmlFile;
49 std::vector<double> StepSize;
50 std::vector<double> ErrorNorm;
51 const int nTimeStepSizes = 7;
54 Teuchos::RCP<const Teuchos::Comm<int> > comm =
55 Teuchos::DefaultComm<int>::getComm();
56 Teuchos::RCP<Teuchos::FancyOStream> my_out =
57 Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
58 my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
59 my_out->setOutputToRootOnly(0);
60 for (
int n=0; n<nTimeStepSizes; n++) {
63 RCP<ParameterList> pList =
64 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos.xml");
67 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
68 RCP<SinCosModel<double> > model =
74 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
75 ParameterList& sens_pl = pl->sublist(
"Sensitivities");
76 sens_pl.set(
"Mass Matrix Is Identity",
false);
77 ParameterList& interp_pl =
78 pl->sublist(
"Default Integrator").sublist(
"Solution History").sublist(
"Interpolator");
79 interp_pl.set(
"Interpolator Type",
"Lagrange");
80 interp_pl.set(
"Order", 0);
83 pl->sublist(
"Default Integrator")
84 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
85 RCP<Tempus::IntegratorAdjointSensitivity<double> > integrator =
86 Tempus::integratorAdjointSensitivity<double>(pl, model);
87 order = integrator->getStepper()->getOrder();
90 double t0 = pl->sublist(
"Default Integrator")
91 .sublist(
"Time Step Control").get<
double>(
"Initial Time");
92 RCP<const Thyra::VectorBase<double> > x0 =
93 model->getExactSolution(t0).get_x();
94 const int num_param = model->get_p_space(0)->dim();
95 RCP<Thyra::MultiVectorBase<double> > DxDp0 =
96 Thyra::createMembers(model->get_x_space(), num_param);
97 for (
int i=0; i<num_param; ++i)
98 Thyra::assign(DxDp0->col(i).ptr(),
99 *(model->getExactSensSolution(i, t0).get_x()));
100 integrator->setInitialState(t0, x0, Teuchos::null, Teuchos::null,
101 DxDp0, Teuchos::null, Teuchos::null);
104 bool integratorStatus = integrator->advanceTime();
105 TEST_ASSERT(integratorStatus)
108 double time = integrator->getTime();
109 double timeFinal =pl->sublist(
"Default Integrator")
110 .sublist(
"Time Step Control").get<
double>(
"Final Time");
111 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
116 RCP<const Thyra::VectorBase<double> > x = integrator->getX();
117 RCP<const Thyra::MultiVectorBase<double> > DgDp = integrator->getDgDp();
118 RCP<Thyra::MultiVectorBase<double> > DxDp =
119 Thyra::createMembers(model->get_x_space(), num_param);
121 Thyra::ConstDetachedMultiVectorView<double> dgdp_view(*DgDp);
122 Thyra::DetachedMultiVectorView<double> dxdp_view(*DxDp);
123 const int num_g = DgDp->domain()->dim();
124 for (
int i=0; i<num_g; ++i)
125 for (
int j=0; j<num_param; ++j)
126 dxdp_view(i,j) = dgdp_view(j,i);
128 RCP<const Thyra::VectorBase<double> > x_exact =
129 model->getExactSolution(time).get_x();
130 RCP<Thyra::MultiVectorBase<double> > DxDp_exact =
131 Thyra::createMembers(model->get_x_space(), num_param);
132 for (
int i=0; i<num_param; ++i)
133 Thyra::assign(DxDp_exact->col(i).ptr(),
134 *(model->getExactSensSolution(i, time).get_x()));
137 if (comm->getRank() == 0 && n == nTimeStepSizes-1) {
138 typedef Thyra::DefaultProductVector<double> DPV;
139 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
141 std::ofstream ftmp(
"Tempus_BackwardEuler_SinCos_AdjSens.dat");
143 integrator->getSolutionHistory();
145 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
146 const double time = solutionState->getTime();
147 RCP<const DPV> x_prod_plot =
148 Teuchos::rcp_dynamic_cast<const DPV>(solutionState->getX());
149 RCP<const Thyra::VectorBase<double> > x_plot =
150 x_prod_plot->getVectorBlock(0);
151 RCP<const DMVPV > adjoint_prod_plot =
152 Teuchos::rcp_dynamic_cast<const DMVPV>(x_prod_plot->getVectorBlock(1));
153 RCP<const Thyra::MultiVectorBase<double> > adjoint_plot =
154 adjoint_prod_plot->getMultiVector();
155 RCP<const Thyra::VectorBase<double> > x_exact_plot =
156 model->getExactSolution(time).get_x();
157 ftmp << std::fixed << std::setprecision(7)
159 << std::setw(11) << get_ele(*(x_plot), 0)
160 << std::setw(11) << get_ele(*(x_plot), 1)
161 << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 0)
162 << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 1)
163 << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 0)
164 << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 1)
165 << std::setw(11) << get_ele(*(x_exact_plot), 0)
166 << std::setw(11) << get_ele(*(x_exact_plot), 1)
173 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
174 RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
175 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
176 Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
177 StepSize.push_back(dt);
178 double L2norm = Thyra::norm_2(*xdiff);
180 Teuchos::Array<double> L2norm_DxDp(num_param);
181 Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
182 for (
int i=0; i<num_param; ++i)
183 L2norm += L2norm_DxDp[i]*L2norm_DxDp[i];
184 L2norm = std::sqrt(L2norm);
185 ErrorNorm.push_back(L2norm);
187 *my_out <<
" n = " << n <<
" dt = " << dt <<
" error = " << L2norm
193 double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
194 *my_out <<
" Stepper = BackwardEuler" << std::endl;
195 *my_out <<
" =========================" << std::endl;
196 *my_out <<
" Expected order: " << order << std::endl;
197 *my_out <<
" Observed order: " << slope << std::endl;
198 *my_out <<
" =========================" << std::endl;
199 TEST_FLOATING_EQUALITY( slope, order, 0.015 );
200 TEST_FLOATING_EQUALITY( ErrorNorm[0], 0.151746, 1.0e-4 );
202 if (comm->getRank() == 0) {
203 std::ofstream ftmp(
"Tempus_BackwardEuler_SinCos_AdjSens-Error.dat");
204 double error0 = 0.8*ErrorNorm[0];
205 for (
int n=0; n<nTimeStepSizes; n++) {
206 ftmp << StepSize[n] <<
" " << ErrorNorm[n] <<
" "
207 << error0*(StepSize[n]/StepSize[0]) << std::endl;