int TriKota::ThyraDirectApplicInterface::derived_map_ac(const Dakota::String& ac_name)
{

  if (App != Teuchos::null) {

    // Test for consistency of problem definition between ModelEval and Dakota
    TEST_FOR_EXCEPTION(numVars > numParameters, std::logic_error,
                       "TriKota_Dakota Adapter Error: ");
    TEST_FOR_EXCEPTION(numFns > numResponses, std::logic_error,
                       "TriKota_Dakota Adapter Error: ");
    TEST_FOR_EXCEPTION(hessFlag, std::logic_error,
                       "TriKota_Dakota Adapter Error: ");

    MEB::InArgs<double> inArgs = App->createInArgs();
    MEB::OutArgs<double> outArgs = App->createOutArgs();

    TEST_FOR_EXCEPTION(gradFlag && !supportsSensitivities, std::logic_error,
                       "TriKota_Dakota Adapter Error: ");

    // Load parameters from Dakota to ModelEval data structure
    {
      Thyra::DetachedVectorView<double> my_p(model_p);
      for (unsigned int i=0; i<numVars; i++) my_p[i]=xC[i];
    }

    // Evaluate model
    inArgs.set_p(0,model_p);
    outArgs.set_g(0,model_g);
    if (gradFlag) outArgs.set_DgDp(0,0,
      MEB::DerivativeMultiVector<double>(model_dgdp,orientation));
    App->evalModel(inArgs, outArgs);

    Thyra::DetachedVectorView<double> my_g(model_g);
    for (unsigned int j=0; j<numFns; j++) fnVals[j]= my_g[j];

    if (gradFlag) {
      if (orientation == MEB::DERIV_MV_BY_COL) {
        for (unsigned int j=0; j<numVars; j++) {
          Thyra::DetachedVectorView<double>
             my_dgdp_j(model_dgdp->col(j));
          for (unsigned int i=0; i<numFns; i++)  fnGrads[i][j]= my_dgdp_j[i];
        }
      }
      else {
        for (unsigned int j=0; j<numFns; j++) {
          Thyra::DetachedVectorView<double>
             my_dgdp_j(model_dgdp->col(j));
          for (unsigned int i=0; i<numVars; i++) fnGrads[j][i]= my_dgdp_j[i]; 
        }
      }
    }
  }
  else {
    TEST_FOR_EXCEPTION(parallelLib.parallel_configuration().ea_parallel_level().server_intra_communicator()
               != MPI_COMM_NULL, std::logic_error,
              "\nTriKota Parallelism Error: ModelEvaluator=null, but analysis_comm != MPI_COMMM_NULL");
  }

  return 0;
}
Ejemplo n.º 2
0
int main(int argc, char *argv[])
{

  using std::endl;
  typedef double Scalar;
  typedef double ScalarMag;
  using Teuchos::describe;
  using Teuchos::RCP;
  using Teuchos::rcp;
  using Teuchos::rcp_implicit_cast;
  using Teuchos::rcp_dynamic_cast;
  using Teuchos::as;
  using Teuchos::ParameterList;
  using Teuchos::CommandLineProcessor;
  typedef Teuchos::ParameterList::PrintOptions PLPrintOptions;
  typedef Thyra::ModelEvaluatorBase MEB;
  typedef Thyra::DefaultMultiVectorProductVectorSpace<Scalar> DMVPVS;
  using Thyra::productVectorBase;

  bool result, success = true;

  Teuchos::GlobalMPISession mpiSession(&argc,&argv);

  RCP<Epetra_Comm> epetra_comm;
#ifdef HAVE_MPI
  epetra_comm = rcp( new Epetra_MpiComm(MPI_COMM_WORLD) );
#else
  epetra_comm = rcp( new Epetra_SerialComm );
#endif // HAVE_MPI

  RCP<Teuchos::FancyOStream>
    out = Teuchos::VerboseObjectBase::getDefaultOStream();

  try {

    //
    // Read commandline options
    //

    CommandLineProcessor clp;
    clp.throwExceptions(false);
    clp.addOutputSetupOptions(true);

    std::string paramsFileName = "";
    clp.setOption( "params-file", &paramsFileName,
      "File name for XML parameters" );

    std::string extraParamsString = "";
    clp.setOption( "extra-params", &extraParamsString,
      "Extra XML parameters" );

    std::string extraParamsFile = "";
    clp.setOption( "extra-params-file", &extraParamsFile, "File containing extra parameters in XML format.");

    double maxStateError = 1e-6;
    clp.setOption( "max-state-error", &maxStateError,
      "The maximum allowed error in the integrated state in relation to the exact state solution" );

    double finalTime = 1e-3;
    clp.setOption( "final-time", &finalTime,
      "Final integration time (initial time is 0.0)" );

    int numTimeSteps = 10;
    clp.setOption( "num-time-steps", &numTimeSteps,
      "Number of (fixed) time steps.  If <= 0.0, then variable time steps are taken" );

    bool useBDF = false;
    clp.setOption( "use-BDF", "use-BE", &useBDF,
      "Use BDF or Backward Euler (BE)" );

    bool useIRK = false;
    clp.setOption( "use-IRK", "use-other", &useIRK,
      "Use IRK or something" );

    bool doFwdSensSolve = false;
    clp.setOption( "fwd-sens-solve", "state-solve", &doFwdSensSolve,
      "Do the forward sensitivity solve or just the state solve" );

    bool doFwdSensErrorControl = false;
    clp.setOption( "fwd-sens-err-cntrl", "no-fwd-sens-err-cntrl", &doFwdSensErrorControl,
      "Do error control on the forward sensitivity solve or not" );

    double maxRestateError = 0.0;
    clp.setOption( "max-restate-error", &maxRestateError,
      "The maximum allowed error between the state integrated by itself verses integrated along with DxDp" );

    double maxSensError = 1e-4;
    clp.setOption( "max-sens-error", &maxSensError,
      "The maximum allowed error in the integrated sensitivity in relation to"
      " the finite-difference sensitivity" );

    Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_DEFAULT;
    setVerbosityLevelOption( "verb-level", &verbLevel,
      "Top-level verbosity level.  By default, this gets deincremented as you go deeper into numerical objects.",
      &clp );

    bool testExactSensitivity = false;
    clp.setOption( "test-exact-sens", "no-test-exact-sens", &testExactSensitivity,
      "Test the exact sensitivity with finite differences or not." );

    bool dumpFinalSolutions = false;
    clp.setOption(
      "dump-final-solutions", "no-dump-final-solutions", &dumpFinalSolutions,
      "Determine if the final solutions are dumpped or not." );

    CommandLineProcessor::EParseCommandLineReturn parse_return = clp.parse(argc,argv);
    if( parse_return != CommandLineProcessor::PARSE_SUCCESSFUL ) return parse_return;

    if ( Teuchos::VERB_DEFAULT == verbLevel )
      verbLevel = Teuchos::VERB_LOW;

    const Teuchos::EVerbosityLevel
      solnVerbLevel = ( dumpFinalSolutions ? Teuchos::VERB_EXTREME : verbLevel );

    //
    // Get the base parameter list that all other parameter lists will be read
    // from.
    //

    RCP<ParameterList>
      paramList = Teuchos::parameterList();
    if (paramsFileName.length())
      updateParametersFromXmlFile( paramsFileName, paramList.ptr() );
    if(extraParamsFile.length())
      Teuchos::updateParametersFromXmlFile( "./"+extraParamsFile, paramList.ptr() );
    if (extraParamsString.length())
      updateParametersFromXmlString( extraParamsString, paramList.ptr() );

    if (testExactSensitivity) {
      paramList->sublist(DiagonalTransientModel_name).set("Exact Solution as Response",true);
    }

    paramList->validateParameters(*getValidParameters(),0); // Only validate top level lists!

    //
    // Create the Stratimikos linear solver factory.
    //
    // This is the linear solve strategy that will be used to solve for the
    // linear system with the W.
    //

    Stratimikos::DefaultLinearSolverBuilder linearSolverBuilder;
    linearSolverBuilder.setParameterList(sublist(paramList,Stratimikos_name));
    RCP<Thyra::LinearOpWithSolveFactoryBase<Scalar> >
      W_factory = createLinearSolveStrategy(linearSolverBuilder);

    //
    // Create the underlying EpetraExt::ModelEvaluator
    //

    RCP<EpetraExt::DiagonalTransientModel>
      epetraStateModel = EpetraExt::diagonalTransientModel(
        epetra_comm,
        sublist(paramList,DiagonalTransientModel_name)
        );

    *out <<"\nepetraStateModel valid options:\n";
    epetraStateModel->getValidParameters()->print(
      *out, PLPrintOptions().indent(2).showTypes(true).showDoc(true)
      );

    //
    // Create the Thyra-wrapped ModelEvaluator
    //

    RCP<Thyra::ModelEvaluator<double> >
      stateModel = epetraModelEvaluator(epetraStateModel,W_factory);

    *out << "\nParameter names = " << *stateModel->get_p_names(0) << "\n";

    //
    // Create the Rythmos stateStepper
    //

    RCP<Rythmos::TimeStepNonlinearSolver<double> >
      nonlinearSolver = Rythmos::timeStepNonlinearSolver<double>();
    RCP<ParameterList>
      nonlinearSolverPL = sublist(paramList,TimeStepNonlinearSolver_name);
    nonlinearSolverPL->get("Default Tol",1e-3*maxStateError); // Set default if not set
    nonlinearSolver->setParameterList(nonlinearSolverPL);

    RCP<Rythmos::StepperBase<Scalar> > stateStepper;

    if (useBDF) {
      stateStepper = rcp(
        new Rythmos::ImplicitBDFStepper<double>(
          stateModel, nonlinearSolver
          )
        );
    }
    else if (useIRK) {
      // We need a separate LOWSFB object for the IRK stepper
      RCP<Thyra::LinearOpWithSolveFactoryBase<Scalar> >
        irk_W_factory = createLinearSolveStrategy(linearSolverBuilder);
      RCP<Rythmos::RKButcherTableauBase<double> > irkbt = Rythmos::createRKBT<double>("Backward Euler");
      stateStepper = Rythmos::implicitRKStepper<double>(
        stateModel, nonlinearSolver, irk_W_factory, irkbt
        );
    }
    else {
      stateStepper = rcp(
        new Rythmos::BackwardEulerStepper<double>(
          stateModel, nonlinearSolver
          )
        );
    }

    *out <<"\nstateStepper:\n" << describe(*stateStepper,verbLevel);
    *out <<"\nstateStepper valid options:\n";
    stateStepper->getValidParameters()->print(
      *out, PLPrintOptions().indent(2).showTypes(true).showDoc(true)
      );

    stateStepper->setParameterList(sublist(paramList,RythmosStepper_name));

    //
    // Setup finite difference objects that will be used for tests
    //

    Thyra::DirectionalFiniteDiffCalculator<Scalar> fdCalc;
    fdCalc.setParameterList(sublist(paramList,FdCalc_name));
    fdCalc.setOStream(out);
    fdCalc.setVerbLevel(verbLevel);

    //
    // Use a StepperAsModelEvaluator to integrate the state
    //

    const MEB::InArgs<Scalar>
      state_ic = stateModel->getNominalValues();
    *out << "\nstate_ic:\n" << describe(state_ic,verbLevel);

    RCP<Rythmos::IntegratorBase<Scalar> > integrator;
    {
      RCP<ParameterList>
        integratorPL = sublist(paramList,RythmosIntegrator_name);
      integratorPL->set( "Take Variable Steps", as<bool>(numTimeSteps < 0) );
      integratorPL->set( "Fixed dt", as<double>((finalTime - state_ic.get_t())/numTimeSteps) );
      RCP<Rythmos::IntegratorBase<Scalar> >
        defaultIntegrator = Rythmos::controlledDefaultIntegrator<Scalar>(
          Rythmos::simpleIntegrationControlStrategy<Scalar>(integratorPL)
          );
      integrator = defaultIntegrator;
    }

    RCP<Rythmos::StepperAsModelEvaluator<Scalar> >
      stateIntegratorAsModel = Rythmos::stepperAsModelEvaluator(
        stateStepper, integrator, state_ic
        );
    stateIntegratorAsModel->setVerbLevel(verbLevel);

    *out << "\nUse the StepperAsModelEvaluator to integrate state x(p,finalTime) ... \n";

    RCP<Thyra::VectorBase<Scalar> > x_final;

    {

      Teuchos::OSTab tab(out);

      x_final = createMember(stateIntegratorAsModel->get_g_space(0));

      eval_g(
        *stateIntegratorAsModel,
        0, *state_ic.get_p(0),
        finalTime,
        0, &*x_final
        );

      *out
        << "\nx_final = x(p,finalTime) evaluated using stateIntegratorAsModel:\n"
        << describe(*x_final,solnVerbLevel);

    }

    //
    // Test the integrated state against the exact analytical state solution
    //

    RCP<const Thyra::VectorBase<Scalar> >
      exact_x_final = create_Vector(
        epetraStateModel->getExactSolution(finalTime),
        stateModel->get_x_space()
        );

    result = Thyra::testRelNormDiffErr(
      "exact_x_final", *exact_x_final, "x_final", *x_final,
      "maxStateError", maxStateError, "warningTol", 1.0, // Don't warn
      &*out, solnVerbLevel
      );
    if (!result) success = false;

    //
    // Solve and test the forward sensitivity computation
    //

    if (doFwdSensSolve) {

      //
      // Create the forward sensitivity stepper
      //

      RCP<Rythmos::ForwardSensitivityStepper<Scalar> > stateAndSensStepper =
        Rythmos::forwardSensitivityStepper<Scalar>();
      if (doFwdSensErrorControl) {
        stateAndSensStepper->initializeDecoupledSteppers(
          stateModel, 0, stateModel->getNominalValues(),
          stateStepper, nonlinearSolver,
          integrator->cloneIntegrator(), finalTime
          );
      }
      else {
        stateAndSensStepper->initializeSyncedSteppers(
          stateModel, 0, stateModel->getNominalValues(),
          stateStepper, nonlinearSolver
          );
        // The above call will result in stateStepper and nonlinearSolver being
        // cloned.  This helps to ensure consistency between the state and
        // sensitivity computations!
      }

      //
      // Set the initial condition for the state and forward sensitivities
      //

      RCP<Thyra::VectorBase<Scalar> > s_bar_init
        = createMember(stateAndSensStepper->getFwdSensModel()->get_x_space());
      assign( s_bar_init.ptr(), 0.0 );
      RCP<Thyra::VectorBase<Scalar> > s_bar_dot_init
        = createMember(stateAndSensStepper->getFwdSensModel()->get_x_space());
      assign( s_bar_dot_init.ptr(), 0.0 );
      // Above, I believe that these are the correct initial conditions for
      // s_bar and s_bar_dot given how the EpetraExt::DiagonalTransientModel
      // is currently implemented!

      RCP<const Rythmos::StateAndForwardSensitivityModelEvaluator<Scalar> >
        stateAndSensModel = stateAndSensStepper->getStateAndFwdSensModel();

      MEB::InArgs<Scalar>
        state_and_sens_ic = stateAndSensStepper->getModel()->createInArgs();

      // Copy time, parameters etc.
      state_and_sens_ic.setArgs(state_ic);
      // Set initial condition for x_bar = [ x; s_bar ]
      state_and_sens_ic.set_x(
        stateAndSensModel->create_x_bar_vec(state_ic.get_x(),s_bar_init)
        );
      // Set initial condition for x_bar_dot = [ x_dot; s_bar_dot ]
      state_and_sens_ic.set_x_dot(
        stateAndSensModel->create_x_bar_vec(state_ic.get_x_dot(),s_bar_dot_init)
        );

      *out << "\nstate_and_sens_ic:\n" << describe(state_and_sens_ic,verbLevel);

      stateAndSensStepper->setInitialCondition(state_and_sens_ic);

      //
      // Use a StepperAsModelEvaluator to integrate the state+sens
      //

      RCP<Rythmos::StepperAsModelEvaluator<Scalar> >
        stateAndSensIntegratorAsModel = Rythmos::stepperAsModelEvaluator(
          rcp_implicit_cast<Rythmos::StepperBase<Scalar> >(stateAndSensStepper),
          integrator, state_and_sens_ic
          );
      stateAndSensIntegratorAsModel->setVerbLevel(verbLevel);

      *out << "\nUse the StepperAsModelEvaluator to integrate state + sens x_bar(p,finalTime) ... \n";

      RCP<Thyra::VectorBase<Scalar> > x_bar_final;

      {

        Teuchos::OSTab tab(out);

        x_bar_final = createMember(stateAndSensIntegratorAsModel->get_g_space(0));

        eval_g(
          *stateAndSensIntegratorAsModel,
          0, *state_ic.get_p(0),
          finalTime,
          0, &*x_bar_final
          );

        *out
          << "\nx_bar_final = x_bar(p,finalTime) evaluated using stateAndSensIntegratorAsModel:\n"
          << describe(*x_bar_final,solnVerbLevel);

      }

      //
      // Test that the state computed above is same as computed initially!
      //

      *out << "\nChecking that x(p,finalTime) computed as part of x_bar above is the same ...\n";

      {

        Teuchos::OSTab tab(out);

        RCP<const Thyra::VectorBase<Scalar> >
          x_in_x_bar_final = productVectorBase<Scalar>(x_bar_final)->getVectorBlock(0);

        result = Thyra::testRelNormDiffErr<Scalar>(
          "x_final", *x_final,
          "x_in_x_bar_final", *x_in_x_bar_final,
          "maxRestateError", maxRestateError,
          "warningTol", 1.0, // Don't warn
          &*out, solnVerbLevel
          );
        if (!result) success = false;

      }

      //
      // Compute DxDp using finite differences
      //

      *out << "\nApproximating DxDp(p,t) using directional finite differences of integrator for x(p,t) ...\n";

      RCP<Thyra::MultiVectorBase<Scalar> > DxDp_fd_final;

      {

        Teuchos::OSTab tab(out);


        MEB::InArgs<Scalar>
          fdBasePoint = stateIntegratorAsModel->createInArgs();

        fdBasePoint.set_t(finalTime);
        fdBasePoint.set_p(0,stateModel->getNominalValues().get_p(0));

        DxDp_fd_final = createMembers(
          stateIntegratorAsModel->get_g_space(0),
          stateIntegratorAsModel->get_p_space(0)->dim()
          );

        typedef Thyra::DirectionalFiniteDiffCalculatorTypes::SelectedDerivatives
          SelectedDerivatives;

        MEB::OutArgs<Scalar> fdOutArgs =
          fdCalc.createOutArgs(
            *stateIntegratorAsModel,
            SelectedDerivatives().supports(MEB::OUT_ARG_DgDp,0,0)
            );
        fdOutArgs.set_DgDp(0,0,DxDp_fd_final);

        // Silence the model evaluators that are called.  The fdCal object
        // will show all of the inputs and outputs for each call.
        stateStepper->setVerbLevel(Teuchos::VERB_NONE);
        stateIntegratorAsModel->setVerbLevel(Teuchos::VERB_NONE);

        fdCalc.calcDerivatives(
          *stateIntegratorAsModel, fdBasePoint,
          stateIntegratorAsModel->createOutArgs(), // Don't bother with function value
          fdOutArgs
          );

        *out
          << "\nFinite difference DxDp_fd_final = DxDp(p,finalTime): "
          << describe(*DxDp_fd_final,solnVerbLevel);

      }

      //
      // Test that the integrated sens and the F.D. sens are similar
      //

      *out << "\nChecking that integrated DxDp(p,finalTime) and finite-diff DxDp(p,finalTime) are similar ...\n";

      {

        Teuchos::OSTab tab(out);

        RCP<const Thyra::VectorBase<Scalar> >
          DxDp_vec_final = Thyra::productVectorBase<Scalar>(x_bar_final)->getVectorBlock(1);

        RCP<const Thyra::VectorBase<Scalar> >
          DxDp_fd_vec_final = Thyra::multiVectorProductVector(
            rcp_dynamic_cast<const Thyra::DefaultMultiVectorProductVectorSpace<Scalar> >(
              DxDp_vec_final->range()
              ),
            DxDp_fd_final
            );

        result = Thyra::testRelNormDiffErr(
          "DxDp_vec_final", *DxDp_vec_final,
          "DxDp_fd_vec_final", *DxDp_fd_vec_final,
          "maxSensError", maxSensError,
          "warningTol", 1.0, // Don't warn
          &*out, solnVerbLevel
          );
        if (!result) success = false;

      }

    }

  }
  TEUCHOS_STANDARD_CATCH_STATEMENTS(true,*out,success);

  if(success)
    *out << "\nEnd Result: TEST PASSED" << endl;
  else
    *out << "\nEnd Result: TEST FAILED" << endl;

  return ( success ? 0 : 1 );

} // end main() [Doxygen looks for this!]
Ejemplo n.º 3
0
void ModelEvaluatorDefaultBase<Scalar>::evalModel(
  const ModelEvaluatorBase::InArgs<Scalar> &inArgs,
  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs
  ) const
{

  using Teuchos::outArg;
  typedef ModelEvaluatorBase MEB;

  lazyInitializeDefaultBase();

  const int l_Np = outArgs.Np();
  const int l_Ng = outArgs.Ng();

  //
  // A) Assert that the inArgs and outArgs object match this class!
  //

#ifdef TEUCHOS_DEBUG
  assertInArgsEvalObjects(*this,inArgs);
  assertOutArgsEvalObjects(*this,outArgs,&inArgs);
#endif  
  
  //
  // B) Setup the OutArgs object for the underlying implementation's
  // evalModelImpl(...) function
  //

  MEB::OutArgs<Scalar> outArgsImpl = this->createOutArgsImpl();
  Array<MultiVectorAdjointPair> DgDp_temp_adjoint_copies;

  {

    outArgsImpl.setArgs(outArgs,true);

    // DfDp(l)
    if (outArgsImpl.supports(MEB::OUT_ARG_f)) {
      for ( int l = 0; l < l_Np; ++l ) {
        const DefaultDerivLinearOpSupport defaultLinearOpSupport =
          DfDp_default_op_support_[l];
        if (defaultLinearOpSupport.provideDefaultLinearOp()) {
          outArgsImpl.set_DfDp( l,
            getOutArgImplForDefaultLinearOpSupport(
              outArgs.get_DfDp(l), defaultLinearOpSupport
              )
            );
        }
        else {
          // DfDp(l) already set by outArgsImpl.setArgs(...)!
        }
      }
    }

    // DgDx_dot(j)
    for ( int j = 0; j < l_Ng; ++j ) {
      const DefaultDerivLinearOpSupport defaultLinearOpSupport =
        DgDx_dot_default_op_support_[j];
      if (defaultLinearOpSupport.provideDefaultLinearOp()) {
        outArgsImpl.set_DgDx_dot( j,
          getOutArgImplForDefaultLinearOpSupport(
            outArgs.get_DgDx_dot(j), defaultLinearOpSupport
            )
          );
      }
      else {
        // DgDx_dot(j) already set by outArgsImpl.setArgs(...)!
      }
    }

    // DgDx(j)
    for ( int j = 0; j < l_Ng; ++j ) {
      const DefaultDerivLinearOpSupport defaultLinearOpSupport =
        DgDx_default_op_support_[j];
      if (defaultLinearOpSupport.provideDefaultLinearOp()) {
        outArgsImpl.set_DgDx( j,
          getOutArgImplForDefaultLinearOpSupport(
            outArgs.get_DgDx(j), defaultLinearOpSupport
            )
          );
      }
      else {
        // DgDx(j) already set by outArgsImpl.setArgs(...)!
      }
    }

    // DgDp(j,l)
    for ( int j = 0; j < l_Ng; ++j ) {
      const Array<DefaultDerivLinearOpSupport> &DgDp_default_op_support_j =
        DgDp_default_op_support_[j];
      const Array<DefaultDerivMvAdjointSupport> &DgDp_default_mv_support_j =
        DgDp_default_mv_support_[j];
      for ( int l = 0; l < l_Np; ++l ) {
        const DefaultDerivLinearOpSupport defaultLinearOpSupport =
          DgDp_default_op_support_j[l];
        const DefaultDerivMvAdjointSupport defaultMvAdjointSupport =
          DgDp_default_mv_support_j[l];
        MEB::Derivative<Scalar> DgDp_j_l;
        if (!outArgs.supports(MEB::OUT_ARG_DgDp,j,l).none())
          DgDp_j_l = outArgs.get_DgDp(j,l);
        if (
          defaultLinearOpSupport.provideDefaultLinearOp()
          && !is_null(DgDp_j_l.getLinearOp())
          )
        {
          outArgsImpl.set_DgDp( j, l,
            getOutArgImplForDefaultLinearOpSupport(
              DgDp_j_l, defaultLinearOpSupport
              )
            );
        }
        else if (
          defaultMvAdjointSupport.provideDefaultAdjoint()
          && !is_null(DgDp_j_l.getMultiVector())
          )
        {
          const RCP<MultiVectorBase<Scalar> > DgDp_j_l_mv = 
            DgDp_j_l.getMultiVector();
          if (
            defaultMvAdjointSupport.mvAdjointCopyOrientation()
            ==
            DgDp_j_l.getMultiVectorOrientation()
            )
          {
            // The orientation of the multi-vector is different so we need to
            // create a temporary copy to pass to evalModelImpl(...) and then
            // copy it back again!
            const RCP<MultiVectorBase<Scalar> > DgDp_j_l_mv_adj =
              createMembers(DgDp_j_l_mv->domain(), DgDp_j_l_mv->range()->dim());
            outArgsImpl.set_DgDp( j, l,
              MEB::Derivative<Scalar>(
                DgDp_j_l_mv_adj,
                getOtherDerivativeMultiVectorOrientation(
                  defaultMvAdjointSupport.mvAdjointCopyOrientation()
                  )
                )
              );
            // Remember these multi-vectors so that we can do the transpose copy
            // back after the evaluation!
            DgDp_temp_adjoint_copies.push_back(
              MultiVectorAdjointPair(DgDp_j_l_mv, DgDp_j_l_mv_adj)
              );
          }
          else {
            // The form of the multi-vector is supported by evalModelImpl(..)
            // and is already set on the outArgsImpl object.
          }
        }
        else {
          // DgDp(j,l) already set by outArgsImpl.setArgs(...)!
        }
      }
    }

    // W
    {
      RCP<LinearOpWithSolveBase<Scalar> > W;
      if ( default_W_support_ && !is_null(W=outArgs.get_W()) ) {
        const RCP<const LinearOpWithSolveFactoryBase<Scalar> >
          W_factory = this->get_W_factory();
        // Extract the underlying W_op object (if it exists)
        RCP<const LinearOpBase<Scalar> > W_op_const;
        uninitializeOp<Scalar>(*W_factory, W.ptr(), outArg(W_op_const));
        RCP<LinearOpBase<Scalar> > W_op;
        if (!is_null(W_op_const)) {
          // Here we remove the const.  This is perfectly safe since
          // w.r.t. this class, we put a non-const object in there and we can
          // expect to change that object after the fact.  That is our
          // prerogative.
          W_op = Teuchos::rcp_const_cast<LinearOpBase<Scalar> >(W_op_const);
        }
        else {
          // The W_op object has not been initialized yet so create it.  The
          // next time through, we should not have to do this!
          W_op = this->create_W_op();
        }
        outArgsImpl.set_W_op(W_op);
      }
    }
    
  }

  //
  // C) Evaluate the underlying model implementation!
  //

  this->evalModelImpl( inArgs, outArgsImpl );

  //
  // D) Post-process the output arguments
  //

  // Do explicit transposes for DgDp(j,l) if needed
  const int numMvAdjointCopies = DgDp_temp_adjoint_copies.size();
  for ( int adj_copy_i = 0; adj_copy_i < numMvAdjointCopies; ++adj_copy_i ) {
    const MultiVectorAdjointPair adjPair =
      DgDp_temp_adjoint_copies[adj_copy_i];
    doExplicitMultiVectorAdjoint( *adjPair.mvImplAdjoint, &*adjPair.mvOuter );
  }
  
  // Update W given W_op and W_factory
  {
    RCP<LinearOpWithSolveBase<Scalar> > W;
    if ( default_W_support_ && !is_null(W=outArgs.get_W()) ) {
      const RCP<const LinearOpWithSolveFactoryBase<Scalar> >
        W_factory = this->get_W_factory();
      W_factory->setOStream(this->getOStream());
      W_factory->setVerbLevel(this->getVerbLevel());
      initializeOp<Scalar>(*W_factory, outArgsImpl.get_W_op().getConst(), W.ptr());
    }
  }
  
}