44#ifndef ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
45#define ROL_COMPOSITECONSTRAINT_SIMOPT_DEF_H
49template<
typename Real>
59 bool isConRedParametrized)
70 stateStore_ = ROL::makePtr<VectorController<Real>>();
73template<
typename Real>
76 bool flag,
int iter) {
81template<
typename Real>
83 bool flag,
int iter) {
85 conVal_->update_1(u, flag, iter);
88template<
typename Real>
90 bool flag,
int iter) {
97template<
typename Real>
105template<
typename Real>
109 conVal_->update_1(u,type,iter);
112template<
typename Real>
115 conRed_->update_2(z,type,iter);
121template<
typename Real>
130template<
typename Real>
139template<
typename Real>
148template<
typename Real>
155 conVal_->applyJacobian_1(jv, v, u, *
Sz_, tol);
158template<
typename Real>
169template<
typename Real>
176 conVal_->applyInverseJacobian_1(ijv, v, u, *
Sz_, tol);
179template<
typename Real>
186 conVal_->applyAdjointJacobian_1(ajv, v, u, *
Sz_, tol);
189template<
typename Real>
200template<
typename Real>
207 conVal_->applyInverseAdjointJacobian_1(ijv, v, u, *
Sz_, tol);
210template<
typename Real>
218 conVal_->applyAdjointHessian_11(ahwv, w, v, u, *
Sz_, tol);
221template<
typename Real>
233template<
typename Real>
245template<
typename Real>
266 dualZ1_->scale(
static_cast<Real
>(-1));
275template<
typename Real>
284template<
typename Real>
290 bool isComputed =
false;
309template<
typename Real>
318 jv.
scale(
static_cast<Real
>(-1));
321template<
typename Real>
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
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
void solve(Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z, Real &tol) override
void update_2(const Vector< Real > &z, bool flag=true, int iter=-1) override
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
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
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
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
ROL::Ptr< VectorController< Real > > stateStore_
void update_1(const Vector< Real > &u, bool flag=true, int iter=-1) override
void value(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
void applyAdjointJacobian_2(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
ROL::Ptr< Vector< Real > > Sz_
void applySens(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
const ROL::Ptr< Constraint_SimOpt< Real > > conVal_
const bool isConRedParametrized_
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
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
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
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
void applyJacobian_1(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
const ROL::Ptr< Constraint_SimOpt< Real > > conRed_
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,...
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.