ROL
ROL_CompositeConstraint_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 #ifndef ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
45 #define ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
46 
47 namespace ROL {
48 
49 template<typename Real>
51  const ROL::Ptr<Constraint_SimOpt<Real>> &conVal,
52  const ROL::Ptr<Constraint_SimOpt<Real>> &conRed,
53  const Vector<Real> &cVal,
54  const Vector<Real> &cRed,
55  const Vector<Real> &u,
56  const Vector<Real> &Sz,
57  const Vector<Real> &z,
58  bool storage,
59  bool isConRedParametrized)
60  : Constraint_SimOpt<Real>(), conVal_(conVal), conRed_(conRed),
61  updateFlag_(true), newUpdate_(false), updateIter_(0), updateType_(UpdateType::Initial),
62  storage_(storage), isConRedParametrized_(isConRedParametrized) {
63  Sz_ = Sz.clone();
64  primRed_ = cRed.clone();
65  dualRed_ = cRed.dual().clone();
66  primZ_ = z.clone();
67  dualZ_ = z.dual().clone();
68  dualZ1_ = z.dual().clone();
69  primU_ = u.clone();
70  stateStore_ = ROL::makePtr<VectorController<Real>>();
71 }
72 
73 template<typename Real>
75  const Vector<Real> &z,
76  bool flag, int iter) {
77  update_1(u,flag,iter);
78  update_2(z,flag,iter);
79 }
80 
81 template<typename Real>
83  bool flag, int iter) {
84  primU_->set(u);
85  conVal_->update_1(u, flag, iter);
86 }
87 
88 template<typename Real>
90  bool flag, int iter) {
91  conRed_->update_2(z,flag,iter);
92  updateFlag_ = flag;
93  updateIter_ = iter;
94  stateStore_->constraintUpdate(true);
95 }
96 
97 template<typename Real>
99  const Vector<Real> &z,
100  UpdateType type, int iter) {
101  update_1(u,type,iter);
102  update_2(z,type,iter);
103 }
104 
105 template<typename Real>
107  UpdateType type, int iter) {
108  primU_->set(u);
109  conVal_->update_1(u,type,iter);
110 }
111 
112 template<typename Real>
114  UpdateType type, int iter) {
115  conRed_->update_2(z,type,iter);
116  updateType_ = type;
117  updateIter_ = iter;
118  stateStore_->constraintUpdate(type);
119 }
120 
121 template<typename Real>
123  const Vector<Real> &z,
124  UpdateType type, int iter) {
125  Real ctol(std::sqrt(ROL_EPSILON<Real>()));
126  solveConRed(*Sz_, z, ctol);
127  conVal_->solve_update(u,*Sz_,type,iter);
128 }
129 
130 template<typename Real>
132  const Vector<Real> &u,
133  const Vector<Real> &z,
134  Real &tol) {
135  solveConRed(*Sz_, z, tol);
136  conVal_->value(c, u, *Sz_, tol);
137 }
138 
139 template<typename Real>
141  Vector<Real> &u,
142  const Vector<Real> &z,
143  Real &tol) {
144  solveConRed(*Sz_, z, tol);
145  conVal_->solve(c, u, *Sz_, tol);
146 }
147 
148 template<typename Real>
150  const Vector<Real> &v,
151  const Vector<Real> &u,
152  const Vector<Real> &z,
153  Real &tol) {
154  solveConRed(*Sz_, z, tol);
155  conVal_->applyJacobian_1(jv, v, u, *Sz_, tol);
156 }
157 
158 template<typename Real>
160  const Vector<Real> &v,
161  const Vector<Real> &u,
162  const Vector<Real> &z,
163  Real &tol) {
164  solveConRed(*Sz_, z, tol);
165  applySens(*primZ_, v, *Sz_, z, tol);
166  conVal_->applyJacobian_2(jv, *primZ_, u, *Sz_, tol);
167 }
168 
169 template<typename Real>
171  const Vector<Real> &v,
172  const Vector<Real> &u,
173  const Vector<Real> &z,
174  Real &tol) {
175  solveConRed(*Sz_, z, tol);
176  conVal_->applyInverseJacobian_1(ijv, v, u, *Sz_, tol);
177 }
178 
179 template<typename Real>
181  const Vector<Real> &v,
182  const Vector<Real> &u,
183  const Vector<Real> &z,
184  Real &tol) {
185  solveConRed(*Sz_, z, tol);
186  conVal_->applyAdjointJacobian_1(ajv, v, u, *Sz_, tol);
187 }
188 
189 template<typename Real>
191  const Vector<Real> &v,
192  const Vector<Real> &u,
193  const Vector<Real> &z,
194  Real &tol) {
195  solveConRed(*Sz_, z, tol);
196  conVal_->applyAdjointJacobian_2(*dualZ_, v, u, *Sz_, tol);
197  applyAdjointSens(ajv, *dualZ_, *Sz_, z, tol);
198 }
199 
200 template<typename Real>
202  const Vector<Real> &v,
203  const Vector<Real> &u,
204  const Vector<Real> &z,
205  Real &tol) {
206  solveConRed(*Sz_, z, tol);
207  conVal_->applyInverseAdjointJacobian_1(ijv, v, u, *Sz_, tol);
208 }
209 
210 template<typename Real>
212  const Vector<Real> &w,
213  const Vector<Real> &v,
214  const Vector<Real> &u,
215  const Vector<Real> &z,
216  Real &tol) {
217  solveConRed(*Sz_, z, tol);
218  conVal_->applyAdjointHessian_11(ahwv, w, v, u, *Sz_, tol);
219 }
220 
221 template<typename Real>
223  const Vector<Real> &w,
224  const Vector<Real> &v,
225  const Vector<Real> &u,
226  const Vector<Real> &z,
227  Real &tol) {
228  solveConRed(*Sz_, z, tol);
229  conVal_->applyAdjointHessian_12(*dualZ_, w, v, u, *Sz_, tol);
230  applyAdjointSens(ahwv, *dualZ_, *Sz_, z, tol);
231 }
232 
233 template<typename Real>
235  const Vector<Real> &w,
236  const Vector<Real> &v,
237  const Vector<Real> &u,
238  const Vector<Real> &z,
239  Real &tol) {
240  solveConRed(*Sz_, z, tol);
241  applySens(*primZ_, v, *Sz_, z, tol);
242  conVal_->applyAdjointHessian_21(ahwv, w, *primZ_, u, *Sz_, tol);
243 }
244 
245 template<typename Real>
247  const Vector<Real> &w,
248  const Vector<Real> &v,
249  const Vector<Real> &u,
250  const Vector<Real> &z,
251  Real &tol) {
252  ahwv.zero();
253  solveConRed(*Sz_, z, tol);
254  applySens(*primZ_, v, *Sz_, z, tol);
255 
256  conVal_->applyAdjointJacobian_2(*dualZ_, w, u, *Sz_, tol);
257  conRed_->applyInverseAdjointJacobian_1(*dualRed_, *dualZ_, *Sz_, z, tol);
258  conRed_->applyAdjointHessian_22(*dualZ_, *dualRed_, v, *Sz_, z, tol);
259  ahwv.axpy(static_cast<Real>(-1), *dualZ_);
260  conRed_->applyAdjointHessian_12(*dualZ_, *dualRed_, *primZ_, *Sz_, z, tol);
261  ahwv.axpy(static_cast<Real>(-1), *dualZ_);
262 
263  conRed_->applyAdjointHessian_11(*dualZ1_, *dualRed_, *primZ_, *Sz_, z, tol);
264  conRed_->applyAdjointHessian_21(*dualZ_, *dualRed_, v, *Sz_, z, tol);
265  dualZ1_->plus(*dualZ_);
266  dualZ1_->scale(static_cast<Real>(-1));
267 
268  conVal_->applyAdjointHessian_22(*dualZ_, w, *primZ_, u, *Sz_, tol);
269  dualZ1_->plus(*dualZ_);
270 
271  applyAdjointSens(*dualZ_, *dualZ1_, *Sz_, z, tol);
272  ahwv.plus(*dualZ_);
273 }
274 
275 template<typename Real>
276 void CompositeConstraint_SimOpt<Real>::setParameter(const std::vector<Real> &param) {
277  conVal_->setParameter(param);
278  if (isConRedParametrized_) {
279  conRed_->setParameter(param);
281  }
282 }
283 
284 template<typename Real>
286  const Vector<Real> &z,
287  Real &tol) {
288  std::vector<Real> param = Constraint_SimOpt<Real>::getParameter();
289  // Check if state has been computed.
290  bool isComputed = false;
291  if (storage_) isComputed = stateStore_->get(Sz,param);
292  // Solve state equation if not done already.
293  if (!isComputed || !storage_) {
294  // Solve state equation.
295  conRed_->solve(*primRed_,Sz,z,tol);
296  // Update equality constraint with new Sim variable.
297  if (newUpdate_) conRed_->update_1(Sz,updateType_,updateIter_);
298  else conRed_->update_1(Sz,updateFlag_,updateIter_);
299  // Update equality constraint.
300  if (newUpdate_) conRed_->update(Sz,z,updateType_,updateIter_);
301  else conRed_->update(Sz,z,updateFlag_,updateIter_);
302  if (newUpdate_) conVal_->update(*primU_,Sz,updateType_,updateIter_);
303  else conVal_->update(*primU_,Sz,updateFlag_,updateIter_);
304  // Store state.
305  if (storage_) stateStore_->set(Sz,param);
306  }
307 }
308 
309 template<typename Real>
311  const Vector<Real> &v,
312  const Vector<Real> &Sz,
313  const Vector<Real> &z,
314  Real &tol) {
315  // Solve linearization of reducible constraint in direction v
316  conRed_->applyJacobian_2(*primRed_, v, Sz, z, tol);
317  conRed_->applyInverseJacobian_1(jv, *primRed_, Sz, z, tol);
318  jv.scale(static_cast<Real>(-1));
319 }
320 
321 template<typename Real>
323  const Vector<Real> &v,
324  const Vector<Real> &Sz,
325  const Vector<Real> &z,
326  Real &tol) {
327  // Solve adjoint of linearized reducible constraint
328  conRed_->applyInverseAdjointJacobian_1(*dualRed_, v, Sz, z, tol);
329  conRed_->applyAdjointJacobian_2(ajv, *dualRed_, Sz, z, tol);
330  ajv.scale(static_cast<Real>(-1));
331 }
332 
333 } // namespace ROL
334 
335 #endif
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...
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).
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:226
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...
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void plus(const Vector &x)=0
Compute , where .
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
void solveConRed(Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
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, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.
const std::vector< Real > getParameter(void) const
void value(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Evaluate the constraint operator at .
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)
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
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 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 setParameter(const std::vector< Real > &param) override
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 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 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 applyAdjointSens(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
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 ...
virtual void setParameter(const std::vector< Real > &param)
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 applySens(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
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, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.
void solve(Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Given , solve for .
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, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.
Defines the constraint operator interface for simulation-based optimization.
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 .
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 is ...
ROL::Ptr< VectorController< Real > > stateStore_