slope 0.28.0
Loading...
Searching...
No Matches
slope Namespace Reference

Namespace containing SLOPE regression implementation. More...

Namespaces

namespace  constants
 Namespace containing constants used in the slope library.
 

Classes

class  Accuracy
 Classification Accuracy scoring metric. More...
 
class  AUC
 Area Under the ROC Curve (AUC-ROC) scoring metric. More...
 
class  Clusters
 Representation of the clusters in SLOPE. More...
 
struct  CvConfig
 Configuration settings for cross-validation. More...
 
struct  CvResult
 Contains overall results from a cross-validation process. More...
 
class  Deviance
 Deviance scoring metric. More...
 
class  Folds
 Manages data partitioning for cross-validation. More...
 
struct  GridResult
 Stores cross-validation results for a specific set of hyperparameters. More...
 
class  Hybrid
 Hybrid CD-PGD solver for SLOPE. More...
 
class  Logistic
 The Logistic class represents a logistic loss function. More...
 
class  Loss
 Loss function interface. More...
 
class  MAE
 Mean Absolute Error (MAE) scoring metric. More...
 
class  MaximizeScore
 Scoring metric that aims to maximize the score value. More...
 
class  MinimizeScore
 Scoring metric that aims to minimize the score value. More...
 
class  MisClass
 Misclassification Rate scoring metric. More...
 
class  MSE
 Mean Squared Error (MSE) scoring metric. More...
 
class  Multinomial
 The Multinomial class represents a multinomial logistic regression loss function. More...
 
class  NoScreening
 No screening rule - uses all variables. More...
 
class  PGD
 Proximal Gradient Descent solver for SLOPE optimization. More...
 
class  Poisson
 The Poisson class represents a Poisson regression loss function. More...
 
class  Quadratic
 Implementation of the Quadratic loss function. More...
 
class  Score
 Base class for scoring metrics used in regularized generalized linear regression. More...
 
class  ScreeningRule
 Base class for screening rules in SLOPE. More...
 
class  Slope
 The SLOPE model. More...
 
class  SlopeFit
 A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting. More...
 
class  SlopePath
 Container class for SLOPE regression solution paths. More...
 
class  SolverBase
 Abstract base class for SLOPE optimization solvers. More...
 
class  SortedL1Norm
 Class representing the Sorted L1 Norm. More...
 
class  StrongScreening
 Implements strong screening rules for SLOPE. More...
 
class  Threads
 Manages OpenMP thread settings across libslope. More...
 
class  Timer
 Timer class for measuring elapsed time with high resolution. More...
 
struct  Warning
 Structure representing a warning with its code and message. More...
 
class  WarningLogger
 Thread-safe warning logger for tracking runtime warnings. More...
 

Typedefs

using ArrayXb = Eigen::Array< bool, Eigen::Dynamic, 1 >
 Dynamic-size column vector of boolean values Wrapper around Eigen::Array<bool, Eigen::Dynamic, 1>
 
using ArrayXXb = Eigen::Array< bool, Eigen::Dynamic, Eigen::Dynamic >
 Dynamic-size matrix of boolean values Wrapper around Eigen::Array<bool, Eigen::Dynamic, Eigen::Dynamic>
 

Enumerations

enum class  JitNormalization { None = 0 , Center = 1 , Scale = 2 , Both = 3 }
 Enums to control predictor standardization behavior. More...
 
enum class  WarningCode { GENERIC_WARNING , DEPRECATED_FEATURE , MAXIT_REACHED , LINE_SEARCH_FAILED }
 Standard warning codes used throughout the slope library. More...
 

Functions

Eigen::SparseMatrix< int > patternMatrix (const Eigen::VectorXd &beta)
 
std::vector< std::map< std::string, double > > createGrid (const std::map< std::string, std::vector< double > > &param_values)
 Creates a grid of parameter combinations from parameter value ranges.
 
void findBestParameters (CvResult &cv_result, const std::unique_ptr< Score > &scorer)
 Identifies the best parameters from cross-validation results.
 
template<typename MatrixType >
CvResult crossValidate (Slope model, MatrixType &x, const Eigen::MatrixXd &y_in, const CvConfig &config=CvConfig())
 Performs cross-validation on a SLOPE model to select optimal hyperparameters.
 
template<typename MatrixType >
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.
 
template<typename MatrixType >
double estimateNoise (MatrixType &x, Eigen::MatrixXd &y, const bool fit_intercept)
 Estimates noise (standard error) in a linear model using OLS residuals.
 
template<typename MatrixType >
SlopePath estimateAlpha (MatrixType &x, Eigen::MatrixXd &y, const Slope &model)
 Estimates the regularization parameter alpha for SLOPE regression.
 
std::vector< int > kktCheck (const Eigen::VectorXd &gradient, const Eigen::VectorXd &beta, const Eigen::ArrayXd &lambda, const std::vector< int > &strong_set)
 Checks KKT conditions for SLOPE optimization.
 
std::string warningCodeToString (WarningCode code)
 Convert a warning code to its string representation.
 
std::unique_ptr< LosssetupLoss (const std::string &loss)
 Factory function to create the appropriate loss function based on the distribution family.
 
Eigen::VectorXd logSumExp (const Eigen::MatrixXd &a)
 
Eigen::MatrixXd softmax (const Eigen::MatrixXd &a)
 
std::vector< int > setUnion (const std::vector< int > &a, const std::vector< int > &b)
 Computes the union of two sorted integer vectors.
 
std::vector< int > setDiff (const std::vector< int > &a, const std::vector< int > &b)
 Computes the set difference of two sorted integer vectors.
 
Eigen::ArrayXd geomSpace (const double start, const double end, const int n)
 Creates an array of n numbers in geometric progression from start to end.
 
Eigen::VectorXd l2Norms (const Eigen::SparseMatrix< double > &x)
 Computes the L2 (Euclidean) norms for each column of a sparse matrix.
 
Eigen::VectorXd l2Norms (const Eigen::MatrixXd &x)
 Computes the L2 (Euclidean) norms for each column of a dense matrix.
 
Eigen::VectorXd means (const Eigen::SparseMatrix< double > &x)
 Computes the arithmetic mean for each column of a sparse matrix.
 
Eigen::VectorXd means (const Eigen::MatrixXd &x)
 Computes the arithmetic mean for each column of a dense matrix.
 
Eigen::VectorXd ranges (const Eigen::SparseMatrix< double > &x)
 Computes the range (max - min) for each column of a matrix.
 
Eigen::VectorXd ranges (const Eigen::MatrixXd &x)
 Computes the range (max - min) for each column of a dense matrix.
 
Eigen::VectorXd maxAbs (const Eigen::SparseMatrix< double > &x)
 Computes the maximum absolute value for each column of a matrix.
 
Eigen::VectorXd maxAbs (const Eigen::MatrixXd &x)
 Computes the maximum absolute value for each column of a dense matrix.
 
