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 2 Stage 2nd order by Runge");
57 RKMethods.push_back(
"RK Explicit Trapezoidal");
58 RKMethods.push_back(
"General ERK");
61 if (method_name !=
"") {
62 auto it = std::find(RKMethods.begin(), RKMethods.end(),
method_name);
63 TEUCHOS_TEST_FOR_EXCEPTION(it == RKMethods.end(), std::logic_error,
67 std::vector<double> RKMethodErrors;
68 if (use_combined_method) {
69 RKMethodErrors.push_back(0.183799);
70 RKMethodErrors.push_back(6.88637e-06);
71 RKMethodErrors.push_back(6.88637e-06);
72 RKMethodErrors.push_back(0.000264154);
73 RKMethodErrors.push_back(5.22798e-05);
74 RKMethodErrors.push_back(0.000261896);
75 RKMethodErrors.push_back(0.000261896);
76 RKMethodErrors.push_back(0.000261896);
77 RKMethodErrors.push_back(0.00934377);
78 RKMethodErrors.push_back(0.00934377);
79 RKMethodErrors.push_back(6.88637e-06);
82 RKMethodErrors.push_back(0.183799);
83 RKMethodErrors.push_back(2.1915e-05);
84 RKMethodErrors.push_back(2.23367e-05);
85 RKMethodErrors.push_back(0.000205051);
86 RKMethodErrors.push_back(2.85141e-05);
87 RKMethodErrors.push_back(0.000126478);
88 RKMethodErrors.push_back(9.64964e-05);
89 RKMethodErrors.push_back(0.000144616);
90 RKMethodErrors.push_back(0.00826159);
91 RKMethodErrors.push_back(0.00710492);
92 RKMethodErrors.push_back(2.1915e-05);
94 Teuchos::RCP<const Teuchos::Comm<int> > comm =
95 Teuchos::DefaultComm<int>::getComm();
96 Teuchos::RCP<Teuchos::FancyOStream> my_out =
97 Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
98 my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
99 my_out->setOutputToRootOnly(0);
101 for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
104 if (method_name !=
"" && RKMethods[m] != method_name)
107 std::string RKMethod_ = RKMethods[m];
108 std::replace(RKMethod_.begin(), RKMethod_.end(),
' ',
'_');
109 std::replace(RKMethod_.begin(), RKMethod_.end(),
'/',
'.');
110 std::vector<double> StepSize;
111 std::vector<double> ErrorNorm;
112 const int nTimeStepSizes = 7;
115 for (
int n=0; n<nTimeStepSizes; n++) {
118 RCP<ParameterList> pList =
119 getParametersFromXmlFile(
"Tempus_ExplicitRK_SinCos.xml");
122 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
123 scm_pl->set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
124 RCP<SinCosModel<double> > model =
125 Teuchos::rcp(
new SinCosModel<double>(scm_pl));
128 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
129 if (RKMethods[m] ==
"General ERK") {
130 pl->sublist(
"Demo Integrator").set(
"Stepper Name",
"Demo Stepper 2");
132 pl->sublist(
"Demo Stepper").set(
"Stepper Type", RKMethods[m]);
139 ParameterList& sens_pl = pl->sublist(
"Sensitivities");
140 if (use_combined_method)
141 sens_pl.set(
"Sensitivity Method",
"Combined");
143 sens_pl.set(
"Sensitivity Method",
"Staggered");
144 sens_pl.set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
145 ParameterList& interp_pl =
146 pl->sublist(
"Demo Integrator").sublist(
"Solution History").sublist(
"Interpolator");
147 interp_pl.set(
"Interpolator Type",
"Lagrange");
148 interp_pl.set(
"Order", 3);
151 pl->sublist(
"Demo Integrator")
152 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
153 RCP<Tempus::IntegratorForwardSensitivity<double> > integrator =
154 Tempus::integratorForwardSensitivity<double>(pl, model);
155 order = integrator->getStepper()->getOrder();
158 double t0 = pl->sublist(
"Demo Integrator")
159 .sublist(
"Time Step Control").get<
double>(
"Initial Time");
162 RCP<Thyra::VectorBase<double> > x0 =
163 model->getNominalValues().get_x()->clone_v();
164 const int num_param = model->get_p_space(0)->dim();
165 RCP<Thyra::MultiVectorBase<double> > DxDp0 =
166 Thyra::createMembers(model->get_x_space(), num_param);
167 for (
int i=0; i<num_param; ++i)
168 Thyra::assign(DxDp0->col(i).ptr(),
169 *(model->getExactSensSolution(i, t0).get_x()));
170 integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
171 DxDp0, Teuchos::null, Teuchos::null);
174 bool integratorStatus = integrator->advanceTime();
175 TEST_ASSERT(integratorStatus)
178 double time = integrator->getTime();
179 double timeFinal = pl->sublist("Demo Integrator")
180 .sublist("Time Step Control").get<
double>("Final Time");
181 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
184 RCP<const Thyra::VectorBase<
double> > x = integrator->getX();
185 RCP<const Thyra::MultiVectorBase<
double> > DxDp = integrator->getDxDp();
186 RCP<const Thyra::VectorBase<
double> > x_exact =
187 model->getExactSolution(time).get_x();
188 RCP<Thyra::MultiVectorBase<
double> > DxDp_exact =
189 Thyra::createMembers(model->get_x_space(), num_param);
190 for (
int i=0; i<num_param; ++i)
191 Thyra::assign(DxDp_exact->col(i).ptr(),
192 *(model->getExactSensSolution(i, time).get_x()));
195 if (comm->getRank() == 0 && n == nTimeStepSizes-1) {
196 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
198 std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_Sens.dat");
200 integrator->getSolutionHistory();
201 RCP< Thyra::MultiVectorBase<double> > DxDp_exact_plot =
202 Thyra::createMembers(model->get_x_space(), num_param);
203 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
204 RCP<const SolutionState<double> > solutionState =
205 (*solutionHistory)[i];
206 double time_i = solutionState->getTime();
207 RCP<const DMVPV> x_prod_plot =
208 Teuchos::rcp_dynamic_cast<
const DMVPV>(solutionState->getX());
209 RCP<const Thyra::VectorBase<double> > x_plot =
210 x_prod_plot->getMultiVector()->col(0);
211 RCP<const Thyra::MultiVectorBase<double> > DxDp_plot =
212 x_prod_plot->getMultiVector()->subView(Teuchos::Range1D(1,num_param));
213 RCP<const Thyra::VectorBase<double> > x_exact_plot =
214 model->getExactSolution(time_i).get_x();
215 for (
int j=0; j<num_param; ++j)
216 Thyra::assign(DxDp_exact_plot->col(j).ptr(),
217 *(model->getExactSensSolution(j, time_i).get_x()));
218 ftmp << std::fixed << std::setprecision(7)
220 << std::setw(11) << get_ele(*(x_plot), 0)
221 << std::setw(11) << get_ele(*(x_plot), 1);
222 for (
int j=0; j<num_param; ++j)
223 ftmp << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 0)
224 << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 1);
225 ftmp << std::setw(11) << get_ele(*(x_exact_plot), 0)
226 << std::setw(11) << get_ele(*(x_exact_plot), 1);
227 for (
int j=0; j<num_param; ++j)
228 ftmp << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 0)
229 << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 1);
236 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
237 RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
238 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
239 Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
240 StepSize.push_back(dt);
241 double L2norm = Thyra::norm_2(*xdiff);
243 Teuchos::Array<double> L2norm_DxDp(num_param);
244 Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
245 for (
int i=0; i<num_param; ++i)
246 L2norm += L2norm_DxDp[i]*L2norm_DxDp[i];
247 L2norm = std::sqrt(L2norm);
248 ErrorNorm.push_back(L2norm);
250 *my_out <<
" n = " << n <<
" dt = " << dt <<
" error = " << L2norm
255 double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
256 *my_out <<
" Stepper = " << RKMethods[m] << std::endl;
257 *my_out <<
" =========================" << std::endl;
258 *my_out <<
" Expected order: " << order << std::endl;
259 *my_out <<
" Observed order: " << slope << std::endl;
260 *my_out <<
" =========================" << std::endl;
261 TEST_FLOATING_EQUALITY( slope, order, 0.04 );
262 TEST_FLOATING_EQUALITY( ErrorNorm[0], RKMethodErrors[m], 1.0e-4 );
264 if (comm->getRank() == 0) {
265 std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_Sens-Error.dat");
266 double error0 = 0.8*ErrorNorm[0];
267 for (
int n=0; n<nTimeStepSizes; n++) {
268 ftmp << StepSize[n] <<
" " << ErrorNorm[n] <<
" "
269 << error0*(pow(StepSize[n]/StepSize[0],order)) << std::endl;
275 Teuchos::TimeMonitor::summarize();
void test_sincos_fsa(const bool use_combined_method, const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out, bool &success)
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
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...