11#ifndef ROL_REDUCED_CONSTRAINT_SIMOPT_DEF_H
12#define ROL_REDUCED_CONSTRAINT_SIMOPT_DEF_H
20 if (newUpdate_) conRed_->update_2(z,updateType_,updateIter_);
21 else conRed_->update_2(z,updateFlag_,updateIter_);
26 if (!isComputed || !storage_) {
28 conRed_->solve(*dualadjoint_,*state_,z,tol);
35 if (newUpdate_) conRed_->update_1(*state_,updateType_,updateIter_);
36 else conRed_->update_1(*state_,updateFlag_,updateIter_);
38 if (newUpdate_) conVal_->update(*state_,z,updateType_,updateIter_);
39 else conVal_->update(*state_,z,updateFlag_,updateIter_);
49 if (!isComputed || !storage_) {
51 conVal_->applyAdjointJacobian_1(*dualstate_,w,*state_,z,tol);
53 conRed_->applyInverseAdjointJacobian_1(*adjoint_,*dualstate_,*state_,z,tol);
54 adjoint_->scale(
static_cast<Real
>(-1));
64 conRed_->applyJacobian_2(*dualadjoint_,v,*state_,z,tol);
65 dualadjoint_->scale(
static_cast<Real
>(-1));
66 conRed_->applyInverseJacobian_1(*state_sens_,*dualadjoint_,*state_,z,tol);
73 conVal_->applyAdjointHessian_11(*dualstate_,w,*state_sens_,*state_,z,tol);
74 conVal_->applyAdjointHessian_21(*dualstate1_,w,v,*state_,z,tol);
75 dualstate_->plus(*dualstate1_);
77 conRed_->applyAdjointHessian_11(*dualstate1_,*adjoint_,*state_sens_,*state_,z,tol);
78 dualstate_->plus(*dualstate1_);
79 conRed_->applyAdjointHessian_21(*dualstate1_,*adjoint_,v,*state_,z,tol);
80 dualstate_->plus(*dualstate1_);
82 dualstate_->scale(
static_cast<Real
>(-1));
83 conRed_->applyInverseAdjointJacobian_1(*adjoint_sens_,*dualstate_,*state_,z,tol);
100 stateStore_( stateStore ),
102 state_( state->clone() ),
103 adjoint_( adjoint->clone() ),
104 residual_( residual->clone() ),
105 state_sens_( state->clone() ),
106 adjoint_sens_( adjoint->clone() ),
107 dualstate_( state->dual().clone() ),
108 dualstate1_( state->dual().clone() ),
109 dualadjoint_( adjoint->dual().clone() ),
110 dualcontrol_( control->dual().clone() ),
111 dualresidual_( residual->dual().clone() ),
112 storage_(storage), useFDhessVec_(useFDhessVec),
113 nupda_(0), nvalu_(0), njaco_(0), najac_(0), nhess_(0),
114 nstat_(0), nadjo_(0), nssen_(0), nasen_(0),
116 newUpdate_(false), isUpdated_(true) {}
135 stateStore_( stateStore ),
137 state_( state->clone() ),
138 adjoint_( adjoint->clone() ),
139 residual_( residual->clone() ),
140 state_sens_( state->clone() ),
141 adjoint_sens_( adjoint->clone() ),
142 dualstate_( dualstate->clone() ),
143 dualstate1_( dualstate->clone() ),
144 dualadjoint_( dualadjoint->clone() ),
145 dualcontrol_( dualcontrol->clone() ),
146 dualresidual_( dualresidual->clone() ),
147 storage_(storage), useFDhessVec_(useFDhessVec),
148 nupda_(0), nvalu_(0), njaco_(0), najac_(0), nhess_(0),
149 nstat_(0), nadjo_(0), nssen_(0), nasen_(0),
151 newUpdate_(false), isUpdated_(true) {}
155 int nupda(0), nvalu(0), njaco(0), najac(0), nhess(0), nstat(0), nadjo(0), nssen(0), nasen(0);
156 if (bman == nullPtr) {
168 auto sumAll = [bman](
int val) {
169 Real global(0), local(val);
170 bman->sumAll(&local,&global,1);
171 return static_cast<int>(global);
173 nupda = sumAll(nupda_);
174 nvalu = sumAll(nvalu_);
175 njaco = sumAll(njaco_);
176 najac = sumAll(najac_);
177 nhess = sumAll(nhess_);
178 nstat = sumAll(nstat_);
179 nadjo = sumAll(nadjo_);
180 nssen = sumAll(nssen_);
181 nasen = sumAll(nasen_);
184 stream << std::string(80,
'=') << std::endl;
185 stream <<
" ROL::Reduced_Objective_SimOpt::summarize" << std::endl;
186 stream <<
" Number of calls to update: " << nupda << std::endl;
187 stream <<
" Number of calls to value: " << nvalu << std::endl;
188 stream <<
" Number of calls to applyJacobian: " << njaco << std::endl;
189 stream <<
" Number of calls to applyAdjointJacobian: " << najac << std::endl;
190 stream <<
" Number of calls to hessvec: " << nhess << std::endl;
191 stream <<
" Number of state solves: " << nstat << std::endl;
192 stream <<
" Number of adjoint solves: " << nadjo << std::endl;
193 stream <<
" Number of state sensitivity solves: " << nssen << std::endl;
194 stream <<
" Number of adjoint sensitivity solves: " << nasen << std::endl;
195 stream << std::string(80,
'=') << std::endl;
201 nupda_ = 0; nvalu_ = 0; njaco_ = 0; najac_ = 0; nhess_ = 0;
202 nstat_ = 0; nadjo_ = 0; nssen_ = 0; nasen_ = 0;
210 stateStore_->constraintUpdate(
true);
211 adjointStore_->constraintUpdate(flag);
220 stateStore_->objectiveUpdate(type);
221 adjointStore_->objectiveUpdate(type);
228 solve_state_equation(z,tol);
230 conVal_->value(c,*state_,z,tol);
238 solve_state_equation(z,tol);
240 solve_state_sensitivity(v,z,tol);
242 conVal_->applyJacobian_1(*residual_,*state_sens_,*state_,z,tol);
244 conVal_->applyJacobian_2(jv,v,*state_,z,tol);
253 solve_state_equation(z,tol);
255 solve_adjoint_equation(w,z,tol);
257 conVal_->applyAdjointJacobian_2(*dualcontrol_,w,*state_,z,tol);
259 conRed_->applyAdjointJacobian_2(ajw,*adjoint_,*state_,z,tol);
260 ajw.
plus(*dualcontrol_);
268 if ( useFDhessVec_ ) {
273 solve_state_equation(z,tol);
275 solve_adjoint_equation(w,z,tol);
277 solve_state_sensitivity(v,z,tol);
279 solve_adjoint_sensitivity(w,v,z,tol);
281 conRed_->applyAdjointJacobian_2(ahwv,*adjoint_sens_,*state_,z,tol);
282 conVal_->applyAdjointHessian_12(*dualcontrol_,w,*state_sens_,*state_,z,tol);
283 ahwv.
plus(*dualcontrol_);
284 conVal_->applyAdjointHessian_22(*dualcontrol_,w,v,*state_,z,tol);
285 ahwv.
plus(*dualcontrol_);
286 conRed_->applyAdjointHessian_12(*dualcontrol_,*adjoint_,*state_sens_,*state_,z,tol);
287 ahwv.
plus(*dualcontrol_);
288 conRed_->applyAdjointHessian_22(*dualcontrol_,*adjoint_,v,*state_,z,tol);
289 ahwv.
plus(*dualcontrol_);
Defines the constraint operator interface for simulation-based optimization.
Defines the general constraint operator interface.
virtual void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the derivative of the adjoint of the constraint Jacobian at to vector in direction ,...
Reduced_Constraint_SimOpt(const ROL::Ptr< Constraint_SimOpt< Real > > &conVal, const ROL::Ptr< Constraint_SimOpt< Real > > &conRed, const ROL::Ptr< VectorController< Real > > &stateStore, const ROL::Ptr< Vector< Real > > &state, const ROL::Ptr< Vector< Real > > &control, const ROL::Ptr< Vector< Real > > &adjoint, const ROL::Ptr< Vector< Real > > &residual, bool storage=true, bool useFDhessVec=false)
Constructor.
void solve_adjoint_sensitivity(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given , the adjoint variable , and a direction , solve the adjoint sensitvity equation for .
void solve_state_sensitivity(const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given which solves the state equation and a direction , solve the state senstivity equation for .
void solve_state_equation(const Vector< Real > &z, Real &tol)
void value(Vector< Real > &c, const Vector< Real > &z, Real &tol)
Given , evaluate the equality constraint where solves .
void summarize(std::ostream &stream, const Ptr< BatchManager< Real > > &bman=nullPtr) const
void applyAdjointHessian(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given , evaluate the Hessian of the objective function in the direction .
void update(const Vector< Real > &z, bool flag=true, int iter=-1)
Update the SimOpt objective function and equality constraint.
void applyAdjointJacobian(Vector< Real > &ajw, const Vector< Real > &w, const Vector< Real > &z, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector .
void solve_adjoint_equation(const Vector< Real > &w, const Vector< Real > &z, Real &tol)
Given which solves the state equation, solve the adjoint equation for .
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
Given , apply the Jacobian to a vector where solves .
Defines the linear algebra or vector space interface.
virtual void plus(const Vector &x)=0
Compute , where .