ROL
Loading...
Searching...
No Matches
ROL_Reduced_Constraint_SimOpt_Def.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
11#ifndef ROL_REDUCED_CONSTRAINT_SIMOPT_DEF_H
12#define ROL_REDUCED_CONSTRAINT_SIMOPT_DEF_H
13
14namespace ROL {
15
16template<class Real>
18 if (!isUpdated_) {
19 // Update equality constraint with new Opt variable.
20 if (newUpdate_) conRed_->update_2(z,updateType_,updateIter_);
21 else conRed_->update_2(z,updateFlag_,updateIter_);
22 }
23 // Check if state has been computed.
24 bool isComputed = storage_ ? stateStore_->get(*state_,Constraint<Real>::getParameter()) : false;
25 // Solve state equation if not done already.
26 if (!isComputed || !storage_) {
27 // Solve state equation.
28 conRed_->solve(*dualadjoint_,*state_,z,tol);
29 nstat_++;
30 // Store state.
31 if (storage_) stateStore_->set(*state_,Constraint<Real>::getParameter());
32 }
33 if (!isUpdated_) {
34 // Update equality constraint with new Sim variable.
35 if (newUpdate_) conRed_->update_1(*state_,updateType_,updateIter_);
36 else conRed_->update_1(*state_,updateFlag_,updateIter_);
37 // Update full objective function.
38 if (newUpdate_) conVal_->update(*state_,z,updateType_,updateIter_);
39 else conVal_->update(*state_,z,updateFlag_,updateIter_);
40 isUpdated_ = true;
41 }
42}
43
44template<class Real>
46 // Check if adjoint has been computed.
47 bool isComputed = storage_ ? adjointStore_->get(*adjoint_,Constraint<Real>::getParameter()) : false;
48 // Solve adjoint equation if not done already.
49 if (!isComputed || !storage_) {
50 // Evaluate the full gradient wrt u
51 conVal_->applyAdjointJacobian_1(*dualstate_,w,*state_,z,tol);
52 // Solve adjoint equation
53 conRed_->applyInverseAdjointJacobian_1(*adjoint_,*dualstate_,*state_,z,tol);
54 adjoint_->scale(static_cast<Real>(-1));
55 nadjo_++;
56 // Store adjoint
57 if (storage_) adjointStore_->set(*adjoint_,Constraint<Real>::getParameter());
58 }
59}
60
61template<class Real>
63 // Solve state sensitivity equation
64 conRed_->applyJacobian_2(*dualadjoint_,v,*state_,z,tol);
65 dualadjoint_->scale(static_cast<Real>(-1));
66 conRed_->applyInverseJacobian_1(*state_sens_,*dualadjoint_,*state_,z,tol);
67 nssen_++;
68}
69
70template<class Real>
72 // Evaluate full hessVec in the direction (s,v)
73 conVal_->applyAdjointHessian_11(*dualstate_,w,*state_sens_,*state_,z,tol);
74 conVal_->applyAdjointHessian_21(*dualstate1_,w,v,*state_,z,tol);
75 dualstate_->plus(*dualstate1_);
76 // Apply adjoint Hessian of constraint
77 conRed_->applyAdjointHessian_11(*dualstate1_,*adjoint_,*state_sens_,*state_,z,tol);
78 dualstate_->plus(*dualstate1_);
79 conRed_->applyAdjointHessian_21(*dualstate1_,*adjoint_,v,*state_,z,tol);
80 dualstate_->plus(*dualstate1_);
81 // Solve adjoint sensitivity equation
82 dualstate_->scale(static_cast<Real>(-1));
83 conRed_->applyInverseAdjointJacobian_1(*adjoint_sens_,*dualstate_,*state_,z,tol);
84 nasen_++;
85}
86
87template<class Real>
89 const ROL::Ptr<Constraint_SimOpt<Real>> &conVal,
90 const ROL::Ptr<Constraint_SimOpt<Real>> &conRed,
91 const ROL::Ptr<VectorController<Real>> &stateStore,
92 const ROL::Ptr<Vector<Real>> &state,
93 const ROL::Ptr<Vector<Real>> &control,
94 const ROL::Ptr<Vector<Real>> &adjoint,
95 const ROL::Ptr<Vector<Real>> &residual,
96 bool storage,
97 bool useFDhessVec)
98 : conVal_( conVal ),
99 conRed_( conRed ),
100 stateStore_( stateStore ),
101 adjointStore_( ROL::makePtr<VectorController<Real>>() ),
102 state_( state->clone() ),
103 adjoint_( adjoint->clone() ),
104 residual_( residual->clone() ),
105 state_sens_( state->clone() ),
106 adjoint_sens_( adjoint->clone() ),
107 dualstate_( state->dual().clone() ),
108 dualstate1_( state->dual().clone() ),
109 dualadjoint_( adjoint->dual().clone() ),
110 dualcontrol_( control->dual().clone() ),
111 dualresidual_( residual->dual().clone() ),
112 storage_(storage), useFDhessVec_(useFDhessVec),
113 nupda_(0), nvalu_(0), njaco_(0), najac_(0), nhess_(0),
114 nstat_(0), nadjo_(0), nssen_(0), nasen_(0),
115 updateFlag_(true), updateIter_(0), updateType_(UpdateType::Initial),
116 newUpdate_(false), isUpdated_(true) {}
117
118template<class Real>
120 const ROL::Ptr<Constraint_SimOpt<Real>> &conVal,
121 const ROL::Ptr<Constraint_SimOpt<Real>> &conRed,
122 const ROL::Ptr<VectorController<Real>> &stateStore,
123 const ROL::Ptr<Vector<Real>> &state,
124 const ROL::Ptr<Vector<Real>> &control,
125 const ROL::Ptr<Vector<Real>> &adjoint,
126 const ROL::Ptr<Vector<Real>> &residual,
127 const ROL::Ptr<Vector<Real>> &dualstate,
128 const ROL::Ptr<Vector<Real>> &dualcontrol,
129 const ROL::Ptr<Vector<Real>> &dualadjoint,
130 const ROL::Ptr<Vector<Real>> &dualresidual,
131 bool storage,
132 bool useFDhessVec)
133 : conVal_( conVal ),
134 conRed_( conRed ),
135 stateStore_( stateStore ),
136 adjointStore_( ROL::makePtr<VectorController<Real>>() ),
137 state_( state->clone() ),
138 adjoint_( adjoint->clone() ),
139 residual_( residual->clone() ),
140 state_sens_( state->clone() ),
141 adjoint_sens_( adjoint->clone() ),
142 dualstate_( dualstate->clone() ),
143 dualstate1_( dualstate->clone() ),
144 dualadjoint_( dualadjoint->clone() ),
145 dualcontrol_( dualcontrol->clone() ),
146 dualresidual_( dualresidual->clone() ),
147 storage_(storage), useFDhessVec_(useFDhessVec),
148 nupda_(0), nvalu_(0), njaco_(0), najac_(0), nhess_(0),
149 nstat_(0), nadjo_(0), nssen_(0), nasen_(0),
150 updateFlag_(true), updateIter_(0), updateType_(UpdateType::Initial),
151 newUpdate_(false), isUpdated_(true) {}
152
153template<class Real>
154void Reduced_Constraint_SimOpt<Real>::summarize(std::ostream &stream, const Ptr<BatchManager<Real>> &bman) const {
155 int nupda(0), nvalu(0), njaco(0), najac(0), nhess(0), nstat(0), nadjo(0), nssen(0), nasen(0);
156 if (bman == nullPtr) {
157 nupda = nupda_;
158 nvalu = nvalu_;
159 njaco = njaco_;
160 najac = najac_;
161 nhess = nhess_;
162 nstat = nstat_;
163 nadjo = nadjo_;
164 nssen = nssen_;
165 nasen = nasen_;
166 }
167 else {
168 auto sumAll = [bman](int val) {
169 Real global(0), local(val);
170 bman->sumAll(&local,&global,1);
171 return static_cast<int>(global);
172 };
173 nupda = sumAll(nupda_);
174 nvalu = sumAll(nvalu_);
175 njaco = sumAll(njaco_);
176 najac = sumAll(najac_);
177 nhess = sumAll(nhess_);
178 nstat = sumAll(nstat_);
179 nadjo = sumAll(nadjo_);
180 nssen = sumAll(nssen_);
181 nasen = sumAll(nasen_);
182 }
183 stream << std::endl;
184 stream << std::string(80,'=') << std::endl;
185 stream << " ROL::Reduced_Objective_SimOpt::summarize" << std::endl;
186 stream << " Number of calls to update: " << nupda << std::endl;
187 stream << " Number of calls to value: " << nvalu << std::endl;
188 stream << " Number of calls to applyJacobian: " << njaco << std::endl;
189 stream << " Number of calls to applyAdjointJacobian: " << najac << std::endl;
190 stream << " Number of calls to hessvec: " << nhess << std::endl;
191 stream << " Number of state solves: " << nstat << std::endl;
192 stream << " Number of adjoint solves: " << nadjo << std::endl;
193 stream << " Number of state sensitivity solves: " << nssen << std::endl;
194 stream << " Number of adjoint sensitivity solves: " << nasen << std::endl;
195 stream << std::string(80,'=') << std::endl;
196 stream << std::endl;
197}
198
199template<class Real>
201 nupda_ = 0; nvalu_ = 0; njaco_ = 0; najac_ = 0; nhess_ = 0;
202 nstat_ = 0; nadjo_ = 0; nssen_ = 0; nasen_ = 0;
203}
204
205template<class Real>
206void Reduced_Constraint_SimOpt<Real>::update( const Vector<Real> &z, bool flag, int iter ) {
207 nupda_++;
208 updateFlag_ = flag;
209 updateIter_ = iter;
210 stateStore_->constraintUpdate(true);
211 adjointStore_->constraintUpdate(flag);
212}
213template<class Real>
215 nupda_++;
216 isUpdated_ = false;
217 newUpdate_ = true;
218 updateType_ = type;
219 updateIter_ = iter;
220 stateStore_->objectiveUpdate(type);
221 adjointStore_->objectiveUpdate(type);
222}
223
224template<class Real>
226 nvalu_++;
227 // Solve state equation
228 solve_state_equation(z,tol);
229 // Get constraint value
230 conVal_->value(c,*state_,z,tol);
231}
232
233template<class Real>
235 const Vector<Real> &z, Real &tol ) {
236 njaco_++;
237 // Solve state equation.
238 solve_state_equation(z,tol);
239 // Solve state sensitivity equation.
240 solve_state_sensitivity(v,z,tol);
241 // Apply Sim Jacobian to state sensitivity.
242 conVal_->applyJacobian_1(*residual_,*state_sens_,*state_,z,tol);
243 // Apply Opt Jacobian to vector.
244 conVal_->applyJacobian_2(jv,v,*state_,z,tol);
245 jv.plus(*residual_);
246}
247
248template<class Real>
250 const Vector<Real> &z, Real &tol ) {
251 najac_++;
252 // Solve state equation
253 solve_state_equation(z,tol);
254 // Solve adjoint equation
255 solve_adjoint_equation(w,z,tol);
256 // Evaluate the full gradient wrt z
257 conVal_->applyAdjointJacobian_2(*dualcontrol_,w,*state_,z,tol);
258 // Build gradient
259 conRed_->applyAdjointJacobian_2(ajw,*adjoint_,*state_,z,tol);
260 ajw.plus(*dualcontrol_);
261}
262
263template<class Real>
265 const Vector<Real> &v, const Vector<Real> &z,
266 Real &tol ) {
267 nhess_++;
268 if ( useFDhessVec_ ) {
270 }
271 else {
272 // Solve state equation
273 solve_state_equation(z,tol);
274 // Solve adjoint equation
275 solve_adjoint_equation(w,z,tol);
276 // Solve state sensitivity equation
277 solve_state_sensitivity(v,z,tol);
278 // Solve adjoint sensitivity equation
279 solve_adjoint_sensitivity(w,v,z,tol);
280 // Build hessVec
281 conRed_->applyAdjointJacobian_2(ahwv,*adjoint_sens_,*state_,z,tol);
282 conVal_->applyAdjointHessian_12(*dualcontrol_,w,*state_sens_,*state_,z,tol);
283 ahwv.plus(*dualcontrol_);
284 conVal_->applyAdjointHessian_22(*dualcontrol_,w,v,*state_,z,tol);
285 ahwv.plus(*dualcontrol_);
286 conRed_->applyAdjointHessian_12(*dualcontrol_,*adjoint_,*state_sens_,*state_,z,tol);
287 ahwv.plus(*dualcontrol_);
288 conRed_->applyAdjointHessian_22(*dualcontrol_,*adjoint_,v,*state_,z,tol);
289 ahwv.plus(*dualcontrol_);
290 }
291}
292
293}
294
295#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 .