ROL
ROL_CompositeConstraint_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#ifndef ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
11#define ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
12
13namespace ROL {
14
15template<typename Real>
17 const ROL::Ptr<Constraint_SimOpt<Real>> &conVal,
18 const ROL::Ptr<Constraint_SimOpt<Real>> &conRed,
19 const Vector<Real> &cVal,
20 const Vector<Real> &cRed,
21 const Vector<Real> &u,
22 const Vector<Real> &Sz,
23 const Vector<Real> &z,
24 bool storage,
25 bool isConRedParametrized)
26 : Constraint_SimOpt<Real>(), conVal_(conVal), conRed_(conRed),
27 updateFlag_(true), newUpdate_(false), updateIter_(0), updateType_(UpdateType::Initial),
28 storage_(storage), isConRedParametrized_(isConRedParametrized) {
29 Sz_ = Sz.clone();
30 primRed_ = cRed.clone();
31 dualRed_ = cRed.dual().clone();
32 primZ_ = z.clone();
33 dualZ_ = z.dual().clone();
34 dualZ1_ = z.dual().clone();
35 primU_ = u.clone();
36 stateStore_ = ROL::makePtr<VectorController<Real>>();
37}
38
39template<typename Real>
41 const Vector<Real> &z,
42 bool flag, int iter) {
43 update_1(u,flag,iter);
44 update_2(z,flag,iter);
45}
46
47template<typename Real>
49 bool flag, int iter) {
50 primU_->set(u);
51 conVal_->update_1(u, flag, iter);
52}
53
54template<typename Real>
56 bool flag, int iter) {
57 conRed_->update_2(z,flag,iter);
58 updateFlag_ = flag;
59 updateIter_ = iter;
60 stateStore_->constraintUpdate(true);
61}
62
63template<typename Real>
65 const Vector<Real> &z,
66 UpdateType type, int iter) {
67 update_1(u,type,iter);
68 update_2(z,type,iter);
69}
70
71template<typename Real>
73 UpdateType type, int iter) {
74 primU_->set(u);
75 conVal_->update_1(u,type,iter);
76}
77
78template<typename Real>
80 UpdateType type, int iter) {
81 conRed_->update_2(z,type,iter);
82 updateType_ = type;
83 updateIter_ = iter;
84 stateStore_->constraintUpdate(type);
85}
86
87template<typename Real>
89 const Vector<Real> &z,
90 UpdateType type, int iter) {
91 Real ctol(std::sqrt(ROL_EPSILON<Real>()));
92 solveConRed(*Sz_, z, ctol);
93 conVal_->solve_update(u,*Sz_,type,iter);
94}
95
96template<typename Real>
98 const Vector<Real> &u,
99 const Vector<Real> &z,
100 Real &tol) {
101 solveConRed(*Sz_, z, tol);
102 conVal_->value(c, u, *Sz_, tol);
103}
104
105template<typename Real>
107 Vector<Real> &u,
108 const Vector<Real> &z,
109 Real &tol) {
110 solveConRed(*Sz_, z, tol);
111 conVal_->solve(c, u, *Sz_, tol);
112}
113
114template<typename Real>
116 const Vector<Real> &v,
117 const Vector<Real> &u,
118 const Vector<Real> &z,
119 Real &tol) {
120 solveConRed(*Sz_, z, tol);
121 conVal_->applyJacobian_1(jv, v, u, *Sz_, tol);
122}
123
124template<typename Real>
126 const Vector<Real> &v,
127 const Vector<Real> &u,
128 const Vector<Real> &z,
129 Real &tol) {
130 solveConRed(*Sz_, z, tol);
131 applySens(*primZ_, v, *Sz_, z, tol);
132 conVal_->applyJacobian_2(jv, *primZ_, u, *Sz_, tol);
133}
134
135template<typename Real>
137 const Vector<Real> &v,
138 const Vector<Real> &u,
139 const Vector<Real> &z,
140 Real &tol) {
141 solveConRed(*Sz_, z, tol);
142 conVal_->applyInverseJacobian_1(ijv, v, u, *Sz_, tol);
143}
144
145template<typename Real>
147 const Vector<Real> &v,
148 const Vector<Real> &u,
149 const Vector<Real> &z,
150 Real &tol) {
151 solveConRed(*Sz_, z, tol);
152 conVal_->applyAdjointJacobian_1(ajv, v, u, *Sz_, tol);
153}
154
155template<typename Real>
157 const Vector<Real> &v,
158 const Vector<Real> &u,
159 const Vector<Real> &z,
160 Real &tol) {
161 solveConRed(*Sz_, z, tol);
162 conVal_->applyAdjointJacobian_2(*dualZ_, v, u, *Sz_, tol);
163 applyAdjointSens(ajv, *dualZ_, *Sz_, z, tol);
164}
165
166template<typename Real>
168 const Vector<Real> &v,
169 const Vector<Real> &u,
170 const Vector<Real> &z,
171 Real &tol) {
172 solveConRed(*Sz_, z, tol);
173 conVal_->applyInverseAdjointJacobian_1(ijv, v, u, *Sz_, tol);
174}
175
176template<typename Real>
178 const Vector<Real> &w,
179 const Vector<Real> &v,
180 const Vector<Real> &u,
181 const Vector<Real> &z,
182 Real &tol) {
183 solveConRed(*Sz_, z, tol);
184 conVal_->applyAdjointHessian_11(ahwv, w, v, u, *Sz_, tol);
185}
186
187template<typename Real>
189 const Vector<Real> &w,
190 const Vector<Real> &v,
191 const Vector<Real> &u,
192 const Vector<Real> &z,
193 Real &tol) {
194 solveConRed(*Sz_, z, tol);
195 conVal_->applyAdjointHessian_12(*dualZ_, w, v, u, *Sz_, tol);
196 applyAdjointSens(ahwv, *dualZ_, *Sz_, z, tol);
197}
198
199template<typename Real>
201 const Vector<Real> &w,
202 const Vector<Real> &v,
203 const Vector<Real> &u,
204 const Vector<Real> &z,
205 Real &tol) {
206 solveConRed(*Sz_, z, tol);
207 applySens(*primZ_, v, *Sz_, z, tol);
208 conVal_->applyAdjointHessian_21(ahwv, w, *primZ_, u, *Sz_, tol);
209}
210
211template<typename Real>
213 const Vector<Real> &w,
214 const Vector<Real> &v,
215 const Vector<Real> &u,
216 const Vector<Real> &z,
217 Real &tol) {
218 ahwv.zero();
219 solveConRed(*Sz_, z, tol);
220 applySens(*primZ_, v, *Sz_, z, tol);
221
222 conVal_->applyAdjointJacobian_2(*dualZ_, w, u, *Sz_, tol);
223 conRed_->applyInverseAdjointJacobian_1(*dualRed_, *dualZ_, *Sz_, z, tol);
224 conRed_->applyAdjointHessian_22(*dualZ_, *dualRed_, v, *Sz_, z, tol);
225 ahwv.axpy(static_cast<Real>(-1), *dualZ_);
226 conRed_->applyAdjointHessian_12(*dualZ_, *dualRed_, *primZ_, *Sz_, z, tol);
227 ahwv.axpy(static_cast<Real>(-1), *dualZ_);
228
229 conRed_->applyAdjointHessian_11(*dualZ1_, *dualRed_, *primZ_, *Sz_, z, tol);
230 conRed_->applyAdjointHessian_21(*dualZ_, *dualRed_, v, *Sz_, z, tol);
231 dualZ1_->plus(*dualZ_);
232 dualZ1_->scale(static_cast<Real>(-1));
233
234 conVal_->applyAdjointHessian_22(*dualZ_, w, *primZ_, u, *Sz_, tol);
235 dualZ1_->plus(*dualZ_);
236
237 applyAdjointSens(*dualZ_, *dualZ1_, *Sz_, z, tol);
238 ahwv.plus(*dualZ_);
239}
240
241template<typename Real>
242void CompositeConstraint_SimOpt<Real>::setParameter(const std::vector<Real> &param) {
243 conVal_->setParameter(param);
244 if (isConRedParametrized_) {
245 conRed_->setParameter(param);
247 }
248}
249
250template<typename Real>
252 const Vector<Real> &z,
253 Real &tol) {
254 std::vector<Real> param = Constraint_SimOpt<Real>::getParameter();
255 // Check if state has been computed.
256 bool isComputed = false;
257 if (storage_) isComputed = stateStore_->get(Sz,param);
258 // Solve state equation if not done already.
259 if (!isComputed || !storage_) {
260 // Solve state equation.
261 conRed_->solve(*primRed_,Sz,z,tol);
262 // Update equality constraint with new Sim variable.
263 if (newUpdate_) conRed_->update_1(Sz,updateType_,updateIter_);
264 else conRed_->update_1(Sz,updateFlag_,updateIter_);
265 // Update equality constraint.
266 if (newUpdate_) conRed_->update(Sz,z,updateType_,updateIter_);
267 else conRed_->update(Sz,z,updateFlag_,updateIter_);
268 if (newUpdate_) conVal_->update(*primU_,Sz,updateType_,updateIter_);
269 else conVal_->update(*primU_,Sz,updateFlag_,updateIter_);
270 // Store state.
271 if (storage_) stateStore_->set(Sz,param);
272 }
273}
274
275template<typename Real>
277 const Vector<Real> &v,
278 const Vector<Real> &Sz,
279 const Vector<Real> &z,
280 Real &tol) {
281 // Solve linearization of reducible constraint in direction v
282 conRed_->applyJacobian_2(*primRed_, v, Sz, z, tol);
283 conRed_->applyInverseJacobian_1(jv, *primRed_, Sz, z, tol);
284 jv.scale(static_cast<Real>(-1));
285}
286
287template<typename Real>
289 const Vector<Real> &v,
290 const Vector<Real> &Sz,
291 const Vector<Real> &z,
292 Real &tol) {
293 // Solve adjoint of linearized reducible constraint
294 conRed_->applyInverseAdjointJacobian_1(*dualRed_, v, Sz, z, tol);
295 conRed_->applyAdjointJacobian_2(ajv, *dualRed_, Sz, z, tol);
296 ajv.scale(static_cast<Real>(-1));
297}
298
299} // namespace ROL
300
301#endif
void applyInverseJacobian_1(Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the inverse partial constraint Jacobian at , , to the vector .
void applyAdjointHessian_11(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at ...
void solve(Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Given , solve for .
void update_2(const Vector< Real > &z, bool flag=true, int iter=-1) override
Update constraint functions with respect to Opt variable. x is the optimization variable,...
void applyAdjointHessian_22(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian ...
void applyJacobian_2(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the partial constraint Jacobian at , , to the vector .
void applyAdjointSens(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
void solve_update(const Vector< Real > &u, const Vector< Real > &z, UpdateType type, int iter=-1) override
Update SimOpt constraint during solve (disconnected from optimization updates).
void solveConRed(Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
void applyAdjointHessian_21(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at...
ROL::Ptr< VectorController< Real > > stateStore_
void update_1(const Vector< Real > &u, bool flag=true, int iter=-1) override
Update constraint functions with respect to Sim variable. x is the optimization variable,...
void value(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Evaluate the constraint operator at .
void applyAdjointJacobian_2(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
void applySens(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
void applyAdjointHessian_12(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at...
void setParameter(const std::vector< Real > &param) override
void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1) override
Update constraint functions. x is the optimization variable, flag = true if optimization variable i...
void applyAdjointJacobian_1(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
void applyInverseAdjointJacobian_1(Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the inverse of the adjoint of the partial constraint Jacobian at , , to the vector .
void applyJacobian_1(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the partial constraint Jacobian at , , to the vector .
CompositeConstraint_SimOpt(const ROL::Ptr< Constraint_SimOpt< Real > > &conVal, const ROL::Ptr< Constraint_SimOpt< Real > > &conRed, const Vector< Real > &cVal, const Vector< Real > &cRed, const Vector< Real > &u, const Vector< Real > &Sz, const Vector< Real > &z, bool storage=true, bool isConRedParametrized=false)
Defines the constraint operator interface for simulation-based optimization.
virtual void setParameter(const std::vector< Real > &param)
const std::vector< Real > getParameter(void) const
Defines the linear algebra or vector space interface.
virtual void scale(const Real alpha)=0
Compute where .
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
virtual void plus(const Vector &x)=0
Compute , where .
virtual void zero()
Set to zero vector.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
virtual void update_2(const Vector< Real > &z, bool flag=true, int iter=-1) override
Update constraint functions with respect to Opt variable. z is the control variable,...