11#include <Eigen/SparseCore>
37 return (T(0) < val) - (val < T(0));
51cumSum(
const T& x,
const bool leading_zero =
false)
53 const size_t start = leading_zero ? 1 : 0;
54 Eigen::ArrayXd out(x.size() + start);
60 std::partial_sum(x.begin(), x.end(), out.begin() + start);
78 return 1.0 / (1.0 + std::exp(-x));
94 assert(x > 0 && x < 1 &&
"Input must be in (0, 1)");
96 return std::log(x) - std::log1p(-x);
110clamp(
const T& x,
const T& lo,
const T& hi)
112 return x < lo ? lo : x > hi ? hi : x;
152 const std::vector<int>& active_set,
153 const Eigen::VectorXd& beta0,
154 const Eigen::VectorXd& beta,
155 const Eigen::VectorXd& x_centers,
156 const Eigen::VectorXd& x_scales,
158 const bool intercept)
162 int m = beta0.size();
164 Eigen::MatrixXd eta = Eigen::MatrixXd::Zero(n, m);
167 bool large_problem = active_set.size() > 100 && n * active_set.size() > 1e7;
168#pragma omp parallel num_threads(Threads::get()) if (large_problem)
171 Eigen::MatrixXd eta_local = Eigen::MatrixXd::Zero(n, m);
174#pragma omp for nowait
176 for (
int i = 0; i < static_cast<int>(active_set.size()); ++i) {
177 int ind = active_set[i];
178 auto [k, j] = std::div(ind, p);
180 switch (jit_normalization) {
182 eta_local.col(k) += x.col(j) * beta(ind) / x_scales(j);
183 eta_local.col(k).array() -= beta(ind) * x_centers(j) / x_scales(j);
187 eta_local.col(k) += x.col(j) * beta(ind);
188 eta_local.col(k).array() -= beta(ind) * x_centers(j);
192 eta_local.col(k) += x.col(j) * beta(ind) / x_scales(j);
196 eta_local.col(k) += x.col(j) * beta(ind);
210 eta.rowwise() += beta0.transpose();
234 const Eigen::MatrixXd& residual,
235 const std::vector<int>& active_set,
236 const Eigen::VectorXd& x_centers,
237 const Eigen::VectorXd& x_scales,
238 const Eigen::VectorXd& w,
241 const int n = x.rows();
242 const int p = x.cols();
243 const int m = residual.cols();
245 assert(gradient.size() == p * m &&
246 "Gradient matrix has incorrect dimensions");
248 Eigen::MatrixXd weighted_residual(n, m);
249 Eigen::ArrayXd wr_sums(m);
252 bool large_problem = active_set.size() > 100 && n * active_set.size() > 1e5;
253#pragma omp parallel for num_threads(Threads::get()) if (large_problem)
255 for (
int k = 0; k < m; ++k) {
256 weighted_residual.col(k) = residual.col(k).cwiseProduct(w);
257 wr_sums(k) = weighted_residual.col(k).sum();
261#pragma omp parallel for num_threads(Threads::get()) if (large_problem)
263 for (
int i = 0; i < static_cast<int>(active_set.size()); ++i) {
264 int ind = active_set[i];
265 auto [k, j] = std::div(ind, p);
267 switch (jit_normalization) {
270 (x.col(j).dot(weighted_residual.col(k)) - x_centers(j) * wr_sums(k)) /
275 (x.col(j).dot(weighted_residual.col(k)) - x_centers(j) * wr_sums(k)) /
280 x.col(j).dot(weighted_residual.col(k)) / (x_scales(j) * n);
283 gradient(ind) = x.col(j).dot(weighted_residual.col(k)) / n;
306 const Eigen::VectorXd& offset,
307 const std::vector<int>& active_set,
308 const Eigen::VectorXd& x_centers,
309 const Eigen::VectorXd& x_scales,
312 const int n = x.rows();
313 const int p = x.cols();
315 for (
size_t i = 0; i < active_set.size(); ++i) {
316 int ind = active_set[i];
317 auto [k, j] = std::div(ind, p);
319 switch (jit_normalization) {
322 offset(k) * (x.col(j).sum() / n - x_centers(j)) / x_scales(j);
325 gradient(ind) -= offset(k) * (x.col(j).sum() / n - x_centers(j));
328 gradient(ind) -= offset(k) * x.col(j).sum() / (n * x_scales(j));
331 gradient(ind) -= offset(k) * x.col(j).sum() / n;
346setUnion(
const std::vector<int>& a,
const std::vector<int>& b);
359setDiff(
const std::vector<int>& a,
const std::vector<int>& b);
376 return std::distance(x.begin(), std::max_element(x.begin(), x.end()));
394 return std::distance(x.begin(), std::min_element(x.begin(), x.end()));
410template<
typename T,
typename Comparator>
417 return std::distance(x.begin(), std::max_element(x.begin(), x.end(), comp));
435geomSpace(
const double start,
const double end,
const int n);
453 const int p = x.cols();
455 Eigen::VectorXd out(p);
457 for (
int j = 0; j < p; ++j) {
458 out(j) = x.col(j).cwiseAbs().sum();
478 const int p = x.cols();
480 Eigen::VectorXd out(p);
482 for (
int j = 0; j < p; ++j) {
483 out(j) = x.col(j).norm();
503 return x.colwise().norm();
519maxAbs(
const Eigen::SparseMatrixBase<T>& x)
521 const int p = x.cols();
523 Eigen::VectorXd out(p);
525 for (
int j = 0; j < p; ++j) {
526 double x_j_maxabs = 0.0;
528 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
529 x_j_maxabs = std::max(x_j_maxabs, std::abs(it.value()));
554 return x.cwiseAbs().colwise().maxCoeff();
570means(
const Eigen::SparseMatrixBase<T>& x)
572 const int n = x.rows();
573 const int p = x.cols();
575 Eigen::VectorXd out(p);
577 for (
int j = 0; j < p; ++j) {
578 out(j) = x.col(j).sum() / n;
599 return x.colwise().mean();
618 const int n = x.rows();
619 const int p = x.cols();
621 Eigen::VectorXd x_means =
means(x);
622 Eigen::VectorXd out(p);
624 for (
int j = 0; j < p; ++j) {
625 double sum_sq_diff = 0.0;
626 const double mean = x_means(j);
629 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
630 double diff = it.value() - mean;
631 sum_sq_diff += diff * diff;
635 int nz_count = x.col(j).nonZeros();
637 sum_sq_diff += (n - nz_count) * mean * mean;
641 out(j) = std::sqrt(sum_sq_diff / n);
666 Eigen::VectorXd x_means =
means(x);
667 Eigen::VectorXd out(p);
669 for (
int j = 0; j < p; ++j) {
670 out(j) = (x.col(j).array() - x_means(j)).matrix().norm();
673 out.array() /= std::sqrt(n);
690ranges(
const Eigen::SparseMatrixBase<T>& x)
692 const int p = x.cols();
694 Eigen::VectorXd out(p);
696 for (
int j = 0; j < p; ++j) {
697 double x_j_max = 0.0;
698 double x_j_min = 0.0;
700 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
701 x_j_max = std::max(x_j_max, it.value());
702 x_j_min = std::min(x_j_min, it.value());
705 out(j) = x_j_max - x_j_min;
728 return x.colwise().maxCoeff() - x.colwise().minCoeff();
743mins(
const Eigen::SparseMatrixBase<T>& x)
745 const int p = x.cols();
747 Eigen::VectorXd out(p);
749 for (
int j = 0; j < p; ++j) {
750 double x_j_min = 0.0;
752 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
753 x_j_min = std::min(x_j_min, it.value());
777mins(
const Eigen::MatrixBase<T>& x)
779 return x.colwise().minCoeff();
803 Eigen::MatrixXd& residual,
806 const Eigen::MatrixXd& w,
807 const Eigen::VectorXd& x_centers,
808 const Eigen::VectorXd& x_scales,
811 using namespace Eigen;
813 const int n = x.rows();
814 const int p = x.cols();
815 const int n_clusters = clusters.
size();
817 Eigen::VectorXd gradient = Eigen::VectorXd::Zero(n_clusters);
819 for (
int j = 0; j < n_clusters; ++j) {
820 double c_old = clusters.
coeff(j);
829 s.reserve(cluster_size);
831 for (
auto c_it = clusters.
cbegin(j); c_it != clusters.
cend(j); ++c_it) {
833 double s_k =
sign(beta(ind));
840 if (cluster_size == 1) {
841 int k = *clusters.
cbegin(j);
843 x, k, w, residual, x_centers, x_scales, s[0], jit_normalization, n);
846 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.
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.
Thread management for parallel computations.