1 #ifndef VIENNACL_LINALG_LANCZOS_HPP_
2 #define VIENNACL_LINALG_LANCZOS_HPP_
36 #include <boost/random.hpp>
37 #include <boost/random/mersenne_twister.hpp>
38 #include <boost/numeric/ublas/matrix.hpp>
39 #include <boost/numeric/ublas/matrix_proxy.hpp>
40 #include <boost/numeric/ublas/matrix_expression.hpp>
41 #include <boost/numeric/ublas/matrix_sparse.hpp>
42 #include <boost/numeric/ublas/vector.hpp>
43 #include <boost/numeric/ublas/operation.hpp>
44 #include <boost/numeric/ublas/vector_expression.hpp>
45 #include <boost/numeric/ublas/io.hpp>
76 vcl_size_t krylov = 100) : factor_(
factor), num_eigenvalues_(numeig), method_(met), krylov_size_(krylov) {}
85 void factor(
double fct) { factor_ = fct; }
88 double factor()
const {
return factor_; }
97 void method(
int met){ method_ = met; }
123 template<
typename MatrixT,
typename VectorT >
135 boost::normal_distribution<CPU_ScalarType> N(0, 1);
136 boost::bernoulli_distribution<CPU_ScalarType> B(0.5);
137 boost::triangle_distribution<CPU_ScalarType> T(-1, 0, 1);
139 boost::variate_generator<boost::mt11213b&, boost::normal_distribution<CPU_ScalarType> > get_N(mt, N);
140 boost::variate_generator<boost::mt11213b&, boost::bernoulli_distribution<CPU_ScalarType> > get_B(mt, B);
141 boost::variate_generator<boost::mt11213b&, boost::triangle_distribution<CPU_ScalarType> > get_T(mt, T);
144 long i, k, retry, reorths;
145 std::vector<long> l_bound(size/2), u_bound(size/2);
147 CPU_ScalarType squ_eps, eta, temp, eps, retry_th;
149 std::vector< std::vector<CPU_ScalarType> > w(2, std::vector<CPU_ScalarType>(size));
150 CPU_ScalarType cpu_beta;
152 boost::numeric::ublas::vector<CPU_ScalarType> s(n);
155 CPU_ScalarType inner_rt;
157 ScalarType vcl_alpha;
158 std::vector<CPU_ScalarType> alphas, betas;
159 boost::numeric::ublas::matrix<CPU_ScalarType> Q(n, size);
162 eps = std::numeric_limits<CPU_ScalarType>::epsilon();
163 squ_eps = std::sqrt(eps);
165 eta = std::exp(std::log(eps) * tag.
factor());
178 alphas.push_back(vcl_alpha);
180 betas.push_back(vcl_beta);
183 for (i = 1;i < static_cast<long>(
size); i++)
185 r = u - vcl_alpha * r;
188 betas.push_back(vcl_beta);
194 w[index][0] = (betas[1] * w[
vcl_size_t(k)][1] + (alphas[0] - vcl_alpha) * w[
vcl_size_t(k)][0] - betas[
vcl_size_t(i) - 1] * w[index][0]) / vcl_beta + eps * 0.3 * get_N() * (betas[1] + vcl_beta);
198 w[index][j] = (betas[j + 1] * w[
vcl_size_t(k)][j + 1] + (alphas[j] - vcl_alpha) * w[
vcl_size_t(k)][j] + betas[j] * w[
vcl_size_t(k)][j - 1] - betas[
vcl_size_t(i) - 1] * w[index][j]) / vcl_beta + eps * 0.3 * get_N() * (betas[j + 1] + vcl_beta);
200 w[index][
vcl_size_t(i) - 1] = 0.6 * eps * CPU_ScalarType(n) * get_N() * betas[1] / vcl_beta;
204 for (vcl_size_t j = 0; j <
vcl_size_t(batches); j++)
209 for (k = l_bound[j];k < u_bound[j];k++)
213 r = r - inner_rt * t;
214 w[index][
vcl_size_t(k)] = 1.5 * eps * get_N();
220 vcl_beta = vcl_beta * temp;
225 for (vcl_size_t j = 0; j <
vcl_size_t(i); j++)
227 if (std::fabs(w[index][j]) >= squ_eps)
231 r = r - inner_rt * t;
232 w[index][j] = 1.5 * eps * get_N();
235 while (k >= 0 && std::fabs(w[index][
vcl_size_t(k)]) > eta)
239 r = r - inner_rt * t;
240 w[index][
vcl_size_t(k)] = 1.5 * eps * get_N();
247 while (k < i && std::fabs(w[index][
vcl_size_t(k)]) > eta)
251 r = r - inner_rt * t;
252 w[index][
vcl_size_t(k)] = 1.5 * eps * get_N();
266 vcl_beta = vcl_beta * temp;
269 while (temp < retry_th)
271 for (vcl_size_t j = 0; j <
vcl_size_t(i); j++)
275 r = r - inner_rt * t;
281 vcl_beta = vcl_beta * temp;
293 alphas.push_back(vcl_alpha);
296 return bisect(alphas, betas);
308 template<
typename MatrixT,
typename VectorT>
318 ScalarType vcl_alpha;
319 std::vector<CPU_ScalarType> alphas, betas;
323 boost::numeric::ublas::vector<CPU_ScalarType> s(r.size()), u_zero(n), q(n);
324 boost::numeric::ublas::matrix<CPU_ScalarType> Q(n, size);
326 u_zero = boost::numeric::ublas::zero_vector<CPU_ScalarType>(n);
340 r = u - vcl_alpha * r;
347 alphas.push_back(vcl_alpha);
348 betas.push_back(vcl_beta);
352 return bisect(alphas, betas);
363 template<
typename MatrixT,
typename VectorT >
372 CPU_NumericType temp;
373 CPU_NumericType norm;
374 NumericType vcl_beta;
375 NumericType vcl_alpha;
376 std::vector<CPU_NumericType> alphas, betas;
379 NumericType inner_rt;
380 boost::numeric::ublas::vector<CPU_NumericType> u_zero(n), s(r.size()), q(n);
381 boost::numeric::ublas::matrix<CPU_NumericType> Q(n, size);
396 r = r - inner_rt * t;
401 vcl_beta = temp * norm;
407 r = u - vcl_alpha * r;
412 alphas.push_back(vcl_alpha);
413 betas.push_back(vcl_beta);
416 return bisect(alphas, betas);
428 template<
typename MatrixT>
429 std::vector< typename viennacl::result_of::cpu_value_type<typename MatrixT::value_type>::type >
434 typedef typename viennacl::result_of::vector_for_matrix<MatrixT>::type VectorT;
437 boost::normal_distribution<CPU_NumericType> N(0, 1);
438 boost::bernoulli_distribution<CPU_NumericType> B(0.5);
439 boost::triangle_distribution<CPU_NumericType> T(-1, 0, 1);
441 boost::variate_generator<boost::mt11213b&, boost::normal_distribution<CPU_NumericType> > get_N(mt, N);
442 boost::variate_generator<boost::mt11213b&, boost::bernoulli_distribution<CPU_NumericType> > get_B(mt, B);
443 boost::variate_generator<boost::mt11213b&, boost::triangle_distribution<CPU_NumericType> > get_T(mt, T);
445 std::vector<CPU_NumericType> eigenvalues;
447 VectorT r(matrix_size);
448 std::vector<CPU_NumericType> s(matrix_size);
451 s[i] = 3.0 * get_B() + get_T() - 1.5;
471 std::vector<CPU_NumericType> largest_eigenvalues;
474 largest_eigenvalues.push_back(eigenvalues[size_krylov-i]);
477 return largest_eigenvalues;
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.
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< typename viennacl::result_of::cpu_value_type< typename MatrixT::value_type >::type > lanczosPRO(MatrixT const &A, VectorT &r, vcl_size_t size, lanczos_tag const &tag)
Implementation of the Lanczos PRO algorithm.
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.)
std::vector< typename viennacl::result_of::cpu_value_type< typename MatrixT::value_type >::type > lanczos(MatrixT const &A, VectorT &r, vcl_size_t size, lanczos_tag)
Implementation of the lanczos algorithm without reorthogonalization.
Implementation of the compressed_matrix class.
void copy_vec_to_vec(viennacl::vector< NumericT > const &src, OtherVectorT &dest)
overloaded function for copying vectors
void num_eigenvalues(vcl_size_t numeig)
Sets the number of eigenvalues.
vcl_size_t num_eigenvalues() const
Returns the number of eigenvalues.
std::vector< typename viennacl::result_of::cpu_value_type< typename MatrixT::value_type >::type > lanczosFRO(MatrixT const &A, VectorT &r, vcl_size_t size, lanczos_tag)
Implementation of the Lanczos FRO algorithm.
T::ERROR_CANNOT_DEDUCE_CPU_SCALAR_TYPE_FOR_T type
std::vector< typename viennacl::result_of::cpu_value_type< typename MatrixT::value_type >::type > eig(MatrixT const &matrix, lanczos_tag const &tag)
Implementation of the calculation of eigenvalues using lanczos.
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
NumericT max(std::vector< NumericT > const &v1)
void krylov_size(vcl_size_t max)
Sets the size of the kylov space.
Implementation of the algorithm for finding eigenvalues of a tridiagonal matrix.
double factor() const
Returns the exponent.
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_column > column(const matrix_base< NumericT, F > &A, unsigned int j)
void factor(double fct)
Sets the exponent of epsilon.
A tag for the lanczos algorithm.