10#include <Eigen/SparseCore>
16#include <unordered_set>
32template<
typename Derived>
34nonZeros(
const Eigen::SparseMatrixBase<Derived>& x)
36 return x.derived().nonZeros();
45template<
typename Derived>
49 return (x.derived().array() != 0).count();
70sort(T& v,
const bool descending =
false)
74 v.data(), v.data() + v.size(), std::greater<typename T::value_type>());
77 v.data(), v.data() + v.size(), std::less<typename T::value_type>());
98 for (
int i = 0; i < x.size(); i++) {
128 vector<int> idx(v.size());
129 iota(idx.begin(), idx.end(), 0);
132 sort(idx.begin(), idx.end(), [&v](
int i,
int j) { return v[i] > v[j]; });
134 sort(idx.begin(), idx.end(), [&v](
int i,
int j) { return v[i] < v[j]; });
153permute(T& values,
const std::vector<int>& ind)
158 T out(values.size());
163 for (
int i = 0; i < values.size(); ++i)
164 out[i] = std::move(values[ind[i]]);
169 values = std::move(out);
189 T out(values.size());
191 for (
int i = 0; i < values.size(); ++i)
192 out[ind[i]] = std::move(values[i]);
194 values = std::move(out);
211move_elements(std::vector<T>& v,
const int from,
const int to,
const int size)
219 assert(from + size <=
static_cast<int>(v.size()));
220 std::rotate(v.begin() + to, v.begin() + from, v.begin() + from + size);
222 assert(to + size <=
static_cast<int>(v.size()));
224 v.begin() + from, v.begin() + from + size, v.begin() + to + size);
242 const std::set<std::string>& valid_options,
243 const std::string& parameter_name);
258subset(
const Eigen::EigenBase<T>& x,
const std::vector<int>& indices)
260 return subset(x.derived(), indices);
275typename Eigen::MatrixBase<T>::PlainObject
276subset(
const Eigen::DenseBase<T>& x,
const std::vector<int>& indices)
278 return x.derived()(indices, all);
296subset(
const Eigen::SparseMatrixBase<T>& x,
const std::vector<int>& indices)
298 std::vector<Eigen::Triplet<double>> triplets;
301 for (
int j = 0; j < x.cols(); ++j) {
302 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
303 auto it_idx = std::find(indices.begin(), indices.end(), it.row());
305 if (it_idx != indices.end()) {
306 int new_row = std::distance(indices.begin(), it_idx);
307 triplets.emplace_back(new_row, j, it.value());
312 T out(indices.size(), x.cols());
313 out.setFromTriplets(triplets.begin(), triplets.end());
332subsetCols(
const Eigen::MatrixBase<T>& x,
const std::vector<int>& indices)
334 return x.derived()(all, indices);
351subsetCols(
const Eigen::SparseMatrixBase<T>& x,
const std::vector<int>& indices)
353 std::vector<Eigen::Triplet<double>> triplets;
356 for (
size_t j_idx = 0; j_idx < indices.size(); ++j_idx) {
357 int j = indices[j_idx];
358 for (
typename T::InnerIterator it(x.derived(), j); it; ++it) {
359 triplets.emplace_back(it.row(), j_idx, it.value());
363 T out(x.rows(), indices.size());
364 out.setFromTriplets(triplets.begin(), triplets.end());
380inline std::unordered_set<double>
383 std::unordered_set<double>
unique;
384 for (Eigen::Index j = 0; j < x.cols(); j++) {
385 for (Eigen::Index i = 0; i < x.rows(); i++) {
403template<
typename Derived>
407 return x.derived().array().isFinite().all();
422template<
typename Derived>
426 for (
int j = 0; j < x.cols(); ++j) {
427 for (
typename Derived::InnerIterator it(x.derived(), j); it; ++it) {
428 if (!std::isfinite(it.value())) {
Eigen compatibility layer for version differences.
Namespace containing SLOPE regression implementation.
std::unordered_set< double > unique(const Eigen::MatrixXd &x)
Create a set of unique values from an Eigen matrix.
void move_elements(std::vector< T > &v, const int from, const int to, const int size)
std::vector< int > which(const T &x)
T subsetCols(const Eigen::MatrixBase< T > &x, const std::vector< int > &indices)
Extract specified columns from a dense matrix.
T subset(const Eigen::EigenBase< T > &x, const std::vector< int > &indices)
Extract a subset of rows from an Eigen matrix.
void validateOption(const std::string &value, const std::set< std::string > &valid_options, const std::string ¶meter_name)
Validates if a given value exists in a set of valid options.
std::vector< int > sortIndex(T &v, const bool descending=false)
void permute(T &values, const std::vector< int > &ind)
void sort(T &v, const bool descending=false)
void inversePermute(T &values, const std::vector< int > &ind)
Eigen::Index nonZeros(const Eigen::SparseMatrixBase< Derived > &x)
bool isFinite(const Eigen::DenseBase< Derived > &x)
Check if all elements in a dense matrix are finite.