44 #ifndef ROL_AUGMENTEDLAGRANGIAN_H 45 #define ROL_AUGMENTEDLAGRANGIAN_H 52 #include "Teuchos_RCP.hpp" 89 const Teuchos::RCP<Objective<Real> >
obj_;
90 Teuchos::RCP<QuadraticPenalty<Real> >
pen_;
126 const Real penaltyParameter,
129 Teuchos::ParameterList &parlist)
130 : obj_(obj), penaltyParameter_(penaltyParameter),
131 fval_(0), nfval_(0), ngval_(0), isValueComputed_(false), isGradientComputed_(false) {
133 gradient_ = optVec.
dual().clone();
134 dualOptVector_ = optVec.
dual().clone();
136 Teuchos::ParameterList& sublist = parlist.sublist(
"Step").sublist(
"Augmented Lagrangian");
137 scaleLagrangian_ = sublist.get(
"Use Scaled Augmented Lagrangian",
false);
138 int HessianApprox = sublist.get(
"Level of Hessian Approximation", 0);
140 pen_ = Teuchos::rcp(
new QuadraticPenalty<Real>(con,multiplier,penaltyParameter,optVec,conVec,scaleLagrangian_,HessianApprox));
149 : obj_(Teuchos::null), pen_(Teuchos::null), dualOptVector_(Teuchos::null),
150 fval_(0), gradient_(Teuchos::null), nfval_(0), ngval_(0),
151 scaleLagrangian_(false), isValueComputed_(false), isGradientComputed_(false) {}
154 obj_->update(x,flag,iter);
155 pen_->update(x,flag,iter);
162 if ( !isValueComputed_ ) {
163 fval_ = obj_->value(x,tol); nfval_++;
164 isValueComputed_ =
true;
167 Real pval = pen_->value(x,tol);
170 if (scaleLagrangian_) {
178 if ( !isGradientComputed_ ) {
179 obj_->gradient(*gradient_,x,tol); ngval_++;
180 isGradientComputed_ =
true;
184 pen_->gradient(*dualOptVector_,x,tol);
186 if ( scaleLagrangian_ ) {
187 g.
scale(static_cast<Real>(1)/penaltyParameter_);
189 g.
plus(*dualOptVector_);
194 obj_->hessVec(hv,v,x,tol);
196 pen_->hessVec(*dualOptVector_,v,x,tol);
198 if ( scaleLagrangian_ ) {
199 hv.
scale(static_cast<Real>(1)/penaltyParameter_);
201 hv.
plus(*dualOptVector_);
206 Real tol = std::sqrt(ROL_EPSILON<Real>());
208 if ( !isValueComputed_ ) {
209 fval_ = obj_->value(x,tol); nfval_++;
210 isValueComputed_ =
true;
217 pen_->getConstraintVec(c,x);
222 return pen_->getNumberConstraintEvaluations();
237 nfval_ = 0; ngval_ = 0;
238 pen_->reset(multiplier,penaltyParameter);
Provides the interface to evaluate objective functions.
Provides the interface to evaluate the augmented Lagrangian.
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 .
Provides the interface to evaluate the quadratic constraint penalty.
Contains definitions of custom data types in ROL.
virtual Real getObjectiveValue(const Vector< Real > &x)
Teuchos::RCP< Vector< Real > > gradient_
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)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Defines the equality constraint operator interface.
AugmentedLagrangian(const Teuchos::RCP< Objective< Real > > &obj, const Teuchos::RCP< EqualityConstraint< Real > > &con, const Vector< Real > &multiplier, const Real penaltyParameter, const Vector< Real > &optVec, const Vector< Real > &conVec, Teuchos::ParameterList &parlist)
Constructor.
virtual int getNumberFunctionEvaluations(void) const
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.
virtual int getNumberConstraintEvaluations(void) const
Teuchos::RCP< Vector< Real > > dualOptVector_
virtual int getNumberGradientEvaluations(void) const
virtual void set(const Vector &x)
Set where .
Teuchos::RCP< QuadraticPenalty< Real > > pen_
virtual void getConstraintVec(Vector< Real > &c, const Vector< Real > &x)
const Teuchos::RCP< Objective< Real > > obj_