slope 6.5.0
Loading...
Searching...
No Matches
utils.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "eigen_compat.h"
9#include <Eigen/Core>
10#include <Eigen/SparseCore>
11#include <algorithm>
12#include <cassert>
13#include <numeric>
14#include <set>
15#include <string>
16#include <unordered_set>
17#include <vector>
18
19namespace slope {
20
21using slope::all;
22
32template<typename Derived>
33Eigen::Index
34nonZeros(const Eigen::SparseMatrixBase<Derived>& x)
35{
36 return x.derived().nonZeros();
37}
38
45template<typename Derived>
46Eigen::Index
47nonZeros(const Eigen::DenseBase<Derived>& x)
48{
49 return (x.derived().array() != 0).count();
50}
51
68template<typename T>
69void
70sort(T& v, const bool descending = false)
71{
72 if (descending) {
73 std::sort(
74 v.data(), v.data() + v.size(), std::greater<typename T::value_type>());
75 } else {
76 std::sort(
77 v.data(), v.data() + v.size(), std::less<typename T::value_type>());
78 }
79}
80
93template<typename T>
94std::vector<int>
95which(const T& x)
96{
97 std::vector<int> out;
98 for (int i = 0; i < x.size(); i++) {
99 if (x[i]) {
100 out.emplace_back(i);
101 }
102 }
103
104 return out;
105}
106
122template<typename T>
123std::vector<int>
124sortIndex(T& v, const bool descending = false)
125{
126 using namespace std;
127
128 vector<int> idx(v.size());
129 iota(idx.begin(), idx.end(), 0);
130
131 if (descending) {
132 sort(idx.begin(), idx.end(), [&v](int i, int j) { return v[i] > v[j]; });
133 } else {
134 sort(idx.begin(), idx.end(), [&v](int i, int j) { return v[i] < v[j]; });
135 }
136
137 return idx;
138}
139
151template<typename T>
152void
153permute(T& values, const std::vector<int>& ind)
154{
158 T out(values.size());
159
163 for (int i = 0; i < values.size(); ++i)
164 out[i] = std::move(values[ind[i]]);
165
169 values = std::move(out);
170}
171
185template<typename T>
186void
187inversePermute(T& values, const std::vector<int>& ind)
188{
189 T out(values.size());
191 for (int i = 0; i < values.size(); ++i)
192 out[ind[i]] = std::move(values[i]);
193
194 values = std::move(out);
195}
196
209template<typename T>
210void
211move_elements(std::vector<T>& v, const int from, const int to, const int size)
212{
213 assert(from >= 0);
214 assert(to >= 0);
215 assert(size >= 0);
216 assert(from != to);
217
218 if (from > to) {
219 assert(from + size <= static_cast<int>(v.size()));
220 std::rotate(v.begin() + to, v.begin() + from, v.begin() + from + size);
221 } else {
222 assert(to + size <= static_cast<int>(v.size()));
223 std::rotate(
224 v.begin() + from, v.begin() + from + size, v.begin() + to + size);
225 }
226}
227
240void
241validateOption(const std::string& value,
242 const std::set<std::string>& valid_options,
243 const std::string& parameter_name);
244
256template<typename T>
257T
258subset(const Eigen::EigenBase<T>& x, const std::vector<int>& indices)
259{
260 return subset(x.derived(), indices);
261}
262
274template<typename T>
275typename Eigen::MatrixBase<T>::PlainObject
276subset(const Eigen::DenseBase<T>& x, const std::vector<int>& indices)
277{
278 return x.derived()(indices, all);
279}
280
294template<typename T>
295T
296subset(const Eigen::SparseMatrixBase<T>& x, const std::vector<int>& indices)
297{
298 std::vector<Eigen::Triplet<double>> triplets;
299 triplets.reserve(slope::nonZeros(x.derived()));
300
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());
304
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());
308 }
309 }
310 }
311
312 T out(indices.size(), x.cols());
313 out.setFromTriplets(triplets.begin(), triplets.end());
314
315 return out;
316}
317
330template<typename T>
331T
332subsetCols(const Eigen::MatrixBase<T>& x, const std::vector<int>& indices)
333{
334 return x.derived()(all, indices);
335}
336
349template<typename T>
350T
351subsetCols(const Eigen::SparseMatrixBase<T>& x, const std::vector<int>& indices)
352{
353 std::vector<Eigen::Triplet<double>> triplets;
354 triplets.reserve(slope::nonZeros(x.derived()));
355
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());
360 }
361 }
362
363 T out(x.rows(), indices.size());
364 out.setFromTriplets(triplets.begin(), triplets.end());
365
366 return out;
367}
368
380inline std::unordered_set<double>
381unique(const Eigen::MatrixXd& x)
382{
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++) {
386 unique.insert(x(i, j));
387 }
388 }
389
390 return unique;
391}
392
403template<typename Derived>
404inline bool
405isFinite(const Eigen::DenseBase<Derived>& x)
406{
407 return x.derived().array().isFinite().all();
408}
409
422template<typename Derived>
423inline bool
424isFinite(const Eigen::SparseMatrixBase<Derived>& x)
425{
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())) {
429 return false;
430 }
431 }
432 }
433 return true;
434}
435
436} // namespace slope
Eigen compatibility layer for version differences.
Namespace containing SLOPE regression implementation.
Definition clusters.h:11
std::unordered_set< double > unique(const Eigen::MatrixXd &x)
Create a set of unique values from an Eigen matrix.
Definition utils.h:381
void move_elements(std::vector< T > &v, const int from, const int to, const int size)
Definition utils.h:211
std::vector< int > which(const T &x)
Definition utils.h:95
T subsetCols(const Eigen::MatrixBase< T > &x, const std::vector< int > &indices)
Extract specified columns from a dense matrix.
Definition utils.h:332
T subset(const Eigen::EigenBase< T > &x, const std::vector< int > &indices)
Extract a subset of rows from an Eigen matrix.
Definition utils.h:258
void validateOption(const std::string &value, const std::set< std::string > &valid_options, const std::string &parameter_name)
Validates if a given value exists in a set of valid options.
std::vector< int > sortIndex(T &v, const bool descending=false)
Definition utils.h:124
void permute(T &values, const std::vector< int > &ind)
Definition utils.h:153
void sort(T &v, const bool descending=false)
Definition utils.h:70
void inversePermute(T &values, const std::vector< int > &ind)
Definition utils.h:187
Eigen::Index nonZeros(const Eigen::SparseMatrixBase< Derived > &x)
Definition utils.h:34
bool isFinite(const Eigen::DenseBase< Derived > &x)
Check if all elements in a dense matrix are finite.
Definition utils.h:405