44 #ifndef ROL_AUGMENTEDLAGRANGIAN_H
45 #define ROL_AUGMENTEDLAGRANGIAN_H
52 #include "ROL_Ptr.hpp"
89 const ROL::Ptr<Objective<Real> >
obj_;
90 ROL::Ptr<QuadraticPenalty<Real> >
pen_;
129 const Real penaltyParameter,
132 ROL::ParameterList &parlist)
141 ROL::ParameterList& sublist = parlist.sublist(
"Step").sublist(
"Augmented Lagrangian");
143 int HessianApprox = sublist.get(
"Level of Hessian Approximation", 0);
145 pen_ = ROL::makePtr<QuadraticPenalty<Real>>(con,multiplier,penaltyParameter,optVec,conVec,
scaleLagrangian_,HessianApprox);
160 obj_->update(x,flag,iter);
161 pen_->update(x,flag,iter);
166 void setScaling(
const Real fscale,
const Real cscale = 1.0) {
168 pen_->setScaling(cscale);
178 Real pval =
pen_->value(x,tol);
206 obj_->hessVec(hv,v,x,tol);
219 Real tol = std::sqrt(ROL_EPSILON<Real>());
229 Real tol = std::sqrt(ROL_EPSILON<Real>());
240 pen_->getConstraintVec(c,x);
245 return pen_->getNumberConstraintEvaluations();
261 pen_->reset(multiplier,penaltyParameter);
Provides the interface to evaluate objective functions.
Provides the interface to evaluate the augmented Lagrangian.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
ROL::Ptr< Vector< Real > > dualOptVector_
virtual void scale(const Real alpha)=0
Compute where .
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual void plus(const Vector &x)=0
Compute , where .
virtual int getNumberConstraintEvaluations(void) const
Contains definitions of custom data types in ROL.
virtual Real getObjectiveValue(const Vector< Real > &x)
virtual Real value(const Vector< Real > &x, Real &tol)
Compute value.
Defines the linear algebra or vector space interface.
virtual void reset(const Vector< Real > &multiplier, const Real penaltyParameter)
AugmentedLagrangian(const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< Constraint< Real > > &con, const Vector< Real > &multiplier, const Real penaltyParameter, const Vector< Real > &optVec, const Vector< Real > &conVec, ROL::ParameterList &parlist)
Constructor.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
AugmentedLagrangian()
Null constructor.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
ROL::Ptr< Vector< Real > > gradient_
virtual int getNumberGradientEvaluations(void) const
const ROL::Ptr< Objective< Real > > obj_
virtual int getNumberFunctionEvaluations(void) const
void setScaling(const Real fscale, const Real cscale=1.0)
const Ptr< const Vector< Real > > getObjectiveGradient(const Vector< Real > &x)
virtual void set(const Vector &x)
Set where .
ROL::Ptr< QuadraticPenalty< Real > > pen_
virtual void getConstraintVec(Vector< Real > &c, const Vector< Real > &x)
Defines the general constraint operator interface.