45 const bool use_combined_method,
46 const bool use_dfdp_as_tangent, Teuchos::FancyOStream& out,
49 std::vector<std::string> stepperTypes;
50 stepperTypes.push_back(
"Partitioned IMEX RK 1st order");
51 stepperTypes.push_back(
"Partitioned IMEX RK SSP2");
52 stepperTypes.push_back(
"Partitioned IMEX RK ARS 233");
53 stepperTypes.push_back(
"General Partitioned IMEX RK");
57 auto it = std::find(stepperTypes.begin(), stepperTypes.end(),
method_name);
58 TEUCHOS_TEST_FOR_EXCEPTION(it == stepperTypes.end(), std::logic_error,
62 std::vector<double> stepperOrders;
63 std::vector<double> stepperErrors;
64 if (use_dfdp_as_tangent) {
65 if (use_combined_method) {
66 stepperOrders.push_back(1.16082);
67 stepperOrders.push_back(1.97231);
68 stepperOrders.push_back(2.5914);
69 stepperOrders.push_back(1.99148);
71 stepperErrors.push_back(0.00820931);
72 stepperErrors.push_back(0.287112);
73 stepperErrors.push_back(0.00646096);
74 stepperErrors.push_back(0.148848);
77 stepperOrders.push_back(1.07932);
78 stepperOrders.push_back(1.97396);
79 stepperOrders.push_back(2.63724);
80 stepperOrders.push_back(1.99133);
82 stepperErrors.push_back(0.055626);
83 stepperErrors.push_back(0.198898);
84 stepperErrors.push_back(0.00614135);
85 stepperErrors.push_back(0.0999881);
89 if (use_combined_method) {
90 stepperOrders.push_back(1.1198);
91 stepperOrders.push_back(1.98931);
92 stepperOrders.push_back(2.60509);
93 stepperOrders.push_back(1.992);
95 stepperErrors.push_back(0.00619674);
96 stepperErrors.push_back(0.294989);
97 stepperErrors.push_back(0.0062125);
98 stepperErrors.push_back(0.142489);
101 stepperOrders.push_back(1.07932);
102 stepperOrders.push_back(1.97396);
103 stepperOrders.push_back(2.63724);
104 stepperOrders.push_back(1.99133);
106 stepperErrors.push_back(0.055626);
107 stepperErrors.push_back(0.198898);
108 stepperErrors.push_back(0.00614135);
109 stepperErrors.push_back(0.0999881);
113 std::vector<double> stepperInitDt;
114 stepperInitDt.push_back(0.0125);
115 stepperInitDt.push_back(0.05);
116 stepperInitDt.push_back(0.05);
117 stepperInitDt.push_back(0.05);
119 Teuchos::RCP<const Teuchos::Comm<int>> comm =
120 Teuchos::DefaultComm<int>::getComm();
122 std::vector<std::string>::size_type m;
123 for (m = 0; m != stepperTypes.size(); m++) {
127 std::string stepperType = stepperTypes[m];
128 std::string stepperName = stepperTypes[m];
129 std::replace(stepperName.begin(), stepperName.end(),
' ',
'_');
130 std::replace(stepperName.begin(), stepperName.end(),
'/',
'.');
132 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
133 std::vector<RCP<Thyra::VectorBase<double>>> sensitivities;
134 std::vector<double> StepSize;
135 std::vector<double> ErrorNorm;
136 const int nTimeStepSizes = 3;
137 double dt = stepperInitDt[m];
139 for (
int n = 0; n < nTimeStepSizes; n++) {
141 RCP<ParameterList> pList =
142 getParametersFromXmlFile(
"Tempus_IMEX_RK_VanDerPol.xml");
145 RCP<ParameterList> vdpmPL = sublist(pList,
"VanDerPolModel",
true);
146 vdpmPL->set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
147 const bool useProductVector =
true;
148 RCP<VanDerPol_IMEX_ExplicitModel<double>> explicitModel = Teuchos::rcp(
152 RCP<VanDerPol_IMEXPart_ImplicitModel<double>> implicitModel =
156 const int numExplicitBlocks = 1;
157 const int parameterIndex = 4;
158 RCP<Tempus::WrapperModelEvaluatorPairPartIMEX_Basic<double>> model =
161 explicitModel, implicitModel, numExplicitBlocks,
165 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
166 ParameterList& sens_pl = pl->sublist(
"Sensitivities");
167 if (use_combined_method)
168 sens_pl.set(
"Sensitivity Method",
"Combined");
170 sens_pl.set(
"Sensitivity Method",
"Staggered");
171 sens_pl.set(
"Reuse State Linear Solver",
true);
173 sens_pl.set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
174 ParameterList& interp_pl = pl->sublist(
"Default Integrator")
175 .sublist(
"Solution History")
176 .sublist(
"Interpolator");
177 interp_pl.set(
"Interpolator Type",
"Lagrange");
178 interp_pl.set(
"Order", 2);
181 if (stepperType ==
"General Partitioned IMEX RK") {
183 pl->sublist(
"Default Integrator")
184 .set(
"Stepper Name",
"General IMEX RK");
187 pl->sublist(
"Default Stepper").set(
"Stepper Type", stepperType);
191 if (n == nTimeStepSizes - 1)
197 pl->sublist(
"Default Integrator")
198 .sublist(
"Time Step Control")
199 .set(
"Initial Time Step", dt);
200 pl->sublist(
"Default Integrator")
201 .sublist(
"Time Step Control")
202 .remove(
"Time Step Control Strategy");
203 RCP<Tempus::IntegratorForwardSensitivity<double>> integrator =
204 Tempus::createIntegratorForwardSensitivity<double>(pl, model);
205 order = integrator->getStepper()->getOrder();
208 bool integratorStatus = integrator->advanceTime();
209 TEST_ASSERT(integratorStatus)
212 double time = integrator->getTime();
213 double timeFinal = pl->sublist(
"Default Integrator")
214 .sublist(
"Time Step Control")
215 .get<
double>(
"Final Time");
216 double tol = 100.0 * std::numeric_limits<double>::epsilon();
217 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
220 auto solution = Thyra::createMember(model->get_x_space());
221 auto sensitivity = Thyra::createMember(model->get_x_space());
222 Thyra::copy(*(integrator->getX()), solution.ptr());
223 Thyra::copy(*(integrator->getDxDp()->col(0)), sensitivity.ptr());
224 solutions.push_back(solution);
225 sensitivities.push_back(sensitivity);
226 StepSize.push_back(dt);
229 if ((n == 0) || (n == nTimeStepSizes - 1)) {
230 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
232 std::string fname =
"Tempus_" + stepperName +
"_VanDerPol_Sens-Ref.dat";
233 if (n == 0) fname =
"Tempus_" + stepperName +
"_VanDerPol_Sens.dat";
234 std::ofstream ftmp(fname);
235 RCP<const SolutionHistory<double>> solutionHistory =
236 integrator->getSolutionHistory();
237 int nStates = solutionHistory->getNumStates();
238 for (
int i = 0; i < nStates; i++) {
239 RCP<const SolutionState<double>> solutionState =
240 (*solutionHistory)[i];
241 RCP<const DMVPV> x_prod =
242 Teuchos::rcp_dynamic_cast<const DMVPV>(solutionState->getX());
243 RCP<const Thyra::VectorBase<double>> x =
244 x_prod->getMultiVector()->col(0);
245 RCP<const Thyra::VectorBase<double>> dxdp =
246 x_prod->getMultiVector()->col(1);
247 double ttime = solutionState->getTime();
248 ftmp << std::fixed << std::setprecision(7) << ttime <<
" "
249 << std::setw(11) << get_ele(*x, 0) <<
" " << std::setw(11)
250 << get_ele(*x, 1) <<
" " << std::setw(11) << get_ele(*dxdp, 0)
251 <<
" " << std::setw(11) << get_ele(*dxdp, 1) << std::endl;
259 auto ref_solution = solutions[solutions.size() - 1];
260 auto ref_sensitivity = sensitivities[solutions.size() - 1];
261 std::vector<double> StepSizeCheck;
262 for (std::size_t i = 0; i < (solutions.size() - 1); ++i) {
263 auto sol = solutions[i];
264 auto sen = sensitivities[i];
265 Thyra::Vp_StV(sol.ptr(), -1.0, *ref_solution);
266 Thyra::Vp_StV(sen.ptr(), -1.0, *ref_sensitivity);
267 const double L2norm_sol = Thyra::norm_2(*sol);
268 const double L2norm_sen = Thyra::norm_2(*sen);
269 const double L2norm =
270 std::sqrt(L2norm_sol * L2norm_sol + L2norm_sen * L2norm_sen);
271 StepSizeCheck.push_back(StepSize[i]);
272 ErrorNorm.push_back(L2norm);
280 computeLinearRegressionLogLog<double>(StepSizeCheck, ErrorNorm);
281 out <<
" Stepper = " << stepperType << std::endl;
282 out <<
" =========================" << std::endl;
283 out <<
" Expected order: " << order << std::endl;
284 out <<
" Observed order: " << slope << std::endl;
285 out <<
" =========================" << std::endl;
286 TEST_FLOATING_EQUALITY(slope, stepperOrders[m], 0.02);
287 TEST_FLOATING_EQUALITY(ErrorNorm[0], stepperErrors[m], 1.0e-4);
291 std::ofstream ftmp(
"Tempus_" + stepperName +
"_VanDerPol_Sens-Error.dat");
292 double error0 = 0.8 * ErrorNorm[0];
293 for (std::size_t n = 0; n < StepSizeCheck.size(); n++) {
294 ftmp << StepSizeCheck[n] <<
" " << ErrorNorm[n] <<
" "
295 << error0 * (pow(StepSize[n] / StepSize[0], order)) << std::endl;
300 Teuchos::TimeMonitor::summarize();