44 #ifndef ROL_INTERIORPOINTPENALTY_H
45 #define ROL_INTERIORPOINTPENALTY_H
49 #include "ROL_ParameterList.hpp"
75 const ROL::Ptr<const V>
lo_;
76 const ROL::Ptr<const V>
up_;
105 return (x>0) ? std::log(x) : Real(0.0);
114 return (x>0) ? 1.0/x : Real(0.0);
123 Real
apply(
const Real &x,
const Real &y )
const {
124 return (x>0) ? y/x : Real(0.0);
132 class Mask :
public Elementwise::BinaryFunction<Real> {
137 Real
apply(
const Real &x,
const Real &y )
const {
149 ROL::ParameterList &parlist ) :
150 obj_(obj),
bnd_(con),
lo_( con->getLowerBound() ),
up_( con->getUpperBound() ) {
156 ValueSet isBoundedBelow( ROL_NINF<Real>(), ValueSet::GREATER_THAN, one, zero );
157 ValueSet isBoundedAbove( ROL_INF<Real>(), ValueSet::LESS_THAN, one, zero );
162 maskL_->applyBinary(isBoundedBelow,*
lo_);
163 maskU_->applyBinary(isBoundedAbove,*
up_);
165 ROL::ParameterList &iplist = parlist.sublist(
"Step").sublist(
"Primal Dual Interior Point");
166 ROL::ParameterList &lblist = iplist.sublist(
"Barrier Objective");
169 kappaD_ = lblist.get(
"Linear Damping Coefficient",1.e-4);
170 mu_ = lblist.get(
"Initial Barrier Parameter",0.1);
175 g_ =
lo_->dual().clone();
214 obj_->update(x,flag,iter);
215 bnd_->update(x,flag,iter);
232 Elementwise::ReductionSum<Real> sum;
233 Elementwise::Multiply<Real> mult;
240 Real linearTerm = 0.0;
251 c_->applyBinary(mult,*
a_);
254 linearTerm +=
c_->reduce(sum);
257 a_->applyUnary(mlog);
270 c_->applyBinary(mult,*
b_);
274 linearTerm +=
c_->reduce(sum);
278 b_->applyUnary(mlog);
283 fval -=
mu_*(aval+bval);
300 obj_->gradient(*
g_,x,tol);
310 a_->applyUnary(mrec);
315 b_->applyUnary(mrec);
347 Elementwise::Multiply<Real> mult;
348 Elementwise::Power<Real> square(2.0);
350 obj_->hessVec(hv,v,x,tol);
354 a_->applyUnary(mrec);
356 a_->applyUnary(square);
358 a_->applyBinary(mult,v);
362 b_->applyUnary(mrec);
364 b_->applyUnary(square);
366 b_->applyBinary(mult,v);
383 #endif // ROL_INTERIORPOINTPENALTY_H
Provides the interface to evaluate objective functions.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update the interior point penalized objective.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
ROL::Ptr< const Vector< Real > > getLowerMask(void) const
ROL::Ptr< const Vector< Real > > getUpperMask(void) const
InteriorPointPenalty(const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< BoundConstraint< Real > > &con, ROL::ParameterList &parlist)
Defines the linear algebra or vector space interface.
Real apply(const Real &x, const Real &y) const
Real apply(const Real &x) const
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
const ROL::Ptr< OBJ > getObjective(void)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Compute action of Hessian on vector.
Provides the interface to evaluate the Interior Pointy log barrier penalty function with upper and lo...
const ROL::Ptr< BND > getBoundConstraint(void)
BoundConstraint< Real > BND
ROL::Ptr< Vector< Real > > getGradient(void) const
Real apply(const Real &x, const Real &y) const
const ROL::Ptr< const V > lo_
const ROL::Ptr< const V > up_
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
const ROL::Ptr< BND > bnd_
Provides the interface to apply upper and lower bound constraints.
int getNumberFunctionEvaluations(void) const
int getNumberGradientEvaluations(void) const
const ROL::Ptr< OBJ > obj_
virtual void set(const Vector &x)
Set where .
Elementwise::ValueSet< Real > ValueSet
Real getObjectiveValue(void) const
Real apply(const Real &x) const