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_IntegratorForwardSensitivity.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"
33 namespace Tempus_Test {
36 using Teuchos::ParameterList;
37 using Teuchos::sublist;
38 using Teuchos::getParametersFromXmlFile;
47 const bool use_dfdp_as_tangent,
48 Teuchos::FancyOStream &out,
bool &success)
50 std::vector<double> StepSize;
51 std::vector<double> ErrorNorm;
52 const int nTimeStepSizes = 7;
55 Teuchos::RCP<const Teuchos::Comm<int> > comm =
56 Teuchos::DefaultComm<int>::getComm();
57 Teuchos::RCP<Teuchos::FancyOStream> my_out =
58 Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
59 my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
60 my_out->setOutputToRootOnly(0);
61 for (
int n=0; n<nTimeStepSizes; n++) {
64 RCP<ParameterList> pList =
65 getParametersFromXmlFile(
"Tempus_BDF2_SinCos_SA.xml");
68 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
69 scm_pl->set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
70 RCP<SinCosModel<double> > model =
71 Teuchos::rcp(
new SinCosModel<double>(scm_pl));
76 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
77 ParameterList& sens_pl = pl->sublist(
"Sensitivities");
78 if (use_combined_method)
79 sens_pl.set(
"Sensitivity Method",
"Combined");
81 sens_pl.set(
"Sensitivity Method",
"Staggered");
84 sens_pl.set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
85 ParameterList& interp_pl =
86 pl->sublist(
"Default Integrator").sublist(
"Solution History").sublist(
"Interpolator");
87 interp_pl.set(
"Interpolator Type",
"Lagrange");
88 interp_pl.set(
"Order", 1);
91 pl->sublist(
"Default Integrator")
92 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
93 RCP<Tempus::IntegratorForwardSensitivity<double> > integrator =
94 Tempus::integratorForwardSensitivity<double>(pl, model);
95 order = integrator->getStepper()->getOrder();
98 double t0 = pl->sublist(
"Default Integrator")
99 .sublist(
"Time Step Control").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").get<
double>("Final Time");
119 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
122 RCP<const Thyra::VectorBase<
double> > x = integrator->getX();
123 RCP<const Thyra::MultiVectorBase<
double> > DxDp = integrator->getDxDp();
124 RCP<const Thyra::VectorBase<
double> > x_exact =
125 model->getExactSolution(time).get_x();
126 RCP<Thyra::MultiVectorBase<
double> > DxDp_exact =
127 Thyra::createMembers(model->get_x_space(), num_param);
128 for (
int i=0; i<num_param; ++i)
129 Thyra::assign(DxDp_exact->col(i).ptr(),
130 *(model->getExactSensSolution(i, time).get_x()));
133 if (comm->getRank() == 0 && n == nTimeStepSizes-1) {
134 typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
136 std::ofstream ftmp(
"Tempus_BDF2_SinCos_Sens.dat");
138 integrator->getSolutionHistory();
139 RCP< Thyra::MultiVectorBase<double> > DxDp_exact_plot =
140 Thyra::createMembers(model->get_x_space(), num_param);
141 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
142 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
143 double time_i = solutionState->getTime();
144 RCP<const DMVPV> x_prod_plot =
145 Teuchos::rcp_dynamic_cast<
const DMVPV>(solutionState->getX());
146 RCP<const Thyra::VectorBase<double> > x_plot =
147 x_prod_plot->getMultiVector()->col(0);
148 RCP<const Thyra::MultiVectorBase<double> > DxDp_plot =
149 x_prod_plot->getMultiVector()->subView(Teuchos::Range1D(1,num_param));
150 RCP<const Thyra::VectorBase<double> > x_exact_plot =
151 model->getExactSolution(time_i).get_x();
152 for (
int j=0; j<num_param; ++j)
153 Thyra::assign(DxDp_exact_plot->col(j).ptr(),
154 *(model->getExactSensSolution(j, time_i).get_x()));
155 ftmp << std::fixed << std::setprecision(7)
157 << std::setw(11) << get_ele(*(x_plot), 0)
158 << std::setw(11) << get_ele(*(x_plot), 1);
159 for (
int j=0; j<num_param; ++j)
160 ftmp << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 0)
161 << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 1);
162 ftmp << std::setw(11) << get_ele(*(x_exact_plot), 0)
163 << std::setw(11) << get_ele(*(x_exact_plot), 1);
164 for (
int j=0; j<num_param; ++j)
165 ftmp << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 0)
166 << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 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 = BDF2" << 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.0344598, 1.0e-4 );
202 if (comm->getRank() == 0) {
203 std::ofstream ftmp(
"Tempus_BDF2_SinCos_Sens-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;
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...