ROL
ROL_NonlinearCG_U.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_NONLINEARCG_U_H
11 #define ROL_NONLINEARCG_U_H
12 
14 #include "ROL_Types.hpp"
15 #include "ROL_NonlinearCG.hpp"
16 
24 
25 namespace ROL {
26 
27 template<typename Real>
28 class NonlinearCG_U : public DescentDirection_U<Real> {
29 private:
30 
31  Ptr<NonlinearCG<Real>> nlcg_;
33  std::string ncgName_;
34 
35 public:
36 
46  NonlinearCG_U( ParameterList &parlist,
47  const Ptr<NonlinearCG<Real>> &nlcg = nullPtr)
49  // Initialize secant object
50  ParameterList& Llist = parlist.sublist("Step").sublist("Line Search");
51  if ( nlcg == nullPtr ) {
52  ncgName_ = Llist.sublist("Descent Method").get("Nonlinear CG Type","Oren-Luenberger");
53  enlcg_
55  nlcg_ = ROL::makePtr<NonlinearCG<Real>>(enlcg_);
56  }
57  else {
58  ncgName_ = Llist.sublist("Descent Method").get("User Defined Nonlinear CG Name",
59  "Unspecified User Define Nonlinear CG Method");
60  }
61  }
62 
63  void compute( Vector<Real> &s, Real &snorm, Real &sdotg, int &iter, int &flag,
64  const Vector<Real> &x, const Vector<Real> &g, Objective<Real> &obj) override {
65  nlcg_->run(s,g,x,obj);
66  //sdotg = -s.dot(g.dual());
67  sdotg = -s.apply(g);
68  if (sdotg >= static_cast<Real>(0)) {
69  s.set(g.dual());
70  //sdotg = -s.dot(g.dual());
71  sdotg = -s.apply(g);
72  }
73  s.scale(static_cast<Real>(-1));
74  snorm = s.norm();
75  iter = 0;
76  flag = 0;
77  }
78 
79  std::string printName(void) const override {
80  std::stringstream name;
81  name << ncgName_ << " Nonlinear CG";
82  return name.str();
83  }
84 }; // class ROL::NonlinearCG_U
85 
86 } // namespace ROL
87 
88 #endif
Provides the interface to evaluate objective functions.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:192
virtual void scale(const Real alpha)=0
Compute where .
virtual Real apply(const Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
Definition: ROL_Vector.hpp:204
std::string printName(void) const override
Contains definitions of custom data types in ROL.
Provides the interface to compute optimization steps with nonlinear CG.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:46
void compute(Vector< Real > &s, Real &snorm, Real &sdotg, int &iter, int &flag, const Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj) override
NonlinearCG_U(ParameterList &parlist, const Ptr< NonlinearCG< Real >> &nlcg=nullPtr)
Constructor.
ENonlinearCG
Enumeration of nonlinear CG algorithms.
Definition: ROL_Types.hpp:536
ENonlinearCG StringToENonlinearCG(std::string s)
Definition: ROL_Types.hpp:608
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:175
virtual Real norm() const =0
Returns where .
Provides the interface to compute unconstrained optimization steps for line search.
Ptr< NonlinearCG< Real > > nlcg_
NonlinearCG object (used for quasi-Newton)