9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
12 #include "Teuchos_DefaultComm.hpp"
14 #include "Thyra_VectorStdOps.hpp"
15 #include "Thyra_MultiVectorStdOps.hpp"
17 #include "Tempus_IntegratorBasic.hpp"
18 #include "Tempus_IntegratorForwardSensitivity.hpp"
20 #include "Thyra_DefaultMultiVectorProductVector.hpp"
21 #include "Thyra_DefaultProductVector.hpp"
23 #include "../TestModels/SinCosModel.hpp"
24 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
29 namespace Tempus_Test {
32 using Teuchos::ParameterList;
33 using Teuchos::sublist;
34 using Teuchos::getParametersFromXmlFile;
43 const bool use_combined_method,
44 const bool use_dfdp_as_tangent,
45 Teuchos::FancyOStream &out,
bool &success)
47 std::vector<std::string> RKMethods;
48 RKMethods.push_back(
"RK Forward Euler");
49 RKMethods.push_back(
"RK Explicit 4 Stage");
50 RKMethods.push_back(
"RK Explicit 3/8 Rule");
51 RKMethods.push_back(
"RK Explicit 4 Stage 3rd order by Runge");
52 RKMethods.push_back(
"RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
53 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order");
54 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order TVD");
55 RKMethods.push_back(
"RK Explicit 3 Stage 3rd order by Heun");
56 RKMethods.push_back(
"RK Explicit Midpoint");
57 RKMethods.push_back(
"RK Explicit Trapezoidal");
58 RKMethods.push_back(
"Heuns Method");
59 RKMethods.push_back(
"General ERK");
62 if (method_name !=
"") {
63 auto it = std::find(RKMethods.begin(), RKMethods.end(),
method_name);
64 TEUCHOS_TEST_FOR_EXCEPTION(it == RKMethods.end(), std::logic_error,
65 "Invalid RK method name '" << method_name <<
"'");
68 std::vector<double> RKMethodErrors;
69 if (use_combined_method) {
70 RKMethodErrors.push_back(0.183799);
71 RKMethodErrors.push_back(6.88637e-06);
72 RKMethodErrors.push_back(6.88637e-06);
73 RKMethodErrors.push_back(0.000264154);
74 RKMethodErrors.push_back(5.22798e-05);
75 RKMethodErrors.push_back(0.000261896);
76 RKMethodErrors.push_back(0.000261896);
77 RKMethodErrors.push_back(0.000261896);
78 RKMethodErrors.push_back(0.00934377);
79 RKMethodErrors.push_back(0.00934377);
80 RKMethodErrors.push_back(0.00934377);
81 RKMethodErrors.push_back(6.88637e-06);
84 RKMethodErrors.push_back(0.183799);
85 RKMethodErrors.push_back(2.1915e-05);
86 RKMethodErrors.push_back(2.23367e-05);
87 RKMethodErrors.push_back(0.000205051);
88 RKMethodErrors.push_back(2.85141e-05);
89 RKMethodErrors.push_back(0.000126478);
90 RKMethodErrors.push_back(9.64964e-05);
91 RKMethodErrors.push_back(0.000144616);
92 RKMethodErrors.push_back(0.00826159);
93 RKMethodErrors.push_back(0.00710492);
94 RKMethodErrors.push_back(0.00710492);
95 RKMethodErrors.push_back(2.1915e-05);
97 Teuchos::RCP<const Teuchos::Comm<int> > comm =
98 Teuchos::DefaultComm<int>::getComm();
99 Teuchos::RCP<Teuchos::FancyOStream> my_out =
100 Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
101 my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
102 my_out->setOutputToRootOnly(0);
104 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
107 if (method_name !=
"" && RKMethods[m] != method_name)
110 std::string RKMethod_ = RKMethods[m];
111 std::replace(RKMethod_.begin(), RKMethod_.end(),
' ',
'_');
112 std::replace(RKMethod_.begin(), RKMethod_.end(),
'/',
'.');
113 std::vector<double> StepSize;
114 std::vector<double> ErrorNorm;
115 const int nTimeStepSizes = 7;
118 for (
int n=0; n<nTimeStepSizes; n++) {
121 RCP<ParameterList> pList =
122 getParametersFromXmlFile(
"Tempus_ExplicitRK_SinCos.xml");
125 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
126 scm_pl->set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
127 RCP<SinCosModel<double> > model =
128 Teuchos::rcp(
new SinCosModel<double>(scm_pl));
131 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
132 if (RKMethods[m] ==
"General ERK") {
133 pl->sublist(
"Demo Integrator").set(
"Stepper Name",
"Demo Stepper 2");
135 pl->sublist(
"Demo Stepper").set(
"Stepper Type", RKMethods[m]);
142 ParameterList& sens_pl = pl->sublist(
"Sensitivities");
143 if (use_combined_method)
144 sens_pl.set(
"Sensitivity Method",
"Combined");
146 sens_pl.set(
"Sensitivity Method",
"Staggered");
147 sens_pl.set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
148 ParameterList& interp_pl =
149 pl->sublist(
"Demo Integrator").sublist(
"Solution History").sublist(
"Interpolator");
150 interp_pl.set(
"Interpolator Type",
"Lagrange");
151 interp_pl.set(
"Order", 3);
154 pl->sublist(
"Demo Integrator")
155 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
156 RCP<Tempus::IntegratorForwardSensitivity<double> > integrator =
157 Tempus::integratorForwardSensitivity<double>(pl, model);
158 order = integrator->getStepper()->getOrder();
161 double t0 = pl->sublist(
"Demo Integrator")
162 .sublist(
"Time Step Control").get<
double>(
"Initial Time");
165 RCP<Thyra::VectorBase<double> > x0 =
166 model->getNominalValues().get_x()->clone_v();
167 const int num_param = model->get_p_space(0)->dim();
168 RCP<Thyra::MultiVectorBase<double> > DxDp0 =
169 Thyra::createMembers(model->get_x_space(), num_param);
170 for (
int i=0; i<num_param; ++i)
171 Thyra::assign(DxDp0->col(i).ptr(),
172 *(model->getExactSensSolution(i, t0).get_x()));
173 integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
174 DxDp0, Teuchos::null, Teuchos::null);
177 bool integratorStatus = integrator->advanceTime();
178 TEST_ASSERT(integratorStatus)
181 double time = integrator->getTime();
182 double timeFinal = pl->sublist("Demo Integrator")
183 .sublist("Time Step Control").get<
double>("Final Time");
184 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
187 RCP<const Thyra::VectorBase<
double> > x = integrator->getX();
188 RCP<const Thyra::MultiVectorBase<
double> > DxDp = integrator->getDxDp();
189 RCP<const Thyra::VectorBase<
double> > x_exact =
190 model->getExactSolution(time).get_x();
191 RCP<Thyra::MultiVectorBase<
double> > DxDp_exact =
192 Thyra::createMembers(model->get_x_space(), num_param);
193 for (
int i=0; i<num_param; ++i)
194 Thyra::assign(DxDp_exact->col(i).ptr(),
195 *(model->getExactSensSolution(i, time).get_x()));
198 if (comm->getRank() == 0 && n == nTimeStepSizes-1) {
199 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
201 std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_Sens.dat");
202 RCP<const SolutionHistory<double> > solutionHistory =
203 integrator->getSolutionHistory();
204 RCP< Thyra::MultiVectorBase<double> > DxDp_exact_plot =
205 Thyra::createMembers(model->get_x_space(), num_param);
206 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
207 RCP<const SolutionState<double> > solutionState =
208 (*solutionHistory)[i];
209 double time_i = solutionState->getTime();
210 RCP<const DMVPV> x_prod_plot =
211 Teuchos::rcp_dynamic_cast<
const DMVPV>(solutionState->getX());
212 RCP<const Thyra::VectorBase<double> > x_plot =
213 x_prod_plot->getMultiVector()->col(0);
214 RCP<const Thyra::MultiVectorBase<double> > DxDp_plot =
215 x_prod_plot->getMultiVector()->subView(Teuchos::Range1D(1,num_param));
216 RCP<const Thyra::VectorBase<double> > x_exact_plot =
217 model->getExactSolution(time_i).get_x();
218 for (
int j=0; j<num_param; ++j)
219 Thyra::assign(DxDp_exact_plot->col(j).ptr(),
220 *(model->getExactSensSolution(j, time_i).get_x()));
221 ftmp << std::fixed << std::setprecision(7)
223 << std::setw(11) << get_ele(*(x_plot), 0)
224 << std::setw(11) << get_ele(*(x_plot), 1);
225 for (
int j=0; j<num_param; ++j)
226 ftmp << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 0)
227 << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 1);
228 ftmp << std::setw(11) << get_ele(*(x_exact_plot), 0)
229 << std::setw(11) << get_ele(*(x_exact_plot), 1);
230 for (
int j=0; j<num_param; ++j)
231 ftmp << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 0)
232 << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 1);
239 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
240 RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
241 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
242 Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
243 StepSize.push_back(dt);
244 double L2norm = Thyra::norm_2(*xdiff);
246 Teuchos::Array<double> L2norm_DxDp(num_param);
247 Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
248 for (
int i=0; i<num_param; ++i)
249 L2norm += L2norm_DxDp[i]*L2norm_DxDp[i];
250 L2norm = std::sqrt(L2norm);
251 ErrorNorm.push_back(L2norm);
253 *my_out <<
" n = " << n <<
" dt = " << dt <<
" error = " << L2norm
258 double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
259 *my_out <<
" Stepper = " << RKMethods[m] << std::endl;
260 *my_out <<
" =========================" << std::endl;
261 *my_out <<
" Expected order: " << order << std::endl;
262 *my_out <<
" Observed order: " << slope << std::endl;
263 *my_out <<
" =========================" << std::endl;
264 TEST_FLOATING_EQUALITY( slope, order, 0.04 );
265 TEST_FLOATING_EQUALITY( ErrorNorm[0], RKMethodErrors[m], 1.0e-4 );
267 if (comm->getRank() == 0) {
268 std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_Sens-Error.dat");
269 double error0 = 0.8*ErrorNorm[0];
270 for (
int n=0; n<nTimeStepSizes; n++) {
271 ftmp << StepSize[n] <<
" " << ErrorNorm[n] <<
" "
272 << error0*(pow(StepSize[n]/StepSize[0],order)) << std::endl;
278 Teuchos::TimeMonitor::summarize();
void test_sincos_fsa(const bool use_combined_method, const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out, bool &success)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...