44 #ifndef ROL_TYPEU_LINESEARCHALGORITHM_DEF_H
45 #define ROL_TYPEU_LINESEARCHALGORITHM_DEF_H
53 template<
typename Real>
58 :
Algorithm<Real>(), desc_(descent), lineSearch_(lineSearch),
66 ParameterList& Llist = parlist.sublist(
"Step").sublist(
"Line Search");
67 ParameterList& Glist = parlist.sublist(
"General");
74 lineSearchName_ = Llist.sublist(
"Line-Search Method").get(
"Type",
"Cubic Interpolation");
79 lineSearchName_ = Llist.sublist(
"Line-Search Method").get(
"User Defined Line Search Name",
"Unspecified User Defined Line Search");
81 if (
desc_ == nullPtr) {
82 ParameterList& dlist = Llist.sublist(
"Descent Method");
85 desc_ = DescentDirectionUFactory<Real>(parlist,secant);
88 descentName_ = Llist.sublist(
"Descent Method").get(
"User Defined Descent Direction Name",
"Unspecified User Defined Descent Direction");
92 template<
typename Real>
96 std::ostream &outStream) {
99 lineSearch_->initialize(x,g);
100 desc_->initialize(x,g);
102 Real ftol = std::sqrt(ROL_EPSILON<Real>());
104 state_->value = obj.
value(x,ftol);
106 obj.
gradient(*state_->gradientVec,x,ftol);
108 state_->gnorm = state_->gradientVec->norm();
109 state_->snorm = ROL_INF<Real>();
112 template<
typename Real>
116 std::ostream &outStream ) {
119 Real ftrial(0), gs(0), tol(std::sqrt(ROL_EPSILON<Real>()));
120 initialize(x,g,obj,outStream);
121 state_->searchSize = one;
122 Ptr<Vector<Real>> gprev = g.
clone();
125 if (verbosity_ > 0) writeOutput(outStream,
true);
127 while (status_->check(*state_)) {
129 desc_->compute(*state_->stepVec,state_->snorm,gs,SPiter_,SPflag_,x,
130 *state_->gradientVec,obj);
133 ftrial = state_->value;
134 ls_nfval_ = 0; ls_ngrad_ = 0;
135 lineSearch_->run(state_->searchSize,ftrial,ls_nfval_,ls_ngrad_,gs,*state_->stepVec,x,obj);
138 if(!acceptLastAlpha_) {
139 lineSearch_->setMaxitUpdate(state_->searchSize,ftrial,state_->value);
143 state_->stepVec->scale(state_->searchSize);
144 state_->snorm *= state_->searchSize;
147 x.
plus(*state_->stepVec);
151 state_->value = obj.
value(x,tol);
153 gprev->set(*state_->gradientVec);
154 obj.
gradient(*state_->gradientVec,x,tol);
156 state_->gnorm = state_->gradientVec->norm();
159 state_->iterateVec->set(x);
160 state_->nfval += ls_nfval_;
161 state_->ngrad += ls_ngrad_;
162 desc_->update(x,*state_->stepVec,*gprev,*state_->gradientVec,state_->snorm,state_->iter);
166 if (verbosity_ > 0) writeOutput(outStream,printHeader_);
171 template<
typename Real>
173 std::ios_base::fmtflags osFlags(os.flags());
174 if (verbosity_ > 1) {
175 os << std::string(109,
'-') << std::endl;
177 os <<
" status output definitions" << std::endl << std::endl;
178 os <<
" iter - Number of iterates (steps taken)" << std::endl;
179 os <<
" value - Objective function value" << std::endl;
180 os <<
" gnorm - Norm of the gradient" << std::endl;
181 os <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
182 os <<
" alpha - Line search step length" << std::endl;
183 os <<
" #fval - Cumulative number of times the objective function was evaluated" << std::endl;
184 os <<
" #grad - Cumulative number of times the gradient was computed" << std::endl;
185 os <<
" ls_#fval - Number of the times the objective function was evaluated during the line search" << std::endl;
186 os <<
" ls_#grad - Number of the times the gradient was evaluated during the line search" << std::endl;
188 os <<
" iterCG - Number of Krylov iterations used to compute search direction" << std::endl;
189 os <<
" flagCG - Krylov solver flag" << std::endl;
191 os << std::string(109,
'-') << std::endl;
195 os << std::setw(6) << std::left <<
"iter";
196 os << std::setw(15) << std::left <<
"value";
197 os << std::setw(15) << std::left <<
"gnorm";
198 os << std::setw(15) << std::left <<
"snorm";
199 os << std::setw(15) << std::left <<
"alpha";
200 os << std::setw(10) << std::left <<
"#fval";
201 os << std::setw(10) << std::left <<
"#grad";
202 os << std::setw(10) << std::left <<
"ls_#fval";
203 os << std::setw(10) << std::left <<
"ls_#grad";
205 os << std::setw(10) << std::left <<
"iterCG";
206 os << std::setw(10) << std::left <<
"flagCG";
212 template<
typename Real>
214 std::ios_base::fmtflags osFlags(os.flags());
215 os << std::endl << desc_->printName();
217 os <<
"Line Search: " << lineSearchName_;
222 template<
typename Real>
224 std::ios_base::fmtflags osFlags(os.flags());
225 os << std::scientific << std::setprecision(6);
226 if ( state_->iter == 0 ) {
229 if ( print_header ) {
232 if ( state_->iter == 0 ) {
234 os << std::setw(6) << std::left << state_->iter;
235 os << std::setw(15) << std::left << state_->value;
236 os << std::setw(15) << std::left << state_->gnorm;
237 os << std::setw(15) << std::left <<
"---";
238 os << std::setw(15) << std::left <<
"---";
239 os << std::setw(10) << std::left << state_->nfval;
240 os << std::setw(10) << std::left << state_->ngrad;
241 os << std::setw(10) << std::left <<
"---";
242 os << std::setw(10) << std::left <<
"---";
244 os << std::setw(10) << std::left <<
"---";
245 os << std::setw(10) << std::left <<
"---";
251 os << std::setw(6) << std::left << state_->iter;
252 os << std::setw(15) << std::left << state_->value;
253 os << std::setw(15) << std::left << state_->gnorm;
254 os << std::setw(15) << std::left << state_->snorm;
255 os << std::setw(15) << std::left << state_->searchSize;
256 os << std::setw(10) << std::left << state_->nfval;
257 os << std::setw(10) << std::left << state_->ngrad;
258 os << std::setw(10) << std::left << ls_nfval_;
259 os << std::setw(10) << std::left << ls_ngrad_;
261 os << std::setw(10) << std::left << SPiter_;
262 os << std::setw(10) << std::left << SPflag_;
Provides interface for and implements line searches.
Provides the interface to evaluate objective functions.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void writeHeader(std::ostream &os) const override
Print iterate header.
bool acceptLastAlpha_
For backwards compatibility. When max function evaluations are reached take last step.
virtual void plus(const Vector &x)=0
Compute , where .
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
std::string lineSearchName_
void writeName(std::ostream &os) const override
Print step name.
void initialize(const Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout)
void initialize(const Vector< Real > &x, const Vector< Real > &g)
Defines the linear algebra or vector space interface.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Ptr< LineSearch_U< Real > > lineSearch_
Line-search object.
EDescentU StringToEDescentU(std::string s)
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
LineSearchAlgorithm(ParameterList &parlist, const Ptr< Secant< Real >> &secant=nullPtr, const Ptr< DescentDirection_U< Real >> &descent=nullPtr, const Ptr< LineSearch_U< Real >> &lineSearch=nullPtr)
Constructor.
std::string ECurvatureConditionUToString(ECurvatureConditionU ls)
Provides an interface to run unconstrained optimization algorithms.
Provides interface for and implements limited-memory secant operators.
Provides an interface to check status of optimization algorithms.
ECurvatureConditionU econd_
enum determines type of curvature condition
void writeOutput(std::ostream &os, const bool print_header=false) const override
Print iterate status.
EDescentU edesc_
enum determines type of descent direction
ELineSearchU StringToELineSearchU(std::string s)
void run(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout) override
Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual opt...
Provides the interface to compute unconstrained optimization steps for line search.
virtual void writeExitStatus(std::ostream &os) const
const Ptr< CombinedStatusTest< Real > > status_
Ptr< DescentDirection_U< Real > > desc_
Unglobalized step object.
ECurvatureConditionU StringToECurvatureConditionU(std::string s)
ELineSearchU els_
enum determines type of line search