46 #ifndef ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
47 #define ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
49 #include "ROL_FunctionBindings.hpp"
56 namespace ph = std::placeholders;
58 template<
typename Real>
74 ts_.t.at(1) = 0.02345;
82 return bind( &
Con::update, &con_, ph::_1, cref(un), cref(z), ts_ );
86 return bind( &
Con::update, &con_, cref(uo), ph::_1, cref(z), ts_ );
90 return bind( &
Con::update, &con_, cref(uo), cref(un), ph::_1, ts_ );
97 ph::_1, ph::_2, cref(un), cref(z), ts_ );
102 ph::_1, cref(uo), ph::_2, cref(z), ts_ );
107 ph::_1, cref(uo), cref(un), ph::_2, ts_ );
112 ph::_1, cref(uo), ph::_2, cref(z), ts_ );
118 return bind( &Con::applyJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
119 cref(un), cref(z), ts_ );
123 return bind( &Con::applyJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
124 ph::_3, cref(z), ts_ );
128 return bind( &Con::applyInverseJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
129 ph::_3, cref(z), ts_ );
133 return bind( &Con::applyJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
134 cref(un), ph::_3, ts_ );
140 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
141 cref(un), cref(z), ts_ );
145 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
146 ph::_3, cref(z), ts_ );
150 return bind( &Con::applyInverseAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
151 ph::_3, cref(z), ts_ );
155 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
156 cref(un), ph::_3, ts_ );
162 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
163 cref(un), cref(z), ts_ );
167 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, cref(uo),
168 ph::_3, cref(z), ts_ );
172 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, cref(uo),
173 cref(un), ph::_3, ts_ );
177 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, ph::_3,
178 cref(un), cref(z), ts_ );
182 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
183 ph::_3, cref(z), ts_ );
187 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
188 cref(un), ph::_3, ts_ );
192 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, ph::_3,
193 cref(un), cref(z), ts_ );
197 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
198 ph::_3, cref(z), ts_ );
202 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
203 cref(un), ph::_3, ts_ );
209 return bind( &Con::applyAdjointHessian_un_un, &con_, ph::_1, cref(l), ph::_2,
210 cref(uo), ph::_3, cref(z), ts_ );
214 return bind( &Con::applyAdjointHessian_un_uo, &con_, ph::_1, cref(l), ph::_2,
215 cref(uo), ph::_3, cref(z), ts_ );
219 return bind( &Con::applyAdjointHessian_un_z, &con_, ph::_1, cref(l), ph::_2,
220 cref(uo), ph::_3, cref(z), ts_ );
226 return bind( &Con::applyAdjointHessian_uo_un, &con_, ph::_1, cref(l), ph::_2,
227 ph::_3, cref(un), cref(z), ts_ );
231 return bind( &Con::applyAdjointHessian_uo_uo, &con_, ph::_1, cref(l), ph::_2,
232 ph::_3, cref(un), cref(z), ts_ );
236 return bind( &Con::applyAdjointHessian_uo_z, &con_, ph::_1, cref(l), ph::_2,
237 ph::_3, cref(un), cref(z), ts_ );
243 return bind( &Con::applyAdjointHessian_z_un, &con_, ph::_1, cref(l), ph::_2,
244 cref(uo), cref(un), ph::_3, ts_ );
248 return bind( &Con::applyAdjointHessian_z_uo, &con_, ph::_1, cref(l), ph::_2,
249 cref(uo), cref(un), ph::_3, ts_ );
253 return bind( &Con::applyAdjointHessian_z_z, &con_, ph::_1, cref(l), ph::_2,
254 cref(uo), cref(un), ph::_3, ts_ );
261 using details::DynamicConstraint_CheckInterface;
263 template<
typename Real>
265 return DynamicConstraint_CheckInterface<Real>(con);
268 template<
typename Real>
271 return DynamicConstraint_CheckInterface<Real>(con,timeStamp);
277 #endif // ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
DynamicConstraint_CheckInterface(Con &con)
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_solve_t< Real > solve_un(const V &uo, const V &z)
f_dderiv_t< Real > adjointHessian_un_un(const V &uo, const V &z, const V &l)
f_update_t< Real > update_un(const V &uo, const V &z)
Defines the linear algebra or vector space interface.
f_dderiv_t< Real > inverseAdjointJacobian_un(const V &uo, const V &z)
f_dderiv_t< Real > jacobian_z(const V &uo, const V &un)
f_vector_t< Real > value_uo(const V &un, const V &z)
f_dderiv_t< Real > adjointJacobian_un_uo(const V &un, const V &z)
f_vector_t< Real > value_un(const V &uo, const V &z)
f_dderiv_t< Real > adjointJacobian_uo_z(const V &uo, const V &un)
f_dderiv_t< Real > jacobian_uo(const V &un, const V &z)
f_dderiv_t< Real > adjointHessian_z_z(const V &uo, const V &un, const V &l)
f_vector_t< Real > value_z(const V &uo, const V &un)
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)
DynamicConstraint_CheckInterface(Con &con, TimeStamp< Real > &ts)
f_dderiv_t< Real > adjointJacobian_uo_uo(const V &un, const V &z)
f_update_t< Real > update_z(const V &uo, const V &un)
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 > jacobian_un(const V &uo, const V &z)
f_dderiv_t< Real > adjointJacobian_un_un(const V &uo, const V &z)
f_update_t< Real > update_uo(const V &un, 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)