slope 6.0.1
Loading...
Searching...
No Matches
pgd.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "../losses/loss.h"
9#include "../math.h"
10#include "../sorted_l1_norm.h"
11#include "solver.h"
12#include <Eigen/Dense>
13#include <Eigen/SparseCore>
14#include <memory>
15
16namespace slope {
17
25class PGD : public SolverBase
26{
27public:
35 bool intercept,
36 const std::string& update_type)
38 , learning_rate(1.0)
39 , learning_rate_decr(0.5)
40 , update_type{ update_type }
41 , t(1.0)
42 {
43 }
44
46 void run(Eigen::VectorXd& beta0,
47 Eigen::VectorXd& beta,
48 Eigen::MatrixXd& eta,
49 const Eigen::ArrayXd& lambda,
50 const std::unique_ptr<Loss>& loss,
51 const SortedL1Norm& penalty,
52 const Eigen::VectorXd& gradient,
53 const std::vector<int>& working_set,
54 const Eigen::MatrixXd& x,
55 const Eigen::VectorXd& x_centers,
56 const Eigen::VectorXd& x_scales,
57 const Eigen::MatrixXd& y) override;
58
60 void run(Eigen::VectorXd& beta0,
61 Eigen::VectorXd& beta,
62 Eigen::MatrixXd& eta,
63 const Eigen::ArrayXd& lambda,
64 const std::unique_ptr<Loss>& loss,
65 const SortedL1Norm& penalty,
66 const Eigen::VectorXd& gradient,
67 const std::vector<int>& working_set,
68 const Eigen::SparseMatrix<double>& x,
69 const Eigen::VectorXd& x_centers,
70 const Eigen::VectorXd& x_scales,
71 const Eigen::MatrixXd& y) override;
72
74 void run(Eigen::VectorXd& beta0,
75 Eigen::VectorXd& beta,
76 Eigen::MatrixXd& eta,
77 const Eigen::ArrayXd& lambda,
78 const std::unique_ptr<Loss>& loss,
79 const SortedL1Norm& penalty,
80 const Eigen::VectorXd& gradient,
81 const std::vector<int>& working_set,
82 const Eigen::Map<Eigen::MatrixXd>& x,
83 const Eigen::VectorXd& x_centers,
84 const Eigen::VectorXd& x_scales,
85 const Eigen::MatrixXd& y) override;
86
88 void run(Eigen::VectorXd& beta0,
89 Eigen::VectorXd& beta,
90 Eigen::MatrixXd& eta,
91 const Eigen::ArrayXd& lambda,
92 const std::unique_ptr<Loss>& loss,
93 const SortedL1Norm& penalty,
94 const Eigen::VectorXd& gradient,
95 const std::vector<int>& working_set,
96 const Eigen::Map<Eigen::SparseMatrix<double>>& x,
97 const Eigen::VectorXd& x_centers,
98 const Eigen::VectorXd& x_scales,
99 const Eigen::MatrixXd& y) override;
100
101private:
102 template<typename MatrixType>
103 void runImpl(Eigen::VectorXd& beta0,
104 Eigen::VectorXd& beta,
105 Eigen::MatrixXd& eta,
106 const Eigen::ArrayXd& lambda,
107 const std::unique_ptr<Loss>& loss,
108 const SortedL1Norm& penalty,
109 const Eigen::VectorXd& gradient,
110 const std::vector<int>& working_set,
111 const MatrixType& x,
112 const Eigen::VectorXd& x_centers,
113 const Eigen::VectorXd& x_scales,
114 const Eigen::MatrixXd& y)
115 {
116 using Eigen::all;
117 using Eigen::MatrixXd;
118 using Eigen::VectorXd;
119
120 int n_working = working_set.size();
121
122 Eigen::VectorXd beta_old = beta(working_set);
123 Eigen::VectorXd beta_diff(n_working);
124 Eigen::VectorXd beta0_old = beta0;
125
126 double g_old = loss->loss(eta, y);
127 double t_old = t;
128
129 Eigen::MatrixXd residual = loss->residual(eta, y);
130 Eigen::VectorXd intercept_grad = residual.colwise().mean();
131 Eigen::VectorXd grad_working = gradient(working_set);
132
133 int line_search_iter = 0;
134 const int max_line_search_iter = 100; // maximum iterations before exit
135
136 while (true) {
137 if (intercept) {
138 beta0 = beta0_old - this->learning_rate * intercept_grad;
139 }
140
141 beta(working_set) =
142 penalty.prox(beta_old - this->learning_rate * grad_working,
143 this->learning_rate * lambda.head(n_working));
144
145 beta_diff = beta(working_set) - beta_old;
146 double beta_diff_norm =
147 beta_diff.dot(grad_working) +
148 (1.0 / (2 * this->learning_rate)) * beta_diff.squaredNorm();
149
150 if (intercept) {
151 Eigen::VectorXd beta0_diff = beta0 - beta0_old;
152
153 beta_diff_norm +=
154 intercept_grad.dot(beta0_diff) +
155 (1.0 / (2 * this->learning_rate)) * beta0_diff.squaredNorm();
156 }
157
158 eta = linearPredictor(x,
159 working_set,
160 beta0,
161 beta,
162 x_centers,
163 x_scales,
165 intercept);
166
167 double g = loss->loss(eta, y);
168 double q = g_old + beta_diff_norm;
169
170 if (q >= g * (1 - 1e-12) || this->learning_rate < 1e-12) {
171 this->learning_rate *= 1.1;
172 break;
173 } else {
174 this->learning_rate *= this->learning_rate_decr;
175 }
176
177 line_search_iter++;
178 if (line_search_iter >= max_line_search_iter) {
179 break;
180 }
181 }
182
183 if (update_type == "fista") {
184 if (beta_prev.size() == 0) {
185 beta_prev = beta;
186 }
187
188 this->t = 0.5 * (1.0 + std::sqrt(1.0 + 4.0 * t_old * t_old));
189 Eigen::VectorXd beta_current = beta(working_set);
190
191 beta(working_set) +=
192 (beta_current - beta_prev(working_set)) * (t_old - 1.0) / this->t;
193 beta_prev(working_set) = beta_current;
194
195 eta = linearPredictor(x,
196 working_set,
197 beta0,
198 beta,
199 x_centers,
200 x_scales,
202 intercept);
203 }
204 }
205
206 double learning_rate;
207 double learning_rate_decr;
208 std::string update_type;
209 double t;
210 Eigen::VectorXd beta_prev;
211};
212
213} // namespace slope
Proximal Gradient Descent solver for SLOPE optimization.
Definition pgd.h:26
void run(Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::MatrixXd &eta, const Eigen::ArrayXd &lambda, const std::unique_ptr< Loss > &loss, const SortedL1Norm &penalty, const Eigen::VectorXd &gradient, const std::vector< int > &working_set, const Eigen::Map< Eigen::MatrixXd > &x, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::MatrixXd &y) override
Pure virtual function defining the solver's optimization routine.
void run(Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::MatrixXd &eta, const Eigen::ArrayXd &lambda, const std::unique_ptr< Loss > &loss, const SortedL1Norm &penalty, const Eigen::VectorXd &gradient, const std::vector< int > &working_set, const Eigen::SparseMatrix< double > &x, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::MatrixXd &y) override
Pure virtual function defining the solver's optimization routine.
void run(Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::MatrixXd &eta, const Eigen::ArrayXd &lambda, const std::unique_ptr< Loss > &loss, const SortedL1Norm &penalty, const Eigen::VectorXd &gradient, const std::vector< int > &working_set, const Eigen::MatrixXd &x, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::MatrixXd &y) override
Pure virtual function defining the solver's optimization routine.
PGD(JitNormalization jit_normalization, bool intercept, const std::string &update_type)
Constructs Proximal Gradient Descent solver for SLOPE optimization.
Definition pgd.h:34
void run(Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::MatrixXd &eta, const Eigen::ArrayXd &lambda, const std::unique_ptr< Loss > &loss, const SortedL1Norm &penalty, const Eigen::VectorXd &gradient, const std::vector< int > &working_set, const Eigen::Map< Eigen::SparseMatrix< double > > &x, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::MatrixXd &y) override
Pure virtual function defining the solver's optimization routine.
Abstract base class for SLOPE optimization solvers.
Definition solver.h:30
JitNormalization jit_normalization
JIT feature normalization strategy.
Definition solver.h:165
bool intercept
If true, fits intercept term.
Definition solver.h:166
Class representing the Sorted L1 Norm.
Eigen::MatrixXd prox(const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda) const
Computes the proximal operator of the Sorted L1 Norm.
Namespace containing SLOPE regression implementation.
Definition clusters.h:11
JitNormalization
Enums to control predictor standardization behavior.
Eigen::MatrixXd linearPredictor(const T &x, const std::vector< int > &active_set, const Eigen::VectorXd &beta0, const Eigen::VectorXd &beta, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization jit_normalization, const bool intercept)
Definition math.h:151
Numerical solver class for SLOPE (Sorted L-One Penalized Estimation)