![]() |
slope 0.29.0
|
Hybrid CD-PGD solver for SLOPE. More...
#include <hybrid.h>
Public Member Functions | |
Hybrid (JitNormalization jit_normalization, bool intercept, bool update_clusters, int cd_iterations) | |
Constructs Hybrid solver for SLOPE optimization. | |
void | run (Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::MatrixXd &eta, const Eigen::ArrayXd &lambda, const std::unique_ptr< Loss > &loss, const SortedL1Norm &penalty, const Eigen::VectorXd &gradient, const std::vector< int > &working_set, const Eigen::MatrixXd &x, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::MatrixXd &y) override |
Pure virtual function defining the solver's optimization routine. | |
void | run (Eigen::VectorXd &beta0, Eigen::VectorXd &beta, Eigen::MatrixXd &eta, const Eigen::ArrayXd &lambda, const std::unique_ptr< Loss > &loss, const SortedL1Norm &penalty, const Eigen::VectorXd &gradient, const std::vector< int > &working_set, const Eigen::SparseMatrix< double > &x, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::MatrixXd &y) override |
Pure virtual function defining the solver's optimization routine. | |
![]() | |
SolverBase (JitNormalization jit_normalization, bool intercept) | |
Constructs a base solver for SLOPE optimization. | |
virtual | ~SolverBase ()=default |
Default desstructor. | |
Additional Inherited Members | |
![]() | |
JitNormalization | jit_normalization |
JIT feature normalization strategy. | |
bool | intercept |
If true, fits intercept term. | |
Hybrid CD-PGD solver for SLOPE.
This solver alternates between coordinate descent (CD) and proximal gradient descent (PGD) steps to solve the SLOPE optimization problem. Because the SLOPE problem is non-separable, CD steps do not work out-of-the-box. In particular, they cannot be used to split he clusters (or would need to be augmented with additional steps). Instead, we use a hybrid method:
The switching between methods is controlled by the cd_iterations parameter, which determines how often PGD steps are taken versus CD steps.
|
inline |
Constructs Hybrid solver for SLOPE optimization.
jit_normalization | Feature normalization strategy |
intercept | If true, fits intercept term |
update_clusters | If true, updates clusters during optimization |
cd_iterations | Frequency of proximal gradient descent updates |
|
overridevirtual |
Pure virtual function defining the solver's optimization routine.
beta0 | Intercept terms for each response |
beta | Coefficients (size p x m) |
eta | Linear predictor matrix (n samples x m responses) |
lambda | Vector of regularization parameters |
loss | Pointer to loss function object |
penalty | Sorted L1 norm object for proximal operations |
gradient | Gradient matrix for loss function |
working_set | Vector of indices for active predictors |
x | Input feature matrix (n samples x p predictors) |
x_centers | Vector of feature means for centering |
x_scales | Vector of feature scales for normalization |
y | Response matrix (n samples x m responses) |
Implements slope::SolverBase.
Definition at line 17 of file hybrid.cpp.
|
overridevirtual |
Pure virtual function defining the solver's optimization routine.
beta0 | Intercept terms for each response |
beta | Coefficients (size p x m) |
eta | Linear predictor matrix (n samples x m responses) |
lambda | Vector of regularization parameters |
loss | Pointer to loss function object |
penalty | Sorted L1 norm object for proximal operations |
gradient | Gradient matrix for loss function |
working_set | Vector of indices for active predictors |
x | Input feature matrix (n samples x p predictors) |
x_centers | Vector of feature means for centering |
x_scales | Vector of feature scales for normalization |
y | Response matrix (n samples x m responses) |
Implements slope::SolverBase.
Definition at line 46 of file hybrid.cpp.