ROL
ROL_NonlinearLeastSquaresObjective_Def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_NONLINEARLEASTSQUARESOBJECTIVE_DEF_H
11#define ROL_NONLINEARLEASTSQUARESOBJECTIVE_DEF_H
12
13namespace ROL {
14
15template<typename Real>
17 const Vector<Real> &optvec,
18 const Vector<Real> &convec,
19 const bool GNH)
20 : con_(con), GaussNewtonHessian_(GNH) {
21 c1_ = convec.clone(); c1dual_ = c1_->dual().clone();
22 c2_ = convec.clone();
23 x_ = optvec.dual().clone();
24}
25
26template<typename Real>
28 Real tol = std::sqrt(ROL_EPSILON<Real>());
29 con_->update(x,type,iter);
30 con_->value(*c1_,x,tol);
31 c1dual_->set(c1_->dual());
32}
33
34template<typename Real>
35void NonlinearLeastSquaresObjective<Real>::update( const Vector<Real> &x, bool flag, int iter ) {
36 Real tol = std::sqrt(ROL_EPSILON<Real>());
37 con_->update(x,flag,iter);
38 con_->value(*c1_,x,tol);
39 c1dual_->set(c1_->dual());
40}
41
42template<typename Real>
44 Real half(0.5);
45 return half*(c1_->dot(*c1_));
46}
47
48template<typename Real>
50 con_->applyAdjointJacobian(g,*c1dual_,x,tol);
51}
52
53template<typename Real>
55 con_->applyJacobian(*c2_,v,x,tol);
56 con_->applyAdjointJacobian(hv,c2_->dual(),x,tol);
57 if ( !GaussNewtonHessian_ ) {
58 con_->applyAdjointHessian(*x_,*c1dual_,v,x,tol);
59 hv.plus(*x_);
60 }
61}
62
63template<typename Real>
65 con_->applyPreconditioner(Pv,v,x,x.dual(),tol);
66}
67
68template<typename Real>
69void NonlinearLeastSquaresObjective<Real>::setParameter(const std::vector<Real> &param) {
71 con_->setParameter(param);
72}
73
74} // namespace ROL
75
76#endif
Defines the general constraint operator interface.
NonlinearLeastSquaresObjective(const Ptr< Constraint< Real > > &con, const Vector< Real > &optvec, const Vector< Real > &convec, const bool GNH=false)
Constructor.
void update(const Vector< Real > &x, UpdateType type, int iter=-1) override
Update objective function.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol) override
Compute gradient.
void setParameter(const std::vector< Real > &param) override
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply Hessian approximation to vector.
void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply preconditioner to vector.
Real value(const Vector< Real > &x, Real &tol) override
Compute value.
virtual void setParameter(const std::vector< Real > &param)
Defines the linear algebra or vector space interface.
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 ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.