int
Piro::PerformDakotaAnalysis(
    Thyra::ModelEvaluatorDefaultBase<double>& piroModel,
    Teuchos::ParameterList& dakotaParams,
    RCP< Thyra::VectorBase<double> >& p)
{
#ifdef HAVE_PIRO_TRIKOTA
  dakotaParams.validateParameters(*Piro::getValidPiroAnalysisDakotaParameters(),0);
  using std::string;

  string dakotaIn  = dakotaParams.get("Input File","dakota.in");
  string dakotaOut = dakotaParams.get("Output File","dakota.out");
  string dakotaErr = dakotaParams.get("Error File","dakota.err");
  string dakotaRes = dakotaParams.get("Restart File","dakota_restart.out");
  string dakotaRestartIn;
  if (dakotaParams.isParameter("Restart File To Read"))
    dakotaRestartIn = dakotaParams.get<string>("Restart File To Read");

  int dakotaRestartEvals= dakotaParams.get("Restart Evals To Read", 0);

  int p_index = dakotaParams.get("Parameter Vector Index", 0);
  int g_index = dakotaParams.get("Response Vector Index", 0);

  TriKota::Driver dakota(dakotaIn, dakotaOut, dakotaErr, dakotaRes,
                         dakotaRestartIn, dakotaRestartEvals);

  RCP<TriKota::ThyraDirectApplicInterface> trikota_interface =
    rcp(new TriKota::ThyraDirectApplicInterface
         (dakota.getProblemDescDB(), rcp(&piroModel,false), p_index, g_index),
	false);

  dakota.run(trikota_interface.get());

  Dakota::RealVector finalValues;
  if (dakota.rankZero())
    finalValues = dakota.getFinalSolution().all_continuous_variables();

  // Copy Dakota parameters into Thyra
  p = Thyra::createMember(piroModel.get_p_space(p_index));
  {
      Thyra::DetachedVectorView<double> global_p(p);
      for (int i = 0; i < finalValues.length(); ++i)
        global_p[i] = finalValues[i];
  }

  return 0;
#else
 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
 *out << "ERROR: Trilinos/Piro was not configured to include Dakota analysis."
      << "\nYou must enable TriKota." << endl;
 return 0;  // should not fail tests
#endif
}
int
Piro::Thyra::PerformDakotaAnalysis(
    ::Thyra::ModelEvaluatorDefaultBase<double>& piroModel,
    Teuchos::ParameterList& dakotaParams,
    RCP< ::Thyra::VectorBase<double> >& p)
{
#ifdef Piro_ENABLE_TriKota
  string dakotaIn = dakotaParams.get("Input File","dakota.in");
  string dakotaOut= dakotaParams.get("Output File","dakota.out");
  string dakotaErr= dakotaParams.get("Error File","dakota.err");
  string dakotaRes= dakotaParams.get("Restart File","dakota_restart.out");

  TriKota::Driver dakota(dakotaIn.c_str(), dakotaOut.c_str(),
                         dakotaErr.c_str(), dakotaRes.c_str());

  RCP<TriKota::ThyraDirectApplicInterface> trikota_interface =
    rcp(new TriKota::ThyraDirectApplicInterface
         (dakota.getProblemDescDB(), rcp(&piroModel,false)), false);

  dakota.run(trikota_interface.get());

  Dakota::RealVector finalValues =
    dakota.getFinalSolution().all_continuous_variables();

  // Copy Dakota parameters into Thyra
  p = ::Thyra::createMember(piroModel.get_p_space(0));
  {
      ::Thyra::DetachedVectorView<double> global_p(p);
      for (int i = 0; i < finalValues.length(); ++i) 
        global_p[i] = finalValues[i];
  }

  return 0;
#else
 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
 *out << "ERROR: Trilinos/Piro was not configured to include Dakota analysis."
      << "\nYou must enable TriKota." << endl;
 return 0;  // should not fail tests
#endif
}
void 
TriKota::MPDirectApplicInterface::
wait_local_evaluations(Dakota::PRPQueue& prp_queue)
{
  using Teuchos::RCP;
  using Teuchos::rcp;
  using Stokhos::ProductEpetraVector;
  using Stokhos::ProductEpetraMultiVector;

  if (model != Teuchos::null) {
    unsigned int queue_size = prp_queue.size();
    unsigned int num_blocks = queue_size / block_size;
    unsigned int remainder = queue_size % block_size;

    *out << "TriKota:: Performing multi-point evaluation of " << queue_size
	 << " responses using a block size of " << block_size
	 << std::endl;
    
    Dakota::PRPQueueIter block_iter = prp_queue.begin();
    unsigned int block = 0;
    while (block_iter != prp_queue.end()) {

      // Load up block p
      Dakota::PRPQueueIter iter = block_iter;
      bool valFlag = true;
      bool gradFlag = true;
      bool hessFlag = true;
      unsigned int blk_sz = block_size;
      if (block == num_blocks && remainder > 0)
	blk_sz = remainder;
      for (unsigned int idx=0; idx<blk_sz; idx++) {
	const Dakota::Variables& vars = iter->variables();
	const Dakota::RealVector& xC  = vars.continuous_variables();
	unsigned int numVars = xC.length();
	for (unsigned int i=0; i<numVars; i++) 
	  (*model_p)[idx][i]=xC[i];

	const Dakota::ActiveSet& set  = iter->active_set();
	short asv = set.request_vector()[0];
	valFlag  = valFlag  && (asv & 1);
	gradFlag = gradFlag && (asv & 2);
	hessFlag = hessFlag && (asv & 4);

	Dakota::Response resp = iter->response();
	Dakota::RealVector fnVals = resp.function_values_view();
	unsigned int numFns = fnVals.length();
	
	// test for consistency of problem definition between ModelEval and 
	// Dakota
	TEUCHOS_TEST_FOR_EXCEPTION(numVars > numParameters, std::logic_error,
				   "TriKota_Dakota Adapter Error: ");
	TEUCHOS_TEST_FOR_EXCEPTION(hessFlag, std::logic_error,
				   "TriKota_Dakota Adapter Error: ");
	TEUCHOS_TEST_FOR_EXCEPTION(numFns > numResponses, std::logic_error,
				   "TriKota_Dakota Adapter Error: ");
	TEUCHOS_TEST_FOR_EXCEPTION(gradFlag && !supportsSensitivities, 
				   std::logic_error,
				   "TriKota_Dakota Adapter Error: ");
	
	++iter;
      }

      // Put in copies of last point for remainder
      if (block == num_blocks && remainder > 0) {
	--iter;
	for (unsigned int idx=remainder; idx<block_size; idx++) {
	  const Dakota::Variables& vars = iter->variables();
	  const Dakota::RealVector& xC  = vars.continuous_variables();
	  unsigned int numVars = xC.length();
	  for (unsigned int i=0; i<numVars; i++) 
	    (*model_p)[idx][i]=xC[i];
	}
      }

      // Evaluate model
      EEME::InArgs inArgs = model->createInArgs();
      EEME::OutArgs outArgs = model->createOutArgs();
      inArgs.set_p_mp(p_index, model_p);
      if (valFlag)
	outArgs.set_g_mp(g_index, model_g);
      if (gradFlag) 
	outArgs.set_DgDp_mp(g_index, p_index,
			    EEME::MPDerivative(model_dgdp, orientation));
      *out << "TriKota:: evaluating block " << block;
      if (gradFlag)
	*out << " with gradients." << std::endl;
      else
	*out << " without gradients." << std::endl;
      model->evalModel(inArgs, outArgs);

      // Copy responses from block g
      iter = block_iter;
      for (unsigned int idx=0; idx<blk_sz; idx++) {
	const Dakota::Variables& vars = iter->variables();
	const Dakota::RealVector& xC  = vars.continuous_variables();
	unsigned int numVars = xC.length();
	Dakota::Response         resp = iter->response(); // shared rep
	Dakota::RealVector fnVals     = resp.function_values_view();
	Dakota::RealMatrix fnGrads    = resp.function_gradients_view();
	unsigned int numFns = fnVals.length();
	
	if (valFlag)
	  for (unsigned int j=0; j<numFns; j++) 
	    fnVals[j]= (*model_g)[idx][j];
	if (gradFlag)  {
	  if (orientation == EEME::DERIV_MV_BY_COL) {
	    for (unsigned int i=0; i<numVars; i++)
	      for (unsigned int j=0; j<numFns; j++)
		fnGrads[j][i]= (*model_dgdp)[idx][i][j];
	  } 
	  else {
	    for (unsigned int i=0; i<numFns; i++)
	      for (unsigned int j=0; j<numVars; j++)
		fnGrads[i][j]= (*model_dgdp)[idx][i][j];
	  }
	}

	// indicate completion of job to ApplicationInterface schedulers
	int fn_eval_id = iter->eval_id();
	completionSet.insert(fn_eval_id);

	++iter;
	++block_iter;
      }
      ++block;
    }
    TEUCHOS_TEST_FOR_EXCEPT((remainder == 0 && block != num_blocks) || 
			    (remainder  > 0 && block != num_blocks+1))
  }