44 #ifndef ROL_COMPOSITEOBJECTIVE_SIMOPT_H 45 #define ROL_COMPOSITEOBJECTIVE_SIMOPT_H 61 const std::vector<Teuchos::RCP<Objective_SimOpt<Real> > >
obj_vec_;
62 const Teuchos::RCP<StdObjective<Real> >
std_obj_;
82 int size = obj_vec_.size();
83 vec_grad1_.clear(); vec_grad1_.resize(size,Teuchos::null);
84 vec_grad2_.clear(); vec_grad2_.resize(size,Teuchos::null);
85 vec_hess1_.clear(); vec_hess1_.resize(size,Teuchos::null);
86 vec_hess2_.clear(); vec_hess2_.resize(size,Teuchos::null);
87 for (
int i = 0; i < size; ++i) {
88 vec_grad1_[i] = u.
dual().clone();
89 vec_grad2_[i] = z.
dual().clone();
90 vec_hess1_[i] = u.
dual().clone();
91 vec_hess2_[i] = z.
dual().clone();
93 isInitialized_ =
true;
99 if (!isValueComputed_) {
100 int size = obj_vec_.size();
101 for (
int i = 0; i < size; ++i) {
102 (*obj_value_)[i] = obj_vec_[i]->value(u,z,tol);
104 isValueComputed_ =
true;
110 if (!isGradientComputed_) {
111 std_obj_->gradient(*(obj_grad_vec_),*(obj_value_vec_),tol);
112 isGradientComputed_ =
true;
118 if (!isGradient1Computed_) {
119 int size = obj_vec_.size();
120 for (
int i = 0; i < size; ++i) {
121 obj_vec_[i]->gradient_1(*(vec_grad1_[i]),u,z,tol);
123 isGradient1Computed_ =
true;
129 if (!isGradient2Computed_) {
130 int size = obj_vec_.size();
131 for (
int i = 0; i < size; ++i) {
132 obj_vec_[i]->gradient_2(*(vec_grad2_[i]),u,z,tol);
134 isGradient2Computed_ =
true;
140 int size = obj_vec_.size();
141 for (
int i = 0; i < size; ++i) {
142 (*obj_gv_)[i] = vec_grad1_[i]->dot(v.
dual());
143 obj_vec_[i]->hessVec_11(*(vec_hess1_[i]),v,u,z,tol);
145 std_obj_->hessVec(*(obj_hess_vec_),*(obj_gv_vec_),*(obj_value_vec_),tol);
151 int size = obj_vec_.size();
152 for (
int i = 0; i < size; ++i) {
153 (*obj_gv_)[i] = vec_grad2_[i]->dot(v.
dual());
154 obj_vec_[i]->hessVec_12(*(vec_hess1_[i]),v,u,z,tol);
156 std_obj_->hessVec(*(obj_hess_vec_),*(obj_gv_vec_),*(obj_value_vec_),tol);
162 int size = obj_vec_.size();
163 for (
int i = 0; i < size; ++i) {
164 (*obj_gv_)[i] = vec_grad1_[i]->dot(v.
dual());
165 obj_vec_[i]->hessVec_21(*(vec_hess2_[i]),v,u,z,tol);
167 std_obj_->hessVec(*(obj_hess_vec_),*(obj_gv_vec_),*(obj_value_vec_),tol);
172 int size = obj_vec_.size();
173 for (
int i = 0; i < size; ++i) {
174 (*obj_gv_)[i] = vec_grad2_[i]->dot(v.
dual());
175 obj_vec_[i]->hessVec_22(*(vec_hess2_[i]),v,u,z,tol);
177 std_obj_->hessVec(*(obj_hess_vec_),*(obj_gv_vec_),*(obj_value_vec_),tol);
183 : obj_vec_(obj_vec), std_obj_(std_obj),
184 isInitialized_(false), isValueComputed_(false),
185 isGradientComputed_(false), isGradient1Computed_(false), isGradient2Computed_(false) {
186 obj_value_ = Teuchos::rcp(
new std::vector<Real>(obj_vec_.size(),0));
188 obj_grad_ = Teuchos::rcp(
new std::vector<Real>(obj_vec_.size(),0));
190 obj_gv_ = Teuchos::rcp(
new std::vector<Real>(obj_vec_.size(),0));
192 obj_hess_ = Teuchos::rcp(
new std::vector<Real>(obj_vec_.size(),0));
197 int size = obj_vec_.size();
198 for (
int i = 0; i < size; ++i) {
199 obj_vec_[i]->update(u,z,flag,iter);
201 isValueComputed_ =
false;
209 return std_obj_->value(*obj_value_vec_,tol);
216 int size = obj_vec_.size();
217 for (
int i = 0; i < size; ++i) {
218 g.
axpy((*obj_grad_)[i],*(vec_grad1_[i]));
225 int size = obj_vec_.size();
226 for (
int i = 0; i < size; ++i) {
227 g.
axpy((*obj_grad_)[i],*(vec_grad2_[i]));
235 int size = obj_vec_.size();
236 for (
int i = 0; i < size; ++i) {
237 hv.
axpy((*obj_grad_)[i],*(vec_hess1_[i]));
238 hv.
axpy((*obj_hess_)[i],*(vec_grad1_[i]));
246 int size = obj_vec_.size();
247 for (
int i = 0; i < size; ++i) {
248 hv.
axpy((*obj_grad_)[i],*(vec_hess1_[i]));
249 hv.
axpy((*obj_hess_)[i],*(vec_grad1_[i]));
257 int size = obj_vec_.size();
258 for (
int i = 0; i < size; ++i) {
259 hv.
axpy((*obj_grad_)[i],*(vec_hess2_[i]));
260 hv.
axpy((*obj_hess_)[i],*(vec_grad2_[i]));
268 int size = obj_vec_.size();
269 for (
int i = 0; i < size; ++i) {
270 hv.
axpy((*obj_grad_)[i],*(vec_hess2_[i]));
271 hv.
axpy((*obj_hess_)[i],*(vec_grad2_[i]));
279 const int size = obj_vec_.size();
280 for (
int i = 0; i < size; ++i) {
281 obj_vec_[i]->setParameter(param);
283 std_obj_->setParameter(param);
284 isValueComputed_ =
false;
285 isGradientComputed_ =
false;
286 isGradient1Computed_ =
false;
287 isGradient2Computed_ =
false;
Provides the interface to evaluate simulation-based objective functions.
bool isGradient2Computed_
Teuchos::RCP< StdVector< Real > > obj_grad_vec_
void computeGradient(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
void hessVec_11(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
Provides the interface to evaluate simulation-based composite objective functions.
const std::vector< Teuchos::RCP< Objective_SimOpt< Real > > > obj_vec_
virtual void zero()
Set to zero vector.
void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
Update objective function. u is an iterate, z is an iterate, flag = true if the iterate has changed...
Defines the linear algebra or vector space interface.
void computeGradient1(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
std::vector< Teuchos::RCP< Vector< Real > > > vec_hess1_
Teuchos::RCP< std::vector< Real > > obj_hess_
void computeGradient2(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Specializes the ROL::Objective interface for objective functions that operate on ROL::StdVector's.
void hessVec_21(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void computeHessVec21(const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
std::vector< Teuchos::RCP< Vector< Real > > > vec_hess2_
void computeHessVec12(const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
std::vector< Teuchos::RCP< Vector< Real > > > vec_grad2_
Teuchos::RCP< StdVector< Real > > obj_value_vec_
Teuchos::RCP< StdVector< Real > > obj_hess_vec_
void setParameter(const std::vector< Real > ¶m)
void gradient_1(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
const Teuchos::RCP< StdObjective< Real > > std_obj_
void gradient_2(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.
void hessVec_12(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
void initialize(const Vector< Real > &u, const Vector< Real > &z)
Real value(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute value.
virtual void setParameter(const std::vector< Real > ¶m)
std::vector< Teuchos::RCP< Vector< Real > > > vec_grad1_
void hessVec_22(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Teuchos::RCP< std::vector< Real > > obj_grad_
void computeValue(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
void computeHessVec11(const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Teuchos::RCP< std::vector< Real > > obj_gv_
Teuchos::RCP< StdVector< Real > > obj_gv_vec_
CompositeObjective_SimOpt(const std::vector< Teuchos::RCP< Objective_SimOpt< Real > > > &obj_vec, const Teuchos::RCP< StdObjective< Real > > &std_obj)
bool isGradient1Computed_
void computeHessVec22(const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Teuchos::RCP< std::vector< Real > > obj_value_