irlba
A C++ library for IRLBA
Loading...
Searching...
No Matches
compute.hpp
Go to the documentation of this file.
1#ifndef IRLBA_COMPUTE_HPP
2#define IRLBA_COMPUTE_HPP
3
4#include <cmath>
5#include <cstdint>
6#include <random>
7#include <stdexcept>
8#include <type_traits>
9
10#include "utils.hpp"
11#include "lanczos.hpp"
12#include "Matrix/simple.hpp"
13
14#include "sanisizer/sanisizer.hpp"
15#include "Eigen/Dense"
16
22namespace irlba {
23
27template<typename EigenMatrix_>
28using JacobiSVD = Eigen::JacobiSVD<EigenMatrix_, Eigen::ComputeThinU | Eigen::ComputeThinV>;
29
30template<class Matrix_, class EigenMatrix_, class EigenVector_>
31void exact(const Matrix_& matrix, const Eigen::Index requested_number, EigenMatrix_& outU, EigenMatrix_& outV, EigenVector_& outD) {
32 JacobiSVD<EigenMatrix_> svd(matrix.rows(), matrix.cols());
33
34 auto realizer = matrix.new_known_realize_workspace();
35 EigenMatrix_ buffer;
36 svd.compute(realizer->realize(buffer));
37
38 outD.resize(requested_number);
39 outD = svd.singularValues().head(requested_number);
40
41 outU.resize(matrix.rows(), requested_number);
42 outU = svd.matrixU().leftCols(requested_number);
43
44 outV.resize(matrix.cols(), requested_number);
45 outV = svd.matrixV().leftCols(requested_number);
46}
47
48// Basically (requested * 2 >= smaller), but avoiding overflow from the product.
49inline bool requested_greater_than_or_equal_to_half_smaller(const Eigen::Index requested, const Eigen::Index smaller) {
50 const Eigen::Index half_smaller = smaller / 2;
51 if (requested == half_smaller) {
52 return smaller % 2 == 0;
53 } else {
54 return requested > half_smaller;
55 }
56}
57
58// Basically min(requested_number + extra_work, smaller), but avoiding overflow from the sum.
59inline Eigen::Index choose_requested_plus_extra_work_or_smaller(const Eigen::Index requested_number, const int extra_work, const Eigen::Index smaller) {
60 if (requested_number >= smaller) {
61 return smaller;
62 } else {
63 // This is guaranteed to fit into an Eigen::Index;
64 // either it's equal to 'smaller' or it's less than 'smaller'.
65 return requested_number + sanisizer::min(smaller - requested_number, extra_work);
66 }
67}
68
69// Setting 'k'. This looks kinda weird, but this is deliberate, see the text below Algorithm 3.1 of Baglama and Reichel.
70//
71// - Our n_converged is their k'.
72// - Our requested_number is their k.
73// - Our work is their m.
74// - Our returned k is their k + k''.
75//
76// They don't mention anything about the value corresponding to our input k, but I guess we always want k to increase.
77// So, if our proposed new value of k is lower than our current value of k, we just pick the latter.
78inline Eigen::Index update_k(Eigen::Index k, const Eigen::Index requested_number, const Eigen::Index n_converged, const Eigen::Index work) {
79 // Here, the goal is to reproduce this code from irlba's convtests() function, but without any risk of wraparound or overflow.
80 //
81 // if (k < requested_number + n_converged) {
82 // k = requested_number + n_converged;
83 // }
84 // if (k > work - 3) {
85 // k = work - 3;
86 // }
87 // if (k < 1) {
88 // k = 1;
89 // }
90
91 if (work <= 3) {
92 return 1;
93 }
94 const Eigen::Index limit = work - 3;
95
96 const auto less_than_requested_plus_converged = [&](const Eigen::Index val) -> bool {
97 return val < requested_number || static_cast<Eigen::Index>(val - requested_number) < n_converged;
98 };
99
100 if (less_than_requested_plus_converged(k)) {
101 if (less_than_requested_plus_converged(limit)) {
102 return limit;
103 } else {
104 const Eigen::Index output = requested_number + n_converged;
105 return std::max(output, static_cast<Eigen::Index>(1));
106 }
107 } else {
108 const Eigen::Index output = std::min(k, limit);
109 return std::max(output, static_cast<Eigen::Index>(1));
110 }
111}
141template<class Matrix_, class EigenMatrix_, class EigenVector_>
142std::pair<bool, int> compute(
143 const Matrix_& matrix,
144 const Eigen::Index number,
145 EigenMatrix_& outU,
146 EigenMatrix_& outV,
147 EigenVector_& outD,
148 const Options<EigenVector_>& options
149) {
150 const Eigen::Index smaller = std::min(matrix.rows(), matrix.cols());
151 Eigen::Index requested_number = sanisizer::cast<Eigen::Index>(number);
152 if (requested_number > smaller) {
153 if (options.cap_number) {
154 requested_number = smaller;
155 } else {
156 throw std::runtime_error("requested number of singular values cannot be greater than the smaller matrix dimension");
157 }
158 } else if (requested_number == smaller && !options.exact_for_large_number) {
159 throw std::runtime_error("requested number of singular values must be less than the smaller matrix dimension for IRLBA iterations");
160 }
161
162 // Optimization if no SVs are requested. Also protect the exact SVD code
163 // from a zero-extent matrix, as this causes Eigen to shit itself.
164 if (requested_number == 0) {
165 outD.resize(requested_number);
166 outU.resize(matrix.rows(), requested_number);
167 outV.resize(matrix.cols(), requested_number);
168 return std::make_pair(true, 0);
169 }
170
171 // Falling back to an exact SVD for small matrices or if the requested number is too large
172 // (not enough of a workspace). Hey, I don't make the rules.
173 if (
174 (options.exact_for_small_matrix && smaller < 6) ||
175 (options.exact_for_large_number && requested_greater_than_or_equal_to_half_smaller(requested_number, smaller))
176 ) {
177 exact(matrix, requested_number, outU, outV, outD);
178 return std::make_pair(true, 0);
179 }
180
181 // We know work must be positive at this point, as we would have returned early if requested_number = 0.
182 // Similar, if smaller = 0, we would have either thrown earlier or set requested_number = 0.
183 const Eigen::Index work = choose_requested_plus_extra_work_or_smaller(requested_number, options.extra_work, smaller);
184
185 // Don't worry about sanitizing dimensions for Eigen constructors,
186 // as the former are stored as Eigen::Index and the latter accepts Eigen::Index inputs.
187 EigenMatrix_ V(matrix.cols(), work);
188 std::mt19937_64 eng(options.seed);
189 if (options.initial.has_value()) {
190 const auto& init = *(options.initial);
191 if (init.size() != V.rows()) {
192 throw std::runtime_error("initialization vector does not have expected number of rows");
193 }
194 V.col(0) = init;
195 } else {
196 fill_with_random_normals(V, 0, eng);
197 }
198 V.col(0) /= V.col(0).norm();
199
200 bool converged = false;
201 int iter = 0;
202 Eigen::Index k = 0;
203 JacobiSVD<EigenMatrix_> svd(work, work);
204
205 LanczosWorkspace<EigenVector_, Matrix_> lpwork(matrix);
206
207 EigenMatrix_ W(matrix.rows(), work);
208 EigenMatrix_ Wtmp(matrix.rows(), work);
209 EigenMatrix_ Vtmp(matrix.cols(), work);
210
211 EigenMatrix_ B(work, work);
212 B.setZero(work, work);
213 EigenVector_ res(work);
214 EigenVector_ F(matrix.cols());
215
216 EigenVector_ prevS(work);
217 typename EigenMatrix_::Scalar svtol = options.singular_value_ratio_tolerance;
218 typename EigenMatrix_::Scalar tol = options.convergence_tolerance;
219 typename EigenMatrix_::Scalar svtol_actual = (svtol >= 0 ? svtol : tol);
220
221 for (; iter < options.max_iterations; ++iter) {
222 // Technically, this is only a 'true' Lanczos bidiagonalization
223 // when k = 0. All other times, we're just recycling the machinery,
224 // see the text below Equation 3.11 in Baglama and Reichel.
225 run_lanczos_bidiagonalization(lpwork, W, V, B, eng, k, options);
226
227// if (iter < 2) {
228// std::cout << "B is currently:\n" << B << std::endl;
229// std::cout << "W is currently:\n" << W << std::endl;
230// std::cout << "V is currently:\n" << V << std::endl;
231// }
232
233 svd.compute(B);
234 const auto& BS = svd.singularValues();
235 const auto& BU = svd.matrixU();
236 const auto& BV = svd.matrixV();
237
238 // Checking for convergence.
239 if (B(work - 1, work - 1) == 0) { // a.k.a. the final value of 'S' from the Lanczos iterations.
240 converged = true;
241 break;
242 }
243
244 const auto& F = lpwork.F; // i.e., the residuals, see Algorithm 2.1 of Baglama and Reichel.
245 auto R_F = F.norm();
246
247 // Computes the convergence criterion defined in on the LHS of
248 // Equation 2.13 of Baglama and Riechel. (e_m is literally the unit
249 // vector along the m-th dimension, so it ends up just being the
250 // m-th row of the U matrix.) We expose it here as we will be
251 // re-using the same values to update B, see below.
252 res = R_F * BU.row(work - 1);
253
254 Eigen::Index n_converged = 0;
255 if (iter > 0) {
256 const auto Smax = *std::max_element(BS.begin(), BS.end());
257 const auto threshold = Smax * tol;
258
259 for (Eigen::Index j = 0; j < work; ++j) {
260 if (std::abs(res[j]) <= threshold) {
261 const auto ratio = std::abs(BS[j] - prevS[j]) / BS[j];
262 if (ratio <= svtol_actual) {
263 ++n_converged;
264 }
265 }
266 }
267
268 if (n_converged >= requested_number) {
269 converged = true;
270 break;
271 }
272 }
273 prevS = BS;
274
275 k = update_k(k, requested_number, n_converged, work);
276 Vtmp.leftCols(k).noalias() = V * BV.leftCols(k);
277 V.leftCols(k) = Vtmp.leftCols(k);
278
279 // See Equation 3.2 of Baglama and Reichel, where our 'V' is their
280 // 'P', and our 'F / R_F' is their 'p_{m+1}' (Equation 2.2). 'F'
281 // was orthogonal to the old 'V' and so it should still be
282 // orthogonal to the new left-most columns of 'V'; the input
283 // expectations of the Lanczos bidiagonalization are still met, so
284 // V is ok to use in the next run_lanczos_bidiagonalization().
285 V.col(k) = F / R_F;
286
287 Wtmp.leftCols(k).noalias() = W * BU.leftCols(k);
288 W.leftCols(k) = Wtmp.leftCols(k);
289
290 B.setZero(work, work);
291 for (I<decltype(k)> l = 0; l < k; ++l) {
292 B(l, l) = BS[l];
293
294 // This assignment looks weird but is deliberate, see Equation
295 // 3.6 of Baglama and Reichel. Happily, this is the same value
296 // used to determine convergence in 2.13, so we just re-use it.
297 // (See the equation just above Equation 3.5; I think they
298 // misplaced a tilde on the final 'u', given no other 'u' has
299 // a B_m superscript as well as a tilde.)
300 B(l, k) = res[l];
301 }
302 }
303
304 // See Equation 2.11 of Baglama and Reichel for how to get from B's
305 // singular triplets to mat's singular triplets.
306 outD.resize(requested_number);
307 outD = svd.singularValues().head(requested_number);
308
309 outU.resize(matrix.rows(), requested_number);
310 outU.noalias() = W * svd.matrixU().leftCols(requested_number);
311
312 outV.resize(matrix.cols(), requested_number);
313 outV.noalias() = V * svd.matrixV().leftCols(requested_number);
314
315 return std::make_pair(converged, (converged ? iter + 1 : iter));
316}
317
341template<class InputEigenMatrix_, class OutputEigenMatrix_, class EigenVector_>
342std::pair<bool, int> compute_simple(
343 const InputEigenMatrix_& matrix,
344 Eigen::Index number,
345 OutputEigenMatrix_& outU,
346 OutputEigenMatrix_& outV,
347 EigenVector_& outD,
348 const Options<EigenVector_>& options
349) {
351
352 // Force the compiler to use virtual dispatch to avoid realizing multiple
353 // instances of this function for each InputEigenMatrix_ type. If you
354 // don't like this, call compute() directly.
356
357 return compute<Interface>(wrapped, number, outU, outV, outD, options);
358}
359
365template<class EigenMatrix_, class EigenVector_>
366struct Results {
372 EigenMatrix_ U;
373
379 EigenMatrix_ V;
380
385 EigenVector_ D;
386
391
396};
397
411template<class EigenMatrix_ = Eigen::MatrixXd, class EigenVector_ = Eigen::VectorXd, class Matrix_>
412Results<EigenMatrix_, EigenVector_> compute(const Matrix_& matrix, Eigen::Index number, const Options<EigenVector_>& options) {
414 const auto stats = compute(matrix, number, output.U, output.V, output.D, options);
415 output.converged = stats.first;
416 output.iterations = stats.second;
417 return output;
418}
419
434template<class OutputEigenMatrix_ = Eigen::MatrixXd, class EigenVector_ = Eigen::VectorXd, class InputEigenMatrix_>
435Results<OutputEigenMatrix_, EigenVector_> compute_simple(const InputEigenMatrix_& matrix, Eigen::Index number, const Options<EigenVector_>& options) {
437 const auto stats = compute_simple(matrix, number, output.U, output.V, output.D, options);
438 output.converged = stats.first;
439 output.iterations = stats.second;
440 return output;
441}
442
443}
444
445#endif
Interface for a matrix to use in compute().
Definition interface.hpp:142
A Matrix-compatible wrapper around a "simple" matrix.
Definition simple.hpp:89
std::pair< bool, int > compute(const Matrix_ &matrix, const Eigen::Index number, EigenMatrix_ &outU, EigenMatrix_ &outV, EigenVector_ &outD, const Options< EigenVector_ > &options)
Definition compute.hpp:142
std::pair< bool, int > compute_simple(const InputEigenMatrix_ &matrix, Eigen::Index number, OutputEigenMatrix_ &outU, OutputEigenMatrix_ &outV, EigenVector_ &outD, const Options< EigenVector_ > &options)
Definition compute.hpp:342
Simple wrapper around an Eigen matrix.
Options for running IRLBA in compute() and pca().
Definition Options.hpp:20
double convergence_tolerance
Definition Options.hpp:34
int extra_work
Definition Options.hpp:47
bool exact_for_large_number
Definition Options.hpp:65
int max_iterations
Definition Options.hpp:53
std::optional< EigenVector_ > initial
Definition Options.hpp:82
bool exact_for_small_matrix
Definition Options.hpp:59
bool cap_number
Definition Options.hpp:71
double singular_value_ratio_tolerance
Definition Options.hpp:41
std::mt19937_64::result_type seed
Definition Options.hpp:76
Results of the IRLBA-based SVD by compute().
Definition compute.hpp:366
EigenVector_ D
Definition compute.hpp:385
int iterations
Definition compute.hpp:390
EigenMatrix_ U
Definition compute.hpp:372
EigenMatrix_ V
Definition compute.hpp:379
bool converged
Definition compute.hpp:395