Eigen::VectorXd mins (const Eigen::SparseMatrix< double > &x)
 Computes the minimum value for each column of a sparse matrix.
 
Eigen::VectorXd mins (const Eigen::MatrixXd &x)
 Computes the minimum value for each column of a dense matrix.
 
Eigen::VectorXd stdDevs (const Eigen::SparseMatrix< double > &x)
 Computes the standard deviation for each column of a matrix.
 
Eigen::VectorXd stdDevs (const Eigen::MatrixXd &x)
 Computes the standard deviation for each column of a matrix.
 
template<typename T >
int sign (T val)
 Returns the sign of a given value.
 
template<typename T >
Eigen::ArrayXd cumSum (const T &x)
 
template<typename T >
sigmoid (const T &x)
 
template<typename T >
logit (const T &x)
 
template<typename T >
clamp (const T &x, const T &lo, const T &hi)
 
template<typename T >
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)
 
template<typename T >
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)
 
template<typename T >
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)
 
template<typename T >
int whichMax (const T &x)
 Returns the index of the maximum element in a container.
 
template<typename T >
int whichMin (const T &x)
 Returns the index of the minimum element in a container.
 
template<typename T , typename Comparator >
int whichBest (const T &x, const Comparator &comp)
 Returns the index of the minimum element in a container.
 
template<typename T >
Eigen::VectorXd l1Norms (const T &x)
 Computes the L1 (Manhattan) norms for each column of a matrix.
 
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)
 
JitNormalization normalize (Eigen::SparseMatrix< double > &x, Eigen::VectorXd &x_centers, Eigen::VectorXd &x_scales, const std::string &centering_type, const std::string &scaling_type, const bool)
 
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.
 
template<typename T >
void computeCenters (Eigen::VectorXd &x_centers, const T &x, const std::string &type)
 
template<typename T >
void computeScales (Eigen::VectorXd &x_scales, const T &x, const std::string &type)
 
double normalQuantile (const double p)
 Computes the quantile of a standard normal distribution using the Beasley-Springer-Moro algorithm.
 
Eigen::ArrayXd lambdaSequence (const int p, const double q, const std::string &type, const int n, const double theta1, const double theta2)
 
Eigen::ArrayXd regularizationPath (const Eigen::ArrayXd &alpha_in, const int path_length, double alpha_min_ratio, double alpha_max)
 
double binaryRocAuc (const Eigen::VectorXd &scores, const Eigen::VectorXd &labels)
 
double rocAuc (const Eigen::MatrixXd &scores, const Eigen::MatrixXd &labels)
 
std::vector< int > activeSet (const Eigen::VectorXd &beta)
 Identifies previously active variables.
 
std::vector< int > strongSet (const Eigen::VectorXd &gradient_prev, const Eigen::ArrayXd &lambda, const Eigen::ArrayXd &lambda_prev)
 Determines the strong set using sequential strong rules.
 
std::unique_ptr< ScreeningRulecreateScreeningRule (const std::string &screening_type)
 Creates a screening rule based on the provided type.
 
std::pair< double, double > computeClusterGradientAndHessian (const Eigen::MatrixXd &x, const int j, const std::vector< int > &s, const Clusters &clusters, const Eigen::VectorXd &w, const Eigen::VectorXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization jit_normalization)
 
std::pair< double, double > computeClusterGradientAndHessian (const Eigen::SparseMatrix< double > &x, const int j, const std::vector< int > &s, const Clusters &clusters, const Eigen::VectorXd &w, const Eigen::VectorXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization jit_normalization)
 
template<typename T >
std::pair< double, double > computeGradientAndHessian (const T &x, const int k, const Eigen::VectorXd &w, const Eigen::VectorXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const double s, const JitNormalization jit_normalization, const int n)
 
template<typename T >
void coordinateDescent (Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::VectorXd &residual, Clusters &clusters, const Eigen::ArrayXd &lambda, const T &x, const Eigen::VectorXd &w, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const bool intercept, const JitNormalization jit_normalization, const bool update_clusters)
 
