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;
 
   44                      const bool use_combined_method,
 
   45                      const bool use_dfdp_as_tangent,
 
   46                      Teuchos::FancyOStream &out, 
bool &success)
 
   48   std::vector<std::string> RKMethods;
 
   49   RKMethods.push_back(
"General DIRK");
 
   50   RKMethods.push_back(
"RK Backward Euler");
 
   51   RKMethods.push_back(
"DIRK 1 Stage Theta Method");
 
   52   RKMethods.push_back(
"RK Implicit 1 Stage 1st order Radau IA");
 
   53   RKMethods.push_back(
"RK Implicit Midpoint");
 
   54   RKMethods.push_back(
"SDIRK 2 Stage 2nd order");
 
   55   RKMethods.push_back(
"RK Implicit 2 Stage 2nd order Lobatto IIIB");
 
   56   RKMethods.push_back(
"SDIRK 2 Stage 3rd order");
 
   57   RKMethods.push_back(
"EDIRK 2 Stage 3rd order");
 
   58   RKMethods.push_back(
"EDIRK 2 Stage Theta Method");
 
   59   RKMethods.push_back(
"SDIRK 3 Stage 4th order");
 
   60   RKMethods.push_back(
"SDIRK 5 Stage 4th order");
 
   61   RKMethods.push_back(
"SDIRK 5 Stage 5th order");
 
   62   RKMethods.push_back(
"SDIRK 2(1) Pair");
 
   63   RKMethods.push_back(
"RK Trapezoidal Rule");
 
   64   RKMethods.push_back(
"RK Crank-Nicolson");
 
   67   if (method_name != 
"") {
 
   68     auto it = std::find(RKMethods.begin(), RKMethods.end(), 
method_name);
 
   69     TEUCHOS_TEST_FOR_EXCEPTION(it == RKMethods.end(), std::logic_error,
 
   70       "Invalid RK method name '" << method_name << 
"'");
 
   73   std::vector<double> RKMethodErrors;
 
   74   if (use_combined_method) {
 
   75     RKMethodErrors.push_back(0.000144507);
 
   76     RKMethodErrors.push_back(0.0428449);
 
   77     RKMethodErrors.push_back(0.000297933);
 
   78     RKMethodErrors.push_back(0.0428449);
 
   79     RKMethodErrors.push_back(0.000297933);
 
   80     RKMethodErrors.push_back(0.000144507);
 
   81     RKMethodErrors.push_back(0.000297933);
 
   82     RKMethodErrors.push_back(8.65434e-06);
 
   83     RKMethodErrors.push_back(1.3468e-06);
 
   84     RKMethodErrors.push_back(0.000297933);
 
   85     RKMethodErrors.push_back(5.44037e-07);
 
   86     RKMethodErrors.push_back(2.77342e-09);
 
   87     RKMethodErrors.push_back(1.21689e-10);
 
   88     RKMethodErrors.push_back(0.000603848);
 
   89     RKMethodErrors.push_back(0.000297933);
 
   90     RKMethodErrors.push_back(0.000297933);
 
   93     RKMethodErrors.push_back(0.000125232);
 
   94     RKMethodErrors.push_back(0.0428449);
 
   95     RKMethodErrors.push_back(0.000221049);
 
   96     RKMethodErrors.push_back(0.0383339);
 
   97     RKMethodErrors.push_back(0.000221049);
 
   98     RKMethodErrors.push_back(0.000125232);
 
   99     RKMethodErrors.push_back(0.000272997);
 
  100     RKMethodErrors.push_back(4.79475e-06);
 
  101     RKMethodErrors.push_back(9.63899e-07);
 
  102     RKMethodErrors.push_back(0.000297933);
 
  103     RKMethodErrors.push_back(2.9362e-07);
 
  104     RKMethodErrors.push_back(9.20081e-08);
 
  105     RKMethodErrors.push_back(9.16252e-08);
 
  106     RKMethodErrors.push_back(0.00043969);
 
  107     RKMethodErrors.push_back(0.000297933);
 
  108     RKMethodErrors.push_back(0.000297933);
 
  111   Teuchos::RCP<const Teuchos::Comm<int> > comm =
 
  112     Teuchos::DefaultComm<int>::getComm();
 
  113   Teuchos::RCP<Teuchos::FancyOStream> my_out =
 
  114     Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
 
  115   my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
 
  116   my_out->setOutputToRootOnly(0);
 
  118   for(std::vector<std::string>::size_type m = 0; m != RKMethods.size(); m++) {
 
  121     if (method_name != 
"" && RKMethods[m] != method_name)
 
  124     std::string RKMethod_ = RKMethods[m];
 
  125     std::replace(RKMethod_.begin(), RKMethod_.end(), 
' ', 
'_');
 
  126     std::replace(RKMethod_.begin(), RKMethod_.end(), 
'/', 
'.');
 
  127     std::vector<double> StepSize;
 
  128     std::vector<double> ErrorNorm;
 
  129     const int nTimeStepSizes = 3; 
 
  132     for (
int n=0; n<nTimeStepSizes; n++) {
 
  135       RCP<ParameterList> pList =
 
  136         getParametersFromXmlFile(
"Tempus_DIRK_SinCos.xml");
 
  139       RCP<ParameterList> scm_pl = sublist(pList, 
"SinCosModel", 
true);
 
  140       scm_pl->set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
 
  141       RCP<SinCosModel<double> > model =
 
  145       RCP<ParameterList> pl = sublist(pList, 
"Tempus", 
true);
 
  146       pl->sublist(
"Default Stepper").set(
"Stepper Type", RKMethods[m]);
 
  147       if (RKMethods[m] == 
"SDIRK 2 Stage 2nd order") {
 
  148         pl->sublist(
"Default Stepper").set(
"gamma", 0.2928932188134524);
 
  149       } 
else if (RKMethods[m] == 
"SDIRK 2 Stage 3rd order") {
 
  150         pl->sublist(
"Default Stepper")
 
  151            .set<std::string>(
"Gamma Type", 
"3rd Order A-stable");
 
  157       ParameterList& sens_pl = pl->sublist(
"Sensitivities");
 
  158       if (use_combined_method)
 
  159         sens_pl.set(
"Sensitivity Method", 
"Combined");
 
  161         sens_pl.set(
"Sensitivity Method", 
"Staggered");
 
  162       sens_pl.set(
"Use DfDp as Tangent", use_dfdp_as_tangent);
 
  163       ParameterList& interp_pl =
 
  164         pl->sublist(
"Default Integrator").sublist(
"Solution History").sublist(
"Interpolator");
 
  165       interp_pl.set(
"Interpolator Type", 
"Lagrange");
 
  166       interp_pl.set(
"Order", 4); 
 
  169       pl->sublist(
"Default Integrator")
 
  170          .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
 
  171       RCP<Tempus::IntegratorForwardSensitivity<double> > integrator =
 
  172         Tempus::integratorForwardSensitivity<double>(pl, model);
 
  173       order = integrator->getStepper()->getOrder();
 
  179       RCP<Thyra::VectorBase<double> > x0 =
 
  180         model->getNominalValues().get_x()->clone_v();
 
  181       const int num_param = model->get_p_space(0)->dim();
 
  182       RCP<Thyra::MultiVectorBase<double> > DxDp0 =
 
  183         Thyra::createMembers(model->get_x_space(), num_param);
 
  184       for (
int i=0; i<num_param; ++i)
 
  185         Thyra::assign(DxDp0->col(i).ptr(),
 
  186                       *(model->getExactSensSolution(i, 0.0).get_x()));
 
  187       integrator->initializeSolutionHistory(0.0, x0, Teuchos::null, Teuchos::null,
 
  188                                   DxDp0, Teuchos::null, Teuchos::null);
 
  191       bool integratorStatus = integrator->advanceTime();
 
  192       TEST_ASSERT(integratorStatus)
 
  195       double time = integrator->getTime();
 
  196       double timeFinal = pl->sublist(
"Default Integrator")
 
  197         .sublist(
"Time Step Control").get<
double>(
"Final Time");
 
  198       double tol = 100.0 * std::numeric_limits<double>::epsilon();
 
  199       TEST_FLOATING_EQUALITY(time, timeFinal, tol);
 
  202       RCP<const Thyra::VectorBase<double> > x = integrator->getX();
 
  203       RCP<const Thyra::MultiVectorBase<double> > DxDp = integrator->getDxDp();
 
  204       RCP<const Thyra::VectorBase<double> > x_exact =
 
  205         model->getExactSolution(time).get_x();
 
  206       RCP<Thyra::MultiVectorBase<double> > DxDp_exact =
 
  207         Thyra::createMembers(model->get_x_space(), num_param);
 
  208       for (
int i=0; i<num_param; ++i)
 
  209         Thyra::assign(DxDp_exact->col(i).ptr(),
 
  210                       *(model->getExactSensSolution(i, time).get_x()));
 
  213       if (comm->getRank() == 0 && n == nTimeStepSizes-1) {
 
  214         typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
 
  216         std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_Sens.dat");
 
  217         RCP<const SolutionHistory<double> > solutionHistory =
 
  218           integrator->getSolutionHistory();
 
  219         RCP< Thyra::MultiVectorBase<double> > DxDp_exact_plot =
 
  220           Thyra::createMembers(model->get_x_space(), num_param);
 
  221         for (
int i=0; i<solutionHistory->getNumStates(); i++) {
 
  222           RCP<const SolutionState<double> > solutionState =
 
  223             (*solutionHistory)[i];
 
  224           double time_i = solutionState->getTime();
 
  225           RCP<const DMVPV> x_prod_plot =
 
  226             Teuchos::rcp_dynamic_cast<
const DMVPV>(solutionState->getX());
 
  227           RCP<const Thyra::VectorBase<double> > x_plot =
 
  228             x_prod_plot->getMultiVector()->col(0);
 
  229           RCP<const Thyra::MultiVectorBase<double> > DxDp_plot =
 
  230             x_prod_plot->getMultiVector()->subView(Teuchos::Range1D(1,num_param));
 
  231           RCP<const Thyra::VectorBase<double> > x_exact_plot =
 
  232             model->getExactSolution(time_i).get_x();
 
  233           for (
int j=0; j<num_param; ++j)
 
  234             Thyra::assign(DxDp_exact_plot->col(j).ptr(),
 
  235                           *(model->getExactSensSolution(j, time_i).get_x()));
 
  236           ftmp << std::fixed << std::setprecision(7)
 
  238                << std::setw(11) << get_ele(*(x_plot), 0)
 
  239                << std::setw(11) << get_ele(*(x_plot), 1);
 
  240           for (
int j=0; j<num_param; ++j)
 
  241             ftmp << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 0)
 
  242                  << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 1);
 
  243           ftmp << std::setw(11) << get_ele(*(x_exact_plot), 0)
 
  244                << std::setw(11) << get_ele(*(x_exact_plot), 1);
 
  245           for (
int j=0; j<num_param; ++j)
 
  246             ftmp << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 0)
 
  247                  << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 1);
 
  254       RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
 
  255       RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
 
  256       Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
 
  257       Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
 
  258       StepSize.push_back(dt);
 
  259       double L2norm = Thyra::norm_2(*xdiff);
 
  261       Teuchos::Array<double> L2norm_DxDp(num_param);
 
  262       Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
 
  263       for (
int i=0; i<num_param; ++i)
 
  264         L2norm += L2norm_DxDp[i]*L2norm_DxDp[i];
 
  265       L2norm = std::sqrt(L2norm);
 
  266       ErrorNorm.push_back(L2norm);
 
  272     if (comm->getRank() == 0) {
 
  273       std::ofstream ftmp(
"Tempus_"+RKMethod_+
"_SinCos_Sens-Error.dat");
 
  274       double error0 = 0.8*ErrorNorm[0];
 
  275       for (
int n=0; n<(int)StepSize.size(); n++) {
 
  276         ftmp << StepSize[n]  << 
"   " << ErrorNorm[n] << 
"   " 
  277              << error0*(pow(StepSize[n]/StepSize[0],order)) << std::endl;
 
  291     double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
 
  292     *my_out << 
"  Stepper = " << RKMethods[m] << std::endl;
 
  293     *my_out << 
"  =========================" << std::endl;
 
  294     *my_out << 
"  Expected order: " << order << std::endl;
 
  295     *my_out << 
"  Observed order: " << slope << std::endl;
 
  296     *my_out << 
"  =========================" << std::endl;
 
  299     double order_expected = use_combined_method ? order : std::min(order,4.0);
 
  300     TEST_FLOATING_EQUALITY( slope, order_expected, 0.03 );
 
  301     TEST_FLOATING_EQUALITY( ErrorNorm[0], RKMethodErrors[m], 5.0e-4 );
 
  304   Teuchos::TimeMonitor::summarize();
 
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation  with a...
 
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...