slope 6.2.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
41inline bool
43{
44 return false;
45}
46
55class Slope
56{
57public:
63 Slope() = default;
64
73 void setSolver(const std::string& solver);
74
80 void setIntercept(bool intercept);
81
89 void setNormalization(const std::string& type);
90
97 void setUpdateClusters(bool update_clusters);
98
105 void setReturnClusters(const bool return_clusters);
106
114 void setAlphaMinRatio(double alpha_min_ratio);
115
130 void setAlphaType(const std::string& alpha_type);
131
138 void setLearningRateDecr(double learning_rate_decr);
139
146 void setQ(double q);
147
154 void setOscarParameters(const double theta1, const double theta2);
155
161 void setTol(double tol);
162
168 void setRelaxTol(double tol);
169
178
186
194 void setMaxIterations(int max_it);
195
201 void setPathLength(int path_length);
202
210 void setHybridCdIterations(int cd_iterations);
211
217 void setHybridCdType(const std::string& cd_type);
218
226 void setLambdaType(const std::string& lambda_type);
227
238 void setLoss(const std::string& loss_type);
239
249 void setScreening(const std::string& screening_type);
250
258 void setModifyX(const bool modify_x);
259
264 void setDevChangeTol(const double dev_change_tol);
265
271 void setDevRatioTol(const double dev_ratio_tol);
272
280 void setMaxClusters(const int max_clusters);
281
286 void setCentering(const std::string& type);
287
293 void setCentering(const Eigen::VectorXd& x_centers);
294
299 void setScaling(const std::string& type);
300
307 void setDiagnostics(const bool collect_diagnostics);
308
314 void setScaling(const Eigen::VectorXd& x_scales);
315
326 void setAlphaEstimationMaxIterations(const int alpha_est_maxit);
327
333 void setRandomSeed(const int seed);
334
340 void setRandomSeed(std::optional<int> seed);
341
345 bool hasRandomSeed() const;
346
350 int getRandomSeed() const;
351
357
361 bool getFitIntercept() const;
362
367 const std::string& getLossType();
368
387 template<typename T>
389 Eigen::EigenBase<T>& x,
390 const Eigen::MatrixXd& y_in,
391 Eigen::ArrayXd alpha = Eigen::ArrayXd::Zero(0),
392 Eigen::ArrayXd lambda = Eigen::ArrayXd::Zero(0),
393 std::function<bool()> check_interrupt = defaultInterruptChecker)
394 {
395 using Eigen::MatrixXd;
396 using Eigen::VectorXd;
397
398 const int n = x.rows();
399 const int p = x.cols();
400
401 const int INTERRUPT_FREQ = 100;
402 bool interrupt = false;
403
404 if (n != y_in.rows()) {
405 throw std::invalid_argument(
406 "x and y_in must have the same number of rows");
407 }
408
409 auto jit_normalization = normalize(x.derived(),
410 this->x_centers,
411 this->x_scales,
412 this->centering_type,
413 this->scaling_type,
414 this->modify_x);
415
416 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
417
418 MatrixXd y = loss->preprocessResponse(y_in);
419
420 const int m = y.cols();
421
422 std::vector<int> full_set(p * m);
423 std::iota(full_set.begin(), full_set.end(), 0);
424
425 VectorXd beta0 = VectorXd::Zero(m);
426 VectorXd beta = VectorXd::Zero(p * m);
427
428 MatrixXd eta = MatrixXd::Zero(n, m); // linear predictor
429
430 if (this->intercept) {
431 beta0 = loss->link(y.colwise().mean()).transpose();
432 eta.rowwise() = beta0.transpose();
433 }
434
435 MatrixXd residual = loss->residual(eta, y);
436 VectorXd gradient(beta.size());
437
438 // Path data
439 bool user_alpha = alpha.size() > 0;
440 bool user_lambda = lambda.size() > 0;
441
442 if (!user_lambda) {
443 lambda = lambdaSequence(
444 p * m, this->q, this->lambda_type, n, this->theta1, this->theta2);
445 } else {
446 if (lambda.size() != beta.size()) {
447 throw std::invalid_argument(
448 "lambda must be the same length as the number of coefficients");
449 }
450 if (lambda.minCoeff() < 0) {
451 throw std::invalid_argument("lambda must be non-negative");
452 }
453 if (!lambda.isFinite().all()) {
454 throw std::invalid_argument("lambda must be finite");
455 }
456 }
457
458 // Setup the regularization sequence and path
459 SortedL1Norm sl1_norm;
460
461 // TODO: Make this part of the slope class
462 auto solver = setupSolver(this->solver_type,
463 this->loss_type,
464 jit_normalization,
465 this->intercept,
466 this->update_clusters,
467 this->cd_iterations,
468 this->cd_type,
469 this->random_seed);
470
471 updateGradient(gradient,
472 x.derived(),
473 residual,
474 full_set,
475 this->x_centers,
476 this->x_scales,
477 Eigen::VectorXd::Ones(n),
478 jit_normalization);
479
480 int alpha_max_ind = whichMax(gradient.cwiseAbs());
481 double alpha_max = sl1_norm.dualNorm(gradient, lambda);
482
483 if (alpha_type == "path" ||
484 (alpha_type == "estimate" && alpha_estimate != 1)) {
485 if (alpha_min_ratio < 0) {
486 alpha_min_ratio = n > gradient.size() ? 1e-4 : 1e-2;
487 }
488
489 alpha =
490 regularizationPath(alpha, path_length, alpha_min_ratio, alpha_max);
491 path_length = alpha.size();
492 } else if (alpha_type == "estimate" && alpha_estimate == -1) {
493 if (loss_type != "quadratic") {
494 throw std::invalid_argument("Automatic alpha estimation is only "
495 "available for the quadratic loss");
496 }
497 }
498
499 // Screening setup
500 std::unique_ptr<ScreeningRule> screening_rule =
501 createScreeningRule(this->screening_type);
502 std::vector<int> working_set =
503 screening_rule->initialize(full_set, alpha_max_ind);
504
505 // Path variables
506 double null_deviance = loss->deviance(eta, y);
507 double dev_prev = null_deviance;
508
509 Timer timer;
510
511 double alpha_prev = std::max(alpha_max, alpha(0));
512
513 std::vector<SlopeFit> fits;
514
515 // Regularization path loop
516 for (int path_step = 0; path_step < this->path_length; ++path_step) {
517 // Check for interrupt at the start of each path step
518 bool local_interrupt = false;
519#ifdef _OPENMP
520#pragma omp critical(check_interrupt)
521#endif
522 {
523 local_interrupt = check_interrupt();
524 }
525 if (local_interrupt) {
526 interrupt = true;
527 break;
528 }
529
530 double alpha_curr = alpha(path_step);
531
532 assert(alpha_curr <= alpha_prev && "Alpha must be decreasing");
533
534 Eigen::ArrayXd lambda_curr = alpha_curr * lambda;
535 Eigen::ArrayXd lambda_prev = alpha_prev * lambda;
536
537 std::vector<double> duals, primals, time;
538 timer.start();
539
540 // Update gradient for the full set
541 // TODO: Only update for non-working set since gradient is updated before
542 // the convergence check in the inner loop for the working set
543 updateGradient(gradient,
544 x.derived(),
545 residual,
546 full_set,
547 x_centers,
548 x_scales,
549 Eigen::VectorXd::Ones(x.rows()),
550 jit_normalization);
551
552 working_set = screening_rule->screen(
553 gradient, lambda_curr, lambda_prev, beta, full_set);
554
555 int it = 0;
556 int total_it = 0;
557 for (; it < this->max_it; ++it, ++total_it) {
558 // Compute primal, dual, and gap
559 residual = loss->residual(eta, y);
560 updateGradient(gradient,
561 x.derived(),
562 residual,
563 working_set,
564 this->x_centers,
565 this->x_scales,
566 Eigen::VectorXd::Ones(n),
567 jit_normalization);
568
569 double primal = loss->loss(eta, y) +
570 sl1_norm.eval(beta(working_set),
571 lambda_curr.head(working_set.size()));
572
573 MatrixXd theta = residual;
574
575 // First compute gradient with potential offset for intercept case
576 VectorXd dual_gradient = gradient;
577
578 // TODO: Can we avoid this copy? Maybe revert offset afterwards or,
579 // alternatively, solve intercept until convergence and then no longer
580 // need the offset at all.
581 if (this->intercept) {
582 VectorXd theta_mean = theta.colwise().mean();
583 theta.rowwise() -= theta_mean.transpose();
584
585 offsetGradient(dual_gradient,
586 x.derived(),
587 theta_mean,
588 working_set,
589 this->x_centers,
590 this->x_scales,
591 jit_normalization);
592 }
593
594 // Common scaling operation
595 double dual_norm = sl1_norm.dualNorm(
596 dual_gradient(working_set), lambda_curr.head(working_set.size()));
597 theta.array() /= std::max(1.0, dual_norm);
598
599 double dual = loss->dual(theta, y, Eigen::VectorXd::Ones(n));
600
601 if (collect_diagnostics) {
602 timer.pause();
603 double true_dual = computeDual(beta,
604 residual,
605 loss,
606 sl1_norm,
607 lambda_curr,
608 x.derived(),
609 y,
610 this->x_centers,
611 this->x_scales,
612 jit_normalization,
613 this->intercept);
614 timer.resume();
615
616 time.emplace_back(timer.elapsed());
617 primals.emplace_back(primal);
618 duals.emplace_back(true_dual);
619 }
620
621 double dual_gap = primal - dual;
622
623 assert(dual_gap > -1e-6 && "Dual gap should be positive");
624
625 double tol_scaled = (std::abs(primal) + constants::EPSILON) * this->tol;
626
627 if (dual_gap <= tol_scaled || it == this->max_it) {
628 bool no_violations =
629 screening_rule->checkKktViolations(gradient,
630 beta,
631 lambda_curr,
632 working_set,
633 x.derived(),
634 residual,
635 this->x_centers,
636 this->x_scales,
637 jit_normalization,
638 full_set);
639 if (no_violations) {
640 break;
641 } else {
642 it = 0; // Restart if there are KKT violations
643 }
644 }
645
646 if (it % INTERRUPT_FREQ == 0) {
647 bool local_interrupt = false;
648#ifdef _OPENMP
649#pragma omp critical(check_interrupt)
650#endif
651 {
652 local_interrupt = check_interrupt();
653 }
654 if (local_interrupt) {
655 interrupt = true;
656 break;
657 }
658 }
659
660 solver->run(beta0,
661 beta,
662 eta,
663 lambda_curr,
664 loss,
665 sl1_norm,
666 gradient,
667 working_set,
668 x.derived(),
669 this->x_centers,
670 this->x_scales,
671 y);
672 }
673
674 if (it == this->max_it) {
677 "Maximum number of iterations reached at step = " +
678 std::to_string(path_step) + ".");
679 }
680
681 alpha_prev = alpha_curr;
682
683 // Compute early stopping criteria
684 double dev = loss->deviance(eta, y);
685 double dev_ratio = 1 - dev / null_deviance;
686 double dev_change = path_step == 0 ? 1.0 : 1 - dev / dev_prev;
687 dev_prev = dev;
688
689 Clusters clusters;
690
691 if (return_clusters) {
692 clusters.update(beta);
693 }
694
695 SlopeFit fit{ beta0,
696 beta.reshaped(p, m).sparseView(),
697 clusters,
698 alpha_curr,
699 lambda,
700 dev,
701 null_deviance,
702 primals,
703 duals,
704 time,
705 total_it,
706 this->centering_type,
707 this->scaling_type,
708 this->intercept,
709 this->x_centers,
710 this->x_scales };
711
712 fits.emplace_back(std::move(fit));
713
714 if (interrupt) {
715 break;
716 }
717
718 if (!user_alpha) {
719 int n_unique = unique(beta.cwiseAbs()).size();
720 if (dev_ratio > dev_ratio_tol || dev_change < dev_change_tol ||
721 n_unique >= this->max_clusters.value_or(n + 1)) {
722 break;
723 }
724 }
725 }
726
727 return fits;
728 }
729
747 template<typename T>
749 Eigen::EigenBase<T>& x,
750 const Eigen::MatrixXd& y_in,
751 const double alpha = 1.0,
752 Eigen::ArrayXd lambda = Eigen::ArrayXd::Zero(0),
753 std::function<bool()> check_interrupt = defaultInterruptChecker)
754 {
755 Eigen::ArrayXd alpha_arr(1);
756 alpha_arr(0) = alpha;
757 SlopePath res = path(x, y_in, alpha_arr, lambda, check_interrupt);
758
759 return { res(0) };
760 };
761
786 template<typename T>
788 Eigen::EigenBase<T>& x,
789 Eigen::MatrixXd& y,
790 std::function<bool()> check_interrupt = defaultInterruptChecker)
791 {
792 int n = x.rows();
793 int p = x.cols();
794
795 // Create a copy with alpha type set to path to avoid recursion
796 Slope model_copy = *this;
797 model_copy.setAlphaType("path");
798
799 std::vector<int> selected;
800 Eigen::ArrayXd alpha(1);
801 SlopePath result;
802
803 // Estimate the noise level, if possible
804 if (n >= p + 30) {
805 alpha(0) = estimateNoise(x, y, this->intercept) / n;
806 this->alpha_estimate = alpha(0);
807 result =
808 model_copy.path(x, y, alpha, Eigen::ArrayXd::Zero(0), check_interrupt);
809 } else {
810 for (int it = 0; it < this->alpha_est_maxit; ++it) {
811 T x_selected = subsetCols(x.derived(), selected);
812
813 std::vector<int> selected_prev = selected;
814 selected.clear();
815
816 alpha(0) = estimateNoise(x_selected, y, this->intercept) / n;
817 this->alpha_estimate = alpha(0);
818
819 result = model_copy.path(
820 x, y, alpha, Eigen::ArrayXd::Zero(0), check_interrupt);
821 auto coefs = result.getCoefs().back();
822
823 for (typename Eigen::SparseMatrix<double>::InnerIterator it(coefs, 0);
824 it;
825 ++it) {
826 selected.emplace_back(it.row());
827 }
828
829 if (selected == selected_prev) {
830 return result;
831 }
832
833 if (static_cast<int>(selected.size()) >= n + this->intercept) {
834 throw std::runtime_error(
835 "selected >= n - 1 variables, cannot estimate variance");
836 }
837 }
838
841 "Maximum iterations reached in alpha estimation");
842 }
843
844 return result;
845 }
846
861 template<typename T>
863 T& x,
864 const Eigen::VectorXd& y_in,
865 const double gamma = 0.0,
866 Eigen::VectorXd beta0 = Eigen::VectorXd(0),
867 Eigen::VectorXd beta = Eigen::VectorXd(0))
868 {
869 using Eigen::MatrixXd;
870 using Eigen::VectorXd;
871
872 int n = x.rows();
873 int p = x.cols();
874
875 if (beta0.size() == 0) {
876 beta0 = fit.getIntercepts(false);
877 }
878
879 if (beta.size() == 0) {
880 beta = fit.getCoefs(false);
881 }
882
883 double alpha = fit.getAlpha();
884
885 Timer timer;
886
887 std::vector<double> primals, duals, time;
888 timer.start();
889
890 auto jit_normalization =
891 normalize(x, x_centers, x_scales, centering_type, scaling_type, modify_x);
892
893 bool update_clusters = false;
894
895 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
896
897 MatrixXd y = loss->preprocessResponse(y_in);
898
899 int m = y.cols();
900
901 Eigen::ArrayXd lambda_cumsum_relax = Eigen::ArrayXd::Zero(p * m + 1);
902
903 auto working_set = activeSet(beta);
904
905 Eigen::MatrixXd eta = linearPredictor(x,
906 working_set,
907 beta0,
908 beta,
909 x_centers,
910 x_scales,
911 jit_normalization,
912 intercept);
913 VectorXd gradient = VectorXd::Zero(p * m);
914 MatrixXd residual(n, m);
915 MatrixXd working_residual(n, m);
916
917 MatrixXd w = MatrixXd::Ones(n, m);
918 MatrixXd w_ones = MatrixXd::Ones(n, m);
919 MatrixXd z = y;
920
921 slope::Clusters clusters(beta);
922
923 std::mt19937 rng;
924
925 if (random_seed.has_value()) {
926 rng.seed(*random_seed);
927 } else {
928 rng.seed(std::random_device{}());
929 }
930
931 int passes = 0;
932
933 for (int irls_it = 0; irls_it < max_it_outer_relax; irls_it++) {
934 residual = loss->residual(eta, y);
935
936 if (collect_diagnostics) {
937 primals.push_back(loss->loss(eta, y));
938 duals.push_back(0.0);
939 time.push_back(timer.elapsed());
940 }
941
942 Eigen::VectorXd cluster_gradient = clusterGradient(beta,
943 residual,
944 clusters,
945 x,
946 w_ones,
947 x_centers,
948 x_scales,
949 jit_normalization);
950
951 double norm_grad = cluster_gradient.lpNorm<Eigen::Infinity>();
952
953 if (norm_grad < tol_relax) {
954 break;
955 }
956
957 loss->updateWeightsAndWorkingResponse(w, z, eta, y);
958 working_residual = eta - z;
959
960 for (int inner_it = 0; inner_it < max_it_inner_relax; ++inner_it) {
961 passes++;
962
963 double max_abs_gradient = coordinateDescent(beta0,
964 beta,
965 working_residual,
966 clusters,
967 lambda_cumsum_relax,
968 x,
969 w,
970 x_centers,
971 x_scales,
972 intercept,
973 jit_normalization,
974 update_clusters,
975 rng,
976 cd_type);
977
978 if (max_abs_gradient < tol_relax) {
979 break;
980 }
981 }
982
983 eta = working_residual + z;
984
985 if (irls_it == max_it_outer_relax) {
987 "Maximum number of IRLS iterations reached.");
988 }
989 }
990
991 double dev = loss->deviance(eta, y);
992
993 if (gamma > 0) {
994 Eigen::VectorXd old_coefs = fit.getCoefs(false);
995 Eigen::VectorXd old_intercept = fit.getIntercepts(false);
996 beta = (1 - gamma) * beta + gamma * old_coefs;
997 }
998
999 SlopeFit fit_out{ beta0,
1000 beta.reshaped(p, m).sparseView(),
1001 clusters,
1002 alpha,
1003 fit.getLambda(),
1004 dev,
1006 primals,
1007 duals,
1008 time,
1009 passes,
1010 centering_type,
1011 scaling_type,
1012 intercept,
1013 x_centers,
1014 x_scales };
1015
1016 return fit_out;
1017 }
1018
1031 template<typename T>
1033 T& x,
1034 const Eigen::VectorXd& y,
1035 const double gamma = 0.0)
1036 {
1037 std::vector<SlopeFit> fits;
1038
1039 for (size_t i = 0; i < path.size(); i++) {
1040 // TODO: Reinstate warm starts. Need to be careful about
1041 // the warm started values though since they have to
1042 // agree with the cluster or we will run into trouble.
1043 // We can probably fix this by using the signs
1044 // of the cluster object rather than the betas though.
1045 auto relaxed_fit = relax(path(i), x, y, gamma);
1046
1047 fits.emplace_back(relaxed_fit);
1048 }
1049
1050 return fits;
1051 }
1052
1053private:
1054 // Parameters
1055 bool collect_diagnostics = false;
1056 bool intercept = true;
1057 bool modify_x = false;
1058 bool return_clusters = true;
1059 bool update_clusters = true;
1060 double alpha_min_ratio = -1;
1061 double dev_change_tol = 1e-5;
1062 double dev_ratio_tol = 0.999;
1063 double learning_rate_decr = 0.5;
1064 double q = 0.1;
1065 double theta1 = 1.0;
1066 double theta2 = 0.5;
1067 double tol = 1e-4;
1068 double tol_relax = 1e-4;
1069 double alpha_estimate = -1;
1070 int alpha_est_maxit = 1000;
1071 int cd_iterations = 10;
1072 int max_it = 1e5;
1073 int max_it_inner_relax = 1e5;
1074 int max_it_outer_relax = 50;
1075 int path_length = 100;
1076 std::optional<int> max_clusters = std::nullopt;
1077 std::optional<int> random_seed = 0;
1078 std::string alpha_type = "path";
1079 std::string cd_type = "permuted";
1080 std::string centering_type = "mean";
1081 std::string lambda_type = "bh";
1082 std::string loss_type = "quadratic";
1083 std::string scaling_type = "sd";
1084 std::string screening_type = "strong";
1085 std::string solver_type = "auto";
1086
1087 // Data
1088 Eigen::VectorXd x_centers;
1089 Eigen::VectorXd x_scales;
1090};
1091
1092} // 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:56
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.
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), std::function< bool()> check_interrupt=defaultInterruptChecker)
Computes SLOPE regression solution path for multiple alpha and lambda values.
Definition slope.h:388
SlopePath estimateAlpha(Eigen::EigenBase< T > &x, Eigen::MatrixXd &y, std::function< bool()> check_interrupt=defaultInterruptChecker)
Estimates the regularization parameter alpha for SLOPE regression.
Definition slope.h:787
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.
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:862
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.
SlopeFit fit(Eigen::EigenBase< T > &x, const Eigen::MatrixXd &y_in, const double alpha=1.0, Eigen::ArrayXd lambda=Eigen::ArrayXd::Zero(0), std::function< bool()> check_interrupt=defaultInterruptChecker)
Fits a single SLOPE regression model for given alpha and lambda values.
Definition slope.h:748
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
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:1032
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:332
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.
bool defaultInterruptChecker()
Default no-op interrupt checker.
Definition slope.h:42
std::unordered_set< double > unique(const Eigen::MatrixXd &x)
Create a set of unique values from an Eigen matrix.
Definition utils.h:380
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:331
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:153
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:234
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:306
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:804
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:376
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.