std::unique_ptr< SolverBasesetupSolver (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.
 
std::tuple< double, int > slopeThreshold (const double x, const int j, const Eigen::ArrayXd lambdas, const Clusters &clusters)
 
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.
 
Eigen::MatrixXd subset (const Eigen::MatrixXd &x, const std::vector< int > &indices)
 Extract a subset of rows from an Eigen matrix.
 
Eigen::SparseMatrix< double > subset (const Eigen::SparseMatrix< double > &x, const std::vector< int > &indices)
 Extract a subset of rows from a sparse Eigen matrix.
 
Eigen::MatrixXd subsetCols (const Eigen::MatrixXd &x, const std::vector< int > &indices)
 Extract specified columns from a dense matrix.
 
Eigen::SparseMatrix< double > subsetCols (const Eigen::SparseMatrix< double > &x, const std::vector< int > &indices)
 Extract specified columns from a sparse matrix.
 
template<typename T >
void sort (T &v, const bool descending=false)
 
template<typename T >
std::vector< int > which (const T &x)
 
template<typename T >
std::vector< int > sortIndex (T &v, const bool descending=false)
 
template<typename T >
void permute (T &values, const std::vector< int > &ind)
 
template<typename T >
void inversePermute (T &values, const std::vector< int > &ind)
 
template<typename T >
void move_elements (std::vector< T > &v, const int from, const int to, const int size)
 
std::unordered_set< double > unique (const Eigen::MatrixXd &x)
 Create a set of unique values from an Eigen matrix.
 

Detailed Description

Namespace containing SLOPE regression implementation.

Typedef Documentation

◆ ArrayXb

typedef Eigen::Array< bool, Eigen::Dynamic, 1 > slope::ArrayXb

Dynamic-size column vector of boolean values Wrapper around Eigen::Array<bool, Eigen::Dynamic, 1>

Definition at line 18 of file kkt_check.h.

◆ ArrayXXb

typedef Eigen::Array< bool, Eigen::Dynamic, Eigen::Dynamic > slope::ArrayXXb

Dynamic-size matrix of boolean values Wrapper around Eigen::Array<bool, Eigen::Dynamic, Eigen::Dynamic>

Definition at line 25 of file kkt_check.h.

Enumeration Type Documentation

◆ JitNormalization

enum class slope::JitNormalization
strong

Enums to control predictor standardization behavior.

Enumerator
None 

No JIT normalization.

Center 

Center JIT.

Scale 

Scale JIT.

Both 

Both.

Definition at line 12 of file jit_normalization.h.

◆ WarningCode

enum class slope::WarningCode
strong

Standard warning codes used throughout the slope library.

These enum values represent the various types of warnings that can be generated during computational operations.

Enumerator
GENERIC_WARNING 

General uncategorized warning.

DEPRECATED_FEATURE 

Feature marked for deprecation.

MAXIT_REACHED 

Maximum iterations reached without convergence.

LINE_SEARCH_FAILED 

Line search algorithm failed to converge.

Definition at line 25 of file logger.h.

Function Documentation

◆ activeSet()

std::vector< int > slope::activeSet ( const Eigen::VectorXd &  beta)

Identifies previously active variables.

Parameters
betaCurrent coefficient matrix
Returns
std::vector<int> Indices of variables with non-zero coefficients

Definition at line 19 of file screening.cpp.

◆ binaryRocAuc()

double slope::binaryRocAuc ( const Eigen::VectorXd &  scores,
const Eigen::VectorXd &  labels 
)

Computes the Area Under the ROC Curve (AUC-ROC) for binary classification.

Parameters
scoresVector of prediction scores/probabilities for each sample
labelsVector of true binary labels (0 or 1)
Returns
The AUC-ROC score between 0 and 1, where 1 indicates perfect prediction and 0.5 indicates random prediction

Both input vectors must have the same length. The function calculates how well the scores discriminate between the positive and negative classes.

Definition at line 10 of file score.cpp.

◆ clamp()

template<typename T >
T slope::clamp ( const T &  x,
const T &  lo,
const T &  hi 
)

Returns the value of x clamped between the specified lower and upper bounds.

Template Parameters
Tthe type of the values being clamped
Parameters
xthe value to be clamped
lothe lower bound
hithe upper bound
Returns
the clamped value of x

Definition at line 101 of file math.h.

◆ computeCenters()

template<typename T >
void slope::computeCenters ( Eigen::VectorXd &  x_centers,
const T &  x,
const std::string &  type 
)

Compute centers.

There are two supported centering types:

  • "none": Do not compute centers and scales.
  • "mean": Use arithmetic means.
Template Parameters
TThe type of the input matrix.
Parameters
x_centersA vector where the computed or provided centers will be stored.
xThe input matrix.
typeA string specifying the normalization type ("none", "manual", or "standardization").
Exceptions
std::invalid_argumentif the provided manual centers or scales have invalid dimensions or contain non-finite values.

Definition at line 33 of file normalize.h.

◆ computeClusterGradientAndHessian() [1/2]

std::pair< double, double > slope::computeClusterGradientAndHessian ( const Eigen::MatrixXd &  x,
const int  j,
const std::vector< int > &  s,
const Clusters clusters,
const Eigen::VectorXd &  w,
const Eigen::VectorXd &  residual,
const Eigen::VectorXd &  x_centers,
const Eigen::VectorXd &  x_scales,
const JitNormalization  jit_normalization 
)

Computes the gradient and Hessian for a cluster of variables in coordinate descent.

This function handles the case when multiple variables are in the same cluster (have the same coefficient magnitude), calculating the combined gradient and Hessian needed for the coordinate descent update.

Parameters
xInput matrix
jCluster index
sVector of signs for each variable in the cluster
clustersThe cluster information object
wVector of weights
residualResidual vector
x_centersVector of feature centers (means)
x_scalesVector of feature scales (standard deviations)
jit_normalizationNormalization strategy (Both, Center, Scale, or None)
Returns
std::pair<double, double> containing:
  • first: Hessian of the loss function for the cluster
  • second: gradient of the loss function for the cluster

Definition at line 6 of file hybrid_cd.cpp.

◆ computeClusterGradientAndHessian() [2/2]

std::pair< double, double > slope::computeClusterGradientAndHessian ( const Eigen::SparseMatrix< double > &  x,
const int  j,
const std::vector< int > &  s,
const Clusters clusters,
const Eigen::VectorXd &  w,
const Eigen::VectorXd &  residual,
const Eigen::VectorXd &  x_centers,
const Eigen::VectorXd &  x_scales,
const JitNormalization  jit_normalization 
)

Computes the gradient and Hessian for a cluster of variables in coordinate descent (sparse matrix version).

This overloaded version handles sparse input matrices, optimizing the computation for this data structure.

Parameters
xInput sparse matrix
jCluster index
sVector of signs for each variable in the cluster
clustersThe cluster information object
wVector of weights
residualResidual vector
x_centersVector of feature centers (means)
x_scalesVector of feature scales (standard deviations)
jit_normalizationNormalization strategy (Both, Center, Scale, or None)
Returns
std::pair<double, double> containing:
  • first: Hessian of the loss function for the cluster
  • second: gradient of the loss function for the cluster

Definition at line 55 of file hybrid_cd.cpp.

◆ computeDual()

template<typename MatrixType >
double slope::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.

Template Parameters
MatrixTypeThe type of the design matrix
Parameters
betaCurrent coefficient vector
residualResidual matrix (y - prediction)
lossPointer to the loss function object
sl1_normSorted L1 norm object
lambdaVector of penalty parameters
xDesign matrix
yResponse matrix
x_centersVector of feature means for centering
x_scalesVector of feature scales for normalization
jit_normalizationJust-in-time normalization settings
interceptBoolean indicating if intercept is included in the model
Returns
double The computed dual objective value

This function computes the dual objective value for the SLOPE optimization problem. It handles both cases with and without intercept terms, applying appropriate normalization and gradient computations.

Definition at line 41 of file diagnostics.h.

◆ computeGradientAndHessian()

template<typename T >
std::pair< double, double > slope::computeGradientAndHessian ( const T &  x,
const int  k,
const Eigen::VectorXd &  w,
const Eigen::VectorXd &  residual,
const Eigen::VectorXd &  x_centers,
const Eigen::VectorXd &  x_scales,
const double  s,
const JitNormalization  jit_normalization,
const int  n 
)

Computes the gradient and Hessian for coordinate descent optimization with different normalization strategies.

Template Parameters
TMatrix type (expected to support col() operations like Eigen matrices)
Parameters
xInput matrix
kColumn index to compute derivatives for
wVector of weights
residualResidual vector
x_centersVector of feature centers (means)
x_scalesVector of feature scales (standard deviations)
sStep size parameter
jit_normalizationNormalization strategy (Both, Center, Scale, or None)
nNumber of samples
Returns
std::pair<double, double> containing:
  • first: gradient of the loss function
  • second: diagonal Hessian element

Definition at line 41 of file hybrid_cd.h.

◆ computeScales()

template<typename T >
void slope::computeScales ( Eigen::VectorXd &  x_scales,
const T &  x,
const std::string &  type 
)

Compute scales

There are two supported scaling types:

  • "none": Do not compute centers and scales.
  • "sd": Compute centers and scales using Welford’s algorithm.
Template Parameters
TThe type of the input matrix.
Parameters
xThe input matrix.
x_scalesA vector where the computed or provided scales will be stored.
typeA string specifying the normalization type ("none", "manual", or "standardization").
Exceptions
std::invalid_argumentif the provided manual centers or scales have invalid dimensions or contain non-finite values.

Definition at line 74 of file normalize.h.

◆ coordinateDescent()

template<typename T >
void slope::coordinateDescent ( Eigen::VectorXd &  beta0,
Eigen::VectorXd &  beta,
Eigen::VectorXd &  residual,
Clusters clusters,
const Eigen::ArrayXd &  lambda,
const T &  x,
const Eigen::VectorXd &  w,
const Eigen::VectorXd &  x_centers,
const Eigen::VectorXd &  x_scales,
const bool  intercept,
const JitNormalization  jit_normalization,
const bool  update_clusters 
)

Coordinate Descent Step

This function takes a coordinate descent step in the hybrid CD/PGD algorithm for SLOPE.

Template Parameters
TThe type of the design matrix. This can be either a dense or sparse.
Parameters
beta0The intercept
betaThe coefficients
residualThe residual vector
clustersThe cluster information, stored in a Cluster object.
lambdaRegularization weights
xThe design matrix
wThe weight vector
x_centersThe center values of the data matrix columns
x_scalesThe scale values of the data matrix columns
interceptShuold an intervept be fit?
jit_normalizationType o fJIT normalization.
update_clustersFlag indicating whether to update the cluster information
See also
Clusters
SortedL1Norm
JitNormalization

Definition at line 187 of file hybrid_cd.h.

◆ createGrid()

std::vector< std::map< std::string, double > > slope::createGrid ( const std::map< std::string, std::vector< double > > &  param_values)

Creates a grid of parameter combinations from parameter value ranges.

Parameters
param_valuesMap of parameter names to vectors of possible values
Returns
std::vector<std::map<std::string, double>> Vector of parameter combinations

This function takes a map where keys are parameter names and values are vectors of possible values for each parameter, then generates all possible combinations.

Definition at line 6 of file cv.cpp.

◆ createScreeningRule()

std::unique_ptr< ScreeningRule > slope::createScreeningRule ( const std::string &  screening_type)

Creates a screening rule based on the provided type.

Parameters
screening_typeType of screening rule to create ("none" or "strong")
Returns
std::unique_ptr<ScreeningRule> A pointer to the created screening rule

Definition at line 257 of file screening.cpp.

◆ crossValidate()

template<typename MatrixType >
CvResult slope::crossValidate ( Slope  model,
MatrixType &  x,
const Eigen::MatrixXd &  y_in,
const CvConfig config = CvConfig() 
)

Performs cross-validation on a SLOPE model to select optimal hyperparameters.

Template Parameters
MatrixTypeType of design matrix (supports both dense and sparse matrices)
Parameters
modelThe SLOPE model to be cross-validated
xThe design matrix containing predictors
y_inThe response matrix
configConfiguration parameters for cross-validation (optional)
Returns
CvResult Object containing cross-validation results, including best parameters, regularization paths, and performance metrics across folds

This function implements k-fold cross-validation for SLOPE models, evaluating model performance across a grid of hyperparameters. For each hyperparameter combination, the function:

  1. Splits the data into training and validation sets according to the fold configuration
  2. Fits the model on each training set and evaluates on the validation set
  3. Computes the specified evaluation metric for each regularization parameter
  4. Averages results across folds to select optimal hyperparameters

The function supports parallel processing with OpenMP when available.

Definition at line 160 of file cv.h.

◆ cumSum()

template<typename T >
Eigen::ArrayXd slope::cumSum ( const T &  x)

Calculates the cumulative sum of the elements in the input array.

Template Parameters
TThe type of the input array.
Parameters
xThe input array.
Returns
An Eigen::ArrayXd containing the cumulative sum of the elements in the input array.

Definition at line 46 of file math.h.

◆ estimateAlpha()

template<typename MatrixType >
SlopePath slope::estimateAlpha ( MatrixType &  x,
Eigen::MatrixXd &  y,
const Slope model 
)

Estimates the regularization parameter alpha for SLOPE regression.

This function implements an algorithm to estimate an appropriate regularization parameter (alpha) for SLOPE, which is a generalization of the lasso. When n >= p + 30, it directly estimates alpha from OLS residuals. Otherwise, it uses an iterative procedure that alternates between estimating alpha and fitting the SLOPE model.

The iterative procedure works by:

  1. Starting with an empty set of selected variables
  2. Estimating alpha based on the selected variables
  3. Fitting a SLOPE model with that alpha
  4. Updating the selected variables based on non-zero coefficients
  5. Repeating until convergence or maximum iterations reached
Template Parameters
MatrixTypeThe type of matrix used to store the design matrix
Parameters
xDesign matrix with n observations and p predictors
yResponse matrix
modelThe SLOPE model object containing method parameters
Returns
A SlopePath object containing the fitted model with estimated alpha
Exceptions
std::runtime_errorIf maximum iterations reached or if too many variables selected

Definition at line 92 of file estimate_alpha.h.

◆ estimateNoise()

template<typename MatrixType >
double slope::estimateNoise ( MatrixType &  x,
Eigen::MatrixXd &  y,
const bool  fit_intercept 
)

Estimates noise (standard error) in a linear model using OLS residuals.

Calculates the standard error of the residuals from an Ordinary Least Squares (OLS) regression. This is computed as sqrt(sum(residuals²) / (n - p + intercept_term)), where n is the number of observations and p is the number of predictors.

Template Parameters
MatrixTypeThe type of matrix used to store the design matrix
Parameters
xDesign matrix with n observations and p predictors
yResponse matrix
fit_interceptWhether to include an intercept term in the model
Returns
The estimated noise level (standard error of residuals)

Definition at line 32 of file estimate_alpha.h.

◆ findBestParameters()

void slope::findBestParameters ( CvResult cv_result,
const std::unique_ptr< Score > &  scorer 
)

Identifies the best parameters from cross-validation results.

Parameters
cv_resultCross-validation results to analyze and update
scorerThe scoring metric used to evaluate performance

This function examines all cross-validation results across the parameter grid and updates the cv_result with information about the best parameter combination, including the best score and corresponding indices.

Definition at line 39 of file cv.cpp.

◆ geomSpace()

Eigen::ArrayXd slope::geomSpace ( const double  start,
const double  end,
const int  n 
)

Creates an array of n numbers in geometric progression from start to end.

Parameters
startThe starting value of the sequence
endThe final value of the sequence
nThe number of points to generate
Returns
Eigen::ArrayXd Array containing n points spaced geometrically between start and end
Note
Similar to numpy.geomspace, generates points that are evenly spaced on a log scale
Exceptions
std::invalid_argumentIf start or end is zero, or if n < 1

Definition at line 46 of file math.cpp.

◆ inversePermute()

template<typename T >
void slope::inversePermute ( T &  values,
const std::vector< int > &  ind 
)

Inverse permutes the elements of a container based on the given indices.

This function takes a container of values and a vector of indices and rearranges the elements of the container according to the indices. The resulting container will have the elements in the order specified by the indices.

Template Parameters
TThe type of the container.
Parameters
valuesThe container of values to be permuted.
indThe vector of indices specifying the new order of the elements.

< The resulting container after permutation.

Definition at line 152 of file utils.h.

◆ kktCheck()

std::vector< int > slope::kktCheck ( const Eigen::VectorXd &  gradient,
const Eigen::VectorXd &  beta,
const Eigen::ArrayXd &  lambda,
const std::vector< int > &  strong_set 
)

Checks KKT conditions for SLOPE optimization.

Parameters
gradientThe gradient of the loss function
betaThe current coefficients
lambdaVector of regularization parameters
strong_setVector of indices in the strong set
Returns
std::vector<int> Indices where KKT conditions are violated

Verifies if the current solution satisfies the KKT optimality conditions for the SLOPE optimization problem. Returns indices where violations occur.

Definition at line 9 of file kkt_check.cpp.

◆ l1Norms()

template<typename T >
Eigen::VectorXd slope::l1Norms ( const T &  x)

Computes the L1 (Manhattan) norms for each column of a matrix.

Template Parameters
TMatrix type (must support cols(), col() operations, compatible with Eigen)
Parameters
xInput matrix whose column L1 norms are to be computed
Returns
Eigen::VectorXd Vector containing L1 norms of each column

For each column j in matrix x, computes the L1 norm:

\[ \|x_j\|_1 = \sum_{i} |x_{ij}| \]

where x_{ij} represents the i-th element of the j-th column.

Definition at line 439 of file math.h.

◆ l2Norms() [1/2]

Eigen::VectorXd slope::l2Norms ( const Eigen::MatrixXd &  x)

Computes the L2 (Euclidean) norms for each column of a dense matrix.

Parameters
xInput dense matrix whose column L2 norms are to be computed
Returns
Eigen::VectorXd Vector containing L2 norms of each column

For each column j in matrix x, computes the L2 norm:

\[ \|x_j\|_2 = \sqrt{\sum_{i} x_{ij}^2} \]

where x_{ij} represents the i-th element of the j-th column.

Definition at line 70 of file math.cpp.

◆ l2Norms() [2/2]

Eigen::VectorXd slope::l2Norms ( const Eigen::SparseMatrix< double > &  x)

Computes the L2 (Euclidean) norms for each column of a sparse matrix.

Parameters
xInput sparse matrix whose column L2 norms are to be computed
Returns
Eigen::VectorXd Vector containing L2 norms of each column

For each column j in matrix x, computes the L2 norm:

\[ \|x_j\|_2 = \sqrt{\sum_{i} x_{ij}^2} \]

where x_{ij} represents the i-th element of the j-th column.

Definition at line 56 of file math.cpp.

◆ lambdaSequence()

Eigen::ArrayXd slope::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 
)

