ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
lanczos.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_LANCZOS_HPP_
2 #define VIENNACL_LINALG_LANCZOS_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2015, Institute for Microelectronics,
6  Institute for Analysis and Scientific Computing,
7  TU Wien.
8  Portions of this software are copyright by UChicago Argonne, LLC.
9 
10  -----------------
11  ViennaCL - The Vienna Computing Library
12  -----------------
13 
14  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
15 
16  (A list of authors and contributors can be found in the manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
27 #include <cmath>
28 #include <vector>
29 #include "viennacl/vector.hpp"
31 #include "viennacl/linalg/prod.hpp"
37 
38 namespace viennacl
39 {
40 namespace linalg
41 {
42 
46 {
47 public:
48 
49  enum
50  {
54  };
55 
64  lanczos_tag(double factor = 0.75,
65  vcl_size_t numeig = 10,
66  int met = 0,
67  vcl_size_t krylov = 100) : factor_(factor), num_eigenvalues_(numeig), method_(met), krylov_size_(krylov) {}
68 
70  void num_eigenvalues(vcl_size_t numeig){ num_eigenvalues_ = numeig; }
71 
73  vcl_size_t num_eigenvalues() const { return num_eigenvalues_; }
74 
76  void factor(double fct) { factor_ = fct; }
77 
79  double factor() const { return factor_; }
80 
82  void krylov_size(vcl_size_t max) { krylov_size_ = max; }
83 
85  vcl_size_t krylov_size() const { return krylov_size_; }
86 
88  void method(int met){ method_ = met; }
89 
91  int method() const { return method_; }
92 
93 
94 private:
95  double factor_;
96  vcl_size_t num_eigenvalues_;
97  int method_; // see enum defined above for possible values
98  vcl_size_t krylov_size_;
99 };
100 
101 
102 namespace detail
103 {
108  template<typename NumericT>
109  void inverse_iteration(std::vector<NumericT> const & alphas, std::vector<NumericT> const & betas,
110  NumericT & eigenvalue, std::vector<NumericT> & eigenvector)
111  {
112  std::vector<NumericT> alpha_sweeped = alphas;
113  for (vcl_size_t i=0; i<alpha_sweeped.size(); ++i)
114  alpha_sweeped[i] -= eigenvalue;
115  for (vcl_size_t row=1; row < alpha_sweeped.size(); ++row)
116  alpha_sweeped[row] -= betas[row] * betas[row] / alpha_sweeped[row-1];
117 
118  // starting guess: ignore last equation
119  eigenvector[alphas.size() - 1] = 1.0;
120 
121  for (vcl_size_t iter=0; iter<1; ++iter)
122  {
123  // solve first n-1 equations (A - \lambda I) y = -beta[n]
124  eigenvector[alphas.size() - 1] /= alpha_sweeped[alphas.size() - 1];
125  for (vcl_size_t row2=1; row2 < alphas.size(); ++row2)
126  {
127  vcl_size_t row = alphas.size() - row2 - 1;
128  eigenvector[row] -= eigenvector[row+1] * betas[row+1];
129  eigenvector[row] /= alpha_sweeped[row];
130  }
131 
132  // normalize eigenvector:
133  NumericT norm_vector = 0;
134  for (vcl_size_t i=0; i<eigenvector.size(); ++i)
135  norm_vector += eigenvector[i] * eigenvector[i];
136  norm_vector = std::sqrt(norm_vector);
137  for (vcl_size_t i=0; i<eigenvector.size(); ++i)
138  eigenvector[i] /= norm_vector;
139  }
140 
141  //eigenvalue = (alphas[0] * eigenvector[0] + betas[1] * eigenvector[1]) / eigenvector[0];
142  }
143 
156  template<typename MatrixT, typename DenseMatrixT, typename NumericT>
157  std::vector<NumericT>
158  lanczosPRO (MatrixT const& A, vector_base<NumericT> & r, DenseMatrixT & eigenvectors_A, vcl_size_t size, lanczos_tag const & tag, bool compute_eigenvectors)
159  {
160  // generation of some random numbers, used for lanczos PRO algorithm
162 
163  std::vector<vcl_size_t> l_bound(size/2), u_bound(size/2);
164  vcl_size_t n = r.size();
165  std::vector<NumericT> w(size), w_old(size); //w_k, w_{k-1}
166 
167  NumericT inner_rt;
168  std::vector<NumericT> alphas, betas;
169  viennacl::matrix<NumericT, viennacl::column_major> Q(n, size); //column-major matrix holding the Krylov basis vectors
170 
171  bool second_step = false;
172  NumericT eps = std::numeric_limits<NumericT>::epsilon();
173  NumericT squ_eps = std::sqrt(eps);
174  NumericT eta = std::exp(std::log(eps) * tag.factor());
175 
177 
178  r /= beta;
179 
180  viennacl::vector_base<NumericT> q_0(Q.handle(), Q.size1(), 0, 1);
181  q_0 = r;
182 
185  alphas.push_back(alpha);
186  w[0] = 1;
187  betas.push_back(beta);
188 
189  vcl_size_t batches = 0;
190  for (vcl_size_t i = 1; i < size; i++) // Main loop for setting up the Krylov space
191  {
192  viennacl::vector_base<NumericT> q_iminus1(Q.handle(), Q.size1(), (i-1) * Q.internal_size1(), 1);
193  r = u - alpha * q_iminus1;
194  beta = viennacl::linalg::norm_2(r);
195 
196  betas.push_back(beta);
197  r = r / beta;
198 
199  //
200  // Update recurrence relation for estimating orthogonality loss
201  //
202  w_old = w;
203  w[0] = (betas[1] * w_old[1] + (alphas[0] - alpha) * w_old[0] - betas[i - 1] * w_old[0]) / beta + eps * 0.3 * get_N() * (betas[1] + beta);
204  for (vcl_size_t j = 1; j < i - 1; j++)
205  w[j] = (betas[j + 1] * w_old[j + 1] + (alphas[j] - alpha) * w_old[j] + betas[j] * w_old[j - 1] - betas[i - 1] * w_old[j]) / beta + eps * 0.3 * get_N() * (betas[j + 1] + beta);
206  w[i-1] = 0.6 * eps * NumericT(n) * get_N() * betas[1] / beta;
207 
208  //
209  // Check whether there has been a need for reorthogonalization detected in the previous iteration.
210  // If so, run the reorthogonalization for each batch
211  //
212  if (second_step)
213  {
214  for (vcl_size_t j = 0; j < batches; j++)
215  {
216  for (vcl_size_t k = l_bound[j] + 1; k < u_bound[j] - 1; k++)
217  {
219  inner_rt = viennacl::linalg::inner_prod(r, q_k);
220  r = r - inner_rt * q_k;
221  w[k] = 1.5 * eps * get_N();
222  }
223  }
225  r = r / temp;
226  beta = beta * temp;
227  second_step = false;
228  }
229  batches = 0;
230 
231  //
232  // Check for semiorthogonality
233  //
234  for (vcl_size_t j = 0; j < i; j++)
235  {
236  if (std::fabs(w[j]) >= squ_eps) // tentative loss of orthonormality, hence reorthonomalize
237  {
239  inner_rt = viennacl::linalg::inner_prod(r, q_j);
240  r = r - inner_rt * q_j;
241  w[j] = 1.5 * eps * get_N();
242  vcl_size_t k = j - 1;
243 
244  // orthogonalization with respect to earlier basis vectors
245  while (std::fabs(w[k]) > eta)
246  {
248  inner_rt = viennacl::linalg::inner_prod(r, q_k);
249  r = r - inner_rt * q_k;
250  w[k] = 1.5 * eps * get_N();
251  if (k == 0) break;
252  k--;
253  }
254  l_bound[batches] = k;
255 
256  // orthogonalization with respect to later basis vectors
257  k = j + 1;
258  while (k < i && std::fabs(w[k]) > eta)
259  {
261  inner_rt = viennacl::linalg::inner_prod(r, q_k);
262  r = r - inner_rt * q_k;
263  w[k] = 1.5 * eps * get_N();
264  k++;
265  }
266  u_bound[batches] = k - 1;
267  batches++;
268 
269  j = k-1; // go to end of current batch
270  }
271  }
272 
273  //
274  // Normalize basis vector and reorthogonalize as needed
275  //
276  if (batches > 0)
277  {
279  r = r / temp;
280  beta = beta * temp;
281  second_step = true;
282  }
283 
284  // store Krylov vector in Q:
286  q_i = r;
287 
288  //
289  // determine and store alpha = <r, u> with u = A q_i - beta q_{i-1}
290  //
291  u = viennacl::linalg::prod(A, r);
292  u += (-beta) * q_iminus1;
293  alpha = viennacl::linalg::inner_prod(u, r);
294  alphas.push_back(alpha);
295  }
296 
297  //
298  // Step 2: Compute eigenvalues of tridiagonal matrix obtained during Lanczos iterations:
299  //
300  std::vector<NumericT> eigenvalues = bisect(alphas, betas);
301 
302  //
303  // Step 3: Compute eigenvectors via inverse iteration. Does not update eigenvalues, so only approximate by nature.
304  //
305  if (compute_eigenvectors)
306  {
307  std::vector<NumericT> eigenvector_tridiag(alphas.size());
308  for (std::size_t i=0; i < tag.num_eigenvalues(); ++i)
309  {
310  // compute eigenvector of tridiagonal matrix via inverse:
311  inverse_iteration(alphas, betas, eigenvalues[eigenvalues.size() - i - 1], eigenvector_tridiag);
312 
313  // eigenvector w of full matrix A. Given as w = Q * u, where u is the eigenvector of the tridiagonal matrix
314  viennacl::vector<NumericT> eigenvector_u(eigenvector_tridiag.size());
315  viennacl::copy(eigenvector_tridiag, eigenvector_u);
316 
317  viennacl::vector_base<NumericT> eigenvector_A(eigenvectors_A.handle(),
318  eigenvectors_A.size1(),
319  eigenvectors_A.row_major() ? i : i * eigenvectors_A.internal_size1(),
320  eigenvectors_A.row_major() ? eigenvectors_A.internal_size2() : 1);
321  eigenvector_A = viennacl::linalg::prod(project(Q,
322  range(0, Q.size1()),
323  range(0, eigenvector_u.size())),
324  eigenvector_u);
325  }
326  }
327 
328  return eigenvalues;
329  }
330 
331 
343  template< typename MatrixT, typename DenseMatrixT, typename NumericT>
344  std::vector<NumericT>
345  lanczos(MatrixT const& A, vector_base<NumericT> & r, DenseMatrixT & eigenvectors_A, vcl_size_t krylov_dim, lanczos_tag const & tag, bool compute_eigenvectors)
346  {
347  std::vector<NumericT> alphas, betas;
349  viennacl::matrix<NumericT, viennacl::column_major> Q(r.size(), krylov_dim + 1); // Krylov basis (each Krylov vector is one column)
350 
351  NumericT norm_r = norm_2(r);
352  NumericT beta = norm_r;
353  r /= norm_r;
354 
355  // first Krylov vector:
356  viennacl::vector_base<NumericT> q0(Q.handle(), Q.size1(), 0, 1);
357  q0 = r;
358 
359  //
360  // Step 1: Run Lanczos' method to obtain tridiagonal matrix
361  //
362  for (vcl_size_t i = 0; i < krylov_dim; i++)
363  {
364  betas.push_back(beta);
365  // last available vector from Krylov basis:
366  viennacl::vector_base<NumericT> q_i(Q.handle(), Q.size1(), i * Q.internal_size1(), 1);
367 
368  // Lanczos algorithm:
369  // - Compute A * q:
370  Aq = viennacl::linalg::prod(A, q_i);
371 
372  // - Form Aq <- Aq - <Aq, q_i> * q_i - beta * q_{i-1}, where beta is ||q_i|| before normalization in previous iteration
373  NumericT alpha = viennacl::linalg::inner_prod(Aq, q_i);
374  Aq -= alpha * q_i;
375 
376  if (i > 0)
377  {
378  viennacl::vector_base<NumericT> q_iminus1(Q.handle(), Q.size1(), (i-1) * Q.internal_size1(), 1);
379  Aq -= beta * q_iminus1;
380 
381  // Extra measures for improved numerical stability?
383  {
384  // Gram-Schmidt (re-)orthogonalization:
385  // TODO: Reuse fast (pipelined) routines from GMRES or GEMV
386  for (vcl_size_t j = 0; j < i; j++)
387  {
388  viennacl::vector_base<NumericT> q_j(Q.handle(), Q.size1(), j * Q.internal_size1(), 1);
389  NumericT inner_rq = viennacl::linalg::inner_prod(Aq, q_j);
390  Aq -= inner_rq * q_j;
391  }
392  }
393  }
394 
395  // normalize Aq and add to Krylov basis at column i+1 in Q:
396  beta = viennacl::linalg::norm_2(Aq);
397  viennacl::vector_base<NumericT> q_iplus1(Q.handle(), Q.size1(), (i+1) * Q.internal_size1(), 1);
398  q_iplus1 = Aq / beta;
399 
400  alphas.push_back(alpha);
401  }
402 
403  //
404  // Step 2: Compute eigenvalues of tridiagonal matrix obtained during Lanczos iterations:
405  //
406  std::vector<NumericT> eigenvalues = bisect(alphas, betas);
407 
408  //
409  // Step 3: Compute eigenvectors via inverse iteration. Does not update eigenvalues, so only approximate by nature.
410  //
411  if (compute_eigenvectors)
412  {
413  std::vector<NumericT> eigenvector_tridiag(alphas.size());
414  for (std::size_t i=0; i < tag.num_eigenvalues(); ++i)
415  {
416  // compute eigenvector of tridiagonal matrix via inverse:
417  inverse_iteration(alphas, betas, eigenvalues[eigenvalues.size() - i - 1], eigenvector_tridiag);
418 
419  // eigenvector w of full matrix A. Given as w = Q * u, where u is the eigenvector of the tridiagonal matrix
420  viennacl::vector<NumericT> eigenvector_u(eigenvector_tridiag.size());
421  viennacl::copy(eigenvector_tridiag, eigenvector_u);
422 
423  viennacl::vector_base<NumericT> eigenvector_A(eigenvectors_A.handle(),
424  eigenvectors_A.size1(),
425  eigenvectors_A.row_major() ? i : i * eigenvectors_A.internal_size1(),
426  eigenvectors_A.row_major() ? eigenvectors_A.internal_size2() : 1);
427  eigenvector_A = viennacl::linalg::prod(project(Q,
428  range(0, Q.size1()),
429  range(0, eigenvector_u.size())),
430  eigenvector_u);
431  }
432  }
433 
434  return eigenvalues;
435  }
436 
437 } // end namespace detail
438 
450 template<typename MatrixT, typename DenseMatrixT>
451 std::vector< typename viennacl::result_of::cpu_value_type<typename MatrixT::value_type>::type >
452 eig(MatrixT const & matrix, DenseMatrixT & eigenvectors_A, lanczos_tag const & tag, bool compute_eigenvectors = true)
453 {
454  typedef typename viennacl::result_of::value_type<MatrixT>::type NumericType;
455  typedef typename viennacl::result_of::cpu_value_type<NumericType>::type CPU_NumericType;
456  typedef typename viennacl::result_of::vector_for_matrix<MatrixT>::type VectorT;
457 
459 
460  std::vector<CPU_NumericType> eigenvalues;
461  vcl_size_t matrix_size = matrix.size1();
462  VectorT r(matrix_size);
463  std::vector<CPU_NumericType> s(matrix_size);
464 
465  for (vcl_size_t i=0; i<s.size(); ++i)
466  s[i] = CPU_NumericType(0.5) + random_gen();
467 
469 
470  vcl_size_t size_krylov = (matrix_size < tag.krylov_size()) ? matrix_size
471  : tag.krylov_size();
472 
473  switch (tag.method())
474  {
476  eigenvalues = detail::lanczosPRO(matrix, r, eigenvectors_A, size_krylov, tag, compute_eigenvectors);
477  break;
480  eigenvalues = detail::lanczos(matrix, r, eigenvectors_A, size_krylov, tag, compute_eigenvectors);
481  break;
482  }
483 
484  std::vector<CPU_NumericType> largest_eigenvalues;
485 
486  for (vcl_size_t i = 1; i<=tag.num_eigenvalues(); i++)
487  largest_eigenvalues.push_back(eigenvalues[size_krylov-i]);
488 
489 
490  return largest_eigenvalues;
491 }
492 
493 
503 template<typename MatrixT>
504 std::vector< typename viennacl::result_of::cpu_value_type<typename MatrixT::value_type>::type >
505 eig(MatrixT const & matrix, lanczos_tag const & tag)
506 {
508 
509  viennacl::matrix<NumericType> eigenvectors(matrix.size1(), tag.num_eigenvalues());
510  return eig(matrix, eigenvectors, tag, false);
511 }
512 
513 } // end namespace linalg
514 } // end namespace viennacl
515 #endif
int method() const
Returns the reorthogonalization method.
Definition: lanczos.hpp:91
T norm_2(std::vector< T, A > const &v1)
Definition: norm_2.hpp:96
A reader and writer for the matrix market format is implemented here.
vcl_size_t krylov_size() const
Returns the size of the kylov space.
Definition: lanczos.hpp:85
std::vector< typename viennacl::result_of::cpu_value_type< typename MatrixT::value_type >::type > eig(MatrixT const &matrix, DenseMatrixT &eigenvectors_A, lanczos_tag const &tag, bool compute_eigenvectors=true)
Implementation of the calculation of eigenvalues using lanczos (with and without reorthogonalization)...
Definition: lanczos.hpp:452
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
void method(int met)
Sets the reorthogonalization method.
Definition: lanczos.hpp:88
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
std::vector< NumericT > lanczosPRO(MatrixT const &A, vector_base< NumericT > &r, DenseMatrixT &eigenvectors_A, vcl_size_t size, lanczos_tag const &tag, bool compute_eigenvectors)
Implementation of the Lanczos PRO algorithm (partial reorthogonalization)
Definition: lanczos.hpp:158
Random number generator for returning normally distributed values.
Definition: random.hpp:57
lanczos_tag(double factor=0.75, vcl_size_t numeig=10, int met=0, vcl_size_t krylov=100)
The constructor.
Definition: lanczos.hpp:64
A dense matrix class.
Definition: forwards.h:375
viennacl::enable_if< viennacl::is_stl< typename viennacl::traits::tag_of< VectorT1 >::type >::value, typename VectorT1::value_type >::type inner_prod(VectorT1 const &v1, VectorT2 const &v2)
Definition: inner_prod.hpp:100
Generic interface for the computation of inner products. See viennacl/linalg/vector_operations.hpp for implementations.
float NumericT
Definition: bisect.cpp:40
std::vector< typename viennacl::result_of::cpu_value_type< typename VectorT::value_type >::type > bisect(VectorT const &alphas, VectorT const &betas)
Implementation of the bisect-algorithm for the calculation of the eigenvalues of a tridiagonal matrix...
Definition: bisect.hpp:78
basic_range range
Definition: forwards.h:424
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
Definition: cpu_ram.hpp:34
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
Definition: size.hpp:235
Random number generator for returning uniformly distributed values in the closed interval [0...
Definition: random.hpp:44
Definition: blas3.hpp:36
Implementation of the compressed_matrix class.
void copy_vec_to_vec(viennacl::vector< NumericT > const &src, OtherVectorT &dest)
overloaded function for copying vectors
Definition: bisect.hpp:44
std::vector< NumericT > lanczos(MatrixT const &A, vector_base< NumericT > &r, DenseMatrixT &eigenvectors_A, vcl_size_t krylov_dim, lanczos_tag const &tag, bool compute_eigenvectors)
Implementation of the Lanczos FRO algorithm.
Definition: lanczos.hpp:345
void num_eigenvalues(vcl_size_t numeig)
Sets the number of eigenvalues.
Definition: lanczos.hpp:70
matrix_range< MatrixType > project(MatrixType const &A, viennacl::range const &r1, viennacl::range const &r2)
vcl_size_t num_eigenvalues() const
Returns the number of eigenvalues.
Definition: lanczos.hpp:73
std::size_t vcl_size_t
Definition: forwards.h:75
handle_type & handle()
Returns the OpenCL handle, non-const-version.
Definition: matrix_def.hpp:244
T::ERROR_CANNOT_DEDUCE_CPU_SCALAR_TYPE_FOR_T type
Definition: result_of.hpp:271
size_type size1() const
Returns the number of rows.
Definition: matrix_def.hpp:224
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_row > row(const matrix_base< NumericT, F > &A, unsigned int i)
Definition: matrix.hpp:900
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
NumericT max(std::vector< NumericT > const &v1)
Definition: maxmin.hpp:47
void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)
Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...
A small collection of sequential random number generators.
void krylov_size(vcl_size_t max)
Sets the size of the kylov space. Must be larger than number of eigenvalues to compute.
Definition: lanczos.hpp:82
size_type size() const
Returns the length of the vector (cf. std::vector)
Definition: vector_def.hpp:118
size_type internal_size1() const
Returns the internal number of rows. Usually required for launching OpenCL kernels only...
Definition: matrix_def.hpp:238
Implementation of the algorithm for finding eigenvalues of a tridiagonal matrix.
double factor() const
Returns the exponent.
Definition: lanczos.hpp:79
void inverse_iteration(std::vector< NumericT > const &alphas, std::vector< NumericT > const &betas, NumericT &eigenvalue, std::vector< NumericT > &eigenvector)
Inverse iteration for finding an eigenvector for an eigenvalue.
Definition: lanczos.hpp:109
void factor(double fct)
Sets the exponent of epsilon. Values between 0.6 and 0.9 usually give best results.
Definition: lanczos.hpp:76
A tag for the lanczos algorithm.
Definition: lanczos.hpp:45