44 #ifndef ROL_STOCHASTIC_CONSTRAINT_H 
   45 #define ROL_STOCHASTIC_CONSTRAINT_H 
   56   Ptr<StochasticObjective<Real>> 
robj_;
 
   57   Ptr<Constraint<Real>>          
con_;
 
   63                ROL::ParameterList             &parlist,
 
   66     robj_ = makePtr<StochasticObjective<Real>>(obj,parlist,sampler,1,index);
 
   67     con_  = makePtr<ConstraintFromObjective<Real>>(
robj_);
 
   72                ROL::ParameterList              &parlist,
 
   76       Ptr<ConstraintFromObjective<Real>> cfo
 
   77         = dynamicPtrCast<ConstraintFromObjective<Real>>(con);
 
   78       robj_ = makePtr<StochasticObjective<Real>>(cfo->getObjective(),
 
   79                 parlist,sampler,1,index);
 
   80       con_  = makePtr<ConstraintFromObjective<Real>>(
robj_);
 
   82     catch (std::exception &e) {
 
   88     return robj_->computeStatistic(x);
 
   92     con_->update(x,flag,iter);
 
  100     con_->applyJacobian(jv,v,x,tol);
 
  104     con_->applyAdjointJacobian(ajv,v,x,tol);
 
  108     con_->applyAdjointHessian(ahuv,u,v,x,tol);
 
Provides the interface to evaluate objective functions. 
 
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector . 
 
StochasticConstraint(const Ptr< Constraint< Real >> &con, const Ptr< SampleGenerator< Real >> &sampler, ROL::ParameterList &parlist, const int index=0)
 
Contains definitions of custom data types in ROL. 
 
Defines the linear algebra or vector space interface. 
 
void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the derivative of the adjoint of the constraint Jacobian at  to vector  in direction ...
 
Ptr< Constraint< Real > > con_
 
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
 
Ptr< SampleGenerator< Real > > sampler_
 
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector . 
 
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator  at . 
 
StochasticConstraint(const Ptr< Objective< Real >> &obj, const Ptr< SampleGenerator< Real >> &sampler, ROL::ParameterList &parlist, const int index=0)
 
Ptr< StochasticObjective< Real > > robj_
 
Real computeStatistic(const Vector< Real > &x) const 
 
Defines the general constraint operator interface.