Generates a sequence of regularization weights for the sorted L1 norm.

Parameters
pThe number of lambda values to generate (number of features)
qThe false discovery rate (FDR) level or quantile value (in (0, 1))
typeThe type of sequence to generate:
  • "bh": Benjamini-Hochberg sequence
  • "gaussian": Gaussian sequence
  • "oscar": Octagonal Shrinkage and Clustering Algorithm for Regression
nNumber of observations (only used for gaussian type)
theta1First parameter for OSCAR weights (default: 1.0)
theta2Second parameter for OSCAR weights (default: 1.0)
Returns
Eigen::ArrayXd containing the generated lambda sequence in decreasing order

Definition at line 10 of file regularization_sequence.cpp.

◆ linearPredictor()

template<typename T >
Eigen::MatrixXd slope::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 
)

Computes the gradient of the loss with respect to \(\beta\).

Template Parameters
TThe type of the input matrix.
Parameters
xThe input matrix.
active_setIndicies for active set
beta0Intercept
betaCoefficients
x_centersThe vector of center values for each column of x.
x_scalesThe vector of scale values for each column of x.
jit_normalizationType of JIT normalization.
interceptWhether to fit an intercept.
Returns
The computed gradient vector.

Definition at line 142 of file math.h.

◆ logit()

template<typename T >
T slope::logit ( const T &  x)

