ROL
ROL_StochasticConstraint.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_STOCHASTIC_CONSTRAINT_H
11#define ROL_STOCHASTIC_CONSTRAINT_H
12
15#include "ROL_Types.hpp"
16
17namespace ROL {
18
19template <class Real>
20class StochasticConstraint : public Constraint<Real> {
21private:
22 Ptr<StochasticObjective<Real>> robj_;
23 Ptr<Constraint<Real>> con_;
24 Ptr<SampleGenerator<Real>> sampler_;
25
26public:
28 const Ptr<SampleGenerator<Real>> &sampler,
29 ParameterList &parlist,
30 const int index = 0)
31 : sampler_(sampler) {
32 robj_ = makePtr<StochasticObjective<Real>>(obj,parlist,sampler,1,index);
33 con_ = makePtr<ConstraintFromObjective<Real>>(robj_);
34 }
35
37 const Ptr<SampleGenerator<Real>> &sampler,
38 ParameterList &parlist,
39 const int index = 0)
40 : sampler_(sampler) {
41 try {
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_);
47 }
48 catch (std::exception &e) {
49 throw Exception::NotImplemented(">>> ROL::StochasticConstraint: Input constraint must be a ConstraintFromObjective!");
50 }
51 }
52
53 Real computeStatistic(const Vector<Real> &x) const {
54 return robj_->computeStatistic(x);
55 }
56
57 void setIndex(int ind) {
58 robj_->setIndex(ind);
59 }
60
61 void update(const Vector<Real> &x, UpdateType type, int iter = -1) {
62 con_->update(x,type,iter);
63 }
64
65 void update(const Vector<Real> &x, bool flag = true, int iter = -1) {
66 con_->update(x,flag,iter);
67 }
68
69 void value(Vector<Real> &c, const Vector<Real> &x, Real &tol) {
70 con_->value(c,x,tol);
71 }
72
73 void applyJacobian(Vector<Real> &jv, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
74 con_->applyJacobian(jv,v,x,tol);
75 }
76
77 void applyAdjointJacobian(Vector<Real> &ajv, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
78 con_->applyAdjointJacobian(ajv,v,x,tol);
79 }
80
81 void applyAdjointHessian(Vector<Real> &ahuv, const Vector<Real> &u, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
82 con_->applyAdjointHessian(ahuv,u,v,x,tol);
83 }
84
85}; // class StochasticConstraint
86
87} // namespace ROL
88
89#endif
Contains definitions of custom data types in ROL.
Defines the general constraint operator interface.
Provides the interface to evaluate objective functions.
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 ,...
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
Real computeStatistic(const Vector< Real > &x) const
Ptr< StochasticObjective< Real > > robj_
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 i...
Ptr< SampleGenerator< Real > > sampler_
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, ParameterList &parlist, const int index=0)
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 .
StochasticConstraint(const Ptr< Objective< Real > > &obj, const Ptr< SampleGenerator< Real > > &sampler, ParameterList &parlist, const int index=0)
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.
Defines the linear algebra or vector space interface.