slope 0.29.0
Loading...
Searching...
No Matches
hybrid_cd.cpp
1#include "hybrid_cd.h"
2
3namespace slope {
4
5std::pair<double, double>
6computeClusterGradientAndHessian(const Eigen::MatrixXd& x,
7 const int j,
8 const std::vector<int>& s,
9 const Clusters& clusters,
10 const Eigen::VectorXd& w,
11 const Eigen::VectorXd& residual,
12 const Eigen::VectorXd& x_centers,
13 const Eigen::VectorXd& x_scales,
14 const JitNormalization jit_normalization)
15{
16 int n = x.rows();
17
18 Eigen::VectorXd x_s = Eigen::VectorXd::Zero(n);
19
20 auto s_it = s.cbegin();
21 auto c_it = clusters.cbegin(j);
22
23 for (; c_it != clusters.cend(j); ++c_it, ++s_it) {
24 int k = *c_it;
25 double s_k = *s_it;
26
27 switch (jit_normalization) {
29 x_s += x.col(k) * (s_k / x_scales(k));
30 x_s.array() -= x_centers(k) * s_k / x_scales(k);
31 break;
32
34 x_s += x.col(k) * s_k;
35 x_s.array() -= x_centers(k) * s_k;
36 break;
37
39 x_s += x.col(k) * (s_k / x_scales(k));
40 break;
41
43 x_s += x.col(k) * s_k;
44 break;
45 }
46 }
47
48 double hessian_j = x_s.cwiseAbs2().dot(w) / n;
49 double gradient_j = x_s.cwiseProduct(w).dot(residual) / n;
50
51 return { hessian_j, gradient_j };
52}
53
54std::pair<double, double>
55computeClusterGradientAndHessian(const Eigen::SparseMatrix<double>& x,
56 const int j,
57 const std::vector<int>& s,
58 const Clusters& clusters,
59 const Eigen::VectorXd& w,
60 const Eigen::VectorXd& residual,
61 const Eigen::VectorXd& x_centers,
62 const Eigen::VectorXd& x_scales,
63 const JitNormalization jit_normalization)
64{
65 int n = x.rows();
66
67 Eigen::VectorXd weighted_residual = w.cwiseProduct(residual);
68
69 std::vector<Eigen::Triplet<double>> triplets;
70
71 Eigen::SparseMatrix<double> x_s(n, 1);
72
73 double offset = 0;
74
75 auto s_it = s.cbegin();
76 auto c_it = clusters.cbegin(j);
77
78 for (; c_it != clusters.cend(j); ++c_it, ++s_it) {
79 int k = *c_it;
80 double s_k = *s_it;
81
82 switch (jit_normalization) {
84 offset += x_centers(k) * s_k;
85 break;
87 offset += x_centers(k) * s_k / x_scales(k);
88 break;
90 break;
92 break;
93 }
94
95 double v = 0;
96
97 for (Eigen::SparseMatrix<double>::InnerIterator it(x, k); it; ++it) {
98 switch (jit_normalization) {
100 [[fallthrough]];
102 v = it.value() * s_k;
103 break;
104
106 [[fallthrough]];
108 v = it.value() * s_k / x_scales(k);
109 break;
110 }
111
112 triplets.emplace_back(it.row(), 0, v);
113 }
114 }
115
116 x_s.setFromTriplets(triplets.begin(), triplets.end());
117
118 double hess = 1;
119 double grad = 0;
120
121 assert(x_s.nonZeros() > 0);
122
123 switch (jit_normalization) {
125 [[fallthrough]];
127 hess = (x_s.col(0).cwiseAbs2().dot(w) - 2 * offset * x_s.col(0).dot(w) +
128 std::pow(offset, 2) * w.sum()) /
129 n;
130 grad = x_s.col(0).dot(weighted_residual) / n -
131 offset * weighted_residual.sum() / n;
132 break;
133
135 [[fallthrough]];
137 hess = x_s.col(0).cwiseAbs2().dot(w) / n;
138 grad = x_s.col(0).dot(weighted_residual) / n;
139 break;
140 }
141
142 return { hess, grad };
143}
144
145} // namespace slope
Representation of the clusters in SLOPE.
Definition clusters.h:18
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.
Definition clusters.cpp:48
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.
Definition clusters.cpp:38
An implementation of the coordinate descent step in the hybrid algorithm for solving SLOPE.
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
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)
Definition hybrid_cd.cpp:6
JitNormalization
Enums to control predictor standardization behavior.
@ None
No JIT normalization.