The logit

The logit function is defined as \(\log(\frac{x}{1 - x})\).

Template Parameters
TThe type of the input.
Parameters
xThe input value.
Returns
The result of the logit function.

Definition at line 83 of file math.h.

◆ logSumExp()

Eigen::VectorXd slope::logSumExp ( const Eigen::MatrixXd &  a)

LogSumExp

Parameters
aA matrix
Returns
\(\log(\sum_i \exp(a_i))\)

Definition at line 7 of file math.cpp.

◆ maxAbs() [1/2]

Eigen::VectorXd slope::maxAbs ( const Eigen::MatrixXd &  x)

Computes the maximum absolute value for each column of a dense matrix.

Parameters
xInput dense matrix whose column-wise maximum absolute values are to be computed
Returns
Eigen::VectorXd Vector containing the maximum absolute value of each column

For each column j in matrix x, computes:

\[ \max_{i} |x_{ij}| \]

where x_{ij} represents the i-th element of the j-th column.

Definition at line 147 of file math.cpp.

◆ maxAbs() [2/2]

Eigen::VectorXd slope::maxAbs ( const Eigen::SparseMatrix< double > &  x)

Computes the maximum absolute value for each column of a matrix.

Parameters
xInput matrix to compute column-wise maximum absolute values
Returns
Eigen::VectorXd Vector containing maximum absolute values for each column

For each column j in matrix x, computes:

\[ \max_{i} |x_{ij}| \]

where x_{ij} represents the i-th element of the j-th column.

Definition at line 126 of file math.cpp.

◆ means() [1/2]

Eigen::VectorXd slope::means ( const Eigen::MatrixXd &  x)

Computes the arithmetic mean for each column of a dense matrix.

Parameters
xInput dense matrix whose column means are to be computed
Returns
Eigen::VectorXd Vector containing the arithmetic mean of each column

For each column j in matrix x, computes the mean:

\[ \bar{x}_j = \frac{1}{n}\sum_{i} x_{ij} \]

where n is the number of rows in the matrix and x_{ij} represents the i-th element of the j-th column.

Definition at line 91 of file math.cpp.

◆ means() [2/2]

Eigen::VectorXd slope::means ( const Eigen::SparseMatrix< double > &  x)

Computes the arithmetic mean for each column of a sparse matrix.

Parameters
xInput sparse matrix whose column means are to be computed
Returns
Eigen::VectorXd Vector containing the arithmetic mean of each column

For each column j in matrix x, computes the mean:

\[ \bar{x}_j = \frac{1}{n}\sum_{i} x_{ij} \]

where n is the number of rows in the matrix and x_{ij} represents the i-th element of the j-th column.

Definition at line 76 of file math.cpp.

◆ mins() [1/2]

Eigen::VectorXd slope::mins ( const Eigen::MatrixXd &  x)

Computes the minimum value for each column of a dense matrix.

