ROL
ROL_Reduced_AugmentedLagrangian_SimOpt.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_REDUCED_AUGMENTEDLAGRANGIAN_SIMOPT_H
11#define ROL_REDUCED_AUGMENTEDLAGRANGIAN_SIMOPT_H
12
18#include "ROL_Vector.hpp"
19#include "ROL_Types.hpp"
20#include "ROL_Ptr.hpp"
21#include <iostream>
22
61namespace ROL {
62
63template <class Real>
65private:
66 ROL::Ptr<AugmentedLagrangian_SimOpt<Real> > augLagSimOpt_;
67 ROL::Ptr<Reduced_Objective_SimOpt<Real> > rAugLagSimOpt_;
68 ROL::Ptr<Vector<Real> > state_;
69
70 // Evaluation counters
71 int ngval_;
72
73public:
75 const ROL::Ptr<Constraint_SimOpt<Real> > &redCon,
76 const ROL::Ptr<Constraint_SimOpt<Real> > &augCon,
77 const ROL::Ptr<Vector<Real> > &state,
78 const ROL::Ptr<Vector<Real> > &control,
79 const ROL::Ptr<Vector<Real> > &adjoint,
80 const ROL::Ptr<Vector<Real> > &augConVec,
81 const ROL::Ptr<Vector<Real> > &multiplier,
82 const Real penaltyParameter,
83 ROL::ParameterList &parlist) : state_(state),
84 ngval_(0) {
85
86 augLagSimOpt_ = ROL::makePtr<AugmentedLagrangian_SimOpt<Real>>(obj,
87 augCon,
88 *multiplier,
89 penaltyParameter,
90 *state,
91 *control,
92 *augConVec,
93 parlist);
94 rAugLagSimOpt_ = ROL::makePtr<Reduced_Objective_SimOpt<Real>>(augLagSimOpt_,redCon,state,control,adjoint);
95 rAugLagSimOpt_->update(*control);
96 Real tol = 1e-8;
97 rAugLagSimOpt_->value(*control,tol);
98 }
99
100 void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
101 rAugLagSimOpt_->update(x,flag,iter);
102 }
103
104 Real value( const Vector<Real> &x, Real &tol ) {
105 return rAugLagSimOpt_->value(x,tol);
106 }
107
108 void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
109 ngval_++;
110 rAugLagSimOpt_->gradient(g,x,tol);
111 }
112
113 void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
114 rAugLagSimOpt_->hessVec(hv,v,x,tol);
115 }
116
117 // Return objective function value
119 return augLagSimOpt_->getObjectiveValue(*state_,x);
120 }
121
122 // Return constraint value
124 augLagSimOpt_->getConstraintVec(c,*state_,x);
125 }
126
127 // Return total number of constraint evaluations
129 return augLagSimOpt_->getNumberConstraintEvaluations();
130 }
131
132 // Return total number of objective evaluations
134 return augLagSimOpt_->getNumberFunctionEvaluations();
135 }
136
137 // Return total number of gradient evaluations
139 return ngval_;
140 //return augLagSimOpt_->getNumberGradientEvaluations();
141 }
142
143 // Reset with upated penalty parameter
144 void reset(const Vector<Real> &multiplier, const Real penaltyParameter) {
145 ngval_ = 0;
146 augLagSimOpt_->reset(multiplier,penaltyParameter);
147 }
148
149// For parametrized (stochastic) objective functions and constraints
150public:
151 void setParameter(const std::vector<Real> &param) {
153 rAugLagSimOpt_->setParameter(param);
154 }
155}; // class AugmentedLagrangian
156
157} // namespace ROL
158
159#endif
Contains definitions of custom data types in ROL.
Provides the interface to evaluate the augmented Lagrangian.
Defines the constraint operator interface for simulation-based optimization.
Provides the interface to evaluate simulation-based objective functions.
virtual void setParameter(const std::vector< Real > &param)
Provides the interface to evaluate the reduced SimOpt augmented Lagrangian.
void reset(const Vector< Real > &multiplier, const Real penaltyParameter)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
ROL::Ptr< Reduced_Objective_SimOpt< Real > > rAugLagSimOpt_
ROL::Ptr< AugmentedLagrangian_SimOpt< Real > > augLagSimOpt_
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
void getConstraintVec(Vector< Real > &c, const Vector< Real > &x)
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Reduced_AugmentedLagrangian_SimOpt(const ROL::Ptr< Objective_SimOpt< Real > > &obj, const ROL::Ptr< Constraint_SimOpt< Real > > &redCon, const ROL::Ptr< Constraint_SimOpt< Real > > &augCon, const ROL::Ptr< Vector< Real > > &state, const ROL::Ptr< Vector< Real > > &control, const ROL::Ptr< Vector< Real > > &adjoint, const ROL::Ptr< Vector< Real > > &augConVec, const ROL::Ptr< Vector< Real > > &multiplier, const Real penaltyParameter, ROL::ParameterList &parlist)
Defines the linear algebra or vector space interface.