slope 6.2.1
Loading...
Searching...
No Matches
screening.h
Go to the documentation of this file.
1
9#pragma once
10
11#include "jit_normalization.h"
12#include <Eigen/SparseCore>
13#include <memory>
14#include <vector>
15
16namespace slope {
17
24std::vector<int>
25activeSet(const Eigen::VectorXd& beta);
26
35std::vector<int>
36strongSet(const Eigen::VectorXd& gradient_prev,
37 const Eigen::ArrayXd& lambda,
38 const Eigen::ArrayXd& lambda_prev);
39
48{
49public:
53 virtual ~ScreeningRule() = default;
54
63 virtual std::vector<int> initialize(const std::vector<int>& full_set,
64 int alpha_max_ind) = 0;
65
76 virtual std::vector<int> screen(Eigen::VectorXd& gradient,
77 const Eigen::ArrayXd& lambda_curr,
78 const Eigen::ArrayXd& lambda_prev,
79 const Eigen::VectorXd& beta,
80 const std::vector<int>& full_set) = 0;
81
98 virtual bool checkKktViolations(Eigen::VectorXd& gradient,
99 const Eigen::VectorXd& beta,
100 const Eigen::ArrayXd& lambda_curr,
101 std::vector<int>& working_set,
102 const Eigen::MatrixXd& x,
103 const Eigen::MatrixXd& residual,
104 const Eigen::VectorXd& x_centers,
105 const Eigen::VectorXd& x_scales,
106 JitNormalization jit_normalization,
107 const std::vector<int>& full_set) = 0;
123 virtual bool checkKktViolations(Eigen::VectorXd& gradient,
124 const Eigen::VectorXd& beta,
125 const Eigen::ArrayXd& lambda_curr,
126 std::vector<int>& working_set,
127 const Eigen::SparseMatrix<double>& x,
128 const Eigen::MatrixXd& residual,
129 const Eigen::VectorXd& x_centers,
130 const Eigen::VectorXd& x_scales,
131 JitNormalization jit_normalization,
132 const std::vector<int>& full_set) = 0;
133
149 virtual bool checkKktViolations(Eigen::VectorXd& gradient,
150 const Eigen::VectorXd& beta,
151 const Eigen::ArrayXd& lambda_curr,
152 std::vector<int>& working_set,
153 const Eigen::Map<Eigen::MatrixXd>& x,
154 const Eigen::MatrixXd& residual,
155 const Eigen::VectorXd& x_centers,
156 const Eigen::VectorXd& x_scales,
157 JitNormalization jit_normalization,
158 const std::vector<int>& full_set) = 0;
159
175 virtual bool checkKktViolations(
176 Eigen::VectorXd& gradient,
177 const Eigen::VectorXd& beta,
178 const Eigen::ArrayXd& lambda_curr,
179 std::vector<int>& working_set,
180 const Eigen::Map<Eigen::SparseMatrix<double>>& x,
181 const Eigen::MatrixXd& residual,
182 const Eigen::VectorXd& x_centers,
183 const Eigen::VectorXd& x_scales,
184 JitNormalization jit_normalization,
185 const std::vector<int>& full_set) = 0;
186
191 virtual std::string toString() const = 0;
192
193protected:
195 std::vector<int> strong_set;
196};
197
203{
204public:
205 std::vector<int> initialize(const std::vector<int>& full_set,
206 int alpha_max_ind) override;
207
208 std::vector<int> screen(Eigen::VectorXd& gradient,
209 const Eigen::ArrayXd& lambda_curr,
210 const Eigen::ArrayXd& lambda_prev,
211 const Eigen::VectorXd& beta,
212 const std::vector<int>& full_set) override;
213
214 bool checkKktViolations(Eigen::VectorXd& gradient,
215 const Eigen::VectorXd& beta,
216 const Eigen::ArrayXd& lambda_curr,
217 std::vector<int>& working_set,
218 const Eigen::MatrixXd& x,
219 const Eigen::MatrixXd& residual,
220 const Eigen::VectorXd& x_centers,
221 const Eigen::VectorXd& x_scales,
222 JitNormalization jit_normalization,
223 const std::vector<int>& full_set) override;
224
225 bool checkKktViolations(Eigen::VectorXd& gradient,
226 const Eigen::VectorXd& beta,
227 const Eigen::ArrayXd& lambda_curr,
228 std::vector<int>& working_set,
229 const Eigen::SparseMatrix<double>& x,
230 const Eigen::MatrixXd& residual,
231 const Eigen::VectorXd& x_centers,
232 const Eigen::VectorXd& x_scales,
233 JitNormalization jit_normalization,
234 const std::vector<int>& full_set) override;
235
236 bool checkKktViolations(Eigen::VectorXd& gradient,
237 const Eigen::VectorXd& beta,
238 const Eigen::ArrayXd& lambda_curr,
239 std::vector<int>& working_set,
240 const Eigen::Map<Eigen::MatrixXd>& x,
241 const Eigen::MatrixXd& residual,
242 const Eigen::VectorXd& x_centers,
243 const Eigen::VectorXd& x_scales,
244 JitNormalization jit_normalization,
245 const std::vector<int>& full_set) override;
246
247 bool checkKktViolations(Eigen::VectorXd& gradient,
248 const Eigen::VectorXd& beta,
249 const Eigen::ArrayXd& lambda_curr,
250 std::vector<int>& working_set,
251 const Eigen::Map<Eigen::SparseMatrix<double>>& x,
252 const Eigen::MatrixXd& residual,
253 const Eigen::VectorXd& x_centers,
254 const Eigen::VectorXd& x_scales,
255 JitNormalization jit_normalization,
256 const std::vector<int>& full_set) override;
257
258 std::string toString() const override;
259};
260
266{
267public:
268 std::vector<int> initialize(const std::vector<int>& full_set,
269 int alpha_max_ind) override;
270
271 std::vector<int> screen(Eigen::VectorXd& gradient,
272 const Eigen::ArrayXd& lambda_curr,
273 const Eigen::ArrayXd& lambda_prev,
274 const Eigen::VectorXd& beta,
275 const std::vector<int>& full_set) override;
276
277 bool checkKktViolations(Eigen::VectorXd& gradient,
278 const Eigen::VectorXd& beta,
279 const Eigen::ArrayXd& lambda_curr,
280 std::vector<int>& working_set,
281 const Eigen::MatrixXd& x,
282 const Eigen::MatrixXd& residual,
283 const Eigen::VectorXd& x_centers,
284 const Eigen::VectorXd& x_scales,
285 JitNormalization jit_normalization,
286 const std::vector<int>& full_set) override;
287
288 bool checkKktViolations(Eigen::VectorXd& gradient,
289 const Eigen::VectorXd& beta,
290 const Eigen::ArrayXd& lambda_curr,
291 std::vector<int>& working_set,
292 const Eigen::Map<Eigen::MatrixXd>& x,
293 const Eigen::MatrixXd& residual,
294 const Eigen::VectorXd& x_centers,
295 const Eigen::VectorXd& x_scales,
296 JitNormalization jit_normalization,
297 const std::vector<int>& full_set) override;
298
299 bool checkKktViolations(Eigen::VectorXd& gradient,
300 const Eigen::VectorXd& beta,
301 const Eigen::ArrayXd& lambda_curr,
302 std::vector<int>& working_set,
303 const Eigen::SparseMatrix<double>& x,
304 const Eigen::MatrixXd& residual,
305 const Eigen::VectorXd& x_centers,
306 const Eigen::VectorXd& x_scales,
307 JitNormalization jit_normalization,
308 const std::vector<int>& full_set) override;
309
310 bool checkKktViolations(Eigen::VectorXd& gradient,
311 const Eigen::VectorXd& beta,
312 const Eigen::ArrayXd& lambda_curr,
313 std::vector<int>& working_set,
314 const Eigen::Map<Eigen::SparseMatrix<double>>& x,
315 const Eigen::MatrixXd& residual,
316 const Eigen::VectorXd& x_centers,
317 const Eigen::VectorXd& x_scales,
318 JitNormalization jit_normalization,
319 const std::vector<int>& full_set) override;
320
321 std::string toString() const override;
322
323private:
324 template<typename MatrixType>
325 bool checkKktViolationsImpl(Eigen::VectorXd& gradient,
326 const Eigen::VectorXd& beta,
327 const Eigen::ArrayXd& lambda_curr,
328 std::vector<int>& working_set,
329 const MatrixType& x,
330 const Eigen::MatrixXd& residual,
331 const Eigen::VectorXd& x_centers,
332 const Eigen::VectorXd& x_scales,
333 JitNormalization jit_normalization,
334 const std::vector<int>& full_set);
335};
336
344std::unique_ptr<ScreeningRule>
345createScreeningRule(const std::string& screening_type);
346
347} // namespace slope
No screening rule - uses all variables.
Definition screening.h:203
std::vector< int > initialize(const std::vector< int > &full_set, int alpha_max_ind) override
Initialize the screening rule at the start of the path algorithm.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::SparseMatrix< double > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations with sparse matrix input.
std::vector< int > screen(Eigen::VectorXd &gradient, const Eigen::ArrayXd &lambda_curr, const Eigen::ArrayXd &lambda_prev, const Eigen::VectorXd &beta, const std::vector< int > &full_set) override
Screen for the next path step.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::Map< Eigen::MatrixXd > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations with sparse matrix input.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::Map< Eigen::SparseMatrix< double > > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations with sparse matrix input.
std::string toString() const override
Get string representation of the screening rule.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::MatrixXd &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations and update working set if necessary.
Base class for screening rules in SLOPE.
Definition screening.h:48
virtual bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::MatrixXd &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set)=0
Check for KKT violations and update working set if necessary.
virtual bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::SparseMatrix< double > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set)=0
Check for KKT violations with sparse matrix input.
std::vector< int > strong_set
Strong set of variables.
Definition screening.h:195
virtual ~ScreeningRule()=default
Virtual destructor.
virtual bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::Map< Eigen::MatrixXd > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set)=0
Check for KKT violations with sparse matrix input.
virtual std::vector< int > initialize(const std::vector< int > &full_set, int alpha_max_ind)=0
Initialize the screening rule at the start of the path algorithm.
virtual std::string toString() const =0
Get string representation of the screening rule.
virtual bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::Map< Eigen::SparseMatrix< double > > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set)=0
Check for KKT violations with sparse matrix input.
virtual std::vector< int > screen(Eigen::VectorXd &gradient, const Eigen::ArrayXd &lambda_curr, const Eigen::ArrayXd &lambda_prev, const Eigen::VectorXd &beta, const std::vector< int > &full_set)=0
Screen for the next path step.
Implements strong screening rules for SLOPE.
Definition screening.h:266
std::vector< int > screen(Eigen::VectorXd &gradient, const Eigen::ArrayXd &lambda_curr, const Eigen::ArrayXd &lambda_prev, const Eigen::VectorXd &beta, const std::vector< int > &full_set) override
Screen for the next path step.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::SparseMatrix< double > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations with sparse matrix input.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::Map< Eigen::MatrixXd > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations with sparse matrix input.
std::vector< int > initialize(const std::vector< int > &full_set, int alpha_max_ind) override
Initialize the screening rule at the start of the path algorithm.
std::string toString() const override
Get string representation of the screening rule.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::Map< Eigen::SparseMatrix< double > > &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations with sparse matrix input.
bool checkKktViolations(Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda_curr, std::vector< int > &working_set, const Eigen::MatrixXd &x, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, JitNormalization jit_normalization, const std::vector< int > &full_set) override
Check for KKT violations and update working set if necessary.
Enums to control predictor standardization behavior.
Namespace containing SLOPE regression implementation.
Definition clusters.h:11
std::unique_ptr< ScreeningRule > createScreeningRule(const std::string &screening_type)
Creates a screening rule based on the provided type.
JitNormalization
Enums to control predictor standardization behavior.
std::vector< int > strongSet(const Eigen::VectorXd &gradient_prev, const Eigen::ArrayXd &lambda, const Eigen::ArrayXd &lambda_prev)
Determines the strong set using sequential strong rules.
std::vector< int > activeSet(const Eigen::VectorXd &beta)
Identifies previously active variables.