Parameters
xInput dense matrix whose column minimums are to be computed
Returns
Eigen::VectorXd Vector containing the minimum value of each column

For each column j in matrix x, computes:

\[ \min_{i}(x_{ij}) \]

where x_{ij} represents the i-th element of the j-th column.

Uses Eigen's built-in column-wise operations for efficient computation on dense matrices.

Definition at line 174 of file math.cpp.

◆ mins() [2/2]

Eigen::VectorXd slope::mins ( const Eigen::SparseMatrix< double > &  x)

Computes the minimum value for each column of a sparse matrix.

Parameters
xInput sparse matrix whose column minimums are to be computed
Returns
Eigen::VectorXd Vector containing the minimum value of each column

For each column j in matrix x, computes:

\[ \min_{i}(x_{ij}) \]

where x_{ij} represents the i-th element of the j-th column.

Definition at line 153 of file math.cpp.

◆ move_elements()

template<typename T >
void slope::move_elements ( std::vector< T > &  v,
const int  from,
const int  to,
const int  size 
)

Moves a range of elements within a vector.

This function moves a range of elements within a vector from one position to another. The elements are moved in a way that the order is preserved.

Template Parameters
TThe type of elements in the vector.
Parameters
vThe vector containing the elements.
fromThe starting index of the range to be moved.
toThe ending index of the range to be moved.
sizeThe size of the range to be moved.

Definition at line 176 of file utils.h.

◆ normalize() [1/2]

JitNormalization slope::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 
)

Normalize a dense matrix by centering and scaling.

The function computes column centers and scaling factors based on the specified normalization type ("none", "manual", or "standardization"). If modify_x is true, the normalization is applied directly to the input matrix.

Parameters
xThe dense input matrix to be normalized.
x_centersA vector that will hold the column centers. It will be resized to match the number of columns.
x_scalesA vector that will hold the column scaling factors. It will be resized to match the number of columns.
centering_typeA string specifying the normalization type ("none", "manual", or "standardization").
scaling_typeA string specifying the normalization type ("none", "manual", or "standardization").
modify_xIf true, modifies x in-place; otherwise, x remains unchanged (centers/scales are still computed).
Returns
true if normalization succeeds, false otherwise.

Definition at line 6 of file normalize.cpp.

◆ normalize() [2/2]

JitNormalization slope::normalize ( Eigen::SparseMatrix< double > &  x,
Eigen::VectorXd &  x_centers,
Eigen::VectorXd &  x_scales,
const std::string &  centering_type,
const std::string &  scaling_type,
const bool  modify_x 
)

Normalize a sparse matrix by scaling only.

To preserve sparsity, centering is not applied. The scaling factors for each column are computed according to the specified normalization type ("none", "manual", or "standardization"). If modify_x is true, the scaling is applied directly to the input matrix.

Parameters
xThe sparse input matrix to be normalized.
x_centersA vector that will hold the column centers. For sparse matrices, centering is typically skipped; this parameter is maintained for consistency.
x_scalesA vector that will hold the column scaling factors. It will be resized to match the number of columns.
centering_typeA string specifying the normalization type ("none", "manual", or "standardization").
scaling_typeA string specifying the normalization type ("none", "manual", or "standardization").
modify_xIf true, performs in-place scaling on x; otherwise, leaves x unchanged.
Returns
true if normalization succeeds, false otherwise.

Definition at line 54 of file normalize.cpp.

◆ normalQuantile()

double slope::normalQuantile ( const double  p)

Computes the quantile of a standard normal distribution using the Beasley-Springer-Moro algorithm.

Parameters
pProbability value
Returns
Quantile

Definition at line 9 of file qnorm.cpp.

◆ offsetGradient()

template<typename T >
void slope::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 
)

Computes the gradient of the loss with respect to \(\beta\).

Template Parameters
TThe type of the input matrix.
Parameters
gradientThe residual vector.
xThe input matrix.
offsetGradient offset
active_setIndices for the active_set
x_centersThe vector of center values for each column of x.
x_scalesThe vector of scale values for each column of x.
jit_normalizationType of JIT normalization just-in-time.

Definition at line 295 of file math.h.

◆ patternMatrix()

Eigen::SparseMatrix< int > slope::Clusters::patternMatrix ( const Eigen::VectorXd &  beta)

Definition at line 257 of file clusters.cpp.

◆ permute()

template<typename T >
void slope::permute ( T &  values,
const std::vector< int > &  ind 
)

Permutes the elements of a container according to the given indices.

This function takes a container of values and permutes its elements according to the given indices. The indices specify the new order of the elements in the container.

Template Parameters
TThe type of the container.
Parameters
valuesThe container of values to be permuted.
indThe vector of indices specifying the new order of the elements.

The container to store the permuted values.

Permute the values according to the given indices.

Assign the permuted values back to the original container.

Definition at line 118 of file utils.h.

◆ ranges() [1/2]

Eigen::VectorXd slope::ranges ( const Eigen::MatrixXd &  x)

Computes the range (max - min) for each column of a dense matrix.

Parameters
xInput dense matrix whose column ranges are to be computed
Returns
Eigen::VectorXd Vector containing the range of each column

For each column j in matrix x, computes:

\[ range_j = \max_{i}(x_{ij}) - \min_{i}(x_{ij}) \]

where x_{ij} represents the i-th element of the j-th column.

Uses Eigen's efficient colwise() operations to compute maximum and minimum values for each column simultaneously.

Definition at line 120 of file math.cpp.

◆ ranges() [2/2]

Eigen::VectorXd slope::ranges ( const Eigen::SparseMatrix< double > &  x)

Computes the range (max - min) for each column of a matrix.

Parameters
xInput matrix whose column ranges are to be computed
Returns
Eigen::VectorXd Vector containing the range of each column

For each column j in matrix x, computes:

\[ range_j = \max_{i}(x_{ij}) - \min_{i}(x_{ij}) \]

where x_{ij} represents the i-th element of the j-th column.

Definition at line 97 of file math.cpp.

◆ regularizationPath()

Eigen::ArrayXd slope::regularizationPath ( const Eigen::ArrayXd &  alpha_in,
const int  path_length,
double  alpha_min_ratio,
const double  alpha_max 
)

Computes a sequence of regularization weights for the SLOPE path.

Parameters
alpha_inAlpha sequence, of length zero if automatic.
path_lengthLength of path.
alpha_min_ratioRatio of minimum to maximum alpha
alpha_maxValue of alpha as which the model is completely sparse
Returns
Eigen::ArrayXd containing the sequence of regularization parameters from strongest (alpha_max) to weakest (alpha_max * alpha_min_ratio)

Definition at line 72 of file regularization_sequence.cpp.

◆ rescaleCoefficients()

std::tuple< Eigen::VectorXd, Eigen::MatrixXd > slope::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.

This function rescales the coefficients by dividing each coefficient by the corresponding scale factor and subtracting the product of the center and the coefficient from the intercept.

