47 std::vector<double> StepSize;
48 std::vector<double> ErrorNorm;
49 const int nTimeStepSizes = 7;
52 Teuchos::RCP<const Teuchos::Comm<int> > comm =
53 Teuchos::DefaultComm<int>::getComm();
54 for (
int n = 0; n < nTimeStepSizes; n++) {
56 RCP<ParameterList> pList =
57 getParametersFromXmlFile(
"Tempus_BackwardEuler_SinCos_ASA.xml");
61 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
62 RCP<SinCosModel<double> > model =
64 RCP<SinCosModelAdjoint<double> > adjoint_model =
70 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
71 ParameterList& sens_pl = pl->sublist(
"Sensitivities");
72 sens_pl.set(
"Mass Matrix Is Identity",
false);
73 ParameterList& interp_pl = pl->sublist(
"Default Integrator")
74 .sublist(
"Solution History")
75 .sublist(
"Interpolator");
76 interp_pl.set(
"Interpolator Type",
"Lagrange");
77 interp_pl.set(
"Order", 0);
80 pl->sublist(
"Default Stepper").set(
"Use FSAL",
false);
84 pl->sublist(
"Default Stepper")
85 .set(
"Initial Condition Consistency Check",
false);
88 pl->sublist(
"Default Integrator")
89 .sublist(
"Time Step Control")
90 .set(
"Initial Time Step", dt);
91 RCP<Tempus::IntegratorAdjointSensitivity<double> > integrator =
92 Tempus::createIntegratorAdjointSensitivity<double>(pl, model,
94 order = integrator->getStepper()->getOrder();
97 double t0 = pl->sublist(
"Default Integrator")
98 .sublist(
"Time Step Control")
99 .get<
double>(
"Initial Time");
100 RCP<const Thyra::VectorBase<double> > x0 =
101 model->getExactSolution(t0).get_x();
102 const int num_param = model->get_p_space(0)->dim();
103 RCP<Thyra::MultiVectorBase<double> > DxDp0 =
104 Thyra::createMembers(model->get_x_space(), num_param);
105 for (
int i = 0; i < num_param; ++i)
106 Thyra::assign(DxDp0->col(i).ptr(),
107 *(model->getExactSensSolution(i, t0).get_x()));
108 integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
109 DxDp0, Teuchos::null, Teuchos::null);
112 bool integratorStatus = integrator->advanceTime();
113 TEST_ASSERT(integratorStatus)
116 double time = integrator->getTime();
117 double timeFinal = pl->sublist(
"Default Integrator")
118 .sublist(
"Time Step Control")
119 .get<
double>(
"Final Time");
120 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
125 RCP<const Thyra::VectorBase<double> > x = integrator->getX();
126 RCP<const Thyra::MultiVectorBase<double> > DgDp = integrator->getDgDp();
127 RCP<Thyra::MultiVectorBase<double> > DxDp =
128 Thyra::createMembers(model->get_x_space(), num_param);
130 Thyra::ConstDetachedMultiVectorView<double> dgdp_view(*DgDp);
131 Thyra::DetachedMultiVectorView<double> dxdp_view(*DxDp);
132 const int num_g = DgDp->domain()->dim();
133 for (
int i = 0; i < num_g; ++i)
134 for (
int j = 0; j < num_param; ++j) dxdp_view(i, j) = dgdp_view(j, i);
136 RCP<const Thyra::VectorBase<double> > x_exact =
137 model->getExactSolution(time).get_x();
138 RCP<Thyra::MultiVectorBase<double> > DxDp_exact =
139 Thyra::createMembers(model->get_x_space(), num_param);
140 for (
int i = 0; i < num_param; ++i)
141 Thyra::assign(DxDp_exact->col(i).ptr(),
142 *(model->getExactSensSolution(i, time).get_x()));
145 if (comm->getRank() == 0 && n == nTimeStepSizes - 1) {
146 typedef Thyra::DefaultProductVector<double> DPV;
147 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
149 std::ofstream ftmp(
"Tempus_BackwardEuler_SinCos_AdjSens.dat");
150 RCP<const SolutionHistory<double> > solutionHistory =
151 integrator->getSolutionHistory();
152 for (
int i = 0; i < solutionHistory->getNumStates(); i++) {
153 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
154 const double time_i = solutionState->getTime();
155 RCP<const DPV> x_prod_plot =
156 Teuchos::rcp_dynamic_cast<const DPV>(solutionState->getX());
157 RCP<const Thyra::VectorBase<double> > x_plot =
158 x_prod_plot->getVectorBlock(0);
159 RCP<const DMVPV> adjoint_prod_plot =
160 Teuchos::rcp_dynamic_cast<const DMVPV>(
161 x_prod_plot->getVectorBlock(1));
162 RCP<const Thyra::MultiVectorBase<double> > adjoint_plot =
163 adjoint_prod_plot->getMultiVector();
164 RCP<const Thyra::VectorBase<double> > x_exact_plot =
165 model->getExactSolution(time_i).get_x();
166 ftmp << std::fixed << std::setprecision(7) << time_i << std::setw(11)
167 << get_ele(*(x_plot), 0) << std::setw(11) << get_ele(*(x_plot), 1)
168 << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 0)
169 << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 1)
170 << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 0)
171 << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 1)
172 << std::setw(11) << get_ele(*(x_exact_plot), 0) << std::setw(11)
173 << get_ele(*(x_exact_plot), 1) << std::endl;
179 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
180 RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
181 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
182 Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
183 StepSize.push_back(dt);
184 double L2norm = Thyra::norm_2(*xdiff);
186 Teuchos::Array<double> L2norm_DxDp(num_param);
187 Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
188 for (
int i = 0; i < num_param; ++i)
189 L2norm += L2norm_DxDp[i] * L2norm_DxDp[i];
190 L2norm = std::sqrt(L2norm);
191 ErrorNorm.push_back(L2norm);
198 double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
199 out <<
" Stepper = BackwardEuler" << std::endl;
200 out <<
" =========================" << std::endl;
201 out <<
" Expected order: " << order << std::endl;
202 out <<
" Observed order: " << slope << std::endl;
203 out <<
" =========================" << std::endl;
204 TEST_FLOATING_EQUALITY(slope, order, 0.015);
205 TEST_FLOATING_EQUALITY(ErrorNorm[0], 0.142525, 1.0e-4);
207 if (comm->getRank() == 0) {
208 std::ofstream ftmp(
"Tempus_BackwardEuler_SinCos_AdjSens-Error.dat");
209 double error0 = 0.8 * ErrorNorm[0];
210 for (
int n = 0; n < nTimeStepSizes; n++) {
211 ftmp << StepSize[n] <<
" " << ErrorNorm[n] <<
" "
212 << error0 * (StepSize[n] / StepSize[0]) << std::endl;