45 #ifndef ROL_REDUCED_OBJECTIVE_SIMOPT_H 46 #define ROL_REDUCED_OBJECTIVE_SIMOPT_H 56 Teuchos::RCP<Objective_SimOpt<Real> >
obj_;
57 Teuchos::RCP<EqualityConstraint_SimOpt<Real> >
con_;
81 if ( state_storage_.count(this->getParameter()) && storage_ ) {
86 con_->update_2(x,flag,iter);
88 con_->solve(*dualadjoint_,*state_,x,tol);
90 con_->update_1(*state_,flag,iter);
92 obj_->update(*state_,x,flag,iter);
95 state_storage_.insert(
96 std::pair<std::vector<Real>,Teuchos::RCP<
Vector<Real> > >(
109 if ( adjoint_storage_.count(this->getParameter()) && storage_ ) {
114 obj_->gradient_1(*dualstate_,*state_,x,tol);
116 con_->applyInverseAdjointJacobian_1(*adjoint_,*dualstate_,*state_,x,tol);
117 adjoint_->scale(-1.0);
120 adjoint_storage_.insert(
121 std::pair<std::vector<Real>,Teuchos::RCP<
Vector<Real> > >(
134 con_->applyJacobian_2(*dualadjoint_,v,*state_,x,tol);
135 dualadjoint_->scale(-1.0);
136 con_->applyInverseJacobian_1(*state_sens_,*dualadjoint_,*state_,x,tol);
148 obj_->hessVec_11(*dualstate_,*state_sens_,*state_,x,tol);
149 obj_->hessVec_12(*dualstate1_,v,*state_,x,tol);
150 dualstate_->plus(*dualstate1_);
152 con_->applyAdjointHessian_11(*dualstate1_,*adjoint_,*state_sens_,*state_,x,tol);
153 dualstate_->plus(*dualstate1_);
154 con_->applyAdjointHessian_21(*dualstate1_,*adjoint_,v,*state_,x,tol);
155 dualstate_->plus(*dualstate1_);
157 dualstate_->scale(-1.0);
158 con_->applyInverseAdjointJacobian_1(*adjoint_sens_,*dualstate_,*state_,x,tol);
174 const bool storage =
true,
175 const bool useFDhessVec =
false)
176 : obj_(obj), con_(con),
177 state_(state), state_sens_(state->clone()),
178 adjoint_(adjoint), adjoint_sens_(adjoint->clone()),
179 dualstate_(state->dual().clone()), dualstate1_(state->dual().clone()),
180 dualadjoint_(adjoint->dual().clone()), dualcontrol_(Teuchos::null),
181 storage_(storage), useFDhessVec_(useFDhessVec), is_initialized_(false) {
182 state_storage_.clear();
183 adjoint_storage_.clear();
204 const bool storage =
true,
205 const bool useFDhessVec =
false)
206 : obj_(obj), con_(con),
207 state_(state), state_sens_(state->clone()),
208 adjoint_(adjoint), adjoint_sens_(adjoint->clone()),
209 dualstate_(dualstate), dualstate1_(dualstate->clone()),
210 dualadjoint_(dualadjoint->clone()), dualcontrol_(Teuchos::null),
211 storage_(storage), useFDhessVec_(useFDhessVec), is_initialized_(false) {
212 state_storage_.clear();
213 adjoint_storage_.clear();
220 state_storage_.clear();
222 adjoint_storage_.clear();
234 return obj_->value(*state_,x,tol);
243 if (!is_initialized_) {
244 dualcontrol_ = g.
clone();
245 is_initialized_ =
true;
252 obj_->gradient_2(*dualcontrol_,*state_,x,tol);
254 con_->applyAdjointJacobian_2(g,*adjoint_,*state_,x,tol);
255 g.
plus(*dualcontrol_);
262 if ( useFDhessVec_ ) {
266 if (!is_initialized_) {
267 dualcontrol_ = hv.
clone();
268 is_initialized_ =
true;
279 con_->applyAdjointJacobian_2(hv,*adjoint_sens_,*state_,x,tol);
280 obj_->hessVec_21(*dualcontrol_,*state_sens_,*state_,x,tol);
281 hv.
plus(*dualcontrol_);
282 obj_->hessVec_22(*dualcontrol_,v,*state_,x,tol);
283 hv.
plus(*dualcontrol_);
284 con_->applyAdjointHessian_12(*dualcontrol_,*adjoint_,*state_sens_,*state_,x,tol);
285 hv.
plus(*dualcontrol_);
286 con_->applyAdjointHessian_22(*dualcontrol_,*adjoint_,v,*state_,x,tol);
287 hv.
plus(*dualcontrol_);
301 con_->setParameter(param);
302 obj_->setParameter(param);
Teuchos::RCP< Objective_SimOpt< Real > > obj_
Provides the interface to evaluate objective functions.
Provides the interface to evaluate simulation-based objective functions.
virtual void plus(const Vector &x)=0
Compute , where .
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply a reduced Hessian preconditioner.
void solve_state_equation(const Vector< Real > &x, Real &tol, bool flag=true, int iter=-1)
Reduced_Objective_SimOpt(const Teuchos::RCP< Objective_SimOpt< Real > > &obj, const Teuchos::RCP< EqualityConstraint_SimOpt< Real > > &con, const Teuchos::RCP< Vector< Real > > &state, const Teuchos::RCP< Vector< Real > > &adjoint, const bool storage=true, const bool useFDhessVec=false)
Constructor.
void setParameter(const std::vector< Real > ¶m)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Teuchos::RCP< Vector< Real > > state_
Teuchos::RCP< Vector< Real > > adjoint_sens_
Teuchos::RCP< Vector< Real > > state_sens_
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Given , evaluate the gradient of the objective function where solves .
Defines the linear algebra or vector space interface.
Defines the equality constraint operator interface for simulation-based optimization.
Teuchos::RCP< Vector< Real > > adjoint_
Teuchos::RCP< Vector< Real > > dualstate_
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Given , evaluate the Hessian of the objective function in the direction .
std::map< std::vector< Real >, Teuchos::RCP< Vector< Real > > > adjoint_storage_
std::map< std::vector< Real >, Teuchos::RCP< Vector< Real > > > state_storage_
Real value(const Vector< Real > &x, Real &tol)
Given , evaluate the objective function where solves .
const std::vector< Real > getParameter(void) const
virtual void setParameter(const std::vector< Real > ¶m)
Teuchos::RCP< Vector< Real > > dualcontrol_
void solve_adjoint_sensitivity(const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Given , the adjoint variable , and a direction , solve the adjoint sensitvity equation for ...
virtual void set(const Vector &x)
Set where .
void solve_state_sensitivity(const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Given which solves the state equation and a direction , solve the state senstivity equation for ...
Teuchos::RCP< EqualityConstraint_SimOpt< Real > > con_
Teuchos::RCP< Vector< Real > > dualstate1_
Teuchos::RCP< Vector< Real > > dualadjoint_
Reduced_Objective_SimOpt(const Teuchos::RCP< Objective_SimOpt< Real > > &obj, const Teuchos::RCP< EqualityConstraint_SimOpt< Real > > &con, const Teuchos::RCP< Vector< Real > > &state, const Teuchos::RCP< Vector< Real > > &adjoint, const Teuchos::RCP< Vector< Real > > &dualstate, const Teuchos::RCP< Vector< Real > > &dualadjoint, const bool storage=true, const bool useFDhessVec=false)
Secondary, general constructor for use with dual optimization vector spaces where the user does not d...
void solve_adjoint_equation(const Vector< Real > &x, Real &tol)
Given which solves the state equation, solve the adjoint equation for .
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update the SimOpt objective function and equality constraint.