9#include "../clusters.h"
40std::pair<double, double>
43 const Eigen::VectorXd& w,
44 const Eigen::VectorXd& residual,
45 const Eigen::VectorXd& x_centers,
46 const Eigen::VectorXd& x_scales,
51 double gradient = 0.0;
54 switch (jit_normalization) {
57 (x.col(k).cwiseProduct(w).dot(residual) -
58 w.dot(residual) * x_centers(k)) /
61 (x.col(k).cwiseAbs2().dot(w) - 2 * x_centers(k) * x.col(k).dot(w) +
62 std::pow(x_centers(k), 2) * w.sum()) /
63 (std::pow(x_scales(k), 2) * n);
68 (x.col(k).cwiseProduct(w).dot(residual) -
69 w.dot(residual) * x_centers(k)) /
72 (x.col(k).cwiseAbs2().dot(w) - 2 * x_centers(k) * x.col(k).dot(w) +
73 std::pow(x_centers(k), 2) * w.sum()) /
79 s * (x.col(k).cwiseProduct(w).dot(residual)) / (n * x_scales(k));
80 hessian = x.col(k).cwiseAbs2().dot(w) / (std::pow(x_scales(k), 2) * n);
84 gradient = s * (x.col(k).cwiseProduct(w).dot(residual)) / n;
85 hessian = x.col(k).cwiseAbs2().dot(w) / n;
89 return { gradient, hessian };
115std::pair<double, double>
118 const std::vector<int>& s,
119 const Clusters& clusters,
120 const Eigen::VectorXd& w,
121 const Eigen::VectorXd& residual,
122 const Eigen::VectorXd& x_centers,
123 const Eigen::VectorXd& x_scales,
148std::pair<double, double>
151 const std::vector<int>& s,
152 const Clusters& clusters,
153 const Eigen::VectorXd& w,
154 const Eigen::VectorXd& residual,
155 const Eigen::VectorXd& x_centers,
156 const Eigen::VectorXd& x_scales,
188 Eigen::VectorXd& beta,
189 Eigen::VectorXd& residual,
191 const Eigen::ArrayXd& lambda,
193 const Eigen::VectorXd& w,
194 const Eigen::VectorXd& x_centers,
195 const Eigen::VectorXd& x_scales,
196 const bool intercept,
198 const bool update_clusters)
200 using namespace Eigen;
202 const int n = x.rows();
204 for (
int j = 0; j < clusters.
n_clusters(); ++j) {
205 double c_old = clusters.
coeff(j);
215 s.reserve(cluster_size);
217 for (
auto c_it = clusters.
cbegin(j); c_it != clusters.
cend(j); ++c_it) {
218 double s_k =
sign(beta(*c_it));
222 double hessian_j = 1;
223 double gradient_j = 0;
226 if (cluster_size == 1) {
227 int k = *clusters.
cbegin(j);
229 x, k, w, residual, x_centers, x_scales, s[0], jit_normalization, n);
232 x, j, s, clusters, w, residual, x_centers, x_scales, jit_normalization);
236 c_old - gradient_j / hessian_j, j, lambda / hessian_j, clusters);
238 double c_diff = c_old - c_tilde;
241 auto s_it = s.cbegin();
242 auto c_it = clusters.
cbegin(j);
243 for (; c_it != clusters.
cend(j); ++c_it, ++s_it) {
248 beta(k) = c_tilde * s_k;
251 switch (jit_normalization) {
253 residual -= x.col(k) * (s_k * c_diff / x_scales(k));
254 residual.array() += x_centers(k) * s_k * c_diff / x_scales(k);
258 residual -= x.col(k) * (s_k * c_diff);
259 residual.array() += x_centers(k) * s_k * c_diff;
263 residual -= x.col(k) * (s_k * c_diff / x_scales(k));
267 residual -= x.col(k) * (s_k * c_diff);
273 if (update_clusters) {
274 clusters.
update(j, new_index, std::abs(c_tilde));
276 clusters.
setCoeff(j, std::abs(c_tilde));
281 double beta0_update = residual.dot(w) / n;
282 residual.array() -= beta0_update;
283 beta0(0) -= beta0_update;
Representation of the clusters in SLOPE.
void update(const int old_index, const int new_index, const double c_new)
Updates the cluster structure when an index is changed.
double coeff(const int i) const
Returns the coefficient of the cluster with the given index.
int cluster_size(const int i) const
Returns the size of the cluster with the given index.
void setCoeff(const int i, const double x)
Sets the coefficient of the cluster with the given index.
int n_clusters() const
Returns the number of clusters.
std::vector< int >::const_iterator cend(const int i) const
Returns a constant iterator pointing to the end of the cluster with the given index.
std::vector< int >::const_iterator cbegin(const int i) const
Returns a constant iterator pointing to the beginning of the cluster with the given index.
Namespace containing SLOPE regression implementation.
int sign(T val)
Returns the sign of a given value.
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)
JitNormalization
Enums to control predictor standardization behavior.
@ None
No JIT normalization.
std::tuple< double, int > slopeThreshold(const double x, const int j, const Eigen::ArrayXd lambdas, const Clusters &clusters)
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::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)
The declaration of the slopeThreshold function.