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_IntegratorAdjointSensitivity.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;
 
   44   std::vector<std::string> RKMethods;
 
   45   RKMethods.push_back(
"RK Forward Euler");
 
   46   RKMethods.push_back(
"RK Explicit 4 Stage");
 
   47   RKMethods.push_back(
"RK Explicit 3/8 Rule");
 
   48   RKMethods.push_back(
"RK Explicit 4 Stage 3rd order by Runge");
 
   49   RKMethods.push_back(
"RK Explicit 5 Stage 3rd order by Kinnmark and Gray");
 
   50   RKMethods.push_back(
"RK Explicit 3 Stage 3rd order");
 
   51   RKMethods.push_back(
"RK Explicit 3 Stage 3rd order TVD");
 
   52   RKMethods.push_back(
"RK Explicit 3 Stage 3rd order by Heun");
 
   53   RKMethods.push_back(
"RK Explicit Midpoint");
 
   54   RKMethods.push_back(
"RK Explicit Trapezoidal");
 
   55   RKMethods.push_back(
"Heuns Method");
 
   56   RKMethods.push_back(
"General ERK");
 
   58   std::vector<double> RKMethodErrors;
 
   59   RKMethodErrors.push_back(0.154904);
 
   60   RKMethodErrors.push_back(4.55982e-06);
 
   61   RKMethodErrors.push_back(4.79132e-06);
 
   62   RKMethodErrors.push_back(0.000113603);
 
   63   RKMethodErrors.push_back(4.98796e-05);
 
   64   RKMethodErrors.push_back(0.00014564);
 
   65   RKMethodErrors.push_back(0.000121968);
 
   66   RKMethodErrors.push_back(0.000109495);
 
   67   RKMethodErrors.push_back(0.00559871);
 
   68   RKMethodErrors.push_back(0.00710492);
 
   69   RKMethodErrors.push_back(0.00710492);
 
   70   RKMethodErrors.push_back(4.55982e-06);
 
   72   Teuchos::RCP<const Teuchos::Comm<int> > comm =
 
   73     Teuchos::DefaultComm<int>::getComm();
 
   74   Teuchos::RCP<Teuchos::FancyOStream> my_out =
 
   75     Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
 
   76   my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
 
   77   my_out->setOutputToRootOnly(0);
 
   79   for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
 
   81     std::string RKMethod_ = RKMethods[m];
 
   82     std::replace(RKMethod_.begin(), RKMethod_.end(), 
' ', 
'_');
 
   83     std::replace(RKMethod_.begin(), RKMethod_.end(), 
'/', 
'.');
 
   84     std::vector<double> StepSize;
 
   85     std::vector<double> ErrorNorm;
 
   86     const int nTimeStepSizes = 7;
 
   89     for (
int n=0; n<nTimeStepSizes; n++) {
 
   92       RCP<ParameterList> pList =
 
   93         getParametersFromXmlFile(
"Tempus_ExplicitRK_SinCos.xml");
 
   96       RCP<ParameterList> scm_pl = sublist(pList, 
"SinCosModel", 
true);
 
   97       RCP<SinCosModel<double> > model =
 
  101       RCP<ParameterList> pl = sublist(pList, 
"Tempus", 
true);
 
  102       if (RKMethods[m] == 
"General ERK") {
 
  103         pl->sublist(
"Demo Integrator").set(
"Stepper Name", 
"Demo Stepper 2");
 
  104         pl->sublist(
"Demo Stepper 2").set(
"Initial Condition Consistency",
 
  106         pl->sublist(
"Demo Stepper 2").set(
"Initial Condition Consistency Check",
 
  109         pl->sublist(
"Demo Stepper").set(
"Stepper Type", RKMethods[m]);
 
  110         pl->sublist(
"Demo Stepper").set(
"Initial Condition Consistency",
 
  112         pl->sublist(
"Demo Stepper").set(
"Initial Condition Consistency Check",
 
  120       ParameterList& sens_pl = pl->sublist(
"Sensitivities");
 
  121       sens_pl.set(
"Mass Matrix Is Identity", 
true); 
 
  122       ParameterList& interp_pl =
 
  123         pl->sublist(
"Demo Integrator").sublist(
"Solution History").sublist(
"Interpolator");
 
  124       interp_pl.set(
"Interpolator Type", 
"Lagrange");
 
  125       interp_pl.set(
"Order", 3); 
 
  128       pl->sublist(
"Demo Integrator")
 
  129         .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
 
  130       RCP<Tempus::IntegratorAdjointSensitivity<double> > integrator =
 
  131         Tempus::integratorAdjointSensitivity<double>(pl, model);
 
  132       order = integrator->getStepper()->getOrder();
 
  135       double t0 = pl->sublist(
"Demo Integrator")
 
  136         .sublist(
"Time Step Control").get<
double>(
"Initial Time");
 
  139       RCP<Thyra::VectorBase<double> > x0 =
 
  140         model->getNominalValues().get_x()->clone_v();
 
  141       const int num_param = model->get_p_space(0)->dim();
 
  142       RCP<Thyra::MultiVectorBase<double> > DxDp0 =
 
  143         Thyra::createMembers(model->get_x_space(), num_param);
 
  144       for (
int i=0; i<num_param; ++i)
 
  145         Thyra::assign(DxDp0->col(i).ptr(),
 
  146                       *(model->getExactSensSolution(i, t0).get_x()));
 
  147       integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
 
  148                                   DxDp0, Teuchos::null, Teuchos::null);
 
  151       bool integratorStatus = integrator->advanceTime();
 
  152       TEST_ASSERT(integratorStatus)
 
  155       double time = integrator->getTime();
 
  156       double timeFinal = pl->sublist(
"Demo Integrator")
 
  157         .sublist(
"Time Step Control").get<
double>(
"Final Time");
 
  158       TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
 
  163       RCP<const Thyra::VectorBase<double> > x = integrator->getX();
 
  164       RCP<const Thyra::MultiVectorBase<double> > DgDp = integrator->getDgDp();
 
  165       RCP<Thyra::MultiVectorBase<double> > DxDp =
 
  166         Thyra::createMembers(model->get_x_space(), num_param);
 
  168         Thyra::ConstDetachedMultiVectorView<double> dgdp_view(*DgDp);
 
  169         Thyra::DetachedMultiVectorView<double> dxdp_view(*DxDp);
 
  170         const int num_g = DgDp->domain()->dim();
 
  171         for (
int i=0; i<num_g; ++i)
 
  172           for (
int j=0; j<num_param; ++j)
 
  173             dxdp_view(i,j) = dgdp_view(j,i);
 
  175       RCP<const Thyra::VectorBase<double> > x_exact =
 
  176         model->getExactSolution(time).get_x();
 
  177       RCP<Thyra::MultiVectorBase<double> > DxDp_exact =
 
  178         Thyra::createMembers(model->get_x_space(), num_param);
 
  179       for (
int i=0; i<num_param; ++i)
 
  180         Thyra::assign(DxDp_exact->col(i).ptr(),
 
  181                       *(model->getExactSensSolution(i, time).get_x()));
 
  184       if (comm->getRank() == 0 && n == nTimeStepSizes-1) {
 
  185         typedef Thyra::DefaultProductVector<double> DPV;
 
  186         typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
 
  188         std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_AdjSens.dat");
 
  189         RCP<const SolutionHistory<double> > solutionHistory =
 
  190           integrator->getSolutionHistory();
 
  191         for (
int i=0; i<solutionHistory->getNumStates(); i++) {
 
  192           RCP<const SolutionState<double> > solutionState =
 
  193             (*solutionHistory)[i];
 
  194           const double time_i = solutionState->getTime();
 
  195           RCP<const DPV> x_prod_plot =
 
  196             Teuchos::rcp_dynamic_cast<
const DPV>(solutionState->getX());
 
  197           RCP<const Thyra::VectorBase<double> > x_plot =
 
  198             x_prod_plot->getVectorBlock(0);
 
  199           RCP<const DMVPV > adjoint_prod_plot =
 
  200             Teuchos::rcp_dynamic_cast<
const DMVPV>(x_prod_plot->getVectorBlock(1));
 
  201           RCP<const Thyra::MultiVectorBase<double> > adjoint_plot =
 
  202             adjoint_prod_plot->getMultiVector();
 
  203           RCP<const Thyra::VectorBase<double> > x_exact_plot =
 
  204             model->getExactSolution(time_i).get_x();
 
  205           ftmp << std::fixed << std::setprecision(7)
 
  207                << std::setw(11) << get_ele(*(x_plot), 0)
 
  208                << std::setw(11) << get_ele(*(x_plot), 1)
 
  209                << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 0)
 
  210                << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 1)
 
  211                << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 0)
 
  212                << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 1)
 
  213                << std::setw(11) << get_ele(*(x_exact_plot), 0)
 
  214                << std::setw(11) << get_ele(*(x_exact_plot), 1)
 
  221       RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
 
  222       RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
 
  223       Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
 
  224       Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
 
  225       StepSize.push_back(dt);
 
  226       double L2norm = Thyra::norm_2(*xdiff);
 
  228       Teuchos::Array<double> L2norm_DxDp(num_param);
 
  229       Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
 
  230       for (
int i=0; i<num_param; ++i)
 
  231         L2norm += L2norm_DxDp[i]*L2norm_DxDp[i];
 
  232       L2norm = std::sqrt(L2norm);
 
  233       ErrorNorm.push_back(L2norm);
 
  240     double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
 
  241     *my_out << 
"  Stepper = " << RKMethods[m] << std::endl;
 
  242     *my_out << 
"  =========================" << std::endl;
 
  243     *my_out << 
"  Expected order: " << order << std::endl;
 
  244     *my_out << 
"  Observed order: " << slope << std::endl;
 
  245     *my_out << 
"  =========================" << std::endl;
 
  246     TEST_FLOATING_EQUALITY( slope, order, 0.015 );
 
  247     TEST_FLOATING_EQUALITY( ErrorNorm[0], RKMethodErrors[m], 1.0e-4 );
 
  249     if (comm->getRank() == 0) {
 
  250       std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_AdjSens-Error.dat");
 
  251       double error0 = 0.8*ErrorNorm[0];
 
  252       for (
int n=0; n<nTimeStepSizes; n++) {
 
  253         ftmp << StepSize[n]  << 
"   " << ErrorNorm[n] << 
"   " 
  254              << error0*(pow(StepSize[n]/StepSize[0],order)) << std::endl;
 
  260   Teuchos::TimeMonitor::summarize();
 
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation  with a...
 
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
 
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...