10#ifndef ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
11#define ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
15template<
typename Real>
25 bool isConRedParametrized)
27 updateFlag_(true), newUpdate_(false), updateIter_(0), updateType_(
UpdateType::
Initial),
28 storage_(storage), isConRedParametrized_(isConRedParametrized) {
36 stateStore_ = ROL::makePtr<VectorController<Real>>();
39template<
typename Real>
42 bool flag,
int iter) {
43 update_1(u,flag,iter);
47template<
typename Real>
49 bool flag,
int iter) {
51 conVal_->update_1(u, flag, iter);
54template<
typename Real>
56 bool flag,
int iter) {
57 conRed_->update_2(z,flag,iter);
60 stateStore_->constraintUpdate(
true);
63template<
typename Real>
67 update_1(u,type,iter);
71template<
typename Real>
75 conVal_->update_1(u,type,iter);
78template<
typename Real>
81 conRed_->update_2(z,type,iter);
84 stateStore_->constraintUpdate(type);
87template<
typename Real>
91 Real ctol(std::sqrt(ROL_EPSILON<Real>()));
92 solveConRed(*Sz_, z, ctol);
93 conVal_->solve_update(u,*Sz_,type,iter);
96template<
typename Real>
101 solveConRed(*Sz_, z, tol);
102 conVal_->value(c, u, *Sz_, tol);
105template<
typename Real>
110 solveConRed(*Sz_, z, tol);
111 conVal_->solve(c, u, *Sz_, tol);
114template<
typename Real>
120 solveConRed(*Sz_, z, tol);
121 conVal_->applyJacobian_1(jv, v, u, *Sz_, tol);
124template<
typename Real>
130 solveConRed(*Sz_, z, tol);
131 applySens(*primZ_, v, *Sz_, z, tol);
132 conVal_->applyJacobian_2(jv, *primZ_, u, *Sz_, tol);
135template<
typename Real>
141 solveConRed(*Sz_, z, tol);
142 conVal_->applyInverseJacobian_1(ijv, v, u, *Sz_, tol);
145template<
typename Real>
151 solveConRed(*Sz_, z, tol);
152 conVal_->applyAdjointJacobian_1(ajv, v, u, *Sz_, tol);
155template<
typename Real>
161 solveConRed(*Sz_, z, tol);
162 conVal_->applyAdjointJacobian_2(*dualZ_, v, u, *Sz_, tol);
163 applyAdjointSens(ajv, *dualZ_, *Sz_, z, tol);
166template<
typename Real>
172 solveConRed(*Sz_, z, tol);
173 conVal_->applyInverseAdjointJacobian_1(ijv, v, u, *Sz_, tol);
176template<
typename Real>
183 solveConRed(*Sz_, z, tol);
184 conVal_->applyAdjointHessian_11(ahwv, w, v, u, *Sz_, tol);
187template<
typename Real>
194 solveConRed(*Sz_, z, tol);
195 conVal_->applyAdjointHessian_12(*dualZ_, w, v, u, *Sz_, tol);
196 applyAdjointSens(ahwv, *dualZ_, *Sz_, z, tol);
199template<
typename Real>
206 solveConRed(*Sz_, z, tol);
207 applySens(*primZ_, v, *Sz_, z, tol);
208 conVal_->applyAdjointHessian_21(ahwv, w, *primZ_, u, *Sz_, tol);
211template<
typename Real>
219 solveConRed(*Sz_, z, tol);
220 applySens(*primZ_, v, *Sz_, z, tol);
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_);
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));
234 conVal_->applyAdjointHessian_22(*dualZ_, w, *primZ_, u, *Sz_, tol);
235 dualZ1_->plus(*dualZ_);
237 applyAdjointSens(*dualZ_, *dualZ1_, *Sz_, z, tol);
241template<
typename Real>
243 conVal_->setParameter(param);
244 if (isConRedParametrized_) {
245 conRed_->setParameter(param);
250template<
typename Real>
256 bool isComputed =
false;
257 if (storage_) isComputed = stateStore_->get(Sz,param);
259 if (!isComputed || !storage_) {
261 conRed_->solve(*primRed_,Sz,z,tol);
263 if (newUpdate_) conRed_->update_1(Sz,updateType_,updateIter_);
264 else conRed_->update_1(Sz,updateFlag_,updateIter_);
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_);
271 if (storage_) stateStore_->set(Sz,param);
275template<
typename Real>
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));
287template<
typename Real>
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));
ROL::Ptr< Vector< Real > > primU_
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 .
ROL::Ptr< Vector< Real > > primRed_
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 ...
ROL::Ptr< Vector< Real > > dualRed_
ROL::Ptr< Vector< Real > > dualZ1_
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...
ROL::Ptr< Vector< Real > > Sz_
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 > ¶m) 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...
ROL::Ptr< Vector< Real > > dualZ_
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...
ROL::Ptr< Vector< Real > > primZ_
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 > ¶m)
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,...