ROL
ROL_DynamicConstraint_CheckInterface.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 #pragma once
46 #ifndef ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
47 #define ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
48 
49 #include "ROL_FunctionBindings.hpp"
51 
52 namespace ROL {
53 namespace details {
54 
55 using namespace std;
56 namespace ph = std::placeholders;
57 
58 template<typename Real>
60 private:
61 
62  using V = Vector<Real>;
64 
66  Real tol_;
68 
69 public:
70 
71  DynamicConstraint_CheckInterface( Con& con ) : con_(con) {
72  ts_.t.resize(2);
73  ts_.t.at(0) = 0.01;
74  ts_.t.at(1) = 0.02345;
75  ts_.k = 0;
76  }
77 
78  DynamicConstraint_CheckInterface( Con& con, TimeStamp<Real>& ts ) : con_(con), ts_(ts) {
79  }
80 
81  f_update_t<Real> update_uo( const V& un, const V& z ) {
82  return bind( &Con::update, &con_, ph::_1, cref(un), cref(z), ts_ );
83  }
84 
85  f_update_t<Real> update_un( const V& uo, const V& z ) {
86  return bind( &Con::update, &con_, cref(uo), ph::_1, cref(z), ts_ );
87  }
88 
89  f_update_t<Real> update_z( const V& uo, const V& un ) {
90  return bind( &Con::update, &con_, cref(uo), cref(un), ph::_1, ts_ );
91  }
92 
93  //----------------------------------------------------------------------------
94 
95  f_vector_t<Real> value_uo( const V& un, const V& z ) {
96  return bind( &Con::value, &con_,
97  ph::_1, ph::_2, cref(un), cref(z), ts_ );
98  }
99 
100  f_vector_t<Real> value_un( const V& uo, const V& z ) {
101  return bind( &Con::value, &con_,
102  ph::_1, cref(uo), ph::_2, cref(z), ts_ );
103  }
104 
105  f_vector_t<Real> value_z( const V& uo, const V& un ) {
106  return bind( &Con::value, &con_,
107  ph::_1, cref(uo), cref(un), ph::_2, ts_ );
108  }
109 
110  f_solve_t<Real> solve_un( const V& uo, const V& z ) {
111  return bind( &Con::solve, &con_,
112  ph::_1, cref(uo), ph::_2, cref(z), ts_ );
113  }
114 
115  //----------------------------------------------------------------------------
116 
117  f_dderiv_t<Real> jacobian_uo( const V& un, const V& z ) {
118  return bind( &Con::applyJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
119  cref(un), cref(z), ts_ );
120  }
121 
122  f_dderiv_t<Real> jacobian_un( const V& uo, const V& z ) {
123  return bind( &Con::applyJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
124  ph::_3, cref(z), ts_ );
125  }
126 
127  f_dderiv_t<Real> inverseJacobian_un( const V& uo, const V& z ) {
128  return bind( &Con::applyInverseJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
129  ph::_3, cref(z), ts_ );
130  }
131 
132  f_dderiv_t<Real> jacobian_z( const V& uo, const V& un ) {
133  return bind( &Con::applyJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
134  cref(un), ph::_3, ts_ );
135  }
136 
137  //----------------------------------------------------------------------------
138 
139  f_dderiv_t<Real> adjointJacobian_uo( const V& un, const V& z ) {
140  return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
141  cref(un), cref(z), ts_ );
142  }
143 
144  f_dderiv_t<Real> adjointJacobian_un( const V& uo, const V& z ) {
145  return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
146  ph::_3, cref(z), ts_ );
147  }
148 
149  f_dderiv_t<Real> inverseAdjointJacobian_un( const V& uo, const V& z ) {
150  return bind( &Con::applyInverseAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
151  ph::_3, cref(z), ts_ );
152  }
153 
154  f_dderiv_t<Real> adjointJacobian_z( const V& uo, const V& un ) {
155  return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
156  cref(un), ph::_3, ts_ );
157  }
158 
159  //----------------------------------------------------------------------------
160 
161  f_dderiv_t<Real> adjointJacobian_uo_uo( const V& un, const V& z ) {
162  return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
163  cref(un), cref(z), ts_ );
164  }
165 
166  f_dderiv_t<Real> adjointJacobian_uo_un( const V& uo, const V& z ) {
167  return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, cref(uo),
168  ph::_3, cref(z), ts_ );
169  }
170 
171  f_dderiv_t<Real> adjointJacobian_uo_z( const V& uo, const V& un ) {
172  return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, cref(uo),
173  cref(un), ph::_3, ts_ );
174  }
175 
176  f_dderiv_t<Real> adjointJacobian_un_uo( const V& un, const V& z ) {
177  return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, ph::_3,
178  cref(un), cref(z), ts_ );
179  }
180 
181  f_dderiv_t<Real> adjointJacobian_un_un( const V& uo, const V& z ) {
182  return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
183  ph::_3, cref(z), ts_ );
184  }
185 
186  f_dderiv_t<Real> adjointJacobian_un_z( const V& uo, const V& un ) {
187  return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
188  cref(un), ph::_3, ts_ );
189  }
190 
191  f_dderiv_t<Real> adjointJacobian_z_uo( const V& un, const V& z ) {
192  return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, ph::_3,
193  cref(un), cref(z), ts_ );
194  }
195 
196  f_dderiv_t<Real> adjointJacobian_z_un( const V& uo, const V& z ) {
197  return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
198  ph::_3, cref(z), ts_ );
199  }
200 
201  f_dderiv_t<Real> adjointJacobian_z_z( const V& uo, const V& un ) {
202  return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
203  cref(un), ph::_3, ts_ );
204  }
205 
206  //----------------------------------------------------------------------------
207 
208  f_dderiv_t<Real> adjointHessian_un_un( const V& uo, const V& z, const V& l ) {
209  return bind( &Con::applyAdjointHessian_un_un, &con_, ph::_1, cref(l), ph::_2,
210  cref(uo), ph::_3, cref(z), ts_ );
211  }
212 
213  f_dderiv_t<Real> adjointHessian_un_uo( const V& uo, const V& z, const V& l ) {
214  return bind( &Con::applyAdjointHessian_un_uo, &con_, ph::_1, cref(l), ph::_2,
215  cref(uo), ph::_3, cref(z), ts_ );
216  }
217 
218  f_dderiv_t<Real> adjointHessian_un_z( const V& uo, const V& z, const V& l ) {
219  return bind( &Con::applyAdjointHessian_un_z, &con_, ph::_1, cref(l), ph::_2,
220  cref(uo), ph::_3, cref(z), ts_ );
221  }
222 
223  //----------------------------------------------------------------------------
224 
225  f_dderiv_t<Real> adjointHessian_uo_un( const V& un, const V& z, const V& l ) {
226  return bind( &Con::applyAdjointHessian_uo_un, &con_, ph::_1, cref(l), ph::_2,
227  ph::_3, cref(un), cref(z), ts_ );
228  }
229 
230  f_dderiv_t<Real> adjointHessian_uo_uo( const V& un, const V& z, const V& l ) {
231  return bind( &Con::applyAdjointHessian_uo_uo, &con_, ph::_1, cref(l), ph::_2,
232  ph::_3, cref(un), cref(z), ts_ );
233  }
234 
235  f_dderiv_t<Real> adjointHessian_uo_z( const V& un, const V& z, const V& l ) {
236  return bind( &Con::applyAdjointHessian_uo_z, &con_, ph::_1, cref(l), ph::_2,
237  ph::_3, cref(un), cref(z), ts_ );
238  }
239 
240  //----------------------------------------------------------------------------
241 
242  f_dderiv_t<Real> adjointHessian_z_un( const V& uo, const V& un, const V& l ) {
243  return bind( &Con::applyAdjointHessian_z_un, &con_, ph::_1, cref(l), ph::_2,
244  cref(uo), cref(un), ph::_3, ts_ );
245  }
246 
247  f_dderiv_t<Real> adjointHessian_z_uo( const V& uo, const V& un, const V& l ) {
248  return bind( &Con::applyAdjointHessian_z_uo, &con_, ph::_1, cref(l), ph::_2,
249  cref(uo), cref(un), ph::_3, ts_ );
250  }
251 
252  f_dderiv_t<Real> adjointHessian_z_z( const V& uo, const V& un, const V& l ) {
253  return bind( &Con::applyAdjointHessian_z_z, &con_, ph::_1, cref(l), ph::_2,
254  cref(uo), cref(un), ph::_3, ts_ );
255  }
256 
257 }; // class DynamicConstraint_CheckInterface
258 
259 } // namespace details
260 
261 using details::DynamicConstraint_CheckInterface;
262 
263 template<typename Real>
264 DynamicConstraint_CheckInterface<Real> make_check( DynamicConstraint<Real>& con ) {
265  return DynamicConstraint_CheckInterface<Real>(con);
266 }
267 
268 template<typename Real>
269 DynamicConstraint_CheckInterface<Real> make_check( DynamicConstraint<Real>& con,
270  TimeStamp<Real>& timeStamp ) {
271  return DynamicConstraint_CheckInterface<Real>(con,timeStamp);
272 }
273 
274 } // namespace ROL
275 
276 
277 #endif // ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
278 
279 
f_dderiv_t< Real > adjointJacobian_un_z(const V &uo, const V &un)
f_dderiv_t< Real > adjointJacobian_un(const V &uo, const V &z)
Defines the time-dependent constraint operator interface for simulation-based optimization.
f_dderiv_t< Real > adjointJacobian_uo(const V &un, const V &z)
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1) override
f_dderiv_t< Real > adjointHessian_uo_z(const V &un, const V &z, const V &l)
f_dderiv_t< Real > adjointJacobian_z_un(const V &uo, const V &z)
f_dderiv_t< Real > adjointHessian_un_uo(const V &uo, const V &z, const V &l)
f_dderiv_t< Real > adjointHessian_un_un(const V &uo, const V &z, const V &l)
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
f_dderiv_t< Real > inverseAdjointJacobian_un(const V &uo, const V &z)
f_dderiv_t< Real > adjointJacobian_un_uo(const V &un, const V &z)
f_dderiv_t< Real > adjointJacobian_uo_z(const V &uo, const V &un)
f_dderiv_t< Real > adjointHessian_z_z(const V &uo, const V &un, const V &l)
f_dderiv_t< Real > inverseJacobian_un(const V &uo, const V &z)
f_dderiv_t< Real > adjointJacobian_z(const V &uo, const V &un)
f_dderiv_t< Real > adjointHessian_uo_uo(const V &un, const V &z, const V &l)
f_dderiv_t< Real > adjointJacobian_uo_un(const V &uo, const V &z)
virtual void solve(Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z) override
f_dderiv_t< Real > adjointHessian_un_z(const V &uo, const V &z, const V &l)
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &sol, const Real &mu)
f_dderiv_t< Real > adjointJacobian_uo_uo(const V &un, const V &z)
f_dderiv_t< Real > adjointHessian_uo_un(const V &un, const V &z, const V &l)
f_dderiv_t< Real > adjointHessian_z_uo(const V &uo, const V &un, const V &l)
f_dderiv_t< Real > adjointJacobian_z_uo(const V &un, const V &z)
f_dderiv_t< Real > adjointJacobian_z_z(const V &uo, const V &un)
f_dderiv_t< Real > adjointJacobian_un_un(const V &uo, const V &z)
f_dderiv_t< Real > adjointHessian_z_un(const V &uo, const V &un, const V &l)
DynamicConstraint_CheckInterface< Real > make_check(DynamicConstraint< Real > &con)