slope 0.29.0
Loading...
Searching...
No Matches
slope.cpp
1#include "slope.h"
2#include "clusters.h"
3#include "constants.h"
4#include "diagnostics.h"
5#include "estimate_alpha.h"
6#include "logger.h"
7#include "losses/loss.h"
8#include "losses/setup_loss.h"
9#include "math.h"
10#include "normalize.h"
12#include "screening.h"
14#include "sorted_l1_norm.h"
15#include "timer.h"
16#include "utils.h"
17#include <Eigen/Core>
18#include <Eigen/SparseCore>
19#include <memory>
20#include <numeric>
21#include <set>
22
23namespace slope {
24
25template<typename T>
26SlopePath
28 const Eigen::MatrixXd& y_in,
29 Eigen::ArrayXd alpha,
30 Eigen::ArrayXd lambda)
31{
32 using Eigen::MatrixXd;
33 using Eigen::VectorXd;
34
35 const int n = x.rows();
36 const int p = x.cols();
37
38 if (n != y_in.rows()) {
39 throw std::invalid_argument("x and y_in must have the same number of rows");
40 }
41
42 auto jit_normalization = normalize(x,
43 this->x_centers,
44 this->x_scales,
45 this->centering_type,
46 this->scaling_type,
47 this->modify_x);
48
49 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
50
51 MatrixXd y = loss->preprocessResponse(y_in);
52
53 const int m = y.cols();
54
55 std::vector<int> full_set(p * m);
56 std::iota(full_set.begin(), full_set.end(), 0);
57
58 VectorXd beta0 = VectorXd::Zero(m);
59 VectorXd beta = VectorXd::Zero(p * m);
60
61 MatrixXd eta = MatrixXd::Zero(n, m); // linear predictor
62
63 if (this->intercept) {
64 beta0 = loss->link(y.colwise().mean()).transpose();
65 eta.rowwise() = beta0.transpose();
66 }
67
68 MatrixXd residual = loss->residual(eta, y);
69 VectorXd gradient(beta.size());
70
71 // Path data
72 bool user_alpha = alpha.size() > 0;
73 bool user_lambda = lambda.size() > 0;
74
75 if (!user_lambda) {
76 lambda = lambdaSequence(
77 p * m, this->q, this->lambda_type, n, this->theta1, this->theta2);
78 } else {
79 if (lambda.size() != beta.size()) {
80 throw std::invalid_argument(
81 "lambda must be the same length as the number of coefficients");
82 }
83 if (lambda.minCoeff() < 0) {
84 throw std::invalid_argument("lambda must be non-negative");
85 }
86 if (!lambda.isFinite().all()) {
87 throw std::invalid_argument("lambda must be finite");
88 }
89 }
90
91 // Setup the regularization sequence and path
92 SortedL1Norm sl1_norm;
93
94 // TODO: Make this part of the slope class
95 auto solver = setupSolver(this->solver_type,
96 this->loss_type,
97 jit_normalization,
98 this->intercept,
99 this->update_clusters,
100 this->cd_iterations);
101
102 updateGradient(gradient,
103 x,
104 residual,
105 full_set,
106 this->x_centers,
107 this->x_scales,
108 Eigen::VectorXd::Ones(n),
109 jit_normalization);
110
111 int alpha_max_ind = whichMax(gradient.cwiseAbs());
112 double alpha_max = sl1_norm.dualNorm(gradient, lambda);
113
114 if (alpha_type == "path") {
115 if (alpha_min_ratio < 0) {
116 alpha_min_ratio = n > gradient.size() ? 1e-4 : 1e-2;
117 }
118
119 alpha = regularizationPath(alpha, path_length, alpha_min_ratio, alpha_max);
120 path_length = alpha.size();
121 } else {
122 if (loss_type != "quadratic") {
123 throw std::invalid_argument(
124 "Automatic alpha estimation is only available for the quadratic loss");
125 }
126 return estimateAlpha(x, y, *this);
127 }
128
129 // Screening setup
130 std::unique_ptr<ScreeningRule> screening_rule =
131 createScreeningRule(this->screening_type);
132 std::vector<int> working_set =
133 screening_rule->initialize(full_set, alpha_max_ind);
134
135 // Path variables
136 double null_deviance = loss->deviance(eta, y);
137 double dev_prev = null_deviance;
138
139 Timer timer;
140
141 double alpha_prev = std::max(alpha_max, alpha(0));
142
143 std::vector<SlopeFit> fits;
144
145 // Regularization path loop
146 for (int path_step = 0; path_step < this->path_length; ++path_step) {
147 double alpha_curr = alpha(path_step);
148
149 assert(alpha_curr <= alpha_prev && "Alpha must be decreasing");
150
151 Eigen::ArrayXd lambda_curr = alpha_curr * lambda;
152 Eigen::ArrayXd lambda_prev = alpha_prev * lambda;
153
154 std::vector<double> duals, primals, time;
155 timer.start();
156
157 // Update gradient for the full set
158 // TODO: Only update for non-working set since gradient is updated before
159 // the convergence check in the inner loop for the working set
160 updateGradient(gradient,
161 x,
162 residual,
163 full_set,
164 x_centers,
165 x_scales,
166 Eigen::VectorXd::Ones(x.rows()),
167 jit_normalization);
168
169 working_set = screening_rule->screen(
170 gradient, lambda_curr, lambda_prev, beta, full_set);
171
172 int it = 0;
173 for (; it < this->max_it; ++it) {
174 // Compute primal, dual, and gap
175 residual = loss->residual(eta, y);
176 updateGradient(gradient,
177 x,
178 residual,
179 working_set,
180 this->x_centers,
181 this->x_scales,
182 Eigen::VectorXd::Ones(n),
183 jit_normalization);
184
185 double primal =
186 loss->loss(eta, y) +
187 sl1_norm.eval(beta(working_set), lambda_curr.head(working_set.size()));
188
189 MatrixXd theta = residual;
190
191 // First compute gradient with potential offset for intercept case
192 VectorXd dual_gradient = gradient;
193
194 // TODO: Can we avoid this copy? Maybe revert offset afterwards or,
195 // alternatively, solve intercept until convergence and then no longer
196 // need the offset at all.
197 if (this->intercept) {
198 VectorXd theta_mean = theta.colwise().mean();
199 theta.rowwise() -= theta_mean.transpose();
200
201 offsetGradient(dual_gradient,
202 x,
203 theta_mean,
204 working_set,
205 this->x_centers,
206 this->x_scales,
207 jit_normalization);
208 }
209
210 // Common scaling operation
211 double dual_norm = sl1_norm.dualNorm(
212 dual_gradient(working_set), lambda_curr.head(working_set.size()));
213 theta.array() /= std::max(1.0, dual_norm);
214
215 double dual = loss->dual(theta, y, Eigen::VectorXd::Ones(n));
216
217 if (collect_diagnostics) {
218 timer.pause();
219 double true_dual = computeDual(beta,
220 residual,
221 loss,
222 sl1_norm,
223 lambda_curr,
224 x,
225 y,
226 this->x_centers,
227 this->x_scales,
228 jit_normalization,
229 this->intercept);
230 timer.resume();
231
232 time.emplace_back(timer.elapsed());
233 primals.emplace_back(primal);
234 duals.emplace_back(true_dual);
235 }
236
237 double dual_gap = primal - dual;
238
239 assert(dual_gap > -1e-6 && "Dual gap should be positive");
240
241 double tol_scaled = (std::abs(primal) + constants::EPSILON) * this->tol;
242
243 if (dual_gap <= tol_scaled) {
244 bool no_violations =
245 screening_rule->checkKktViolations(gradient,
246 beta,
247 lambda_curr,
248 working_set,
249 x,
250 residual,
251 this->x_centers,
252 this->x_scales,
253 jit_normalization,
254 full_set);
255 if (no_violations) {
256 break;
257 }
258 }
259
260 solver->run(beta0,
261 beta,
262 eta,
263 lambda_curr,
264 loss,
265 sl1_norm,
266 gradient,
267 working_set,
268 x,
269 this->x_centers,
270 this->x_scales,
271 y);
272 }
273
274 if (it == this->max_it) {
277 "Maximum number of iterations reached at step = " +
278 std::to_string(path_step) + ".");
279 }
280
281 // Store everything for this step of the path
282 auto [beta0_out, beta_out] = rescaleCoefficients(
283 beta0, beta, this->x_centers, this->x_scales, this->intercept);
284
285 alpha_prev = alpha_curr;
286
287 // Compute early stopping criteria
288 double dev = loss->deviance(eta, y);
289 double dev_ratio = 1 - dev / null_deviance;
290 double dev_change = path_step == 0 ? 1.0 : 1 - dev / dev_prev;
291 dev_prev = dev;
292
293 Clusters clusters;
294
295 if (return_clusters) {
296 clusters.update(beta);
297 }
298
300 beta0_out, beta_out.sparseView(), clusters, alpha_curr, lambda,
301 dev, null_deviance, primals, duals, time,
302 it, this->centering_type, this->scaling_type
303 };
304
305 fits.emplace_back(std::move(fit));
306
307 if (!user_alpha) {
308 int n_unique = unique(beta.cwiseAbs()).size();
309 if (dev_ratio > dev_ratio_tol || dev_change < dev_change_tol ||
310 n_unique >= this->max_clusters.value_or(n + 1)) {
311 break;
312 }
313 }
314 }
315
316 return fits;
317}
318
319template<typename T>
322 const Eigen::MatrixXd& y_in,
323 const double alpha,
324 Eigen::ArrayXd lambda)
325{
326 Eigen::ArrayXd alpha_arr(1);
327 alpha_arr(0) = alpha;
328 SlopePath res = path(x, y_in, alpha_arr, lambda);
329
330 return { res(0) };
331};
332
333void
334Slope::setSolver(const std::string& solver)
335{
336 validateOption(solver, { "auto", "pgd", "hybrid", "fista" }, "solver");
337 this->solver_type = solver;
338}
339
340void
341Slope::setIntercept(bool intercept)
342{
343 this->intercept = intercept;
344}
345
346void
347Slope::setNormalization(const std::string& type)
348{
350 type, { "standardization", "min_max", "max_abs", "none" }, "type");
351
352 if (type == "standardization") {
353 this->centering_type = "mean";
354 this->scaling_type = "sd";
355 } else if (type == "none") {
356 this->centering_type = "none";
357 this->scaling_type = "none";
358 } else if (type == "min_max") {
359 this->centering_type = "min";
360 this->scaling_type = "range";
361 } else if (type == "max_abs") {
362 this->centering_type = "none";
363 this->scaling_type = "max_abs";
364 }
365}
366
367void
368Slope::setCentering(const std::string& type)
369{
370 validateOption(type, { "mean", "min", "none" }, "type");
371 this->centering_type = type;
372}
373
374void
375Slope::setCentering(const Eigen::VectorXd& x_centers)
376{
377 this->x_centers = x_centers;
378 this->centering_type = "manual";
379}
380
381void
382Slope::setScaling(const std::string& type)
383{
385 type, { "sd", "l1", "l2", "range", "max_abs", "none" }, "type");
386 this->scaling_type = type;
387}
388
389void
390Slope::setScaling(const Eigen::VectorXd& x_scales)
391{
392 this->x_scales = x_scales;
393 this->scaling_type = "manual";
394}
395
396void
397Slope::setUpdateClusters(bool update_clusters)
398{
399 this->update_clusters = update_clusters;
400}
401
402void
403Slope::setReturnClusters(const bool return_clusters)
404{
405 this->return_clusters = return_clusters;
406}
407
408void
409Slope::setAlphaType(const std::string& alpha_type)
410{
411 validateOption(alpha_type, { "path", "estimate" }, "alpha_type");
412 this->alpha_type = alpha_type;
413}
414
415void
416Slope::setAlphaMinRatio(double alpha_min_ratio)
417{
418 if (alpha_min_ratio <= 0 || alpha_min_ratio >= 1) {
419 throw std::invalid_argument("alpha_min_ratio must be in (0, 1)");
420 }
421 this->alpha_min_ratio = alpha_min_ratio;
422}
423
424void
425Slope::setLearningRateDecr(double learning_rate_decr)
426{
427 if (learning_rate_decr <= 0 || learning_rate_decr >= 1) {
428 throw std::invalid_argument("learning_rate_decr must be in (0, 1)");
429 }
430 this->learning_rate_decr = learning_rate_decr;
431}
432
433void
434Slope::setQ(double q)
435{
436 if (q < 0 || q > 1) {
437 throw std::invalid_argument("q must be between 0 and 1");
438 }
439 this->q = q;
440}
441
442void
443Slope::setOscarParameters(const double theta1, const double theta2)
444{
445 if (theta1 < 0) {
446 throw std::invalid_argument("theta1 must be between 0 and 1");
447 }
448
449 if (theta2 < 0) {
450 throw std::invalid_argument("theta2 must be between 0 and 1");
451 }
452
453 this->theta1 = theta1;
454 this->theta2 = theta2;
455}
456
457void
458Slope::setTol(double tol)
459{
460 if (tol < 0) {
461 throw std::invalid_argument("tol must be non-negative");
462 }
463 this->tol = tol;
464}
465
466void
468{
469 if (max_it < 1) {
470 throw std::invalid_argument("max_it_outer must be >= 1");
471 }
472 this->max_it = max_it;
473}
474
475void
476Slope::setPathLength(int path_length)
477{
478 if (path_length < 1) {
479 throw std::invalid_argument("path_length must be >= 1");
480 }
481 this->path_length = path_length;
482}
483
484void
486{
487 if (cd_iterations < 0) {
488 throw std::invalid_argument("cd_iterations must be >= 0");
489 }
490 this->cd_iterations = cd_iterations;
491}
492
493void
494Slope::setLambdaType(const std::string& lambda_type)
495{
497 lambda_type, { "bh", "gaussian", "oscar", "lasso" }, "lambda_type");
498
499 this->lambda_type = lambda_type;
500}
501
502void
503Slope::setLoss(const std::string& loss_type)
504{
505 validateOption(loss_type,
506 { "quadratic", "logistic", "poisson", "multinomial" },
507 "loss_type");
508 this->loss_type = loss_type;
509}
510
511void
512Slope::setScreening(const std::string& screening_type)
513{
514 validateOption(screening_type, { "strong", "none" }, "screening_type");
515 this->screening_type = screening_type;
516}
517
518void
519Slope::setModifyX(const bool modify_x)
520{
521 this->modify_x = modify_x;
522}
523
524void
525Slope::setDevChangeTol(const double dev_change_tol)
526{
527 if (dev_change_tol < 0 || dev_change_tol > 1) {
528 throw std::invalid_argument("dev_change_tol must be in [0, 1]");
529 }
530
531 this->dev_change_tol = dev_change_tol;
532}
533
534void
535Slope::setDevRatioTol(const double dev_ratio_tol)
536{
537 if (dev_ratio_tol < 0 || dev_ratio_tol > 1) {
538 throw std::invalid_argument("dev_ratio_tol must be in [0, 1]");
539 }
540
541 this->dev_ratio_tol = dev_ratio_tol;
542}
543
544void
545Slope::setMaxClusters(const int max_clusters)
546{
547 if (max_clusters < 1) {
548 throw std::invalid_argument("max_clusters must be >= 1");
549 }
550
551 this->max_clusters = max_clusters;
552}
553
554void
555Slope::setDiagnostics(const bool collect_diagnostics)
556{
557 this->collect_diagnostics = collect_diagnostics;
558}
559
560void
562{
563 this->alpha_est_maxit = alpha_est_maxit;
564}
565
566int
568{
569 return alpha_est_maxit;
570}
571
572bool
574{
575 return intercept;
576}
577
578const std::string&
580{
581 return loss_type;
582}
583
585// Explicit template instantiations
586template SlopePath
587Slope::path<Eigen::MatrixXd>(Eigen::MatrixXd&,
588 const Eigen::MatrixXd&,
589 Eigen::ArrayXd,
590 Eigen::ArrayXd);
591
592template SlopePath
593Slope::path<Eigen::SparseMatrix<double>>(Eigen::SparseMatrix<double>&,
594 const Eigen::MatrixXd&,
595 Eigen::ArrayXd,
596 Eigen::ArrayXd);
597
598template SlopeFit
599Slope::fit<Eigen::MatrixXd>(Eigen::MatrixXd&,
600 const Eigen::MatrixXd&,
601 const double,
602 Eigen::ArrayXd);
603
604template SlopeFit
605Slope::fit<Eigen::SparseMatrix<double>>(Eigen::SparseMatrix<double>&,
606 const Eigen::MatrixXd&,
607 const double,
608 Eigen::ArrayXd);
610
611} // namespace slope
Representation of the clusters in SLOPE.
Definition clusters.h:18
void update(const int old_index, const int new_index, const double c_new)
Updates the cluster structure when an index is changed.
Definition clusters.cpp:164
A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting.
Definition slope_fit.h:27
Container class for SLOPE regression solution paths.
Definition slope_path.h:31
void setSolver(const std::string &solver)
Sets the numerical solver used to fit the model.
Definition slope.cpp:334
void setAlphaMinRatio(double alpha_min_ratio)
Sets the alpha min ratio.
Definition slope.cpp:416
void setMaxIterations(int max_it)
Sets the maximum number of iterations.
Definition slope.cpp:467
void setAlphaEstimationMaxIterations(const int alpha_est_maxit)
Sets the maximum number of iterations for the alpha estimation procedure.
Definition slope.cpp:561
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:535
void setScaling(const std::string &type)
Sets the scaling type.
Definition slope.cpp:382
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:321
void setDevChangeTol(const double dev_change_tol)
Sets tolerance in deviance change for early stopping.
Definition slope.cpp:525
const std::string & getLossType()
Get currently defined loss type.
Definition slope.cpp:579
void setReturnClusters(const bool return_clusters)
Sets the update clusters flag.
Definition slope.cpp:403
void setMaxClusters(const int max_clusters)
Sets tolerance in deviance change for early stopping.
Definition slope.cpp:545
void setIntercept(bool intercept)
Sets the intercept flag.
Definition slope.cpp:341
bool getFitIntercept() const
Returns the intercept flag.
Definition slope.cpp:573
void setDiagnostics(const bool collect_diagnostics)
Toggles collection of diagnostics.
Definition slope.cpp:555
void setNormalization(const std::string &type)
Sets normalization type for the design matrix.
Definition slope.cpp:347
int getAlphaEstimationMaxIterations() const
Gets the maximum number of iterations allowed for the alpha estimation procedure.
Definition slope.cpp:567
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:512
void setOscarParameters(const double theta1, const double theta2)
Sets OSCAR parameters.
Definition slope.cpp:443
void setLambdaType(const std::string &lambda_type)
Sets the lambda type for regularization weights.
Definition slope.cpp:494
void setAlphaType(const std::string &alpha_type)
Sets the alpha type.
Definition slope.cpp:409
void setModifyX(const bool modify_x)
Controls if x should be modified-in-place.
Definition slope.cpp:519
void setCentering(const std::string &type)
Sets the center points for feature normalization.
Definition slope.cpp:368
void setLearningRateDecr(double learning_rate_decr)
Sets the learning rate decrement.
Definition slope.cpp:425
void setLoss(const std::string &loss_type)
Sets the loss function type.
Definition slope.cpp:503
void setPathLength(int path_length)
Sets the path length.
Definition slope.cpp:476
void setHybridCdIterations(int cd_iterations)
Sets the frequence of proximal gradient descent steps.
Definition slope.cpp:485
void setTol(double tol)
Sets the tolerance value.
Definition slope.cpp:458
void setQ(double q)
Sets the q value.
Definition slope.cpp:434
void setUpdateClusters(bool update_clusters)
Sets the update clusters flag.
Definition slope.cpp:397
Class representing the Sorted L1 Norm.
double eval(const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda) const
Evaluates the Sorted L1 Norm.
double dualNorm(const Eigen::VectorXd &a, const Eigen::ArrayXd &lambda) const
Computes the dual norm of a vector.
Timer class for measuring elapsed time with high resolution.
Definition timer.h:19
void start()
Starts the timer by recording the current time point.
Definition timer.cpp:6
void resume()
Resumes the timer after a pause.
Definition timer.cpp:25
double elapsed() const
Returns the elapsed time in seconds since start() was called.
Definition timer.cpp:34
void pause()
Pauses the timer.
Definition timer.cpp:14
static void addWarning(WarningCode code, const std::string &message)
Log a new warning.
Definition logger.cpp:31
The declaration of the Clusters class.
Definitions of constants used in libslope.
Diagnostics for SLOPE optimization.
Functions for estimating noise level and regularization parameter alpha.
Thread-safe warning logging facility for the slope library.
The declartion of the Objctive class and its subclasses, which represent the data-fitting part of the...
Mathematical support functions for the slope package.
constexpr double EPSILON
Small value used for floating-point comparisons to handle precision issues.
Definition constants.h:27
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
std::unique_ptr< Loss > setupLoss(const std::string &loss)
Factory function to create the appropriate loss function based on the distribution family.
std::unordered_set< double > unique(const Eigen::MatrixXd &x)
Create a set of unique values from an Eigen matrix.
Definition utils.h:272
void validateOption(const std::string &value, const std::set< std::string > &valid_options, const std::string &parameter_name)
Validates if a given value exists in a set of valid options.
Definition utils.cpp:7
std::unique_ptr< ScreeningRule > createScreeningRule(const std::string &screening_type)
Creates a screening rule based on the provided type.
Eigen::ArrayXd lambdaSequence(const int p, const double q, const std::string &type, const int n, const double theta1, const double theta2)
JitNormalization normalize(Eigen::MatrixXd &x, Eigen::VectorXd &x_centers, Eigen::VectorXd &x_scales, const std::string &centering_type, const std::string &scaling_type, const bool modify_x)
Definition normalize.cpp:6
double computeDual(const Eigen::VectorXd &beta, const Eigen::MatrixXd &residual, const std::unique_ptr< Loss > &loss, const SortedL1Norm &sl1_norm, const Eigen::ArrayXd &lambda, const MatrixType &x, const Eigen::MatrixXd &y, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization &jit_normalization, const bool intercept)
Computes the dual objective function value for SLOPE optimization.
Definition diagnostics.h:41
@ MAXIT_REACHED
Maximum iterations reached without convergence.
void updateGradient(Eigen::VectorXd &gradient, const T &x, const Eigen::MatrixXd &residual, const std::vector< int > &active_set, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::VectorXd &w, const JitNormalization jit_normalization)
Definition math.h:223
void offsetGradient(Eigen::VectorXd &gradient, const T &x, const Eigen::VectorXd &offset, const std::vector< int > &active_set, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization jit_normalization)
Definition math.h:295
Eigen::ArrayXd regularizationPath(const Eigen::ArrayXd &alpha_in, const int path_length, double alpha_min_ratio, double alpha_max)
SlopePath estimateAlpha(MatrixType &x, Eigen::MatrixXd &y, const Slope &model)
Estimates the regularization parameter alpha for SLOPE regression.
std::unique_ptr< SolverBase > setupSolver(const std::string &solver_type, const std::string &loss, JitNormalization jit_normalization, bool intercept, bool update_clusters, int cd_iterations)
Factory function to create and configure a SLOPE solver.
int whichMax(const T &x)
Returns the index of the maximum element in a container.
Definition math.h:365
std::tuple< Eigen::VectorXd, Eigen::MatrixXd > rescaleCoefficients(const Eigen::VectorXd &beta0, const Eigen::VectorXd &beta, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const bool intercept)
Rescales the coefficients using the given parameters.
Definition normalize.cpp:94
Functions to normalize the design matrix and rescale coefficients in case the design was normalized.
Functions for generating regularization sequences for SLOPE.
Screening rules for SLOPE regression optimization.
Factory function to create the appropriate loss function based on.
Factory function to create and configure a SLOPE solver.
SLOPE (Sorted L-One Penalized Estimation) optimization.
The declaration of the SortedL1Norm class.
Simple high-resolution timer class for performance measurements.
Various utility functions.