10 #ifndef ROL_STOCHASTIC_CONSTRAINT_H
11 #define ROL_STOCHASTIC_CONSTRAINT_H
22 Ptr<StochasticObjective<Real>>
robj_;
23 Ptr<Constraint<Real>>
con_;
29 ParameterList &parlist,
32 robj_ = makePtr<StochasticObjective<Real>>(obj,parlist,sampler,1,index);
33 con_ = makePtr<ConstraintFromObjective<Real>>(
robj_);
38 ParameterList &parlist,
42 Ptr<ConstraintFromObjective<Real>> cfo
43 = dynamicPtrCast<ConstraintFromObjective<Real>>(con);
44 robj_ = makePtr<StochasticObjective<Real>>(cfo->getObjective(),
45 parlist,sampler,1,index);
46 con_ = makePtr<ConstraintFromObjective<Real>>(
robj_);
48 catch (std::exception &e) {
54 return robj_->computeStatistic(x);
62 con_->update(x,type,iter);
66 con_->update(x,flag,iter);
74 con_->applyJacobian(jv,v,x,tol);
78 con_->applyAdjointJacobian(ajv,v,x,tol);
82 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< Objective< Real >> &obj, const Ptr< SampleGenerator< Real >> &sampler, ParameterList &parlist, const int index=0)
Contains definitions of custom data types in ROL.
StochasticConstraint(const Ptr< Constraint< Real >> &con, const Ptr< SampleGenerator< Real >> &sampler, ParameterList &parlist, const int index=0)
Defines the linear algebra or vector space interface.
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.
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 .
Ptr< StochasticObjective< Real > > robj_
Real computeStatistic(const Vector< Real > &x) const
Defines the general constraint operator interface.