slope 6.0.1
Loading...
Searching...
No Matches
slope.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "clusters.h"
9#include "constants.h"
10#include "diagnostics.h"
11#include "estimate_alpha.h"
12#include "logger.h"
13#include "losses/loss.h"
14#include "losses/setup_loss.h"
15#include "math.h"
16#include "normalize.h"
18#include "screening.h"
19#include "slope_fit.h"
20#include "slope_path.h"
21#include "solvers/hybrid_cd.h"
23#include "sorted_l1_norm.h"
24#include "timer.h"
25#include <Eigen/Core>
26#include <Eigen/SparseCore>
27#include <cassert>
28#include <memory>
29#include <numeric>
30#include <optional>
31
35namespace slope {
36
45class Slope
46{
47public:
53 Slope() = default;
54
63 void setSolver(const std::string& solver);
64
70 void setIntercept(bool intercept);
71
79 void setNormalization(const std::string& type);
80
87 void setUpdateClusters(bool update_clusters);
88
95 void setReturnClusters(const bool return_clusters);
96
104 void setAlphaMinRatio(double alpha_min_ratio);
105
120 void setAlphaType(const std::string& alpha_type);
121
128 void setLearningRateDecr(double learning_rate_decr);
129
136 void setQ(double q);
137
144 void setOscarParameters(const double theta1, const double theta2);
145
151 void setTol(double tol);
152
158 void setRelaxTol(double tol);
159
168
176
184 void setMaxIterations(int max_it);
185
191 void setPathLength(int path_length);
192
200 void setHybridCdIterations(int cd_iterations);
201
207 void setHybridCdType(const std::string& cd_type);
208
216 void setLambdaType(const std::string& lambda_type);
217
228 void setLoss(const std::string& loss_type);
229
239 void setScreening(const std::string& screening_type);
240
248 void setModifyX(const bool modify_x);
249
254 void setDevChangeTol(const double dev_change_tol);
255
261 void setDevRatioTol(const double dev_ratio_tol);
262
270 void setMaxClusters(const int max_clusters);
271
276 void setCentering(const std::string& type);
277
283 void setCentering(const Eigen::VectorXd& x_centers);
284
289 void setScaling(const std::string& type);
290
297 void setDiagnostics(const bool collect_diagnostics);
298
304 void setScaling(const Eigen::VectorXd& x_scales);
305
316 void setAlphaEstimationMaxIterations(const int alpha_est_maxit);
317
323 void setRandomSeed(const int seed);
324
330 void setRandomSeed(std::optional<int> seed);
331
335 bool hasRandomSeed() const;
336
340 int getRandomSeed() const;
341
347
351 bool getFitIntercept() const;
352
357 const std::string& getLossType();
358
375 template<typename T>
376 SlopePath path(Eigen::EigenBase<T>& x,
377 const Eigen::MatrixXd& y_in,
378 Eigen::ArrayXd alpha = Eigen::ArrayXd::Zero(0),
379 Eigen::ArrayXd lambda = Eigen::ArrayXd::Zero(0))
380 {
381 using Eigen::MatrixXd;
382 using Eigen::VectorXd;
383
384 const int n = x.rows();
385 const int p = x.cols();
386
387 if (n != y_in.rows()) {
388 throw std::invalid_argument(
389 "x and y_in must have the same number of rows");
390 }
391
392 auto jit_normalization = normalize(x.derived(),
393 this->x_centers,
394 this->x_scales,
395 this->centering_type,
396 this->scaling_type,
397 this->modify_x);
398
399 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
400
401 MatrixXd y = loss->preprocessResponse(y_in);
402
403 const int m = y.cols();
404
405 std::vector<int> full_set(p * m);
406 std::iota(full_set.begin(), full_set.end(), 0);
407
408 VectorXd beta0 = VectorXd::Zero(m);
409 VectorXd beta = VectorXd::Zero(p * m);
410
411 MatrixXd eta = MatrixXd::Zero(n, m); // linear predictor
412
413 if (this->intercept) {
414 beta0 = loss->link(y.colwise().mean()).transpose();
415 eta.rowwise() = beta0.transpose();
416 }
417
418 MatrixXd residual = loss->residual(eta, y);
419 VectorXd gradient(beta.size());
420
421 // Path data
422 bool user_alpha = alpha.size() > 0;
423 bool user_lambda = lambda.size() > 0;
424
425 if (!user_lambda) {
426 lambda = lambdaSequence(
427 p * m, this->q, this->lambda_type, n, this->theta1, this->theta2);
428 } else {
429 if (lambda.size() != beta.size()) {
430 throw std::invalid_argument(
431 "lambda must be the same length as the number of coefficients");
432 }
433 if (lambda.minCoeff() < 0) {
434 throw std::invalid_argument("lambda must be non-negative");
435 }
436 if (!lambda.isFinite().all()) {
437 throw std::invalid_argument("lambda must be finite");
438 }
439 }
440
441 // Setup the regularization sequence and path
442 SortedL1Norm sl1_norm;
443
444 // TODO: Make this part of the slope class
445 auto solver = setupSolver(this->solver_type,
446 this->loss_type,
447 jit_normalization,
448 this->intercept,
449 this->update_clusters,
450 this->cd_iterations,
451 this->cd_type,
452 this->random_seed);
453
454 updateGradient(gradient,
455 x.derived(),
456 residual,
457 full_set,
458 this->x_centers,
459 this->x_scales,
460 Eigen::VectorXd::Ones(n),
461 jit_normalization);
462
463 int alpha_max_ind = whichMax(gradient.cwiseAbs());
464 double alpha_max = sl1_norm.dualNorm(gradient, lambda);
465
466 if (alpha_type == "path" ||
467 (alpha_type == "estimate" && alpha_estimate != 1)) {
468 if (alpha_min_ratio < 0) {
469 alpha_min_ratio = n > gradient.size() ? 1e-4 : 1e-2;
470 }
471
472 alpha =
473 regularizationPath(alpha, path_length, alpha_min_ratio, alpha_max);
474 path_length = alpha.size();
475 } else if (alpha_type == "estimate" && alpha_estimate == -1) {
476 if (loss_type != "quadratic") {
477 throw std::invalid_argument("Automatic alpha estimation is only "
478 "available for the quadratic loss");
479 }
480 }
481
482 // Screening setup
483 std::unique_ptr<ScreeningRule> screening_rule =
484 createScreeningRule(this->screening_type);
485 std::vector<int> working_set =
486 screening_rule->initialize(full_set, alpha_max_ind);
487
488 // Path variables
489 double null_deviance = loss->deviance(eta, y);
490 double dev_prev = null_deviance;
491
492 Timer timer;
493
494 double alpha_prev = std::max(alpha_max, alpha(0));
495
496 std::vector<SlopeFit> fits;
497
498 // Regularization path loop
499 for (int path_step = 0; path_step < this->path_length; ++path_step) {
500 double alpha_curr = alpha(path_step);
501
502 assert(alpha_curr <= alpha_prev && "Alpha must be decreasing");
503
504 Eigen::ArrayXd lambda_curr = alpha_curr * lambda;
505 Eigen::ArrayXd lambda_prev = alpha_prev * lambda;
506
507 std::vector<double> duals, primals, time;
508 timer.start();
509
510 // Update gradient for the full set
511 // TODO: Only update for non-working set since gradient is updated before
512 // the convergence check in the inner loop for the working set
513 updateGradient(gradient,
514 x.derived(),
515 residual,
516 full_set,
517 x_centers,
518 x_scales,
519 Eigen::VectorXd::Ones(x.rows()),
520 jit_normalization);
521
522 working_set = screening_rule->screen(
523 gradient, lambda_curr, lambda_prev, beta, full_set);
524
525 int it = 0;
526 int total_it = 0;
527 for (; it < this->max_it; ++it, ++total_it) {
528 // Compute primal, dual, and gap
529 residual = loss->residual(eta, y);
530 updateGradient(gradient,
531 x.derived(),
532 residual,
533 working_set,
534 this->x_centers,
535 this->x_scales,
536 Eigen::VectorXd::Ones(n),
537 jit_normalization);
538
539 double primal = loss->loss(eta, y) +
540 sl1_norm.eval(beta(working_set),
541 lambda_curr.head(working_set.size()));
542
543 MatrixXd theta = residual;
544
545 // First compute gradient with potential offset for intercept case
546 VectorXd dual_gradient = gradient;
547
548 // TODO: Can we avoid this copy? Maybe revert offset afterwards or,
549 // alternatively, solve intercept until convergence and then no longer
550 // need the offset at all.
551 if (this->intercept) {
552 VectorXd theta_mean = theta.colwise().mean();
553 theta.rowwise() -= theta_mean.transpose();
554
555 offsetGradient(dual_gradient,
556 x.derived(),
557 theta_mean,
558 working_set,
559 this->x_centers,
560 this->x_scales,
561 jit_normalization);
562 }
563
564 // Common scaling operation
565 double dual_norm = sl1_norm.dualNorm(
566 dual_gradient(working_set), lambda_curr.head(working_set.size()));
567 theta.array() /= std::max(1.0, dual_norm);
568
569 double dual = loss->dual(theta, y, Eigen::VectorXd::Ones(n));
570
571 if (collect_diagnostics) {
572 timer.pause();
573 double true_dual = computeDual(beta,
574 residual,
575 loss,
576 sl1_norm,
577 lambda_curr,
578 x.derived(),
579 y,
580 this->x_centers,
581 this->x_scales,
582 jit_normalization,
583 this->intercept);
584 timer.resume();
585
586 time.emplace_back(timer.elapsed());
587 primals.emplace_back(primal);
588 duals.emplace_back(true_dual);
589 }
590
591 double dual_gap = primal - dual;
592
593 assert(dual_gap > -1e-6 && "Dual gap should be positive");
594
595 double tol_scaled = (std::abs(primal) + constants::EPSILON) * this->tol;
596
597 if (dual_gap <= tol_scaled || it == this->max_it) {
598 bool no_violations =
599 screening_rule->checkKktViolations(gradient,
600 beta,
601 lambda_curr,
602 working_set,
603 x.derived(),
604 residual,
605 this->x_centers,
606 this->x_scales,
607 jit_normalization,
608 full_set);
609 if (no_violations) {
610 break;
611 } else {
612 it = 0; // Restart if there are KKT violations
613 }
614 }
615
616 solver->run(beta0,
617 beta,
618 eta,
619 lambda_curr,
620 loss,
621 sl1_norm,
622 gradient,
623 working_set,
624 x.derived(),
625 this->x_centers,
626 this->x_scales,
627 y);
628 }
629
630 if (it == this->max_it) {
633 "Maximum number of iterations reached at step = " +
634 std::to_string(path_step) + ".");
635 }
636
637 alpha_prev = alpha_curr;
638
639 // Compute early stopping criteria
640 double dev = loss->deviance(eta, y);
641 double dev_ratio = 1 - dev / null_deviance;
642 double dev_change = path_step == 0 ? 1.0 : 1 - dev / dev_prev;
643 dev_prev = dev;
644
645 Clusters clusters;
646
647 if (return_clusters) {
648 clusters.update(beta);
649 }
650
651 SlopeFit fit{ beta0,
652 beta.reshaped(p, m).sparseView(),
653 clusters,
654 alpha_curr,
655 lambda,
656 dev,
657 null_deviance,
658 primals,
659 duals,
660 time,
661 total_it,
662 this->centering_type,
663 this->scaling_type,
664 this->intercept,
665 this->x_centers,
666 this->x_scales };
667
668 fits.emplace_back(std::move(fit));
669
670 if (!user_alpha) {
671 int n_unique = unique(beta.cwiseAbs()).size();
672 if (dev_ratio > dev_ratio_tol || dev_change < dev_change_tol ||
673 n_unique >= this->max_clusters.value_or(n + 1)) {
674 break;
675 }
676 }
677 }
678
679 return fits;
680 }
681
697 template<typename T>
698 SlopeFit fit(Eigen::EigenBase<T>& x,
699 const Eigen::MatrixXd& y_in,
700 const double alpha = 1.0,
701 Eigen::ArrayXd lambda = Eigen::ArrayXd::Zero(0))
702 {
703 Eigen::ArrayXd alpha_arr(1);
704 alpha_arr(0) = alpha;
705 SlopePath res = path(x, y_in, alpha_arr, lambda);
706
707 return { res(0) };
708 };
709
733 template<typename T>
734 SlopePath estimateAlpha(Eigen::EigenBase<T>& x, Eigen::MatrixXd& y)
735 {
736 int n = x.rows();
737 int p = x.cols();
738
739 // Create a copy with alpha type set to path to avoid recursion
740 Slope model_copy = *this;
741 model_copy.setAlphaType("path");
742
743 std::vector<int> selected;
744 Eigen::ArrayXd alpha(1);
745 SlopePath result;
746
747 // Estimate the noise level, if possible
748 if (n >= p + 30) {
749 alpha(0) = estimateNoise(x, y, this->intercept) / n;
750 this->alpha_estimate = alpha(0);
751 result = model_copy.path(x, y, alpha);
752 } else {
753 for (int it = 0; it < this->alpha_est_maxit; ++it) {
754 T x_selected = subsetCols(x.derived(), selected);
755
756 std::vector<int> selected_prev = selected;
757 selected.clear();
758
759 alpha(0) = estimateNoise(x_selected, y, this->intercept) / n;
760 this->alpha_estimate = alpha(0);
761
762 result = model_copy.path(x, y, alpha);
763 auto coefs = result.getCoefs().back();
764
765 for (typename Eigen::SparseMatrix<double>::InnerIterator it(coefs, 0);
766 it;
767 ++it) {
768 selected.emplace_back(it.row());
769 }
770
771 if (selected == selected_prev) {
772 return result;
773 }
774
775 if (static_cast<int>(selected.size()) >= n + this->intercept) {
776 throw std::runtime_error(
777 "selected >= n - 1 variables, cannot estimate variance");
778 }
779 }
780
783 "Maximum iterations reached in alpha estimation");
784 }
785
786 return result;
787 }
788
803 template<typename T>
805 T& x,
806 const Eigen::VectorXd& y_in,
807 const double gamma = 0.0,
808 Eigen::VectorXd beta0 = Eigen::VectorXd(0),
809 Eigen::VectorXd beta = Eigen::VectorXd(0))
810 {
811 using Eigen::MatrixXd;
812 using Eigen::VectorXd;
813
814 int n = x.rows();
815 int p = x.cols();
816
817 if (beta0.size() == 0) {
818 beta0 = fit.getIntercepts(false);
819 }
820
821 if (beta.size() == 0) {
822 beta = fit.getCoefs(false);
823 }
824
825 double alpha = fit.getAlpha();
826
827 Timer timer;
828
829 std::vector<double> primals, duals, time;
830 timer.start();
831
832 auto jit_normalization =
833 normalize(x, x_centers, x_scales, centering_type, scaling_type, modify_x);
834
835 bool update_clusters = false;
836
837 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
838
839 MatrixXd y = loss->preprocessResponse(y_in);
840
841 int m = y.cols();
842
843 Eigen::ArrayXd lambda_cumsum_relax = Eigen::ArrayXd::Zero(p * m + 1);
844
845 auto working_set = activeSet(beta);
846
847 Eigen::MatrixXd eta = linearPredictor(x,
848 working_set,
849 beta0,
850 beta,
851 x_centers,
852 x_scales,
853 jit_normalization,
854 intercept);
855 VectorXd gradient = VectorXd::Zero(p * m);
856 MatrixXd residual(n, m);
857 MatrixXd working_residual(n, m);
858
859 MatrixXd w = MatrixXd::Ones(n, m);
860 MatrixXd w_ones = MatrixXd::Ones(n, m);
861 MatrixXd z = y;
862
863 slope::Clusters clusters(beta);
864
865 std::mt19937 rng;
866
867 if (random_seed.has_value()) {
868 rng.seed(*random_seed);
869 } else {
870 rng.seed(std::random_device{}());
871 }
872
873 int passes = 0;
874
875 for (int irls_it = 0; irls_it < max_it_outer_relax; irls_it++) {
876 residual = loss->residual(eta, y);
877
878 if (collect_diagnostics) {
879 primals.push_back(loss->loss(eta, y));
880 duals.push_back(0.0);
881 time.push_back(timer.elapsed());
882 }
883
884 Eigen::VectorXd cluster_gradient = clusterGradient(beta,
885 residual,
886 clusters,
887 x,
888 w_ones,
889 x_centers,
890 x_scales,
891 jit_normalization);
892
893 double norm_grad = cluster_gradient.lpNorm<Eigen::Infinity>();
894
895 if (norm_grad < tol_relax) {
896 break;
897 }
898
899 loss->updateWeightsAndWorkingResponse(w, z, eta, y);
900 working_residual = eta - z;
901
902 for (int inner_it = 0; inner_it < max_it_inner_relax; ++inner_it) {
903 passes++;
904
905 double max_abs_gradient = coordinateDescent(beta0,
906 beta,
907 working_residual,
908 clusters,
909 lambda_cumsum_relax,
910 x,
911 w,
912 x_centers,
913 x_scales,
914 intercept,
915 jit_normalization,
916 update_clusters,
917 rng,
918 cd_type);
919
920 if (max_abs_gradient < tol_relax) {
921 break;
922 }
923 }
924
925 eta = working_residual + z;
926
927 if (irls_it == max_it_outer_relax) {
929 "Maximum number of IRLS iterations reached.");
930 }
931 }
932
933 double dev = loss->deviance(eta, y);
934
935 if (gamma > 0) {
936 Eigen::VectorXd old_coefs = fit.getCoefs(false);
937 Eigen::VectorXd old_intercept = fit.getIntercepts(false);
938 beta = (1 - gamma) * beta + gamma * old_coefs;
939 }
940
941 SlopeFit fit_out{ beta0,
942 beta.reshaped(p, m).sparseView(),
943 clusters,
944 alpha,
945 fit.getLambda(),
946 dev,
948 primals,
949 duals,
950 time,
951 passes,
952 centering_type,
953 scaling_type,
954 intercept,
955 x_centers,
956 x_scales };
957
958 return fit_out;
959 }
960
973 template<typename T>
975 T& x,
976 const Eigen::VectorXd& y,
977 const double gamma = 0.0)
978 {
979 std::vector<SlopeFit> fits;
980
981 for (size_t i = 0; i < path.size(); i++) {
982 // TODO: Reinstate warm starts. Need to be careful about
983 // the warm started values though since they have to
984 // agree with the cluster or we will run into trouble.
985 // We can probably fix this by using the signs
986 // of the cluster object rather than the betas though.
987 auto relaxed_fit = relax(path(i), x, y, gamma);
988
989 fits.emplace_back(relaxed_fit);
990 }
991
992 return fits;
993 }
994
995private:
996 // Parameters
997 bool collect_diagnostics = false;
998 bool intercept = true;
999 bool modify_x = false;
1000 bool return_clusters = true;
1001 bool update_clusters = true;
1002 double alpha_min_ratio = -1;
1003 double dev_change_tol = 1e-5;
1004 double dev_ratio_tol = 0.999;
1005 double learning_rate_decr = 0.5;
1006 double q = 0.1;
1007 double theta1 = 1.0;
1008 double theta2 = 0.5;
1009 double tol = 1e-4;
1010 double tol_relax = 1e-4;
1011 double alpha_estimate = -1;
1012 int alpha_est_maxit = 1000;
1013 int cd_iterations = 10;
1014 int max_it = 1e5;
1015 int max_it_inner_relax = 1e5;
1016 int max_it_outer_relax = 50;
1017 int path_length = 100;
1018 std::optional<int> max_clusters = std::nullopt;
1019 std::optional<int> random_seed = 0;
1020 std::string alpha_type = "path";
1021 std::string cd_type = "permuted";
1022 std::string centering_type = "mean";
1023 std::string lambda_type = "bh";
1024 std::string loss_type = "quadratic";
1025 std::string scaling_type = "sd";
1026 std::string screening_type = "strong";
1027 std::string solver_type = "auto";
1028
1029 // Data
1030 Eigen::VectorXd x_centers;
1031 Eigen::VectorXd x_scales;
1032};
1033
1034} // 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.
A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting.
Definition slope_fit.h:27
Eigen::VectorXd getIntercepts(const bool original_scale=true) const
Gets the intercept terms for this SLOPE fit.
Definition slope_fit.h:113
const Eigen::ArrayXd & getLambda() const
Gets the lambda (regularization) parameter used.
Definition slope_fit.h:156
double getAlpha() const
Gets the alpha (mixing) parameter used.
Definition slope_fit.h:161
Eigen::SparseMatrix< double > getCoefs(const bool original_scale=true) const
Gets the sparse coefficient matrix for this fit.
Definition slope_fit.h:133
double getNullDeviance() const
Gets the null model deviance.
Definition slope_fit.h:176
Container class for SLOPE regression solution paths.
Definition slope_path.h:31
std::size_t size() const
Gets the number of solutions in the path.
Definition slope_path.h:255
std::vector< Eigen::SparseMatrix< double > > getCoefs(const bool original_scale=true) const
Returns the vector of coefficient matrices for each solution in the path.
Definition slope_path.h:96
The SLOPE model.
Definition slope.h:46
void setSolver(const std::string &solver)
Sets the numerical solver used to fit the model.
void setAlphaMinRatio(double alpha_min_ratio)
Sets the alpha min ratio.
void setMaxIterations(int max_it)
Sets the maximum number of iterations.
void setAlphaEstimationMaxIterations(const int alpha_est_maxit)
Sets the maximum number of iterations for the alpha estimation procedure.
void setRelaxMaxInnerIterations(int max_it)
Sets the maximum number of inner iterations for the relaxed solver.
void setRandomSeed(std::optional< int > seed)
Sets the random seed.
int getRandomSeed() const
Gets the random seed.
void setDevRatioTol(const double dev_ratio_tol)
Sets tolerance in deviance change for early stopping.
void setScaling(const std::string &type)
Sets the scaling type.
void setRandomSeed(const int seed)
Sets the random seed.
void setDevChangeTol(const double dev_change_tol)
Sets tolerance in deviance change for early stopping.
const std::string & getLossType()
Get currently defined loss type.
SlopePath path(Eigen::EigenBase< 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.h:376
SlopePath estimateAlpha(Eigen::EigenBase< T > &x, Eigen::MatrixXd &y)
Estimates the regularization parameter alpha for SLOPE regression.
Definition slope.h:734
void setReturnClusters(const bool return_clusters)
Sets the return clusters flag.
void setRelaxMaxOuterIterations(int max_it)
Sets the maximum number of outer (IRLS) iterations for the relaxed solver.
void setMaxClusters(const int max_clusters)
Sets the maximum number of clusters.
void setIntercept(bool intercept)
Sets the intercept flag.
bool getFitIntercept() const
Returns the intercept flag.
void setDiagnostics(const bool collect_diagnostics)
Toggles collection of diagnostics.
void setNormalization(const std::string &type)
Sets normalization type for the design matrix.
SlopeFit relax(const SlopeFit &fit, T &x, const Eigen::VectorXd &y_in, const double gamma=0.0, Eigen::VectorXd beta0=Eigen::VectorXd(0), Eigen::VectorXd beta=Eigen::VectorXd(0))
Relaxes a fitted SLOPE model.
Definition slope.h:804
int getAlphaEstimationMaxIterations() const
Gets the maximum number of iterations allowed for the alpha estimation procedure.
void setScreening(const std::string &screening_type)
Sets the type of feature screening used, which discards predictors that are unlikely to be active.
void setOscarParameters(const double theta1, const double theta2)
Sets OSCAR parameters.
void setLambdaType(const std::string &lambda_type)
Sets the lambda type for regularization weights.
void setAlphaType(const std::string &alpha_type)
Sets the alpha type.
void setModifyX(const bool modify_x)
Controls if x should be modified-in-place.
void setCentering(const std::string &type)
Sets the center points for feature normalization.
Slope()=default
SlopeFit fit(Eigen::EigenBase< 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.h:698
void setHybridCdType(const std::string &cd_type)
Sets the frequence of proximal gradient descent steps.
bool hasRandomSeed() const
Checks if a random seed is set.
void setLearningRateDecr(double learning_rate_decr)
Sets the learning rate decrement.
void setLoss(const std::string &loss_type)
Sets the loss function type.
SlopePath relax(const SlopePath &path, T &x, const Eigen::VectorXd &y, const double gamma=0.0)
Relaxes a fitted SLOPE path.
Definition slope.h:974
void setScaling(const Eigen::VectorXd &x_scales)
Sets the scaling factors for feature normalization.
void setRelaxTol(double tol)
Sets the tolerance value for the relaxed SLOPE solver.
void setPathLength(int path_length)
Sets the path length.
void setCentering(const Eigen::VectorXd &x_centers)
Sets the center points for feature normalization.
void setHybridCdIterations(int cd_iterations)
Sets the frequence of proximal gradient descent steps.
void setTol(double tol)
Sets the tolerance value.
void setQ(double q)
Sets the q value.
void setUpdateClusters(bool update_clusters)
Sets the update clusters flag.
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.
void resume()
Resumes the timer after a pause.
double elapsed() const
Returns the elapsed time in seconds since start() was called.
void pause()
Pauses the timer.
static void addWarning(WarningCode code, const std::string &message)
Log a new warning.
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.
An implementation of the coordinate descent step in the hybrid algorithm for solving SLOPE.
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.h:11
Eigen::ArrayXd regularizationPath(const Eigen::ArrayXd &alpha_in, const int path_length, double alpha_min_ratio, const double alpha_max)
double coordinateDescent(Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::MatrixXd &residual, Clusters &clusters, const Eigen::ArrayXd &lambda_cumsum, const T &x, const Eigen::MatrixXd &w, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const bool intercept, const JitNormalization jit_normalization, const bool update_clusters, std::mt19937 &rng, const std::string &cd_type="cyclical")
Definition hybrid_cd.h:331
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, const std::string &cd_type, std::optional< int > random_seed=std::nullopt)
Factory function to create and configure a SLOPE solver.
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:348
double estimateNoise(Eigen::EigenBase< T > &x, Eigen::MatrixXd &y, const bool fit_intercept)
Estimates noise (standard error) in a linear model using OLS residuals.
T subsetCols(const Eigen::MatrixBase< T > &x, const std::vector< int > &indices)
Extract specified columns from a dense matrix.
Definition utils.h:299
std::unique_ptr< ScreeningRule > createScreeningRule(const std::string &screening_type)
Creates a screening rule based on the provided type.
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
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.
Eigen::ArrayXd lambdaSequence(const int p, const double q, const std::string &type, const int n=-1, const double theta1=1.0, const double theta2=1.0)
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:232
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:304
Eigen::VectorXd clusterGradient(Eigen::VectorXd &beta, Eigen::MatrixXd &residual, Clusters &clusters, const T &x, const Eigen::MatrixXd &w, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization jit_normalization)
Definition math.h:802
std::vector< int > activeSet(const Eigen::VectorXd &beta)
Identifies previously active variables.
JitNormalization normalize(Eigen::MatrixBase< T > &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.h:123
int whichMax(const T &x)
Returns the index of the maximum element in a container.
Definition math.h:374
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) fitting results.
Defines the SlopePath class for storing and accessing SLOPE regression solution paths.
The declaration of the SortedL1Norm class.
Simple high-resolution timer class for performance measurements.