1 #ifndef VIENNACL_LINALG_LANCZOS_HPP_ 
    2 #define VIENNACL_LINALG_LANCZOS_HPP_ 
   67               vcl_size_t krylov = 100) : factor_(
factor), num_eigenvalues_(numeig), method_(met), krylov_size_(krylov) {}
 
   76   void factor(
double fct) { factor_ = fct; }
 
   79   double factor()
 const { 
return factor_; }
 
   88   void method(
int met){ method_ = met; }
 
   91   int method()
 const { 
return method_; }
 
  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)
 
  112     std::vector<NumericT> alpha_sweeped = alphas;
 
  113     for (
vcl_size_t i=0; i<alpha_sweeped.size(); ++i)
 
  114       alpha_sweeped[i] -= eigenvalue;
 
  116       alpha_sweeped[
row] -= betas[
row] * betas[
row] / alpha_sweeped[
row-1];
 
  119     eigenvector[alphas.size() - 1] = 1.0;
 
  124       eigenvector[alphas.size() - 1] /= alpha_sweeped[alphas.size() - 1];
 
  125       for (
vcl_size_t row2=1; row2 < alphas.size(); ++row2)
 
  128         eigenvector[
row] -= eigenvector[row+1] * betas[row+1];
 
  129         eigenvector[
row] /= alpha_sweeped[
row];
 
  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;
 
  156   template<
typename MatrixT, 
typename DenseMatrixT, 
typename NumericT>
 
  157   std::vector<NumericT>
 
  163     std::vector<vcl_size_t> l_bound(size/2), u_bound(size/2);
 
  165     std::vector<NumericT> w(size), w_old(size); 
 
  168     std::vector<NumericT> alphas, betas;
 
  171     bool second_step = 
false;
 
  172     NumericT eps = std::numeric_limits<NumericT>::epsilon();
 
  185     alphas.push_back(alpha);
 
  187     betas.push_back(beta);
 
  193       r = u - alpha * q_iminus1;
 
  196       betas.push_back(beta);
 
  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);
 
  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;
 
  216           for (
vcl_size_t k = l_bound[j] + 1; k < u_bound[j] - 1; k++)
 
  220             r = r - inner_rt * q_k;
 
  221             w[k] = 1.5 * eps * get_N();
 
  236         if (std::fabs(w[j]) >= squ_eps) 
 
  240           r = r - inner_rt * q_j;
 
  241           w[j] = 1.5 * eps * get_N();
 
  245           while (std::fabs(w[k]) > eta)
 
  249             r = r - inner_rt * q_k;
 
  250             w[k] = 1.5 * eps * get_N();
 
  254           l_bound[batches] = k;
 
  258           while (k < i && std::fabs(w[k]) > eta)
 
  262             r = r - inner_rt * q_k;
 
  263             w[k] = 1.5 * eps * get_N();
 
  266           u_bound[batches] = k - 1;
 
  292       u += (-beta) * q_iminus1;
 
  294       alphas.push_back(alpha);
 
  300     std::vector<NumericT> eigenvalues = 
bisect(alphas, betas);
 
  305     if (compute_eigenvectors)
 
  307       std::vector<NumericT> eigenvector_tridiag(alphas.size());
 
  311         inverse_iteration(alphas, betas, eigenvalues[eigenvalues.size() - i - 1], eigenvector_tridiag);
 
  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);
 
  323                                                        range(0, eigenvector_u.size())),
 
  343   template< 
typename MatrixT, 
typename DenseMatrixT, 
typename NumericT>
 
  344   std::vector<NumericT>
 
  347     std::vector<NumericT> alphas, betas;
 
  364       betas.push_back(beta);
 
  379         Aq -= beta * q_iminus1;
 
  390             Aq -= inner_rq * q_j;
 
  398       q_iplus1 = Aq / beta;
 
  400       alphas.push_back(alpha);
 
  406     std::vector<NumericT> eigenvalues = 
bisect(alphas, betas);
 
  411     if (compute_eigenvectors)
 
  413       std::vector<NumericT> eigenvector_tridiag(alphas.size());
 
  417         inverse_iteration(alphas, betas, eigenvalues[eigenvalues.size() - i - 1], eigenvector_tridiag);
 
  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);
 
  429                                                        range(0, eigenvector_u.size())),
 
  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)
 
  456   typedef typename viennacl::result_of::vector_for_matrix<MatrixT>::type    VectorT;
 
  460   std::vector<CPU_NumericType> eigenvalues;
 
  462   VectorT r(matrix_size);
 
  463   std::vector<CPU_NumericType> s(matrix_size);
 
  466     s[i] = CPU_NumericType(0.5) + random_gen();
 
  476     eigenvalues = 
detail::lanczosPRO(matrix, r, eigenvectors_A, size_krylov, tag, compute_eigenvectors);
 
  480     eigenvalues = 
detail::lanczos(matrix, r, eigenvectors_A, size_krylov, tag, compute_eigenvectors);
 
  484   std::vector<CPU_NumericType> largest_eigenvalues;
 
  487     largest_eigenvalues.push_back(eigenvalues[size_krylov-i]);
 
  490   return largest_eigenvalues;
 
  503 template<
typename MatrixT>
 
  504 std::vector< typename viennacl::result_of::cpu_value_type<typename MatrixT::value_type>::type >
 
  510   return eig(matrix, eigenvectors, tag, 
false);
 
int method() const 
Returns the reorthogonalization method. 
T norm_2(std::vector< T, A > const &v1)
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. 
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)...
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
void method(int met)
Sets the reorthogonalization method. 
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) 
lanczos_tag(double factor=0.75, vcl_size_t numeig=10, int met=0, vcl_size_t krylov=100)
The constructor. 
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)
Generic interface for the computation of inner products. See viennacl/linalg/vector_operations.hpp for implementations. 
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...
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.) 
Implementation of the compressed_matrix class. 
void copy_vec_to_vec(viennacl::vector< NumericT > const &src, OtherVectorT &dest)
overloaded function for copying vectors 
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. 
void num_eigenvalues(vcl_size_t numeig)
Sets the number of eigenvalues. 
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. 
handle_type & handle()
Returns the OpenCL handle, non-const-version. 
T::ERROR_CANNOT_DEDUCE_CPU_SCALAR_TYPE_FOR_T type
size_type size1() const
Returns the number of rows. 
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_row > row(const matrix_base< NumericT, F > &A, unsigned int i)
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
NumericT max(std::vector< NumericT > const &v1)
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. 
size_type size() const 
Returns the length of the vector (cf. std::vector) 
size_type internal_size1() const
Returns the internal number of rows. Usually required for launching OpenCL kernels only...
Implementation of the algorithm for finding eigenvalues of a tridiagonal matrix. 
double factor() const 
Returns the exponent. 
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. 
void factor(double fct)
Sets the exponent of epsilon. Values between 0.6 and 0.9 usually give best results. 
A tag for the lanczos algorithm.