41 #ifndef ROL_INTERIORPOINT_H 42 #define ROL_INTERIORPOINT_H 64 const static size_type
OPT = 0;
65 const static size_type
SLACK = 1;
67 Teuchos::RCP<Objective<Real> >
obj_;
76 using Teuchos::dyn_cast;
77 return dyn_cast<
const PV>(x);
82 using Teuchos::dyn_cast;
83 return dyn_cast<PV>(x);
90 obj_(obj),barrier_(barrier),mu_(mu),nfval_(0),ngval_(0) {
108 Teuchos::RCP<const V> xopt = xpv.
get(OPT);
109 Teuchos::RCP<const V> s = xpv.
get(SLACK);
112 return obj_->value(*xopt,tol) + mu_*barrier_->value(*s,tol);
121 Teuchos::RCP<const V> xopt = xpv.
get(OPT);
122 Teuchos::RCP<const V> s = xpv.
get(SLACK);
124 Teuchos::RCP<V> gopt = gpv.
get(OPT);
125 Teuchos::RCP<V> gs = gpv.
get(SLACK);
127 obj_->gradient(*gopt, *xopt, tol);
128 barrier_->gradient(*gs, *s, tol);
141 Teuchos::RCP<const V> xopt = xpv.
get(OPT);
142 Teuchos::RCP<const V> s = xpv.
get(SLACK);
144 Teuchos::RCP<const V> vopt = vpv.
get(OPT);
145 Teuchos::RCP<const V> vs = vpv.
get(SLACK);
147 Teuchos::RCP<V> hvopt = hvpv.
get(OPT);
148 Teuchos::RCP<V> hvs = hvpv.
get(SLACK);
150 obj_->hessVec(*hvopt, *vopt, *xopt, tol);
151 barrier_->hessVec(*hvs, *vs, *s, tol);
173 const static size_type
OPT = 0;
176 const static size_type INEQ = 0;
177 const static size_type EQUAL = 1;
179 Teuchos::RCP<EqualityConstraint<Real> >
incon_;
180 Teuchos::RCP<EqualityConstraint<Real> >
eqcon_;
187 using Teuchos::dyn_cast;
188 return dyn_cast<
const PV>(x);
193 using Teuchos::dyn_cast;
194 return dyn_cast<PV>(x);
202 incon_(incon), eqcon_(eqcon),hasEquality_(true), ncval_(0) {}
206 incon_(incon), hasEquality_(false), ncval_(0) {}
221 RCP<const V> xopt = xpv.
get(OPT);
222 RCP<const V> s = xpv.
get(SLACK);
226 RCP<V> ci = cpv.
get(INEQ);
227 incon_->value(*ci, *xopt, tol);
231 RCP<V> ce = cpv.
get(EQUAL);
232 eqcon_->value(*ce, *xopt, tol);
249 RCP<const V> xopt = xpv.
get(OPT);
250 RCP<const V> s = xpv.
get(SLACK);
252 RCP<const V> vopt = vpv.
get(OPT);
253 RCP<const V> vs = vpv.
get(SLACK);
257 RCP<V> jvi = jvpv.
get(INEQ);
258 incon_->applyJacobian(*jvi, *vopt, *xopt, tol);
262 RCP<V> jve = jvpv.
get(EQUAL);
263 eqcon_->applyJacobian(*jve, *vopt, *xopt, tol);
279 RCP<const V> xopt = xpv.
get(OPT);
280 RCP<const V> s = xpv.
get(SLACK);
282 RCP<V> ajvopt = ajvpv.
get(OPT);
283 RCP<V> ajvs = ajvpv.
get(SLACK);
287 RCP<const V> vi = vpv.
get(INEQ);
289 incon_->applyAdjointJacobian(*ajvopt,*vi,*xopt,tol);
297 RCP<const V> ve = vpv.
get(EQUAL);
298 RCP<V> temp = ajvopt->clone();
299 eqcon_->applyAdjointJacobian(*temp,*ve,*xopt,tol);
318 RCP<const V> xopt = xpv.
get(OPT);
319 RCP<const V> s = xpv.
get(SLACK);
321 RCP<const V> vopt = vpv.
get(OPT);
323 RCP<V> ahuvopt = ahuvpv.
get(OPT);
324 RCP<V> ahuvs = ahuvpv.
get(SLACK);
326 RCP<V> temp = ahuvopt->clone();
330 RCP<const V> ui = upv.
get(INEQ);
332 incon_->applyAdjointHessian(*ahuvopt,*ui,*vopt,*xopt,tol);
336 RCP<const V> ue = upv.
get(EQUAL);
337 eqcon_->applyAdjointHessian(*temp,*ue,*vopt,*xopt,tol);
338 ahuvopt->plus(*temp);
358 const static size_type
OPT = 0;
361 Teuchos::RCP<BoundConstraint<Real> >
bc_;
373 PV lowerpv = Teuchos::dyn_cast<PV>(lower_);
374 PV upperpv = Teuchos::dyn_cast<PV>(upper_);
376 Teuchos::RCP<V> lower_opt = lowerpv.
get(OPT);
377 Teuchos::RCP<V> lower_slack = lowerpv.
get(SLACK);
379 Teuchos::RCP<V> upper_opt = upperpv.get(OPT);
380 Teuchos::RCP<V> upper_slack = upperpv.get(SLACK);
382 Elementwise::Fill<Real> setToMax(std::numeric_limits<Real>::max());
383 Elementwise::Fill<Real> setToMin(std::numeric_limits<Real>::lowest());
385 lower_opt->applyUnary( setToMin );
388 upper_opt->applyUnary( setToMax );
389 upper_slack->applyUnary( setToMax );
400 bc_->pruneUpperActive( v, x, eps );
405 bc_->pruneUpperActive( v, g, x, eps );
409 bc_->pruneLowerActive( v, x, eps );
414 bc_->pruneLowerActive( v, g, x, eps );
426 return bc_->isFeasible(v);
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 .
Provides the interface to evaluate objective functions.
PartitionedVector< Real > PV
Teuchos::RCP< EqualityConstraint< Real > > eqcon_
int getNumberFunctionEvaluations(void)
Has both inequality and equality constraints. Treat inequality constraint as equality with slack vari...
Defines the linear algebra of vector space on a generic partitioned vector.
int getNumberGradientEvaluations(void)
Teuchos::RCP< BoundConstraint< Real > > bc_
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
void pruneUpperActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -binding set.
InteriorPointEqualityConstraint(const Teuchos::RCP< EqualityConstraint< Real > > &incon)
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Teuchos::RCP< Vector< Real > > upper_
InteriorPointEqualityConstraint(const Teuchos::RCP< EqualityConstraint< Real > > &incon, const Teuchos::RCP< EqualityConstraint< Real > > &eqcon)
Require positivity of slack variables.
Defines the linear algebra or vector space interface.
const PV & partition(const V &x)
void pruneLowerActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -binding set.
Teuchos::RCP< EqualityConstraint< Real > > incon_
InteriorPointBoundConstraint(const Vector< Real > &x)
void updatePenalty(Real mu)
static const size_type SLACK
Defines the equality constraint operator interface.
PartitionedVector< Real > PV
void pruneUpperActive(Vector< Real > &v, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -active set.
Teuchos::RCP< Vector< Real > > lower_
Teuchos::RCP< const Vector< Real > > get(size_type i) const
void pruneLowerActive(Vector< Real > &v, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -active set.
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 .
void project(Vector< Real > &x)
Project optimization variables onto the bounds.
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
Adds barrier term to generic objective.
Teuchos::RCP< Objective< Real > > barrier_
Provides the interface to apply upper and lower bound constraints.
bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
Teuchos::RCP< Objective< Real > > obj_
Real value(const Vector< Real > &x, Real &tol)
Compute value.
std::vector< PV >::size_type size_type
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual void set(const Vector &x)
Set where .
void setVectorToUpperBound(Vector< Real > &u)
Set the input vector to the upper bound.
PartitionedVector< Real > PV
int getNumberConstraintEvaluations(void)
InteriorPointObjective(const Teuchos::RCP< Objective< Real > > &obj, const Teuchos::RCP< Objective< Real > > &barrier, Real mu)
const PV & partition(const V &x)
void setVectorToLowerBound(Vector< Real > &l)
Set the input vector to the lower bound.
static const size_type OPT