11 #ifndef ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
12 #define ROL_DYNAMICCONSTRAINT_CHECKINTERFACE_HPP
14 #include "ROL_FunctionBindings.hpp"
21 namespace ph = std::placeholders;
23 template<
typename Real>
39 ts_.t.at(1) = 0.02345;
47 return bind( &
Con::update, &con_, ph::_1, cref(un), cref(z), ts_ );
51 return bind( &
Con::update, &con_, cref(uo), ph::_1, cref(z), ts_ );
55 return bind( &
Con::update, &con_, cref(uo), cref(un), ph::_1, ts_ );
62 ph::_1, ph::_2, cref(un), cref(z), ts_ );
67 ph::_1, cref(uo), ph::_2, cref(z), ts_ );
70 f_vector_t<Real>
value_z(
const V& uo,
const V& un ) {
72 ph::_1, cref(uo), cref(un), ph::_2, ts_ );
77 ph::_1, cref(uo), ph::_2, cref(z), ts_ );
83 return bind( &Con::applyJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
84 cref(un), cref(z), ts_ );
88 return bind( &Con::applyJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
89 ph::_3, cref(z), ts_ );
93 return bind( &Con::applyInverseJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
94 ph::_3, cref(z), ts_ );
98 return bind( &Con::applyJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
99 cref(un), ph::_3, ts_ );
105 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
106 cref(un), cref(z), ts_ );
110 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
111 ph::_3, cref(z), ts_ );
115 return bind( &Con::applyInverseAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
116 ph::_3, cref(z), ts_ );
120 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
121 cref(un), ph::_3, ts_ );
127 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, ph::_3,
128 cref(un), cref(z), ts_ );
132 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, cref(uo),
133 ph::_3, cref(z), ts_ );
137 return bind( &Con::applyAdjointJacobian_uo, &con_, ph::_1, ph::_2, cref(uo),
138 cref(un), ph::_3, ts_ );
142 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, ph::_3,
143 cref(un), cref(z), ts_ );
147 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
148 ph::_3, cref(z), ts_ );
152 return bind( &Con::applyAdjointJacobian_un, &con_, ph::_1, ph::_2, cref(uo),
153 cref(un), ph::_3, ts_ );
157 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, ph::_3,
158 cref(un), cref(z), ts_ );
162 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
163 ph::_3, cref(z), ts_ );
167 return bind( &Con::applyAdjointJacobian_z, &con_, ph::_1, ph::_2, cref(uo),
168 cref(un), ph::_3, ts_ );
174 return bind( &Con::applyAdjointHessian_un_un, &con_, ph::_1, cref(l), ph::_2,
175 cref(uo), ph::_3, cref(z), ts_ );
179 return bind( &Con::applyAdjointHessian_un_uo, &con_, ph::_1, cref(l), ph::_2,
180 cref(uo), ph::_3, cref(z), ts_ );
184 return bind( &Con::applyAdjointHessian_un_z, &con_, ph::_1, cref(l), ph::_2,
185 cref(uo), ph::_3, cref(z), ts_ );
191 return bind( &Con::applyAdjointHessian_uo_un, &con_, ph::_1, cref(l), ph::_2,
192 ph::_3, cref(un), cref(z), ts_ );
196 return bind( &Con::applyAdjointHessian_uo_uo, &con_, ph::_1, cref(l), ph::_2,
197 ph::_3, cref(un), cref(z), ts_ );
201 return bind( &Con::applyAdjointHessian_uo_z, &con_, ph::_1, cref(l), ph::_2,
202 ph::_3, cref(un), cref(z), ts_ );
208 return bind( &Con::applyAdjointHessian_z_un, &con_, ph::_1, cref(l), ph::_2,
209 cref(uo), cref(un), ph::_3, ts_ );
213 return bind( &Con::applyAdjointHessian_z_uo, &con_, ph::_1, cref(l), ph::_2,
214 cref(uo), cref(un), ph::_3, ts_ );
218 return bind( &Con::applyAdjointHessian_z_z, &con_, ph::_1, cref(l), ph::_2,
219 cref(uo), cref(un), ph::_3, ts_ );
226 using details::DynamicConstraint_CheckInterface;
228 template<
typename Real>
230 return DynamicConstraint_CheckInterface<Real>(con);
233 template<
typename Real>
236 return DynamicConstraint_CheckInterface<Real>(con,timeStamp);
242 #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)