ROL
ROL_TypeU_TrustRegionAlgorithm.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_TYPEU_TRUSTREGIONALGORITHM_H
11#define ROL_TYPEU_TRUSTREGIONALGORITHM_H
12
15#include "ROL_TrustRegion_U.hpp"
17#include "ROL_Secant.hpp"
18
24namespace ROL {
25namespace TypeU {
26
27template<typename Real>
28class TrustRegionAlgorithm : public Algorithm<Real> {
29private:
30 // TRUST REGION INFORMATION
31 Ptr<TrustRegion_U<Real>> solver_;
32 Ptr<TrustRegionModel_U<Real>> model_;
34 Real delMax_;
35 Real eta0_;
36 Real eta1_;
37 Real eta2_;
38 Real gamma0_;
39 Real gamma1_;
40 Real gamma2_;
41 Real TRsafe_;
42 Real eps_;
44 int SPflag_;
45 int SPiter_;
46
47 // NONMONOTONE INFORMATION
48 bool useNM_;
50
51 // SECANT INFORMATION
55
56 // INEXACT COMPUTATION PARAMETERS
57 std::vector<bool> useInexact_;
58 Real scale0_;
59 Real scale1_;
62 Real gtol_;
63
64 // VERBOSITY SETTING
67
68 using Algorithm<Real>::state_;
69 using Algorithm<Real>::status_;
70 using Algorithm<Real>::initialize;
71
72public:
73
74 TrustRegionAlgorithm( ParameterList &parlist,
75 const Ptr<Secant<Real>> &secant = nullPtr );
76
77 using Algorithm<Real>::run;
78 void run( Vector<Real> &x,
79 const Vector<Real> &g,
80 Objective<Real> &obj,
81 std::ostream &outStream = std::cout) override;
82
83 void writeHeader( std::ostream& os ) const override;
84
85 void writeName( std::ostream& os ) const override;
86
87 void writeOutput( std::ostream& os, const bool print_header = false ) const override;
88
89private:
90
91 void initialize(const Vector<Real> &x, const Vector<Real> &g, Vector<Real> &Bg,
92 Objective<Real> &obj, std::ostream &outStream = std::cout);
93
94 Real computeValue(const Vector<Real> &x, Objective<Real> &obj, Real pRed);
95
108 void computeGradient(const Vector<Real> &x, Objective<Real> &obj, bool accept);
109
110}; // class ROL::TrustRegionAlgorithm
111} // namespace TypeU
112} // namespace ROL
113
115
116#endif
Contains definitions of enums for trust region algorithms.
Provides the interface to evaluate objective functions.
Provides interface for and implements limited-memory secant operators.
Provides an interface to run unconstrained optimization algorithms.
const Ptr< CombinedStatusTest< Real > > status_
const Ptr< AlgorithmState< Real > > state_
Provides an interface to run trust-region methods for unconstrained optimization algorithms.
Real scale1_
Scale for inexact gradient computation.
bool printHeader_
Print header at every iteration.
std::vector< bool > useInexact_
Flags for inexact (0) objective function, (1) gradient, (2) Hessian.
int SPiter_
Subproblem solver iteration count.
Real gamma0_
Radius decrease rate (negative rho).
Real TRsafe_
Safeguard size for numerically evaluating ratio.
void writeHeader(std::ostream &os) const override
Print iterate header.
Real scale0_
Scale for inexact gradient computation.
void writeOutput(std::ostream &os, const bool print_header=false) const override
Print iterate status.
Real eps_
Safeguard for numerically evaluating ratio.
int verbosity_
Print additional information to screen if > 0.
void writeName(std::ostream &os) const override
Print step name.
void computeGradient(const Vector< Real > &x, Objective< Real > &obj, bool accept)
Compute gradient to iteratively satisfy inexactness condition.
Real computeValue(const Vector< Real > &x, Objective< Real > &obj, Real pRed)
Real gamma1_
Radius decrease rate (positive rho).
int SPflag_
Subproblem solver termination flag.
TRUtils::ETRFlag TRflag_
Trust-region exit flag.
Ptr< TrustRegionModel_U< Real > > model_
Container for trust-region model.
void initialize(const Vector< Real > &x, const Vector< Real > &g, Vector< Real > &Bg, Objective< Real > &obj, std::ostream &outStream=std::cout)
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...
ETrustRegionU etr_
Trust-region subproblem solver type.
Ptr< TrustRegion_U< Real > > solver_
Container for trust-region solver object.
Defines the linear algebra or vector space interface.