9 Eigen::VectorXd max_vals = a.rowwise().maxCoeff();
10 Eigen::ArrayXd sum_exp =
11 (a.colwise() - max_vals).array().exp().rowwise().sum();
19 Eigen::VectorXd shift = a.rowwise().maxCoeff();
20 Eigen::MatrixXd exp_a = (a.colwise() - shift).array().exp();
21 Eigen::ArrayXd row_sums = exp_a.rowwise().sum();
22 return exp_a.array().colwise() / row_sums;
26setUnion(
const std::vector<int>& a,
const std::vector<int>& b)
30 a.begin(), a.end(), b.begin(), b.end(), std::back_inserter(out));
36setDiff(
const std::vector<int>& a,
const std::vector<int>& b)
40 a.begin(), a.end(), b.begin(), b.end(), std::back_inserter(out));
46geomSpace(
const double start,
const double end,
const int n)
49 return Eigen::ArrayXd::Constant(1, start);
51 return Eigen::ArrayXd::LinSpaced(n, std::log(start), std::log(end)).exp();
56l2Norms(
const Eigen::SparseMatrix<double>& x)
58 const int p = x.cols();
60 Eigen::VectorXd out(p);
62 for (
int j = 0; j < p; ++j) {
63 out(j) = x.col(j).norm();
72 return x.colwise().norm();
76means(
const Eigen::SparseMatrix<double>& x)
78 const int n = x.rows();
79 const int p = x.cols();
81 Eigen::VectorXd out(p);
83 for (
int j = 0; j < p; ++j) {
84 out(j) = x.col(j).sum() / n;
93 return x.colwise().mean();
97ranges(
const Eigen::SparseMatrix<double>& x)
99 const int p = x.cols();
101 Eigen::VectorXd out(p);
103 for (
int j = 0; j < p; ++j) {
104 double x_j_max = 0.0;
105 double x_j_min = 0.0;
107 for (
typename Eigen::SparseMatrix<double>::InnerIterator it(x, j); it;
109 x_j_max = std::max(x_j_max, it.value());
110 x_j_min = std::min(x_j_min, it.value());
113 out(j) = x_j_max - x_j_min;
122 return x.colwise().maxCoeff() - x.colwise().minCoeff();
126maxAbs(
const Eigen::SparseMatrix<double>& x)
128 const int p = x.cols();
130 Eigen::VectorXd out(p);
132 for (
int j = 0; j < p; ++j) {
133 double x_j_maxabs = 0.0;
135 for (
typename Eigen::SparseMatrix<double>::InnerIterator it(x, j); it;
137 x_j_maxabs = std::max(x_j_maxabs, std::abs(it.value()));
149 return x.cwiseAbs().colwise().maxCoeff();
153mins(
const Eigen::SparseMatrix<double>& x)
155 const int p = x.cols();
157 Eigen::VectorXd out(p);
159 for (
int j = 0; j < p; ++j) {
160 double x_j_min = 0.0;
162 for (
typename Eigen::SparseMatrix<double>::InnerIterator it(x, j); it;
164 x_j_min = std::min(x_j_min, it.value());
176 return x.colwise().minCoeff();
182 const int n = x.rows();
183 const int p = x.cols();
185 Eigen::VectorXd x_means =
means(x);
186 Eigen::VectorXd out(p);
188 for (
int j = 0; j < p; ++j) {
189 double sum_sq_diff = 0.0;
190 const double mean = x_means(j);
193 for (Eigen::SparseMatrix<double>::InnerIterator it(x, j); it; ++it) {
194 double diff = it.value() - mean;
195 sum_sq_diff += diff * diff;
199 int nz_count = x.col(j).nonZeros();
201 sum_sq_diff += (n - nz_count) * mean * mean;
205 out(j) = std::sqrt(sum_sq_diff / n);
217 Eigen::VectorXd x_means =
means(x);
218 Eigen::VectorXd out(p);
220 for (
int j = 0; j < p; ++j) {
221 out(j) = (x.col(j).array() - x_means(j)).matrix().norm();
224 out.array() /= std::sqrt(n);
Definitions of constants used in libslope.
Mathematical support functions for the slope package.
constexpr double P_MIN
Minimum allowed probability value to avoid numerical underflow.
Namespace containing SLOPE regression implementation.
Eigen::VectorXd maxAbs(const Eigen::SparseMatrix< double > &x)
Computes the maximum absolute value for each column of a matrix.
Eigen::VectorXd means(const Eigen::SparseMatrix< double > &x)
Computes the arithmetic mean for each column of a sparse matrix.
Eigen::VectorXd logSumExp(const Eigen::MatrixXd &a)
Eigen::MatrixXd softmax(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::VectorXd stdDevs(const Eigen::SparseMatrix< double > &x)
Computes the standard deviation for each column of a matrix.
std::vector< int > setUnion(const std::vector< int > &a, const std::vector< int > &b)
Computes the union of two sorted integer vectors.
Eigen::VectorXd ranges(const Eigen::SparseMatrix< double > &x)
Computes the range (max - min) for each column of a matrix.
Eigen::VectorXd l2Norms(const Eigen::SparseMatrix< double > &x)
Computes the L2 (Euclidean) norms for each column of a sparse matrix.
Eigen::VectorXd mins(const Eigen::SparseMatrix< double > &x)
Computes the minimum value for each column of a sparse 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.