void
Albany::SamplingBasedScalarResponseFunction::
evaluateSGResponse(
  const double curr_time,
  const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
  const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
  const Stokhos::EpetraVectorOrthogPoly& sg_x,
  const Teuchos::Array<ParamVec>& p,
  const Teuchos::Array<int>& sg_p_index,
  const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
  Stokhos::EpetraVectorOrthogPoly& sg_g)
{
  RCP<const Epetra_BlockMap> x_map = sg_x.coefficientMap();
  RCP<Epetra_Vector> xdot;
  if (sg_xdot != NULL)
    xdot = rcp(new Epetra_Vector(*x_map));
  RCP<Epetra_Vector> xdotdot;
  if (sg_xdotdot != NULL)
    xdotdot = rcp(new Epetra_Vector(*x_map));
  Epetra_Vector x(*x_map);
  Teuchos::Array<ParamVec> pp = p;
  
  RCP<const Epetra_BlockMap> g_map = sg_g.coefficientMap();
  Epetra_Vector g(*g_map);

  // Get quadrature data
  const Teuchos::Array<double>& norms = basis->norm_squared();
  const Teuchos::Array< Teuchos::Array<double> >& points = 
    quad->getQuadPoints();
  const Teuchos::Array<double>& weights = quad->getQuadWeights();
  const Teuchos::Array< Teuchos::Array<double> >& vals = 
    quad->getBasisAtQuadPoints();
  int nqp = points.size();

  // Compute sg_g via quadrature
  sg_g.init(0.0);
  SGConverter c(this, commT);
  for (int qp=0; qp<nqp; qp++) {

    // Evaluate sg_x, sg_xdot at quadrature point
    sg_x.evaluate(vals[qp], x);
    if (sg_xdot != NULL)
      sg_xdot->evaluate(vals[qp], *xdot);
    if (sg_xdotdot != NULL)
      sg_xdotdot->evaluate(vals[qp], *xdotdot);

    // Evaluate parameters at quadrature point
    for (int i=0; i<sg_p_index.size(); i++) {
      int ii = sg_p_index[i];
      for (unsigned int j=0; j<pp[ii].size(); j++)
	pp[ii][j].baseValue = sg_p_vals[ii][j].evaluate(points[qp], vals[qp]);
    }

    // Compute response at quadrature point
    c.evaluateResponse(curr_time, xdot.get(), xdotdot.get(), x, pp, g);

    // Add result into integral
    sg_g.sumIntoAllTerms(weights[qp], vals[qp], norms, g);
  }
}
void
Albany::SamplingBasedScalarResponseFunction::
evaluateSGGradient(
  const double current_time,
  const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
  const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
  const Stokhos::EpetraVectorOrthogPoly& sg_x,
  const Teuchos::Array<ParamVec>& p,
  const Teuchos::Array<int>& sg_p_index,
  const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
  ParamVec* deriv_p,
  Stokhos::EpetraVectorOrthogPoly* sg_g,
  Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dx,
  Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdot,
  Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdotdot,
  Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dp)
{
  RCP<const Epetra_BlockMap> x_map = sg_x.coefficientMap();
  RCP<Epetra_Vector> xdot;
  if (sg_xdot != NULL)
    xdot = rcp(new Epetra_Vector(*x_map));
  RCP<Epetra_Vector> xdotdot;
  if (sg_xdotdot != NULL)
    xdotdot = rcp(new Epetra_Vector(*x_map));
  Epetra_Vector x(*x_map);
  Teuchos::Array<ParamVec> pp = p;

  RCP<Epetra_Vector> g;
  if (sg_g != NULL) {
    sg_g->init(0.0);
    g = rcp(new Epetra_Vector(*(sg_g->coefficientMap())));
  }

  RCP<Epetra_MultiVector> dg_dx;
  if (sg_dg_dx != NULL) {
    sg_dg_dx->init(0.0);
    dg_dx = rcp(new Epetra_MultiVector(*(sg_dg_dx->coefficientMap()), 
				       sg_dg_dx->numVectors()));
  }

  RCP<Epetra_MultiVector> dg_dxdot;
  if (sg_dg_dxdot != NULL) {
    sg_dg_dxdot->init(0.0);
    dg_dxdot = rcp(new Epetra_MultiVector(*(sg_dg_dxdot->coefficientMap()), 
					  sg_dg_dxdot->numVectors()));
  }

  RCP<Epetra_MultiVector> dg_dxdotdot;
  if (sg_dg_dxdotdot != NULL) {
    sg_dg_dxdotdot->init(0.0);
    dg_dxdotdot = rcp(new Epetra_MultiVector(*(sg_dg_dxdotdot->coefficientMap()), 
					  sg_dg_dxdotdot->numVectors()));
  }

  RCP<Epetra_MultiVector> dg_dp;
  if (sg_dg_dp != NULL) {
    sg_dg_dp->init(0.0);
    dg_dp = rcp(new Epetra_MultiVector(*(sg_dg_dp->coefficientMap()), 
				       sg_dg_dp->numVectors()));
  }

  // Get quadrature data
  const Teuchos::Array<double>& norms = basis->norm_squared();
  const Teuchos::Array< Teuchos::Array<double> >& points = 
    quad->getQuadPoints();
  const Teuchos::Array<double>& weights = quad->getQuadWeights();
  const Teuchos::Array< Teuchos::Array<double> >& vals = 
    quad->getBasisAtQuadPoints();
  int nqp = points.size();

  // Compute sg_g via quadrature
  for (int qp=0; qp<nqp; qp++) {

    // Evaluate sg_x, sg_xdot at quadrature point
    sg_x.evaluate(vals[qp], x);
    if (sg_xdot != NULL)
      sg_xdot->evaluate(vals[qp], *xdot);
    if (sg_xdotdot != NULL)
      sg_xdotdot->evaluate(vals[qp], *xdotdot);

    // Evaluate parameters at quadrature point
    for (int i=0; i<sg_p_index.size(); i++) {
      int ii = sg_p_index[i];
      for (unsigned int j=0; j<pp[ii].size(); j++) {
	pp[ii][j].baseValue = sg_p_vals[ii][j].evaluate(points[qp], vals[qp]);
	if (deriv_p != NULL) {
	  for (unsigned int k=0; k<deriv_p->size(); k++)
	    if ((*deriv_p)[k].family->getName() == pp[ii][j].family->getName())
	      (*deriv_p)[k].baseValue = pp[ii][j].baseValue;
	}
      }
    }

    // Compute response at quadrature point
    evaluateGradient(current_time, xdot.get(), xdotdot.get(), x, pp, deriv_p,
		     g.get(), dg_dx.get(), dg_dxdot.get(), dg_dxdotdot.get(), dg_dp.get());

    // Add result into integral
    if (sg_g != NULL)
      sg_g->sumIntoAllTerms(weights[qp], vals[qp], norms, *g);
    if (sg_dg_dx != NULL)
      sg_dg_dx->sumIntoAllTerms(weights[qp], vals[qp], norms, *dg_dx);
    if (sg_dg_dxdot != NULL)
      sg_dg_dxdot->sumIntoAllTerms(weights[qp], vals[qp], norms, *dg_dxdot);
    if (sg_dg_dxdotdot != NULL)
      sg_dg_dxdotdot->sumIntoAllTerms(weights[qp], vals[qp], norms, *dg_dxdotdot);
    if (sg_dg_dp != NULL)
      sg_dg_dp->sumIntoAllTerms(weights[qp], vals[qp], norms, *dg_dp);
  }
}