13#include <Eigen/SparseCore>
39 return (T(0) < val) - (val < T(0));
53cumSum(
const T& x,
const bool leading_zero =
false)
55 const size_t start = leading_zero ? 1 : 0;
56 Eigen::ArrayXd out(x.size() + start);
62 std::partial_sum(x.begin(), x.end(), out.begin() + start);
80 return 1.0 / (1.0 + std::exp(-x));
96 assert(x > 0 && x < 1 &&
"Input must be in (0, 1)");
98 return std::log(x) - std::log1p(-x);
112clamp(
const T& x,
const T& lo,
const T& hi)
114 return x < lo ? lo : x > hi ? hi : x;
154 const std::vector<int>& active_set,
155 const Eigen::VectorXd& beta0,
156 const Eigen::VectorXd& beta,
157 const Eigen::VectorXd& x_centers,
158 const Eigen::VectorXd& x_scales,
160 const bool intercept)
164 int m = beta0.size();
166 Eigen::MatrixXd eta = Eigen::MatrixXd::Zero(n, m);
169 bool large_problem = active_set.size() > 100 && n * active_set.size() > 1e7;
170#pragma omp parallel num_threads(Threads::get()) if (large_problem)
173 Eigen::MatrixXd eta_local = Eigen::MatrixXd::Zero(n, m);
176#pragma omp for nowait
178 for (
int i = 0; i < static_cast<int>(active_set.size()); ++i) {
179 int ind = active_set[i];
180 auto [k, j] = std::div(ind, p);
182 switch (jit_normalization) {
184 eta_local.col(k) += x.col(j) * beta(ind) / x_scales(j);
185 eta_local.col(k).array() -= beta(ind) * x_centers(j) / x_scales(j);
189 eta_local.col(k) += x.col(j) * beta(ind);
190 eta_local.col(k).array() -= beta(ind) * x_centers(j);
194 eta_local.col(k) += x.col(j) * beta(ind) / x_scales(j);
198 eta_local.col(k) += x.col(j) * beta(ind);
212 eta.rowwise() += beta0.transpose();
236 const Eigen::MatrixXd& residual,
237 const std::vector<int>& active_set,
238 const Eigen::VectorXd& x_centers,
239 const Eigen::VectorXd& x_scales,
240 const Eigen::VectorXd& w,
243 const int n = x.rows();
244 const int p = x.cols();
245 const int m = residual.cols();
247 assert(gradient.size() == p * m &&
248 "Gradient matrix has incorrect dimensions");
250 Eigen::MatrixXd weighted_residual(n, m);
251 Eigen::ArrayXd wr_sums(m);
254 bool large_problem = active_set.size() > 100 && n * active_set.size() > 1e5;
255#pragma omp parallel for num_threads(Threads::get()) if (large_problem)
257 for (
int k = 0; k < m; ++k) {
258 weighted_residual.col(k) = residual.col(k).cwiseProduct(w);
259 wr_sums(k) = weighted_residual.col(k).sum();
263#pragma omp parallel for num_threads(Threads::get()) if (large_problem)
265 for (
int i = 0; i < static_cast<int>(active_set.size()); ++i) {
266 int ind = active_set[i];
267 auto [k, j] = std::div(ind, p);
269 switch (jit_normalization) {
272 (x.col(j).dot(weighted_residual.col(k)) - x_centers(j) * wr_sums(k)) /
277 (x.col(j).dot(weighted_residual.col(k)) - x_centers(j) * wr_sums(k)) /
282 x.col(j).dot(weighted_residual.col(k)) / (x_scales(j) * n);
285 gradient(ind) = x.col(j).dot(weighted_residual.col(k)) / n;
308 const Eigen::VectorXd& offset,
309 const std::vector<int>& active_set,
310 const Eigen::VectorXd& x_centers,
311 const Eigen::VectorXd& x_scales,
314 const int n = x.rows();
315 const int p = x.cols();
317 for (
size_t i = 0; i < active_set.size(); ++i) {
318 int ind = active_set[i];
319 auto [k, j] = std::div(ind, p);
321 switch (jit_normalization) {
324 offset(k) * (x.col(j).sum() / n - x_centers(j)) / x_scales(j);
327 gradient(ind) -= offset(k) * (x.col(j).sum() / n - x_centers(j));
330 gradient(ind) -= offset(k) * x.col(j).sum() / (n * x_scales(j));
333 gradient(ind) -= offset(k) * x.col(j).sum() / n;
348setUnion(
const std::vector<int>& a,
const std::vector<int>& b);
361setDiff(
const std::vector<int>& a,
const std::vector<int>& b);
378 return std::distance(x.begin(), std::max_element(x.begin(), x.end()));
396 return std::distance(x.begin(), std::min_element(x.begin(), x.end()));
412template<
typename T,
typename Comparator>
419 return std::distance(x.begin(), std::max_element(x.begin(), x.end(), comp));
437geomSpace(
const double start,
const double end,
const int n);
455 const int p = x.cols();
457 Eigen::VectorXd out(p);
459 for (
int j = 0; j < p; ++j) {
460 out(j) = x.col(j).cwiseAbs().sum();
480 const int p = x.cols();
482 Eigen::VectorXd out(p);
484 for (
int j = 0; j < p; ++j) {
485 out(j) = x.col(j).norm();
505 return x.colwise().norm();
521maxAbs(
const Eigen::SparseMatrixBase<T>& x)
523 const int p = x.cols();
525 Eigen::VectorXd out(p);
527 for (
int j = 0; j < p; ++j) {
528 double x_j_maxabs = 0.0;
530 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
531 x_j_maxabs = std::max(x_j_maxabs, std::abs(it.value()));
556 return x.cwiseAbs().colwise().maxCoeff();
572means(
const Eigen::SparseMatrixBase<T>& x)
574 const int n = x.rows();
575 const int p = x.cols();
577 Eigen::VectorXd out(p);
579 for (
int j = 0; j < p; ++j) {
580 out(j) = x.col(j).sum() / n;
601 return x.colwise().mean();
620 const int n = x.rows();
621 const int p = x.cols();
623 Eigen::VectorXd x_means =
means(x);
624 Eigen::VectorXd out(p);
626 for (
int j = 0; j < p; ++j) {
627 double sum_sq_diff = 0.0;
628 const double mean = x_means(j);
631 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
632 double diff = it.value() - mean;
633 sum_sq_diff += diff * diff;
639 sum_sq_diff += (n - nz_count) * mean * mean;
643 out(j) = std::sqrt(sum_sq_diff / n);
668 Eigen::VectorXd x_means =
means(x);
669 Eigen::VectorXd out(p);
671 for (
int j = 0; j < p; ++j) {
672 out(j) = (x.col(j).array() - x_means(j)).matrix().norm();
675 out.array() /= std::sqrt(n);
692ranges(
const Eigen::SparseMatrixBase<T>& x)
694 const int p = x.cols();
696 Eigen::VectorXd out(p);
698 for (
int j = 0; j < p; ++j) {
699 double x_j_max = 0.0;
700 double x_j_min = 0.0;
702 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
703 x_j_max = std::max(x_j_max, it.value());
704 x_j_min = std::min(x_j_min, it.value());
707 out(j) = x_j_max - x_j_min;
730 return x.colwise().maxCoeff() - x.colwise().minCoeff();
745mins(
const Eigen::SparseMatrixBase<T>& x)
747 const int p = x.cols();
749 Eigen::VectorXd out(p);
751 for (
int j = 0; j < p; ++j) {
752 double x_j_min = 0.0;
754 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
755 x_j_min = std::min(x_j_min, it.value());
779mins(
const Eigen::MatrixBase<T>& x)
781 return x.colwise().minCoeff();
805 Eigen::MatrixXd& residual,
808 const Eigen::MatrixXd& w,
809 const Eigen::VectorXd& x_centers,
810 const Eigen::VectorXd& x_scales,
813 using namespace Eigen;
815 const int n = x.rows();
816 const int p = x.cols();
817 const int n_clusters = clusters.
size();
819 Eigen::VectorXd gradient = Eigen::VectorXd::Zero(n_clusters);
821 for (
int j = 0; j < n_clusters; ++j) {
822 double c_old = clusters.
coeff(j);
831 s.reserve(cluster_size);
833 for (
auto c_it = clusters.
cbegin(j); c_it != clusters.
cend(j); ++c_it) {
835 double s_k =
sign(beta(ind));
842 if (cluster_size == 1) {
843 int k = *clusters.
cbegin(j);
845 x, k, w, residual, x_centers, x_scales, s[0], jit_normalization, n);
848 x, j, s, clusters, w, residual, x_centers, x_scales, jit_normalization);
Representation of the clusters in SLOPE.
double coeff(const int i) const
Returns the coefficient of the cluster with the given index.
std::size_t size()
Returns the number of clusters.
int cluster_size(const int i) const
Returns the size of the cluster with the given index.
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.
The declaration of the Clusters class.
Eigen compatibility layer for version differences.
Enums to control predictor standardization behavior.
Namespace containing SLOPE regression implementation.
Eigen::VectorXd means(const Eigen::SparseMatrixBase< T > &x)
Computes the arithmetic mean for each column of a sparse matrix.
T clamp(const T &x, const T &lo, const T &hi)
Eigen::VectorXd l1Norms(const T &x)
Computes the L1 (Manhattan) norms for each column of a matrix.
int sign(T val)
Returns the sign of a given value.
int whichMin(const T &x)
Returns the index of the minimum element in a container.
std::pair< double, double > computeGradientAndHessian(const T &x, const int ind, const Eigen::MatrixXd &w, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const double s, const JitNormalization jit_normalization, const int n)
JitNormalization
Enums to control predictor standardization behavior.
@ None
No JIT normalization.
Eigen::VectorXd logSumExp(const Eigen::MatrixXd &a)
std::vector< int > setDiff(const std::vector< int > &a, const std::vector< int > &b)
Computes the set difference of two sorted integer vectors.
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)
Eigen::VectorXd ranges(const Eigen::SparseMatrixBase< T > &x)
Computes the range (max - min) for each column of a matrix.
Eigen::MatrixXd softmax(const Eigen::MatrixXd &x)
std::vector< int > setUnion(const std::vector< int > &a, const std::vector< int > &b)
Computes the union of two sorted integer vectors.
Eigen::VectorXd mins(const Eigen::SparseMatrixBase< T > &x)
Computes the minimum value for each column of a sparse matrix.
Eigen::ArrayXd cumSum(const T &x, const bool leading_zero=false)
int whichBest(const T &x, const Comparator &comp)
Returns the index of the minimum element in a container.
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)
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)
Eigen::VectorXd stdDevs(const Eigen::SparseMatrixBase< T > &x)
Computes the standard deviation for each column of a matrix.
Eigen::VectorXd l2Norms(const Eigen::SparseMatrixBase< T > &x)
Computes the L2 (Euclidean) norms for each column of a sparse matrix.
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)
Eigen::VectorXd maxAbs(const Eigen::SparseMatrixBase< T > &x)
Computes the maximum absolute value for each column of a matrix.
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.
std::pair< double, double > computeClusterGradientAndHessian(const Eigen::MatrixBase< T > &x, const int c_ind, const std::vector< int > &s, const Clusters &clusters, const Eigen::MatrixXd &w, const Eigen::MatrixXd &residual, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization jit_normalization)
int whichMax(const T &x)
Returns the index of the maximum element in a container.
Eigen::Index nonZeros(const Eigen::SparseMatrixBase< Derived > &x)
Thread management for parallel computations.
Various utility functions.