ROL
ROL_Reduced_Constraint_SimOpt_Def.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Rapid Optimization Library (ROL) Package
5// Copyright (2014) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact lead developers:
38// Drew Kouri (dpkouri@sandia.gov) and
39// Denis Ridzal (dridzal@sandia.gov)
40//
41// ************************************************************************
42// @HEADER
43
44
45#ifndef ROL_REDUCED_CONSTRAINT_SIMOPT_DEF_H
46#define ROL_REDUCED_CONSTRAINT_SIMOPT_DEF_H
47
48namespace ROL {
49
50template<class Real>
52 if (!isUpdated_) {
53 // Update equality constraint with new Opt variable.
54 if (newUpdate_) conRed_->update_2(z,updateType_,updateIter_);
55 else conRed_->update_2(z,updateFlag_,updateIter_);
56 }
57 // Check if state has been computed.
58 bool isComputed = storage_ ? stateStore_->get(*state_,Constraint<Real>::getParameter()) : false;
59 // Solve state equation if not done already.
60 if (!isComputed || !storage_) {
61 // Solve state equation.
62 conRed_->solve(*dualadjoint_,*state_,z,tol);
63 nstat_++;
64 // Store state.
65 if (storage_) stateStore_->set(*state_,Constraint<Real>::getParameter());
66 }
67 if (!isUpdated_) {
68 // Update equality constraint with new Sim variable.
69 if (newUpdate_) conRed_->update_1(*state_,updateType_,updateIter_);
70 else conRed_->update_1(*state_,updateFlag_,updateIter_);
71 // Update full objective function.
72 if (newUpdate_) conVal_->update(*state_,z,updateType_,updateIter_);
73 else conVal_->update(*state_,z,updateFlag_,updateIter_);
74 isUpdated_ = true;
75 }
76}
77
78template<class Real>
80 // Check if adjoint has been computed.
81 bool isComputed = storage_ ? adjointStore_->get(*adjoint_,Constraint<Real>::getParameter()) : false;
82 // Solve adjoint equation if not done already.
83 if (!isComputed || !storage_) {
84 // Evaluate the full gradient wrt u
85 conVal_->applyAdjointJacobian_1(*dualstate_,w,*state_,z,tol);
86 // Solve adjoint equation
87 conRed_->applyInverseAdjointJacobian_1(*adjoint_,*dualstate_,*state_,z,tol);
88 adjoint_->scale(static_cast<Real>(-1));
89 nadjo_++;
90 // Store adjoint
91 if (storage_) adjointStore_->set(*adjoint_,Constraint<Real>::getParameter());
92 }
93}
94
95template<class Real>
97 // Solve state sensitivity equation
98 conRed_->applyJacobian_2(*dualadjoint_,v,*state_,z,tol);
99 dualadjoint_->scale(static_cast<Real>(-1));
100 conRed_->applyInverseJacobian_1(*state_sens_,*dualadjoint_,*state_,z,tol);
101 nssen_++;
102}
103
104template<class Real>
106 // Evaluate full hessVec in the direction (s,v)
107 conVal_->applyAdjointHessian_11(*dualstate_,w,*state_sens_,*state_,z,tol);
108 conVal_->applyAdjointHessian_21(*dualstate1_,w,v,*state_,z,tol);
109 dualstate_->plus(*dualstate1_);
110 // Apply adjoint Hessian of constraint
111 conRed_->applyAdjointHessian_11(*dualstate1_,*adjoint_,*state_sens_,*state_,z,tol);
112 dualstate_->plus(*dualstate1_);
113 conRed_->applyAdjointHessian_21(*dualstate1_,*adjoint_,v,*state_,z,tol);
114 dualstate_->plus(*dualstate1_);
115 // Solve adjoint sensitivity equation
116 dualstate_->scale(static_cast<Real>(-1));
117 conRed_->applyInverseAdjointJacobian_1(*adjoint_sens_,*dualstate_,*state_,z,tol);
118 nasen_++;
119}
120
121template<class Real>
123 const ROL::Ptr<Constraint_SimOpt<Real>> &conVal,
124 const ROL::Ptr<Constraint_SimOpt<Real>> &conRed,
125 const ROL::Ptr<VectorController<Real>> &stateStore,
126 const ROL::Ptr<Vector<Real>> &state,
127 const ROL::Ptr<Vector<Real>> &control,
128 const ROL::Ptr<Vector<Real>> &adjoint,
129 const ROL::Ptr<Vector<Real>> &residual,
130 bool storage,
131 bool useFDhessVec)
132 : conVal_( conVal ),
133 conRed_( conRed ),
134 stateStore_( stateStore ),
135 adjointStore_( ROL::makePtr<VectorController<Real>>() ),
136 state_( state->clone() ),
137 adjoint_( adjoint->clone() ),
138 residual_( residual->clone() ),
139 state_sens_( state->clone() ),
140 adjoint_sens_( adjoint->clone() ),
141 dualstate_( state->dual().clone() ),
142 dualstate1_( state->dual().clone() ),
143 dualadjoint_( adjoint->dual().clone() ),
144 dualcontrol_( control->dual().clone() ),
145 dualresidual_( residual->dual().clone() ),
146 storage_(storage), useFDhessVec_(useFDhessVec),
147 nupda_(0), nvalu_(0), njaco_(0), najac_(0), nhess_(0),
148 nstat_(0), nadjo_(0), nssen_(0), nasen_(0),
149 updateFlag_(true), updateIter_(0), updateType_(UpdateType::Initial),
150 newUpdate_(false), isUpdated_(true) {}
151
152template<class Real>
154 const ROL::Ptr<Constraint_SimOpt<Real>> &conVal,
155 const ROL::Ptr<Constraint_SimOpt<Real>> &conRed,
156 const ROL::Ptr<VectorController<Real>> &stateStore,
157 const ROL::Ptr<Vector<Real>> &state,
158 const ROL::Ptr<Vector<Real>> &control,
159 const ROL::Ptr<Vector<Real>> &adjoint,
160 const ROL::Ptr<Vector<Real>> &residual,
161 const ROL::Ptr<Vector<Real>> &dualstate,
162 const ROL::Ptr<Vector<Real>> &dualcontrol,
163 const ROL::Ptr<Vector<Real>> &dualadjoint,
164 const ROL::Ptr<Vector<Real>> &dualresidual,
165 bool storage,
166 bool useFDhessVec)
167 : conVal_( conVal ),
168 conRed_( conRed ),
169 stateStore_( stateStore ),
170 adjointStore_( ROL::makePtr<VectorController<Real>>() ),
171 state_( state->clone() ),
172 adjoint_( adjoint->clone() ),
173 residual_( residual->clone() ),
174 state_sens_( state->clone() ),
175 adjoint_sens_( adjoint->clone() ),
176 dualstate_( dualstate->clone() ),
177 dualstate1_( dualstate->clone() ),
178 dualadjoint_( dualadjoint->clone() ),
179 dualcontrol_( dualcontrol->clone() ),
180 dualresidual_( dualresidual->clone() ),
181 storage_(storage), useFDhessVec_(useFDhessVec),
182 nupda_(0), nvalu_(0), njaco_(0), najac_(0), nhess_(0),
183 nstat_(0), nadjo_(0), nssen_(0), nasen_(0),
184 updateFlag_(true), updateIter_(0), updateType_(UpdateType::Initial),
185 newUpdate_(false), isUpdated_(true) {}
186
187template<class Real>
188void Reduced_Constraint_SimOpt<Real>::summarize(std::ostream &stream, const Ptr<BatchManager<Real>> &bman) const {
189 int nupda(0), nvalu(0), njaco(0), najac(0), nhess(0), nstat(0), nadjo(0), nssen(0), nasen(0);
190 if (bman == nullPtr) {
191 nupda = nupda_;
192 nvalu = nvalu_;
193 njaco = njaco_;
194 najac = najac_;
195 nhess = nhess_;
196 nstat = nstat_;
197 nadjo = nadjo_;
198 nssen = nssen_;
199 nasen = nasen_;
200 }
201 else {
202 auto sumAll = [bman](int val) {
203 Real global(0), local(val);
204 bman->sumAll(&local,&global,1);
205 return static_cast<int>(global);
206 };
207 nupda = sumAll(nupda_);
208 nvalu = sumAll(nvalu_);
209 njaco = sumAll(njaco_);
210 najac = sumAll(najac_);
211 nhess = sumAll(nhess_);
212 nstat = sumAll(nstat_);
213 nadjo = sumAll(nadjo_);
214 nssen = sumAll(nssen_);
215 nasen = sumAll(nasen_);
216 }
217 stream << std::endl;
218 stream << std::string(80,'=') << std::endl;
219 stream << " ROL::Reduced_Objective_SimOpt::summarize" << std::endl;
220 stream << " Number of calls to update: " << nupda << std::endl;
221 stream << " Number of calls to value: " << nvalu << std::endl;
222 stream << " Number of calls to applyJacobian: " << njaco << std::endl;
223 stream << " Number of calls to applyAdjointJacobian: " << najac << std::endl;
224 stream << " Number of calls to hessvec: " << nhess << std::endl;
225 stream << " Number of state solves: " << nstat << std::endl;
226 stream << " Number of adjoint solves: " << nadjo << std::endl;
227 stream << " Number of state sensitivity solves: " << nssen << std::endl;
228 stream << " Number of adjoint sensitivity solves: " << nasen << std::endl;
229 stream << std::string(80,'=') << std::endl;
230 stream << std::endl;
231}
232
233template<class Real>
235 nupda_ = 0; nvalu_ = 0; njaco_ = 0; najac_ = 0; nhess_ = 0;
236 nstat_ = 0; nadjo_ = 0; nssen_ = 0; nasen_ = 0;
237}
238
239template<class Real>
240void Reduced_Constraint_SimOpt<Real>::update( const Vector<Real> &z, bool flag, int iter ) {
241 nupda_++;
242 updateFlag_ = flag;
243 updateIter_ = iter;
244 stateStore_->constraintUpdate(true);
245 adjointStore_->constraintUpdate(flag);
246}
247template<class Real>
249 nupda_++;
250 isUpdated_ = false;
251 newUpdate_ = true;
252 updateType_ = type;
253 updateIter_ = iter;
254 stateStore_->objectiveUpdate(type);
255 adjointStore_->objectiveUpdate(type);
256}
257
258template<class Real>
260 nvalu_++;
261 // Solve state equation
262 solve_state_equation(z,tol);
263 // Get constraint value
264 conVal_->value(c,*state_,z,tol);
265}
266
267template<class Real>
269 const Vector<Real> &z, Real &tol ) {
270 njaco_++;
271 // Solve state equation.
272 solve_state_equation(z,tol);
273 // Solve state sensitivity equation.
274 solve_state_sensitivity(v,z,tol);
275 // Apply Sim Jacobian to state sensitivity.
276 conVal_->applyJacobian_1(*residual_,*state_sens_,*state_,z,tol);
277 // Apply Opt Jacobian to vector.
278 conVal_->applyJacobian_2(jv,v,*state_,z,tol);
279 jv.plus(*residual_);
280}
281
282template<class Real>
284 const Vector<Real> &z, Real &tol ) {
285 najac_++;
286 // Solve state equation
287 solve_state_equation(z,tol);
288 // Solve adjoint equation
289 solve_adjoint_equation(w,z,tol);
290 // Evaluate the full gradient wrt z
291 conVal_->applyAdjointJacobian_2(*dualcontrol_,w,*state_,z,tol);
292 // Build gradient
293 conRed_->applyAdjointJacobian_2(ajw,*adjoint_,*state_,z,tol);
294 ajw.plus(*dualcontrol_);
295}
296
297template<class Real>
299 const Vector<Real> &v, const Vector<Real> &z,
300 Real &tol ) {
301 nhess_++;
302 if ( useFDhessVec_ ) {
304 }
305 else {
306 // Solve state equation
307 solve_state_equation(z,tol);
308 // Solve adjoint equation
309 solve_adjoint_equation(w,z,tol);
310 // Solve state sensitivity equation
311 solve_state_sensitivity(v,z,tol);
312 // Solve adjoint sensitivity equation
313 solve_adjoint_sensitivity(w,v,z,tol);
314 // Build hessVec
315 conRed_->applyAdjointJacobian_2(ahwv,*adjoint_sens_,*state_,z,tol);
316 conVal_->applyAdjointHessian_12(*dualcontrol_,w,*state_sens_,*state_,z,tol);
317 ahwv.plus(*dualcontrol_);
318 conVal_->applyAdjointHessian_22(*dualcontrol_,w,v,*state_,z,tol);
319 ahwv.plus(*dualcontrol_);
320 conRed_->applyAdjointHessian_12(*dualcontrol_,*adjoint_,*state_sens_,*state_,z,tol);
321 ahwv.plus(*dualcontrol_);
322 conRed_->applyAdjointHessian_22(*dualcontrol_,*adjoint_,v,*state_,z,tol);
323 ahwv.plus(*dualcontrol_);
324 }
325}
326
327}
328
329#endif
Defines the constraint operator interface for simulation-based optimization.
Defines the general constraint operator interface.
virtual 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 ,...
Reduced_Constraint_SimOpt(const ROL::Ptr< Constraint_SimOpt< Real > > &conVal, const ROL::Ptr< Constraint_SimOpt< Real > > &conRed, const ROL::Ptr< VectorController< Real > > &stateStore, const ROL::Ptr< Vector< Real > > &state, const ROL::Ptr< Vector< Real > > &control, const ROL::Ptr< Vector< Real > > &adjoint, const ROL::Ptr< Vector< Real > > &residual, bool storage=true, bool useFDhessVec=false)
Constructor.
void solve_adjoint_sensitivity(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given , the adjoint variable , and a direction , solve the adjoint sensitvity equation for .
void solve_state_sensitivity(const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given which solves the state equation and a direction , solve the state senstivity equation for .
void solve_state_equation(const Vector< Real > &z, Real &tol)
void value(Vector< Real > &c, const Vector< Real > &z, Real &tol)
Given , evaluate the equality constraint where solves .
void summarize(std::ostream &stream, const Ptr< BatchManager< Real > > &bman=nullPtr) const
void applyAdjointHessian(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given , evaluate the Hessian of the objective function in the direction .
void update(const Vector< Real > &z, bool flag=true, int iter=-1)
Update the SimOpt objective function and equality constraint.
void applyAdjointJacobian(Vector< Real > &ajw, const Vector< Real > &w, const Vector< Real > &z, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector .
void solve_adjoint_equation(const Vector< Real > &w, const Vector< Real > &z, Real &tol)
Given which solves the state equation, solve the adjoint equation for .
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given , apply the Jacobian to a vector where solves .
Defines the linear algebra or vector space interface.
virtual void plus(const Vector &x)=0
Compute , where .