Parameters
beta0The intercept coefficient.
betaThe vector of coefficients.
x_centersThe vector of center values.
x_scalesThe vector of scale factors.
interceptShould an intercept be fit?
Returns
A tuple containing the rescaled intercept and coefficients.
Note
The input vectors beta, x_centers, and x_scales must have the same size.
The output vector beta will be modified in-place.
See also
SlopeParameters

Definition at line 94 of file normalize.cpp.

◆ rocAuc()

double slope::rocAuc ( const Eigen::MatrixXd &  scores,
const Eigen::MatrixXd &  labels 
)

Computes the Area Under the ROC Curve (AUC-ROC) for binary and multi-class classification.

Parameters
scoresMatrix of prediction scores/probabilities (samples x classes)
labelsMatrix of true labels in one-hot encoded format (samples x classes)
Returns
The average AUC-ROC score across all classes, between 0 and 1

For binary classification, this reduces to standard binary AUC-ROC. For multi-class, computes one-vs-rest AUC-ROC for each class and averages. Both input matrices must have the same dimensions.

Definition at line 54 of file score.cpp.

◆ setDiff()

std::vector< int > slope::setDiff ( const std::vector< int > &  a,
const std::vector< int > &  b 
)

Computes the set difference of two sorted integer vectors.

Parameters
aFirst sorted vector of integers (set to subtract from)
bSecond sorted vector of integers (set to subtract)
Returns
std::vector<int> Vector containing elements in a that do not appear in b, maintaining sorted order

Returns A \ B = {x ∈ A | x ∉ B}

Definition at line 36 of file math.cpp.

◆ setUnion()

std::vector< int > slope::setUnion ( const std::vector< int > &  a,
const std::vector< int > &  b 
)

Computes the union of two sorted integer vectors.

Parameters
aFirst sorted vector of integers
bSecond sorted vector of integers
Returns
std::vector<int> Vector containing all elements that appear in either a or b, without duplicates and in sorted order

Definition at line 26 of file math.cpp.

◆ setupLoss()

std::unique_ptr< Loss > slope::setupLoss ( const std::string &  loss)

Factory function to create the appropriate loss function based on the distribution family.

This function creates and returns an loss function object based on the specified statistical distribution family. The supported families are:

  • "logistic": For binary classification problems (logistic regression)
  • "poisson": For count data modeling (Poisson regression)
  • "multinomial": For multi-class classification problems
  • "quadratic": For continuous response variables (linear regression, default if unspecified)
Parameters
lossA string specifying the loss type ("logistic", "poisson", "multinomial", or "quadratic")
Returns
std::unique_ptr<Loss> A unique pointer to the appropriate loss function object

Definition at line 10 of file setup_loss.cpp.

◆ setupSolver()

std::unique_ptr< SolverBase > slope::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.

Creates a solver object based on the specified type and parameters. The solver implements the Sorted L1 Penalized Estimation (SLOPE) algorithm with various configurations possible.

Parameters
solver_typeType of solver to use (e.g., "pgd", "admm")
lossLoss type
jit_normalizationType of JIT normalization
interceptWhether to fit an intercept term
update_clustersWhether to update cluster assignments during optimization (Hybrid solver)
cd_iterationsFrequency of proximal gradient descent updates (Hybrid solver)
Returns
std::unique_ptr<SolverBase> A unique pointer to the configured solver
See also
Loss
JitNormalization
SolverBase
PGD
Hybrid

Definition at line 10 of file setup_solver.cpp.

◆ sigmoid()

template<typename T >
T slope::sigmoid ( const T &  x)

Calculates the sigmoid function for the given input.

The sigmoid function is defined as 1 / (1 + exp(-x)).

Template Parameters
TThe type of the input.
Parameters
xThe input value.
Returns
The result of the sigmoid function.

Definition at line 67 of file math.h.

◆ sign()

template<typename T >
int slope::sign ( val)

Returns the sign of a given value.

This function determines the sign of the input value by comparing it to zero. It returns -1 if the value is negative, 0 if the value is zero, and 1 if the value is positive.

Template Parameters
TThe type of the input value.
Parameters
valThe value for which the sign needs to be determined.
Returns
-1 if the value is negative, 0 if the value is zero, and 1 if the value is positive.

Definition at line 31 of file math.h.

◆ slopeThreshold()

std::tuple< double, int > slope::slopeThreshold ( const double  x,
const int  j,
const Eigen::ArrayXd  lambdas,
const Clusters clusters 
)

Calculates slope thresholding for a given input

This function calculates the slope threshold for a given x and j, using the provided lambdas and clusters. This is used in the coordinate descent update step.

Parameters
xThe value of x.
jThe value of j.
lambdasThe array of lambdas.
clustersThe clusters object.
Returns
A tuple containing the slope threshold and the index.

Definition at line 8 of file slope_threshold.cpp.

◆ softmax()

Eigen::MatrixXd slope::softmax ( const Eigen::MatrixXd &  x)

Softmax

Computes the softmax function for the given input matrix.

Parameters
xA matrix
Returns
\(\exp(a) / \sum_i \exp(a_i)\)

Definition at line 17 of file math.cpp.

◆ sort()

template<typename T >
void slope::sort ( T &  v,
const bool  descending = false 
)

Sorts the elements in a container in ascending or descending order.

This function sorts the elements in the container v in either ascending or descending order, depending on the value of the descending parameter.

Template Parameters
TThe type of the container.
Parameters
vThe container to be sorted.
descendingFlag indicating whether to sort in descending order (default is ascending order).
Note
The elements in the container must be comparable using the < and > operators.
See also
std::sort

Definition at line 37 of file utils.h.

◆ sortIndex()

template<typename T >
std::vector< int > slope::sortIndex ( T &  v,
const bool  descending = false 
)

Sorts the elements of a vector and returns the indices of the sorted elements.

This function sorts the elements of a vector in ascending or descending order and returns the indices of the sorted elements. The sorting is done using the std::sort function from the C++ standard library.

Template Parameters
TThe type of the vector elements.
Parameters
vThe vector to be sorted.
descendingFlag indicating whether to sort in descending order. Default is false.
Returns
A vector of indices representing the sorted order of the elements in the input vector.

Definition at line 89 of file utils.h.

◆ stdDevs() [1/2]

Eigen::VectorXd slope::stdDevs ( const Eigen::MatrixXd &  x)

Computes the standard deviation for each column of a matrix.

Parameters
xInput matrix whose column standard deviations are to be computed
Returns
Eigen::VectorXd Vector containing the standard deviation of each column

For each column j in matrix x, computes the standard deviation:

\[ \sigma_j = \sqrt{\frac{1}{n}\sum_{i} (x_{ij} - \bar{x}_j)^2} \]

where n is the number of rows, x_{ij} represents the i-th element of the j-th column, and \(\bar{x}_j\) is the mean of column j.

