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)
70 pcon_->update(x,flag,iter);
73 pcon_->update(x,type,iter);
82 ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
83 ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
88 catch (
const std::bad_cast &e) {}
91 std::vector<Real> param;
92 Real weight(0), one(1);
94 param =
sampler_->getMyPoint(static_cast<int>(i));
95 weight =
sampler_->getMyWeight(static_cast<int>(i));
96 pcon_->setParameter(param);
98 pcon_->value(*(pc.
get(i)), *(pu.
get(i)), *zptr, tol);
100 pc.
get(i)->scale(weight);
115 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
116 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
121 catch (
const std::bad_cast &e) {}
125 ROL::Ptr<const Vector<Real> > vuptr = vuz.
get_1();
126 ROL::Ptr<const Vector<Real> > vzptr = vuz.
get_2();
131 catch (
const std::bad_cast &e) {}
134 std::vector<Real> param;
135 Real weight(0), one(1);
137 param =
sampler_->getMyPoint(static_cast<int>(i));
138 weight =
sampler_->getMyWeight(static_cast<int>(i));
139 pcon_->setParameter(param);
143 pcon_->applyJacobian(*(pjv.
get(i)), vi, xi, tol);
145 pjv.
get(i)->scale(weight);
157 ROL::Ptr<Vector<Real> > ajvuptr = ajvuz.
get_1();
158 ROL::Ptr<Vector<Real> > ajvzptr = ajvuz.
get_2();
163 catch (
const std::bad_cast &e) {}
167 ROL::Ptr<const Vector<Real> > xuptr = xuz.
get_1();
168 ROL::Ptr<const Vector<Real> > xzptr = xuz.
get_2();
173 catch (
const std::bad_cast &e) {}
178 std::vector<Real> param;
179 Real weight(0), one(1);
180 ROL::Ptr<Vector<Real> > tmp1 = ajvzptr->clone();
181 ROL::Ptr<Vector<Real> > tmp2 = ajvzptr->clone();
183 param =
sampler_->getMyPoint(static_cast<int>(i));
184 weight =
sampler_->getMyWeight(static_cast<int>(i));
185 pcon_->setParameter(param);
189 pcon_->applyAdjointJacobian(ajvi, *(pv.
get(i)), xi, tol);
207 ROL::Ptr<Vector<Real> > ahuvuptr = ahuvuz.
get_1();
208 ROL::Ptr<Vector<Real> > ahuvzptr = ahuvuz.
get_2();
213 catch (
const std::bad_cast &e) {}
219 ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
220 ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
225 catch (
const std::bad_cast &e) {}
229 ROL::Ptr<const Vector<Real> > xuptr = xuz.
get_1();
230 ROL::Ptr<const Vector<Real> > xzptr = xuz.
get_2();
235 catch (
const std::bad_cast &e) {}
238 std::vector<Real> param;
239 Real weight(0), one(1);
240 ROL::Ptr<Vector<Real> > tmp1 = ahuvzptr->clone();
241 ROL::Ptr<Vector<Real> > tmp2 = ahuvzptr->clone();
243 param =
sampler_->getMyPoint(static_cast<int>(i));
244 weight =
sampler_->getMyWeight(static_cast<int>(i));
245 pcon_->setParameter(param);
250 pcon_->applyAdjointHessian(ahuvi, *(pu.
get(i)), vi, xi, tol);
269 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
270 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
275 catch (
const std::bad_cast &e) {}
279 ROL::Ptr<const Vector<Real> > guptr = guz.
get_1();
280 ROL::Ptr<const Vector<Real> > gzptr = guz.
get_2();
285 catch (
const std::bad_cast &e) {}
290 std::vector<Real> param;
291 Real weight(0), one(1);
293 param =
sampler_->getMyPoint(static_cast<int>(i));
294 weight =
sampler_->getMyWeight(static_cast<int>(i));
295 pcon_->setParameter(param);
299 pcon_->applyPreconditioner(*(ppv.
get(i)), *(pv.
get(i)), xi, gi, tol);
301 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.
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 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 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< const Vector< Real > > getVector(void) const
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