slope 6.5.0
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 if (!isFinite(x.derived())) {
410 throw std::invalid_argument("x must not contain NA, NaN, or Inf values");
411 }
412
413 if (!y_in.array().isFinite().all()) {
414 throw std::invalid_argument("y must not contain NA, NaN, or Inf values");
415 }
416
417 auto jit_normalization = normalize(x.derived(),
418 this->x_centers,
419 this->x_scales,
420 this->centering_type,
421 this->scaling_type,
422 this->modify_x);
423
424 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
425
426 MatrixXd y = loss->preprocessResponse(y_in);
427
428 const int m = y.cols();
429
430 std::vector<int> full_set(p * m);
431 std::iota(full_set.begin(), full_set.end(), 0);
432
433 VectorXd beta0 = VectorXd::Zero(m);
434 VectorXd beta = VectorXd::Zero(p * m);
435
436 MatrixXd eta = MatrixXd::Zero(n, m); // linear predictor
437
438 if (this->intercept) {
439 beta0 = loss->link(y.colwise().mean()).transpose();
440 eta.rowwise() = beta0.transpose();
441 }
442
443 MatrixXd residual = loss->residual(eta, y);
444 VectorXd gradient(beta.size());
445
446 // Path data
447 bool user_alpha = alpha.size() > 0;
448 bool user_lambda = lambda.size() > 0;
449
450 if (!user_lambda) {
451 lambda = lambdaSequence(
452 p * m, this->q, this->lambda_type, n, this->theta1, this->theta2);
453 } else {
454 if (lambda.size() != beta.size()) {
455 throw std::invalid_argument(
456 "lambda must be the same length as the number of coefficients");
457 }
458 if (lambda.minCoeff() < 0) {
459 throw std::invalid_argument("lambda must be non-negative");
460 }
461 if (!lambda.isFinite().all()) {
462 throw std::invalid_argument("lambda must be finite");
463 }
464 // Check that lambda is in decreasing order
465 for (int i = 1; i < lambda.size(); ++i) {
466 if (lambda(i) > lambda(i - 1)) {
467 throw std::invalid_argument(
468 "lambda must be in decreasing order");
469 }
470 }
471 }
472
473 // Setup the regularization sequence and path
474 SortedL1Norm sl1_norm;
475
476 // TODO: Make this part of the slope class
477 auto solver = setupSolver(this->solver_type,
478 this->loss_type,
479 jit_normalization,
480 this->intercept,
481 this->update_clusters,
482 this->cd_iterations,
483 this->cd_type,
484 this->random_seed);
485
486 updateGradient(gradient,
487 x.derived(),
488 residual,
489 full_set,
490 this->x_centers,
491 this->x_scales,
492 Eigen::VectorXd::Ones(n),
493 jit_normalization);
494
495 int alpha_max_ind = whichMax(gradient.cwiseAbs());
496 double alpha_max = sl1_norm.dualNorm(gradient, lambda);
497
498 if (alpha_type == "path" ||
499 (alpha_type == "estimate" && alpha_estimate != 1)) {
500 if (alpha_min_ratio < 0) {
501 alpha_min_ratio = n > gradient.size() ? 1e-4 : 1e-2;
502 }
503
504 alpha =
505 regularizationPath(alpha, path_length, alpha_min_ratio, alpha_max);
506 path_length = alpha.size();
507 } else if (alpha_type == "estimate" && alpha_estimate == -1) {
508 if (loss_type != "quadratic") {
509 throw std::invalid_argument("Automatic alpha estimation is only "
510 "available for the quadratic loss");
511 }
512 }
513
514 // Screening setup
515 std::unique_ptr<ScreeningRule> screening_rule =
516 createScreeningRule(this->screening_type);
517 std::vector<int> working_set =
518 screening_rule->initialize(full_set, alpha_max_ind);
519
520 // Path variables
521 double null_deviance = loss->deviance(eta, y);
522 double dev_prev = null_deviance;
523
524 Timer timer;
525
526 double alpha_prev = std::max(alpha_max, alpha(0));
527
528 std::vector<SlopeFit> fits;
529
530 // Regularization path loop
531 for (int path_step = 0; path_step < this->path_length; ++path_step) {
532 // Check for interrupt at the start of each path step
533 bool local_interrupt = false;
534#ifdef _OPENMP
535#pragma omp critical(check_interrupt)
536#endif
537 {
538 local_interrupt = check_interrupt();
539 }
540 if (local_interrupt) {
541 interrupt = true;
542 break;
543 }
544
545 double alpha_curr = alpha(path_step);
546
547 assert(alpha_curr <= alpha_prev && "Alpha must be decreasing");
548
549 Eigen::ArrayXd lambda_curr = alpha_curr * lambda;
550 Eigen::ArrayXd lambda_prev = alpha_prev * lambda;
551
552 std::vector<double> duals, primals, time;
553 timer.start();
554
555 // Update gradient for the full set
556 // TODO: Only update for non-working set since gradient is updated before
557 // the convergence check in the inner loop for the working set
558 updateGradient(gradient,
559 x.derived(),
560 residual,
561 full_set,
562 x_centers,
563 x_scales,
564 Eigen::VectorXd::Ones(x.rows()),
565 jit_normalization);
566
567 working_set = screening_rule->screen(
568 gradient, lambda_curr, lambda_prev, beta, full_set);
569
570 int it = 0;
571 int total_it = 0;
572 for (; it < this->max_it; ++it, ++total_it) {
573 // Compute primal, dual, and gap
574 residual = loss->residual(eta, y);
575 updateGradient(gradient,
576 x.derived(),
577 residual,
578 working_set,
579 this->x_centers,
580 this->x_scales,
581 Eigen::VectorXd::Ones(n),
582 jit_normalization);
583
584 double primal = loss->loss(eta, y) +
585 sl1_norm.eval(beta(working_set),
586 lambda_curr.head(working_set.size()));
587
588 MatrixXd theta = residual;
589
590 // First compute gradient with potential offset for intercept case
591 VectorXd dual_gradient = gradient;
592
593 // TODO: Can we avoid this copy? Maybe revert offset afterwards or,
594 // alternatively, solve intercept until convergence and then no longer
595 // need the offset at all.
596 if (this->intercept) {
597 VectorXd theta_mean = theta.colwise().mean();
598 theta.rowwise() -= theta_mean.transpose();
599
600 offsetGradient(dual_gradient,
601 x.derived(),
602 theta_mean,
603 working_set,
604 this->x_centers,
605 this->x_scales,
606 jit_normalization);
607 }
608
609 // Common scaling operation
610 double dual_norm = sl1_norm.dualNorm(
611 dual_gradient(working_set), lambda_curr.head(working_set.size()));
612 theta.array() /= std::max(1.0, dual_norm);
613
614 double dual = loss->dual(theta, y, Eigen::VectorXd::Ones(n));
615
616 if (collect_diagnostics) {
617 timer.pause();
618 double true_dual = computeDual(beta,
619 residual,
620 loss,
621 sl1_norm,
622 lambda_curr,
623 x.derived(),
624 y,
625 this->x_centers,
626 this->x_scales,
627 jit_normalization,
628 this->intercept);
629 timer.resume();
630
631 time.emplace_back(timer.elapsed());
632 primals.emplace_back(primal);
633 duals.emplace_back(true_dual);
634 }
635
636 double dual_gap = primal - dual;
637
638 assert(dual_gap > -1e-6 && "Dual gap should be positive");
639
640 double tol_scaled = (std::abs(primal) + constants::EPSILON) * this->tol;
641
642 if (dual_gap <= tol_scaled || it == this->max_it) {
643 bool no_violations =
644 screening_rule->checkKktViolations(gradient,
645 beta,
646 lambda_curr,
647 working_set,
648 x.derived(),
649 residual,
650 this->x_centers,
651 this->x_scales,
652 jit_normalization,
653 full_set);
654 if (no_violations) {
655 break;
656 } else {
657 it = 0; // Restart if there are KKT violations
658 }
659 }
660
661 if (it % INTERRUPT_FREQ == 0) {
662 bool local_interrupt = false;
663#ifdef _OPENMP
664#pragma omp critical(check_interrupt)
665#endif
666 {
667 local_interrupt = check_interrupt();
668 }
669 if (local_interrupt) {
670 interrupt = true;
671 break;
672 }
673 }
674
675 solver->run(beta0,
676 beta,
677 eta,
678 lambda_curr,
679 loss,
680 sl1_norm,
681 gradient,
682 working_set,
683 x.derived(),
684 this->x_centers,
685 this->x_scales,
686 y);
687 }
688
689 if (it == this->max_it) {
692 "Maximum number of iterations reached at step = " +
693 std::to_string(path_step) + ".");
694 }
695
696 alpha_prev = alpha_curr;
697
698 // Compute early stopping criteria
699 double dev = loss->deviance(eta, y);
700 double dev_ratio = 1 - dev / null_deviance;
701 double dev_change = path_step == 0 ? 1.0 : 1 - dev / dev_prev;
702 dev_prev = dev;
703
704 Clusters clusters;
705
706 if (return_clusters) {
707 clusters.update(beta);
708 }
709
710 SlopeFit fit{ beta0,
711 beta.reshaped(p, m).sparseView(),
712 clusters,
713 alpha_curr,
714 lambda,
715 dev,
716 null_deviance,
717 primals,
718 duals,
719 time,
720 total_it,
721 this->centering_type,
722 this->scaling_type,
723 this->intercept,
724 this->x_centers,
725 this->x_scales };
726
727 fits.emplace_back(std::move(fit));
728
729 if (interrupt) {
730 break;
731 }
732
733 if (!user_alpha) {
734 int n_unique = unique(beta.cwiseAbs()).size();
735 if (dev_ratio > dev_ratio_tol || dev_change < dev_change_tol ||
736 n_unique >= this->max_clusters.value_or(n + 1)) {
737 break;
738 }
739 }
740 }
741
742 return fits;
743 }
744
762 template<typename T>
764 Eigen::EigenBase<T>& x,
765 const Eigen::MatrixXd& y_in,
766 const double alpha = 1.0,
767 Eigen::ArrayXd lambda = Eigen::ArrayXd::Zero(0),
768 std::function<bool()> check_interrupt = defaultInterruptChecker)
769 {
770 Eigen::ArrayXd alpha_arr(1);
771 alpha_arr(0) = alpha;
772 SlopePath res = path(x, y_in, alpha_arr, lambda, check_interrupt);
773
774 return { res(0) };
775 };
776
801 template<typename T>
803 Eigen::EigenBase<T>& x,
804 Eigen::MatrixXd& y,
805 std::function<bool()> check_interrupt = defaultInterruptChecker)
806 {
807 int n = x.rows();
808 int p = x.cols();
809
810 // Create a copy with alpha type set to path to avoid recursion
811 Slope model_copy = *this;
812 model_copy.setAlphaType("path");
813
814 std::vector<int> selected;
815 Eigen::ArrayXd alpha(1);
816 SlopePath result;
817
818 // Estimate the noise level, if possible
819 if (n >= p + 30) {
820 alpha(0) = estimateNoise(x, y, this->intercept) / n;
821 this->alpha_estimate = alpha(0);
822 result =
823 model_copy.path(x, y, alpha, Eigen::ArrayXd::Zero(0), check_interrupt);
824 } else {
825 for (int it = 0; it < this->alpha_est_maxit; ++it) {
826 T x_selected = subsetCols(x.derived(), selected);
827
828 std::vector<int> selected_prev = selected;
829 selected.clear();
830
831 alpha(0) = estimateNoise(x_selected, y, this->intercept) / n;
832 this->alpha_estimate = alpha(0);
833
834 result = model_copy.path(
835 x, y, alpha, Eigen::ArrayXd::Zero(0), check_interrupt);
836 auto coefs = result.getCoefs().back();
837
838 for (typename Eigen::SparseMatrix<double>::InnerIterator it(coefs, 0);
839 it;
840 ++it) {
841 selected.emplace_back(it.row());
842 }
843
844 if (selected == selected_prev) {
845 return result;
846 }
847
848 if (static_cast<int>(selected.size()) >= n + this->intercept) {
849 throw std::runtime_error(
850 "selected >= n - 1 variables, cannot estimate variance");
851 }
852 }
853
856 "Maximum iterations reached in alpha estimation");
857 }
858
859 return result;
860 }
861
876 template<typename T>
878 T& x,
879 const Eigen::VectorXd& y_in,
880 const double gamma = 0.0,
881 Eigen::VectorXd beta0 = Eigen::VectorXd(0),
882 Eigen::VectorXd beta = Eigen::VectorXd(0))
883 {
884 using Eigen::MatrixXd;
885 using Eigen::VectorXd;
886
887 int n = x.rows();
888 int p = x.cols();
889
890 if (beta0.size() == 0) {
891 beta0 = fit.getIntercepts(false);
892 }
893
894 if (beta.size() == 0) {
895 beta = fit.getCoefs(false);
896 }
897
898 double alpha = fit.getAlpha();
899
900 Timer timer;
901
902 std::vector<double> primals, duals, time;
903 timer.start();
904
905 auto jit_normalization =
906 normalize(x, x_centers, x_scales, centering_type, scaling_type, modify_x);
907
908 bool update_clusters = false;
909
910 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
911
912 MatrixXd y = loss->preprocessResponse(y_in);
913
914 int m = y.cols();
915
916 Eigen::ArrayXd lambda_cumsum_relax = Eigen::ArrayXd::Zero(p * m + 1);
917
918 auto working_set = activeSet(beta);
919
920 Eigen::MatrixXd eta = linearPredictor(x,
921 working_set,
922 beta0,
923 beta,
924 x_centers,
925 x_scales,
926 jit_normalization,
927 intercept);
928 VectorXd gradient = VectorXd::Zero(p * m);
929 MatrixXd residual(n, m);
930 MatrixXd working_residual(n, m);
931
932 MatrixXd w = MatrixXd::Ones(n, m);
933 MatrixXd w_ones = MatrixXd::Ones(n, m);
934 MatrixXd z = y;
935
936 slope::Clusters clusters(beta);
937
938 std::mt19937 rng;
939
940 if (random_seed.has_value()) {
941 rng.seed(*random_seed);
942 } else {
943 rng.seed(std::random_device{}());
944 }
945
946 int passes = 0;
947
948 for (int irls_it = 0; irls_it < max_it_outer_relax; irls_it++) {
949 residual = loss->residual(eta, y);
950
951 if (collect_diagnostics) {
952 primals.push_back(loss->loss(eta, y));
953 duals.push_back(0.0);
954 time.push_back(timer.elapsed());
955 }
956
957 Eigen::VectorXd cluster_gradient = clusterGradient(beta,
958 residual,
959 clusters,
960 x,
961 w_ones,
962 x_centers,
963 x_scales,
964 jit_normalization);
965
966 double norm_grad = cluster_gradient.lpNorm<Eigen::Infinity>();
967
968 if (norm_grad < tol_relax) {
969 break;
970 }
971
972 loss->updateWeightsAndWorkingResponse(w, z, eta, y);
973 working_residual = eta - z;
974
975 for (int inner_it = 0; inner_it < max_it_inner_relax; ++inner_it) {
976 passes++;
977
978 double max_abs_gradient = coordinateDescent(beta0,
979 beta,
980 working_residual,
981 clusters,
982 lambda_cumsum_relax,
983 x,
984 w,
985 x_centers,
986 x_scales,
987 intercept,
988 jit_normalization,
989 update_clusters,
990 rng,
991 cd_type);
992
993 if (max_abs_gradient < tol_relax) {
994 break;
995 }
996 }
997
998 eta = working_residual + z;
999
1000 if (irls_it == max_it_outer_relax) {
1002 "Maximum number of IRLS iterations reached.");
1003 }
1004 }
1005
1006 double dev = loss->deviance(eta, y);
1007
1008 if (gamma > 0) {
1009 Eigen::VectorXd old_coefs = fit.getCoefs(false);
1010 Eigen::VectorXd old_intercept = fit.getIntercepts(false);
1011 beta = (1 - gamma) * beta + gamma * old_coefs;
1012 }
1013
1014 SlopeFit fit_out{ beta0,
1015 beta.reshaped(p, m).sparseView(),
1016 clusters,
1017 alpha,
1018 fit.getLambda(),
1019 dev,
1021 primals,
1022 duals,
1023 time,
1024 passes,
1025 centering_type,
1026 scaling_type,
1027 intercept,
1028 x_centers,
1029 x_scales };
1030
1031 return fit_out;
1032 }
1033
1046 template<typename T>
1048 T& x,
1049 const Eigen::VectorXd& y,
1050 const double gamma = 0.0)
1051 {
1052 std::vector<SlopeFit> fits;
1053
1054 for (size_t i = 0; i < path.size(); i++) {
1055 // TODO: Reinstate warm starts. Need to be careful about
1056 // the warm started values though since they have to
1057 // agree with the cluster or we will run into trouble.
1058 // We can probably fix this by using the signs
1059 // of the cluster object rather than the betas though.
1060 auto relaxed_fit = relax(path(i), x, y, gamma);
1061
1062 fits.emplace_back(relaxed_fit);
1063 }
1064
1065 return fits;
1066 }
1067
1068private:
1069 // Parameters
1070 bool collect_diagnostics = false;
1071 bool intercept = true;
1072 bool modify_x = false;
1073 bool return_clusters = true;
1074 bool update_clusters = true;
1075 double alpha_min_ratio = -1;
1076 double dev_change_tol = 1e-5;
1077 double dev_ratio_tol = 0.999;
1078 double learning_rate_decr = 0.5;
1079 double q = 0.1;
1080 double theta1 = 1.0;
1081 double theta2 = 0.5;
1082 double tol = 1e-4;
1083 double tol_relax = 1e-4;
1084 double alpha_estimate = -1;
1085 int alpha_est_maxit = 1000;
1086 int cd_iterations = 10;
1087 int max_it = 1e5;
1088 int max_it_inner_relax = 1e5;
1089 int max_it_outer_relax = 50;
1090 int path_length = 100;
1091 std::optional<int> max_clusters = std::nullopt;
1092 std::optional<int> random_seed = 0;
1093 std::string alpha_type = "path";
1094 std::string cd_type = "permuted";
1095 std::string centering_type = "mean";
1096 std::string lambda_type = "bh";
1097 std::string loss_type = "quadratic";
1098 std::string scaling_type = "sd";
1099 std::string screening_type = "strong";
1100 std::string solver_type = "auto";
1101
1102 // Data
1103 Eigen::VectorXd x_centers;
1104 Eigen::VectorXd x_scales;
1105};
1106
1107} // 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:32
std::size_t size() const
Gets the number of solutions in the path.
Definition slope_path.h:256
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:97
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:802
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:877
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:763
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:1047
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:333
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:381
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:332
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:154
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:235
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:307
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:805
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:377
bool isFinite(const Eigen::DenseBase< Derived > &x)
Check if all elements in a dense matrix are finite.
Definition utils.h:405
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.