44#ifndef ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
45#define ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
49template<
typename Real>
59 bool isConRedParametrized)
61 updateFlag_(true), newUpdate_(false), updateIter_(0), updateType_(
UpdateType::
Initial),
62 storage_(storage), isConRedParametrized_(isConRedParametrized) {
70 stateStore_ = ROL::makePtr<VectorController<Real>>();
73template<
typename Real>
76 bool flag,
int iter) {
77 update_1(u,flag,iter);
81template<
typename Real>
83 bool flag,
int iter) {
85 conVal_->update_1(u, flag, iter);
88template<
typename Real>
90 bool flag,
int iter) {
91 conRed_->update_2(z,flag,iter);
94 stateStore_->constraintUpdate(
true);
97template<
typename Real>
101 update_1(u,type,iter);
105template<
typename Real>
109 conVal_->update_1(u,type,iter);
112template<
typename Real>
115 conRed_->update_2(z,type,iter);
118 stateStore_->constraintUpdate(type);
121template<
typename Real>
125 Real ctol(std::sqrt(ROL_EPSILON<Real>()));
126 solveConRed(*Sz_, z, ctol);
127 conVal_->solve_update(u,*Sz_,type,iter);
130template<
typename Real>
135 solveConRed(*Sz_, z, tol);
136 conVal_->value(c, u, *Sz_, tol);
139template<
typename Real>
144 solveConRed(*Sz_, z, tol);
145 conVal_->solve(c, u, *Sz_, tol);
148template<
typename Real>
154 solveConRed(*Sz_, z, tol);
155 conVal_->applyJacobian_1(jv, v, u, *Sz_, tol);
158template<
typename Real>
164 solveConRed(*Sz_, z, tol);
165 applySens(*primZ_, v, *Sz_, z, tol);
166 conVal_->applyJacobian_2(jv, *primZ_, u, *Sz_, tol);
169template<
typename Real>
175 solveConRed(*Sz_, z, tol);
176 conVal_->applyInverseJacobian_1(ijv, v, u, *Sz_, tol);
179template<
typename Real>
185 solveConRed(*Sz_, z, tol);
186 conVal_->applyAdjointJacobian_1(ajv, v, u, *Sz_, tol);
189template<
typename Real>
195 solveConRed(*Sz_, z, tol);
196 conVal_->applyAdjointJacobian_2(*dualZ_, v, u, *Sz_, tol);
197 applyAdjointSens(ajv, *dualZ_, *Sz_, z, tol);
200template<
typename Real>
206 solveConRed(*Sz_, z, tol);
207 conVal_->applyInverseAdjointJacobian_1(ijv, v, u, *Sz_, tol);
210template<
typename Real>
217 solveConRed(*Sz_, z, tol);
218 conVal_->applyAdjointHessian_11(ahwv, w, v, u, *Sz_, tol);
221template<
typename Real>
228 solveConRed(*Sz_, z, tol);
229 conVal_->applyAdjointHessian_12(*dualZ_, w, v, u, *Sz_, tol);
230 applyAdjointSens(ahwv, *dualZ_, *Sz_, z, tol);
233template<
typename Real>
240 solveConRed(*Sz_, z, tol);
241 applySens(*primZ_, v, *Sz_, z, tol);
242 conVal_->applyAdjointHessian_21(ahwv, w, *primZ_, u, *Sz_, tol);
245template<
typename Real>
253 solveConRed(*Sz_, z, tol);
254 applySens(*primZ_, v, *Sz_, z, tol);
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_);
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));
268 conVal_->applyAdjointHessian_22(*dualZ_, w, *primZ_, u, *Sz_, tol);
269 dualZ1_->plus(*dualZ_);
271 applyAdjointSens(*dualZ_, *dualZ1_, *Sz_, z, tol);
275template<
typename Real>
277 conVal_->setParameter(param);
278 if (isConRedParametrized_) {
279 conRed_->setParameter(param);
284template<
typename Real>
290 bool isComputed =
false;
291 if (storage_) isComputed = stateStore_->get(Sz,param);
293 if (!isComputed || !storage_) {
295 conRed_->solve(*primRed_,Sz,z,tol);
297 if (newUpdate_) conRed_->update_1(Sz,updateType_,updateIter_);
298 else conRed_->update_1(Sz,updateFlag_,updateIter_);
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_);
305 if (storage_) stateStore_->set(Sz,param);
309template<
typename Real>
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));
321template<
typename Real>
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));
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,...