Definition at line 212 of file math.cpp.

◆ stdDevs() [2/2]

Eigen::VectorXd slope::stdDevs ( const Eigen::SparseMatrix< double > &  x)

Computes the standard deviation for each column of a matrix.

Parameters
xInput matrix whose column standard deviations are to be computed
Returns
Eigen::VectorXd Vector containing the standard deviation of each column

For each column j in matrix x, computes the standard deviation:

\[ \sigma_j = \sqrt{\frac{1}{n}\sum_{i} (x_{ij} - \bar{x}_j)^2} \]

where n is the number of rows, x_{ij} represents the i-th element of the j-th column, and \(\bar{x}_j\) is the mean of column j.

Definition at line 180 of file math.cpp.

◆ strongSet()

std::vector< int > slope::strongSet ( const Eigen::VectorXd &  gradient_prev,
const Eigen::ArrayXd &  lambda,
const Eigen::ArrayXd &  lambda_prev 
)

Determines the strong set using sequential strong rules.

Parameters
gradient_prevGradient from previous solution
lambdaCurrent lambda sequence
lambda_prevPrevious lambda sequence
Returns
std::vector<int> Indices of variables in the strong set

Definition at line 27 of file screening.cpp.

◆ subset() [1/2]

Eigen::MatrixXd slope::subset ( const Eigen::MatrixXd &  x,
const std::vector< int > &  indices 
)

Extract a subset of rows from an Eigen matrix.

Parameters
xThe input matrix to extract rows from
indicesA vector of row indices to extract
Returns
Eigen::MatrixXd A new matrix containing only the specified rows

This function creates a new matrix containing only the rows specified in the indices vector, preserving their order. The number of columns remains the same.

Definition at line 26 of file utils.cpp.

◆ subset() [2/2]

Eigen::SparseMatrix< double > slope::subset ( const Eigen::SparseMatrix< double > &  x,
const std::vector< int > &  indices 
)

Extract a subset of rows from a sparse Eigen matrix.

Parameters
xThe input sparse matrix to extract rows from
indicesA vector of row indices to extract
Returns
Eigen::SparseMatrix<double> A new sparse matrix containing only the specified rows

This function creates a new sparse matrix containing only the rows specified in the indices vector, preserving their order. The number of columns remains the same. The sparsity structure is maintained in the extracted rows.

Definition at line 32 of file utils.cpp.

◆ subsetCols() [1/2]

Eigen::MatrixXd slope::subsetCols ( const Eigen::MatrixXd &  x,
const std::vector< int > &  indices 
)

Extract specified columns from a dense matrix.

Creates a new matrix containing only the columns specified in the indices vector.

Parameters
xInput matrix
indicesVector of column indices to extract (0-based)
Returns
Eigen::MatrixXd New matrix with only the selected columns

Definition at line 54 of file utils.cpp.

◆ subsetCols() [2/2]

Eigen::SparseMatrix< double > slope::subsetCols ( const Eigen::SparseMatrix< double > &  x,
const std::vector< int > &  indices 
)

Extract specified columns from a sparse matrix.

Creates a new sparse matrix containing only the columns specified in the indices vector. The column ordering in the output matches the order of indices.

Parameters
xInput sparse matrix
indicesVector of column indices to extract (0-based)
Returns
Eigen::SparseMatrix<double> New sparse matrix with only the selected columns

Definition at line 60 of file utils.cpp.

◆ unique()

std::unordered_set< double > slope::unique ( const Eigen::MatrixXd &  x)
inline

Create a set of unique values from an Eigen matrix.

Parameters
xThe input matrix to extract unique values from
Returns
std::unordered_set<double> A set containing all unique values found in the matrix

This function iterates through all elements of the input matrix and collects all unique values into an unordered set. The order of elements in the returned set is not guaranteed.

Definition at line 272 of file utils.h.

◆ updateGradient()

template<typename T >
void slope::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 
)

Computes the gradient of the loss with respect to \(\beta\).

Template Parameters
TThe type of the input matrix.
Parameters
gradientThe gradient vector.
xThe input matrix.
residualThe residual matrix.
active_setThe indices for the active set.
x_centersThe vector of center values for each column of x.
x_scalesThe vector of scale values for each column of x.
wWorking weights
jit_normalizationType of JIT normalization just-in-time.

Definition at line 223 of file math.h.

◆ validateOption()

void slope::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.

Throws an informative error if the value is not found in the valid options, listing all valid possibilities in the error message.

Parameters
valueThe value to validate
valid_optionsSet of valid options
parameter_nameName of the parameter being validated (used in error message).
Exceptions
std::invalid_argumentIf value is not in valid_options

Definition at line 7 of file utils.cpp.

◆ warningCodeToString()

std::string slope::warningCodeToString ( WarningCode  code)

Convert a warning code to its string representation.

Parameters
codeThe warning code to convert
Returns
std::string String representation of the warning code

Definition at line 14 of file logger.cpp.

◆ which()

template<typename T >
std::vector< int > slope::which ( const T &  x)

Returns indices of true values in a boolean container.

Template Parameters
TContainer type supporting size() and operator[] (e.g., std::vector<bool>, std::array<bool>)
Parameters
xInput container with boolean-convertible values
Returns
std::vector<int> containing indices where x[i] evaluates to true

Example: std::vector<bool> v = {true, false, true, false, true}; auto indices = which(v); // returns {0, 2, 4}

Definition at line 60 of file utils.h.

◆ whichBest()

template<typename T , typename Comparator >
int slope::whichBest ( const T &  x,
const Comparator &  comp 
)

Returns the index of the minimum element in a container.

Template Parameters
TContainer type that supports iterators and std::min_element
Parameters
xContainer whose maximum element's index is to be found
compComparator that returns true if the first argument is worse than the second argument
Returns
int Zero-based index position of the maximum element

Uses std::max_element to find the iterator to the maximum element, then converts to index position using std::distance. For containers with multiple maximum elements, returns the first occurrence.

Definition at line 403 of file math.h.

◆ whichMax()

template<typename T >
int slope::whichMax ( const T &  x)

Returns the index of the maximum element in a container.

Template Parameters
TContainer type that supports iterators and std::max_element
Parameters
xContainer whose maximum element's index is to be found
Returns
int Zero-based index position of the maximum element

Uses std::max_element to find the iterator to the maximum element, then converts to index position using std::distance. For containers with multiple maximum elements, returns the first occurrence.

Definition at line 365 of file math.h.

◆ whichMin()

template<typename T >
int slope::whichMin ( const T &  x)

Returns the index of the minimum element in a container.

Template Parameters
TContainer type that supports iterators and std::min_element
Parameters
xContainer whose maximum element's index is to be found
Returns
int Zero-based index position of the maximum element

Uses std::max_element to find the iterator to the maximum element, then converts to index position using std::distance. For containers with multiple maximum elements, returns the first occurrence.

Definition at line 383 of file math.h.