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& 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 // Falling back to an exact SVD for small matrices or if the requested number is too large
163 // (not enough of a workspace). Hey, I don't make the rules.
164 if (
165 (options.exact_for_small_matrix && smaller < 6) ||
166 (options.exact_for_large_number && requested_greater_than_or_equal_to_half_smaller(requested_number, smaller))
167 ) {
168 exact(matrix, requested_number, outU, outV, outD);
169 return std::make_pair(true, 0);
170 }
171
172 Eigen::Index work = choose_requested_plus_extra_work_or_smaller(requested_number, options.extra_work, smaller);
173 if (work == 0) {
174 throw std::runtime_error("number of requested dimensions must be positive");
175 }
176
177 // Don't worry about sanitizing dimensions for Eigen constructors,
178 // as the former are stored as Eigen::Index and the latter accepts Eigen::Index inputs.
179 EigenMatrix_ V(matrix.cols(), work);
180 std::mt19937_64 eng(options.seed);
181 if (options.initial) {
182 auto init = reinterpret_cast<EigenVector_*>(options.initial);
183 if (init->size() != V.rows()) {
184 throw std::runtime_error("initialization vector does not have expected number of rows");
185 }
186 V.col(0) = *init;
187 } else {
188 fill_with_random_normals(V, 0, eng);
189 }
190 V.col(0) /= V.col(0).norm();
191
192 bool converged = false;
193 int iter = 0;
194 Eigen::Index k = 0;
195 JacobiSVD<EigenMatrix_> svd(work, work);
196
197 LanczosWorkspace<EigenVector_, Matrix_> lpwork(matrix);
198
199 EigenMatrix_ W(matrix.rows(), work);
200 EigenMatrix_ Wtmp(matrix.rows(), work);
201 EigenMatrix_ Vtmp(matrix.cols(), work);
202
203 EigenMatrix_ B(work, work);
204 B.setZero(work, work);
205 EigenVector_ res(work);
206 EigenVector_ F(matrix.cols());
207
208 EigenVector_ prevS(work);
209 typename EigenMatrix_::Scalar svtol = options.singular_value_ratio_tolerance;
210 typename EigenMatrix_::Scalar tol = options.convergence_tolerance;
211 typename EigenMatrix_::Scalar svtol_actual = (svtol >= 0 ? svtol : tol);
212
213 for (; iter < options.max_iterations; ++iter) {
214 // Technically, this is only a 'true' Lanczos bidiagonalization
215 // when k = 0. All other times, we're just recycling the machinery,
216 // see the text below Equation 3.11 in Baglama and Reichel.
217 run_lanczos_bidiagonalization(lpwork, W, V, B, eng, k, options);
218
219// if (iter < 2) {
220// std::cout << "B is currently:\n" << B << std::endl;
221// std::cout << "W is currently:\n" << W << std::endl;
222// std::cout << "V is currently:\n" << V << std::endl;
223// }
224
225 svd.compute(B);
226 const auto& BS = svd.singularValues();
227 const auto& BU = svd.matrixU();
228 const auto& BV = svd.matrixV();
229
230 // Checking for convergence.
231 if (B(work - 1, work - 1) == 0) { // a.k.a. the final value of 'S' from the Lanczos iterations.
232 converged = true;
233 break;
234 }
235
236 const auto& F = lpwork.F; // i.e., the residuals, see Algorithm 2.1 of Baglama and Reichel.
237 auto R_F = F.norm();
238
239 // Computes the convergence criterion defined in on the LHS of
240 // Equation 2.13 of Baglama and Riechel. (e_m is literally the unit
241 // vector along the m-th dimension, so it ends up just being the
242 // m-th row of the U matrix.) We expose it here as we will be
243 // re-using the same values to update B, see below.
244 res = R_F * BU.row(work - 1);
245
246 Eigen::Index n_converged = 0;
247 if (iter > 0) {
248 const auto Smax = *std::max_element(BS.begin(), BS.end());
249 const auto threshold = Smax * tol;
250
251 for (Eigen::Index j = 0; j < work; ++j) {
252 if (std::abs(res[j]) <= threshold) {
253 const auto ratio = std::abs(BS[j] - prevS[j]) / BS[j];
254 if (ratio <= svtol_actual) {
255 ++n_converged;
256 }
257 }
258 }
259
260 if (n_converged >= requested_number) {
261 converged = true;
262 break;
263 }
264 }
265 prevS = BS;
266
267 k = update_k(k, requested_number, n_converged, work);
268 Vtmp.leftCols(k).noalias() = V * BV.leftCols(k);
269 V.leftCols(k) = Vtmp.leftCols(k);
270
271 // See Equation 3.2 of Baglama and Reichel, where our 'V' is their
272 // 'P', and our 'F / R_F' is their 'p_{m+1}' (Equation 2.2). 'F'
273 // was orthogonal to the old 'V' and so it should still be
274 // orthogonal to the new left-most columns of 'V'; the input
275 // expectations of the Lanczos bidiagonalization are still met, so
276 // V is ok to use in the next run_lanczos_bidiagonalization().
277 V.col(k) = F / R_F;
278
279 Wtmp.leftCols(k).noalias() = W * BU.leftCols(k);
280 W.leftCols(k) = Wtmp.leftCols(k);
281
282 B.setZero(work, work);
283 for (I<decltype(k)> l = 0; l < k; ++l) {
284 B(l, l) = BS[l];
285
286 // This assignment looks weird but is deliberate, see Equation
287 // 3.6 of Baglama and Reichel. Happily, this is the same value
288 // used to determine convergence in 2.13, so we just re-use it.
289 // (See the equation just above Equation 3.5; I think they
290 // misplaced a tilde on the final 'u', given no other 'u' has
291 // a B_m superscript as well as a tilde.)
292 B(l, k) = res[l];
293 }
294 }
295
296 // See Equation 2.11 of Baglama and Reichel for how to get from B's
297 // singular triplets to mat's singular triplets.
298 outD.resize(requested_number);
299 outD = svd.singularValues().head(requested_number);
300
301 outU.resize(matrix.rows(), requested_number);
302 outU.noalias() = W * svd.matrixU().leftCols(requested_number);
303
304 outV.resize(matrix.cols(), requested_number);
305 outV.noalias() = V * svd.matrixV().leftCols(requested_number);
306
307 return std::make_pair(converged, (converged ? iter + 1 : iter));
308}
309
333template<class InputEigenMatrix_, class OutputEigenMatrix_, class EigenVector_>
334std::pair<bool, int> compute_simple(
335 const InputEigenMatrix_& matrix,
336 Eigen::Index number,
337 OutputEigenMatrix_& outU,
338 OutputEigenMatrix_& outV,
339 EigenVector_& outD,
340 const Options& options
341) {
343
344 // Force the compiler to use virtual dispatch to avoid realizing multiple
345 // instances of this function for each InputEigenMatrix_ type. If you
346 // don't like this, call compute() directly.
348
349 return compute<Interface>(wrapped, number, outU, outV, outD, options);
350}
351
357template<class EigenMatrix_, class EigenVector_>
358struct Results {
364 EigenMatrix_ U;
365
371 EigenMatrix_ V;
372
377 EigenVector_ D;
378
383
388};
389
403template<class EigenMatrix_ = Eigen::MatrixXd, class EigenVector_ = Eigen::VectorXd, class Matrix_>
404Results<EigenMatrix_, EigenVector_> compute(const Matrix_& matrix, Eigen::Index number, const Options& options) {
406 const auto stats = compute(matrix, number, output.U, output.V, output.D, options);
407 output.converged = stats.first;
408 output.iterations = stats.second;
409 return output;
410}
411
426template<class OutputEigenMatrix_ = Eigen::MatrixXd, class EigenVector_ = Eigen::VectorXd, class InputEigenMatrix_>
427Results<OutputEigenMatrix_, EigenVector_> compute_simple(const InputEigenMatrix_& matrix, Eigen::Index number, const Options& options) {
429 const auto stats = compute_simple(matrix, number, output.U, output.V, output.D, options);
430 output.converged = stats.first;
431 output.iterations = stats.second;
432 return output;
433}
434
435}
436
437#endif
Interface for a matrix to use in compute().
Definition interface.hpp:142
A Matrix-compatible wrapper around a "simple" matrix.
Definition simple.hpp:125
Implements IRLBA for approximate SVD.
Definition compute.hpp:22
std::pair< bool, int > compute(const Matrix_ &matrix, const Eigen::Index number, EigenMatrix_ &outU, EigenMatrix_ &outV, EigenVector_ &outD, const Options &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 &options)
Definition compute.hpp:334
Simple wrapper around an Eigen matrix.
Options for the IRLBA algorithm.
Definition Options.hpp:16
double singular_value_ratio_tolerance
Definition Options.hpp:37
double convergence_tolerance
Definition Options.hpp:30
bool cap_number
Definition Options.hpp:67
bool exact_for_small_matrix
Definition Options.hpp:55
void * initial
Definition Options.hpp:78
int extra_work
Definition Options.hpp:43
bool exact_for_large_number
Definition Options.hpp:61
std::mt19937_64::result_type seed
Definition Options.hpp:72
int max_iterations
Definition Options.hpp:49
Results of the IRLBA-based SVD by compute().
Definition compute.hpp:358
EigenVector_ D
Definition compute.hpp:377
int iterations
Definition compute.hpp:382
EigenMatrix_ U
Definition compute.hpp:364
EigenMatrix_ V
Definition compute.hpp:371
bool converged
Definition compute.hpp:387