slope 0.28.0
Loading...
Searching...
No Matches
slope.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "slope_fit.h"
9#include "slope_path.h"
10#include <Eigen/Core>
11#include <Eigen/SparseCore>
12#include <cassert>
13#include <optional>
14
18namespace slope {
19
28class Slope
29{
30public:
37 : intercept(true)
38 , modify_x(false)
39 , update_clusters(false)
40 , collect_diagnostics(false)
41 , return_clusters(true)
42 , alpha_min_ratio(-1) // TODO: Use std::optional for alpha_min_ratio
43 , dev_change_tol(1e-5)
44 , dev_ratio_tol(0.999)
45 , learning_rate_decr(0.5)
46 , q(0.1)
47 , tol(1e-4)
48 , alpha_est_maxit(1000)
49 , max_it(1e4)
50 , path_length(100)
51 , cd_iterations(10)
52 , max_clusters(std::optional<int>())
53 , alpha_type("path")
54 , lambda_type("bh")
55 , centering_type("mean")
56 , scaling_type("sd")
57 , loss_type("quadratic")
58 , screening_type("strong")
59 , solver_type("auto")
60 {
61 }
62
71 void setSolver(const std::string& solver);
72
78 void setIntercept(bool intercept);
79
87 void setNormalization(const std::string& type);
88
95 void setUpdateClusters(bool update_clusters);
96
103 void setReturnClusters(const bool return_clusters);
104
112 void setAlphaMinRatio(double alpha_min_ratio);
113
128 void setAlphaType(const std::string& alpha_type);
129
136 void setLearningRateDecr(double learning_rate_decr);
137
144 void setQ(double q);
145
152 void setOscarParameters(const double theta1, const double theta2);
153
159 void setTol(double tol);
160
168 void setMaxIterations(int max_it);
169
175 void setPathLength(int path_length);
176
184 void setHybridCdIterations(int cd_iterations);
185
193 void setLambdaType(const std::string& lambda_type);
194
205 void setLoss(const std::string& loss_type);
206
216 void setScreening(const std::string& screening_type);
217
225 void setModifyX(const bool modify_x);
226
231 void setDevChangeTol(const double dev_change_tol);
232
238 void setDevRatioTol(const double dev_ratio_tol);
239
247 void setMaxClusters(const int max_clusters);
248
253 void setCentering(const std::string& type);
254
260 void setCentering(const Eigen::VectorXd& x_centers);
261
266 void setScaling(const std::string& type);
267
274 void setDiagnostics(const bool collect_diagnostics);
275
281 void setScaling(const Eigen::VectorXd& x_scales);
282
293 void setAlphaEstimationMaxIterations(const int alpha_est_maxit);
294
300
304 bool getFitIntercept() const;
305
310 const std::string& getLossType();
311
328 template<typename T>
329 SlopePath path(T& x,
330 const Eigen::MatrixXd& y_in,
331 Eigen::ArrayXd alpha = Eigen::ArrayXd::Zero(0),
332 Eigen::ArrayXd lambda = Eigen::ArrayXd::Zero(0));
333
349 template<typename T>
350 SlopeFit fit(T& x,
351 const Eigen::MatrixXd& y_in,
352 const double alpha = 1.0,
353 Eigen::ArrayXd lambda = Eigen::ArrayXd::Zero(0));
354
355private:
356 // Parameters
357 bool intercept;
358 bool modify_x;
359 bool update_clusters;
360 bool collect_diagnostics;
361 bool return_clusters;
362 double alpha_min_ratio;
363 double dev_change_tol;
364 double dev_ratio_tol;
365 double learning_rate_decr;
366 double q;
367 double theta1;
368 double theta2;
369 double tol;
370 int alpha_est_maxit;
371 int max_it;
372 int path_length;
373 int cd_iterations;
374 std::optional<int> max_clusters;
375 std::string alpha_type;
376 std::string lambda_type;
377 std::string centering_type;
378 std::string scaling_type;
379 std::string loss_type;
380 std::string screening_type;
381 std::string solver_type;
382
383 // Data
384 Eigen::VectorXd x_centers;
385 Eigen::VectorXd x_scales;
386};
387
388} // namespace slope
A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting.
Definition slope_fit.h:26
Container class for SLOPE regression solution paths.
Definition slope_path.h:31
The SLOPE model.
Definition slope.h:29
void setSolver(const std::string &solver)
Sets the numerical solver used to fit the model.
Definition slope.cpp:335
void setAlphaMinRatio(double alpha_min_ratio)
Sets the alpha min ratio.
Definition slope.cpp:417
void setMaxIterations(int max_it)
Sets the maximum number of iterations.
Definition slope.cpp:468
void setAlphaEstimationMaxIterations(const int alpha_est_maxit)
Sets the maximum number of iterations for the alpha estimation procedure.
Definition slope.cpp:562
SlopePath path(T &x, const Eigen::MatrixXd &y_in, Eigen::ArrayXd alpha=Eigen::ArrayXd::Zero(0), Eigen::ArrayXd lambda=Eigen::ArrayXd::Zero(0))
Computes SLOPE regression solution path for multiple alpha and lambda values.
Definition slope.cpp:27
void setDevRatioTol(const double dev_ratio_tol)
Sets tolerance in deviance change for early stopping.
Definition slope.cpp:536
void setScaling(const std::string &type)
Sets the scaling type.
Definition slope.cpp:383
SlopeFit fit(T &x, const Eigen::MatrixXd &y_in, const double alpha=1.0, Eigen::ArrayXd lambda=Eigen::ArrayXd::Zero(0))
Fits a single SLOPE regression model for given alpha and lambda values.
Definition slope.cpp:322
void setDevChangeTol(const double dev_change_tol)
Sets tolerance in deviance change for early stopping.
Definition slope.cpp:526
const std::string & getLossType()
Get currently defined loss type.
Definition slope.cpp:580
void setReturnClusters(const bool return_clusters)
Sets the update clusters flag.
Definition slope.cpp:404
void setMaxClusters(const int max_clusters)
Sets tolerance in deviance change for early stopping.
Definition slope.cpp:546
void setIntercept(bool intercept)
Sets the intercept flag.
Definition slope.cpp:342
bool getFitIntercept() const
Returns the intercept flag.
Definition slope.cpp:574
void setDiagnostics(const bool collect_diagnostics)
Toggles collection of diagnostics.
Definition slope.cpp:556
void setNormalization(const std::string &type)
Sets normalization type for the design matrix.
Definition slope.cpp:348
int getAlphaEstimationMaxIterations() const
Gets the maximum number of iterations allowed for the alpha estimation procedure.
Definition slope.cpp:568
void setScreening(const std::string &screening_type)
Sets the type of feature screening used, which discards predictors that are unlikely to be active.
Definition slope.cpp:513
void setOscarParameters(const double theta1, const double theta2)
Sets OSCAR parameters.
Definition slope.cpp:444
void setLambdaType(const std::string &lambda_type)
Sets the lambda type for regularization weights.
Definition slope.cpp:495
void setAlphaType(const std::string &alpha_type)
Sets the alpha type.
Definition slope.cpp:410
void setModifyX(const bool modify_x)
Controls if x should be modified-in-place.
Definition slope.cpp:520
void setCentering(const std::string &type)
Sets the center points for feature normalization.
Definition slope.cpp:369
void setLearningRateDecr(double learning_rate_decr)
Sets the learning rate decrement.
Definition slope.cpp:426
void setLoss(const std::string &loss_type)
Sets the loss function type.
Definition slope.cpp:504
void setPathLength(int path_length)
Sets the path length.
Definition slope.cpp:477
void setHybridCdIterations(int cd_iterations)
Sets the frequence of proximal gradient descent steps.
Definition slope.cpp:486
void setTol(double tol)
Sets the tolerance value.
Definition slope.cpp:459
void setQ(double q)
Sets the q value.
Definition slope.cpp:435
void setUpdateClusters(bool update_clusters)
Sets the update clusters flag.
Definition slope.cpp:398
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
SLOPE (Sorted L-One Penalized Estimation) fitting results.
Defines the SlopePath class for storing and accessing SLOPE regression solution paths.