ROL
ROL_AlmostSureConstraint.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_ALMOST_SURE_CONSTRAINT_H
11#define ROL_ALMOST_SURE_CONSTRAINT_H
12
13#include "ROL_RiskVector.hpp"
14#include "ROL_Constraint.hpp"
17
18namespace ROL {
19
20template <class Real>
21class AlmostSureConstraint : public Constraint<Real> {
22private:
23 const Ptr<SampleGenerator<Real>> sampler_;
24 const Ptr<Constraint<Real>> con_;
25
26 Ptr<Vector<Real>> scratch1_;
27 Ptr<Vector<Real>> scratch2_;
29
30public:
32
34 const Ptr<Constraint<Real>> &con)
35 : sampler_(sampler), con_(con), isInitialized_(false) {}
36
37 void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
38 con_->update(x,flag,iter);
39 }
40 void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
41 con_->update(x,type,iter);
42 }
43
45 const Vector<Real> &x,
46 Real &tol) {
47 c.zero();
48 SimulatedVector<Real> &pc = dynamic_cast<SimulatedVector<Real>&>(c);
49
50 std::vector<Real> param;
51 for (int i = 0; i < sampler_->numMySamples(); ++i) {
52 param = sampler_->getMyPoint(i);
53 con_->setParameter(param);
54 con_->update(x);
55 con_->value(*(pc.get(i)), x, tol);
56 }
57 }
58
60 const Vector<Real> &v,
61 const Vector<Real> &x,
62 Real &tol) {
63 jv.zero();
64 SimulatedVector<Real> &pjv = dynamic_cast<SimulatedVector<Real>&>(jv);
65
66 std::vector<Real> param;
67 for (int i = 0; i < sampler_->numMySamples(); ++i) {
68 param = sampler_->getMyPoint(i);
69 con_->setParameter(param);
70 con_->update(x);
71 con_->applyJacobian(*(pjv.get(i)), v, x, tol);
72 }
73 }
74
76 const Vector<Real> &v,
77 const Vector<Real> &x,
78 Real &tol) {
79 ajv.zero();
80 if (!isInitialized_) {
81 scratch1_ = ajv.clone();
82 scratch2_ = ajv.clone();
83 isInitialized_ = true;
84 }
85 const SimulatedVector<Real> &pv = dynamic_cast<const SimulatedVector<Real>&>(v);
86
87 std::vector<Real> param;
88 scratch1_->zero(); scratch2_->zero();
89 for (int i = 0; i < sampler_->numMySamples(); ++i) {
90 param = sampler_->getMyPoint(i);
91 con_->setParameter(param);
92 con_->update(x);
93 con_->applyAdjointJacobian(*scratch1_, *(pv.get(i)), x, tol);
95 }
96 sampler_->sumAll(*scratch2_, ajv);
97 }
98
100 const Vector<Real> &u,
101 const Vector<Real> &v,
102 const Vector<Real> &x,
103 Real &tol) {
104 ahuv.zero();
105 if (!isInitialized_) {
106 scratch1_ = ahuv.clone();
107 scratch2_ = ahuv.clone();
108 isInitialized_ = true;
109 }
110 const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(u);
111
112 std::vector<Real> param;
113 scratch1_->zero(); scratch2_->zero();
114 for (int i = 0; i < sampler_->numMySamples(); ++i) {
115 param = sampler_->getMyPoint(i);
116 con_->setParameter(param);
117 con_->update(x);
118 con_->applyAdjointHessian(*scratch1_, *(pu.get(i)), v, x, tol);
120 }
121 sampler_->sumAll(*scratch2_, ahuv);
122 }
123
125 const Vector<Real> &v,
126 const Vector<Real> &x,
127 const Vector<Real> &g,
128 Real &tol) {
129 Pv.zero();
130 SimulatedVector<Real> &ppv = dynamic_cast<SimulatedVector<Real>&>(Pv);
131 const SimulatedVector<Real> &pv = dynamic_cast<const SimulatedVector<Real>&>(v);
132
133 std::vector<Real> param;
134 scratch1_->zero(); scratch2_->zero();
135 for (int i = 0; i < sampler_->numMySamples(); ++i) {
136 param = sampler_->getMyPoint(i);
137 con_->setParameter(param);
138 con_->update(x);
139 con_->applyPreconditioner(*(ppv.get(i)), *(pv.get(i)), x, g, tol);
140 }
141 }
142
143}; // class AlmostSureConstraint
144
145} // namespace ROL
146
147#endif
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...
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 ,...
AlmostSureConstraint(const Ptr< SampleGenerator< Real > > &sampler, const Ptr< Constraint< Real > > &con)
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
void applyPreconditioner(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
Apply a constraint preconditioner at , , to vector . Ideally, this preconditioner satisfies the follo...
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
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 .
const Ptr< SampleGenerator< Real > > sampler_
const Ptr< Constraint< Real > > con_
Defines the general constraint operator interface.
Defines the linear algebra of a vector space on a generic partitioned vector where the individual vec...
void zero()
Set to zero vector.
ROL::Ptr< const Vector< Real > > get(size_type i) const
void plus(const V &x)
Compute , where .
Defines the linear algebra or vector space interface.
virtual void zero()
Set to zero vector.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.