44 #ifndef ROL_SIMULATED_CONSTRAINT_H
45 #define ROL_SIMULATED_CONSTRAINT_H
56 const ROL::Ptr<SampleGenerator<Real> >
sampler_;
57 const ROL::Ptr<Constraint_SimOpt<Real> >
pcon_;
66 const bool useWeights =
true)
77 ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
78 ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
83 catch (
const std::bad_cast &e) {}
86 std::vector<Real> param;
87 Real weight(0), one(1);
89 param =
sampler_->getMyPoint(static_cast<int>(i));
90 weight =
sampler_->getMyWeight(static_cast<int>(i));
91 pcon_->setParameter(param);
93 pcon_->value(*(pc.
get(i)), *(pu.
get(i)), *zptr, tol);
95 pc.
get(i)->scale(weight);
110 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
111 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
116 catch (
const std::bad_cast &e) {}
120 ROL::Ptr<const Vector<Real> > vuptr = vuz.
get_1();
121 ROL::Ptr<const Vector<Real> > vzptr = vuz.
get_2();
126 catch (
const std::bad_cast &e) {}
129 std::vector<Real> param;
130 Real weight(0), one(1);
132 param =
sampler_->getMyPoint(static_cast<int>(i));
133 weight =
sampler_->getMyWeight(static_cast<int>(i));
134 pcon_->setParameter(param);
138 pcon_->applyJacobian(*(pjv.
get(i)), vi, xi, tol);
140 pjv.
get(i)->scale(weight);
152 ROL::Ptr<Vector<Real> > ajvuptr = ajvuz.
get_1();
153 ROL::Ptr<Vector<Real> > ajvzptr = ajvuz.
get_2();
158 catch (
const std::bad_cast &e) {}
162 ROL::Ptr<const Vector<Real> > xuptr = xuz.
get_1();
163 ROL::Ptr<const Vector<Real> > xzptr = xuz.
get_2();
168 catch (
const std::bad_cast &e) {}
173 std::vector<Real> param;
174 Real weight(0), one(1);
175 ROL::Ptr<Vector<Real> > tmp1 = ajvzptr->clone();
176 ROL::Ptr<Vector<Real> > tmp2 = ajvzptr->clone();
178 param =
sampler_->getMyPoint(static_cast<int>(i));
179 weight =
sampler_->getMyWeight(static_cast<int>(i));
180 pcon_->setParameter(param);
184 pcon_->applyAdjointJacobian(ajvi, *(pv.
get(i)), xi, tol);
202 ROL::Ptr<Vector<Real> > ahuvuptr = ahuvuz.
get_1();
203 ROL::Ptr<Vector<Real> > ahuvzptr = ahuvuz.
get_2();
208 catch (
const std::bad_cast &e) {}
214 ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
215 ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
220 catch (
const std::bad_cast &e) {}
224 ROL::Ptr<const Vector<Real> > xuptr = xuz.
get_1();
225 ROL::Ptr<const Vector<Real> > xzptr = xuz.
get_2();
230 catch (
const std::bad_cast &e) {}
233 std::vector<Real> param;
234 Real weight(0), one(1);
235 ROL::Ptr<Vector<Real> > tmp1 = ahuvzptr->clone();
236 ROL::Ptr<Vector<Real> > tmp2 = ahuvzptr->clone();
238 param =
sampler_->getMyPoint(static_cast<int>(i));
239 weight =
sampler_->getMyWeight(static_cast<int>(i));
240 pcon_->setParameter(param);
245 pcon_->applyAdjointHessian(ahuvi, *(pu.
get(i)), vi, xi, tol);
264 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
265 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
270 catch (
const std::bad_cast &e) {}
274 ROL::Ptr<const Vector<Real> > guptr = guz.
get_1();
275 ROL::Ptr<const Vector<Real> > gzptr = guz.
get_2();
280 catch (
const std::bad_cast &e) {}
285 std::vector<Real> param;
286 Real weight(0), one(1);
288 param =
sampler_->getMyPoint(static_cast<int>(i));
289 weight =
sampler_->getMyWeight(static_cast<int>(i));
290 pcon_->setParameter(param);
294 pcon_->applyPreconditioner(*(ppv.
get(i)), *(pv.
get(i)), xi, gi, tol);
296 ppv.
get(i)->scale(one/(weight*weight));
typename PV< Real >::size_type size_type
virtual void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
ROL::Ptr< const Vector< Real > > get_2() const
Defines the linear algebra or vector space interface for simulation-based optimization.
ROL::Ptr< const Vector< Real > > getVector(void) const
virtual 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 ...
const ROL::Ptr< Constraint_SimOpt< Real > > pcon_
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
const ROL::Ptr< SampleGenerator< Real > > sampler_
Defines the linear algebra of a vector space on a generic partitioned vector where the individual vec...
ROL::Ptr< const Vector< Real > > get(size_type i) const
virtual 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...
size_type numVectors() const
virtual 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 .
SimulatedConstraint(const ROL::Ptr< SampleGenerator< Real > > &sampler, const ROL::Ptr< Constraint_SimOpt< Real > > &pcon, const bool useWeights=true)
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
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 ...
virtual ~SimulatedConstraint()
Defines the constraint operator interface for simulation-based optimization.
Defines the general constraint operator interface.
ROL::Ptr< const Vector< Real > > get_1() const