ViennaCL - The Vienna Computing Library  1.7.1
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
qr_method.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2016, Institute for Microelectronics,
3  Institute for Analysis and Scientific Computing,
4  TU Wien.
5  Portions of this software are copyright by UChicago Argonne, LLC.
6 
7  -----------------
8  ViennaCL - The Vienna Computing Library
9  -----------------
10 
11  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
12 
13  (A list of authors and contributors can be found in the PDF manual)
14 
15  License: MIT (X11), see file LICENSE in the base directory
16 ============================================================================= */
17 
18 
23 /*
24 Solutions for testdata were generated with Scilab line:
25 
26 M=fscanfMat('nsm1.example');e=spec(M);e=gsort(e);rr=real(e);ii=imag(e);e=cat(1, rr, ii); s=strcat(string(e), ' ');write('tmp', s);
27 */
28 
29 
30 //#define VIENNACL_DEBUG_ALL
31 #include <iostream>
32 #include <fstream>
33 #include <stdexcept>
34 #include <vector>
35 
36 #include "viennacl/linalg/prod.hpp"
38 
39 #include "viennacl/tools/timer.hpp"
40 
41 #include <boost/numeric/ublas/vector.hpp>
42 #include <boost/numeric/ublas/matrix.hpp>
43 
44 namespace ublas = boost::numeric::ublas;
45 
46 void read_matrix_size(std::fstream& f, std::size_t& sz);
47 
48 void read_matrix_size(std::fstream& f, std::size_t& sz)
49 {
50  if(!f.is_open())
51  {
52  throw std::invalid_argument("File is not opened");
53  }
54 
55  f >> sz;
56 }
57 
58 template <typename NumericT, typename MatrixLayout>
60 {
61  if(!f.is_open())
62  {
63  throw std::invalid_argument("File is not opened");
64  }
65 
66  boost::numeric::ublas::matrix<NumericT> h_A(A.size1(), A.size2());
67 
68  for(std::size_t i = 0; i < h_A.size1(); i++) {
69  for(std::size_t j = 0; j < h_A.size2(); j++) {
70  NumericT val = 0.0;
71  f >> val;
72  h_A(i, j) = val;
73  }
74  }
75 
76  viennacl::copy(h_A, A);
77 }
78 
79 template<typename NumericT>
80 void read_vector_body(std::fstream& f, std::vector<NumericT>& v) {
81  if(!f.is_open())
82  throw std::invalid_argument("File is not opened");
83 
84  for(std::size_t i = 0; i < v.size(); i++)
85  {
86  NumericT val = 0.0;
87  f >> val;
88  v[i] = val;
89  }
90 }
91 
92 template<typename NumericT, typename MatrixLayout>
94 {
95  ublas::matrix<NumericT> A(A_orig.size1(), A_orig.size2());
96  viennacl::copy(A_orig, A);
97 
98  for (unsigned int i = 0; i < A.size1(); i++) {
99  for (unsigned int j = 0; j < A.size2(); j++) {
100  if ((std::abs(A(i, j)) > EPS) && ((i - 1) != j) && (i != j) && ((i + 1) != j))
101  {
102  // std::cout << "Failed at " << i << " " << j << " " << A(i, j) << "\n";
103  return false;
104  }
105  }
106  }
107  return true;
108 }
109 
110 template <typename NumericT, typename MatrixLayout>
112 {
113  ublas::matrix<NumericT> A(A_orig.size1(), A_orig.size2());
114  viennacl::copy(A_orig, A);
115 
116  for (std::size_t i = 0; i < A.size1(); i++) {
117  for (std::size_t j = 0; j < A.size2(); j++) {
118  if ((std::abs(A(i, j)) > EPS) && (i > (j + 1)))
119  {
120  // std::cout << "Failed at " << i << " " << j << " " << A(i, j) << "\n";
121  return false;
122  }
123  }
124  }
125  return true;
126 }
127 
128 template<typename NumericT>
129 NumericT matrix_compare(ublas::matrix<NumericT>& res,
130  ublas::matrix<NumericT>& ref)
131 {
132  NumericT diff = 0.0;
133  NumericT mx = 0.0;
134 
135  for(std::size_t i = 0; i < res.size1(); i++)
136  {
137  for(std::size_t j = 0; j < res.size2(); j++)
138  {
139  diff = std::max(diff, std::abs(res(i, j) - ref(i, j)));
140  mx = std::max(mx, res(i, j));
141  }
142  }
143 
144  return diff / mx;
145 }
146 
147 template<typename NumericT>
148 NumericT vector_compare(std::vector<NumericT> & res,
149  std::vector<NumericT> & ref)
150 {
151  std::sort(ref.begin(), ref.end());
152  std::sort(res.begin(), res.end());
153 
154  NumericT diff = 0.0;
155  NumericT mx = 0.0;
156  for(size_t i = 0; i < res.size(); i++)
157  {
158  diff = std::max(diff, std::abs(res[i] - ref[i]));
159  mx = std::max(mx, res[i]);
160  }
161 
162  return diff / mx;
163 }
164 
165 template <typename NumericT, typename MatrixLayout>
167 {
168  for (unsigned int i = 0; i < A.size1(); i++) {
169  for (unsigned int j = 0; j < A.size2(); j++)
170  std::cout << std::fixed << A(i, j) << "\t";
171  std::cout << "\n";
172  }
173 }
174 
175 template <typename NumericT, typename MatrixLayout>
176 void test_eigen(const std::string& fn, bool is_symm, NumericT EPS)
177 {
178  std::cout << "Reading..." << "\n";
179  std::size_t sz;
180  // read file
181  std::fstream f(fn.c_str(), std::fstream::in);
182  //read size of input matrix
183  read_matrix_size(f, sz);
184 
186  if (is_row)
187  std::cout << "Testing row-major matrix of size " << sz << "-by-" << sz << std::endl;
188  else
189  std::cout << "Testing column-major matrix of size " << sz << "-by-" << sz << std::endl;
190 
191  viennacl::matrix<NumericT> A_input(sz, sz), A_ref(sz, sz), Q(sz, sz);
192  // reference vector with reference values from file
193  std::vector<NumericT> eigen_ref_re(sz);
194  // calculated real eigenvalues
195  std::vector<NumericT> eigen_re(sz);
196  // calculated im. eigenvalues
197  std::vector<NumericT> eigen_im(sz);
198 
199  // read input matrix from file
200  read_matrix_body(f, A_input);
201  // read reference eigenvalues from file
202  read_vector_body(f, eigen_ref_re);
203 
204 
205  f.close();
206 
207  A_ref = A_input;
208 
209  std::cout << "Calculation..." << "\n";
210 
212  timer.start();
213  // Start the calculation
214  if(is_symm)
215  viennacl::linalg::qr_method_sym(A_input, Q, eigen_re);
216  else
217  viennacl::linalg::qr_method_nsm(A_input, Q, eigen_re, eigen_im);
218 /*
219 
220  std::cout << "\n\n Matrix A: \n\n";
221  matrix_print(A_input);
222  std::cout << "\n\n";
223 
224  std::cout << "\n\n Matrix Q: \n\n";
225  matrix_print(Q);
226  std::cout << "\n\n";
227 */
228 
229  double time_spend = timer.get();
230 
231  std::cout << "Verification..." << "\n";
232 
233  bool is_hessenberg = check_hessenberg(A_input, EPS);
234  bool is_tridiag = check_tridiag(A_input, EPS);
235 
236  ublas::matrix<NumericT> A_ref_ublas(sz, sz), A_input_ublas(sz, sz), Q_ublas(sz, sz), result1(sz, sz), result2(sz, sz);
237  viennacl::copy(A_ref, A_ref_ublas);
238  viennacl::copy(A_input, A_input_ublas);
239  viennacl::copy(Q, Q_ublas);
240 
241  // compute result1 = ublas::prod(Q_ublas, A_input_ublas); (terribly slow when using ublas directly)
242  for (std::size_t i=0; i<result1.size1(); ++i)
243  for (std::size_t j=0; j<result1.size2(); ++j)
244  {
245  NumericT value = 0;
246  for (std::size_t k=0; k<Q_ublas.size2(); ++k)
247  value += Q_ublas(i, k) * A_input_ublas(k, j);
248  result1(i,j) = value;
249  }
250  // compute result2 = ublas::prod(A_ref_ublas, Q_ublas); (terribly slow when using ublas directly)
251  for (std::size_t i=0; i<result2.size1(); ++i)
252  for (std::size_t j=0; j<result2.size2(); ++j)
253  {
254  NumericT value = 0;
255  for (std::size_t k=0; k<A_ref_ublas.size2(); ++k)
256  value += A_ref_ublas(i, k) * Q_ublas(k, j);
257  result2(i,j) = value;
258  }
259 
260 
261  NumericT prods_diff = matrix_compare(result1, result2);
262  NumericT eigen_diff = vector_compare(eigen_re, eigen_ref_re);
263 
264 
265  bool is_ok = is_hessenberg;
266 
267  if(is_symm)
268  is_ok = is_ok && is_tridiag;
269 
270  is_ok = is_ok && (eigen_diff < EPS);
271  is_ok = is_ok && (prods_diff < EPS);
272 
273  // std::cout << A_ref << "\n";
274  // std::cout << A_input << "\n";
275  // std::cout << Q << "\n";
276  // std::cout << eigen_re << "\n";
277  // std::cout << eigen_im << "\n";
278  // std::cout << eigen_ref_re << "\n";
279  // std::cout << eigen_ref_im << "\n";
280 
281  // std::cout << result1 << "\n";
282  // std::cout << result2 << "\n";
283  // std::cout << eigen_ref << "\n";
284  // std::cout << eigen << "\n";
285 
286  printf("%6s [%dx%d] %40s time = %.4f\n", is_ok?"[[OK]]":"[FAIL]", (int)A_ref.size1(), (int)A_ref.size2(), fn.c_str(), time_spend);
287  printf("tridiagonal = %d, hessenberg = %d prod-diff = %f eigen-diff = %f\n", is_tridiag, is_hessenberg, prods_diff, eigen_diff);
288  std::cout << std::endl << std::endl;
289 
290  if (!is_ok)
291  exit(EXIT_FAILURE);
292 
293 }
294 
295 int main()
296 {
297  float epsilon1 = 0.0001f;
298 
299  std::cout << "# Testing setup:" << std::endl;
300  std::cout << " eps: " << epsilon1 << std::endl;
301  std::cout << " numeric: double" << std::endl;
302  std::cout << std::endl;
303  test_eigen<float, viennacl::row_major >("../examples/testdata/eigen/symm5.example", true, epsilon1);
304  test_eigen<float, viennacl::column_major>("../examples/testdata/eigen/symm5.example", true, epsilon1);
305 
306  #ifdef VIENNACL_WITH_OPENCL
308  #endif
309  {
310  double epsilon2 = 1e-5;
311 
312  std::cout << "# Testing setup:" << std::endl;
313  std::cout << " eps: " << epsilon2 << std::endl;
314  std::cout << " numeric: double" << std::endl;
315  std::cout << std::endl;
316  test_eigen<double, viennacl::row_major >("../examples/testdata/eigen/symm5.example", true, epsilon2);
317  test_eigen<double, viennacl::column_major>("../examples/testdata/eigen/symm5.example", true, epsilon2);
318  }
319 
320  //test_eigen<viennacl::row_major>("../../examples/testdata/eigen/symm3.example", true); // Computation of this matrix takes very long
321  //test_eigen<viennacl::column_major>("../../examples/testdata/eigen/symm3.example", true);
322 
323  //test_eigen<viennacl::row_major>("../examples/testdata/eigen/nsm2.example", false);
324  //test_eigen<viennacl::row_major>("../../examples/testdata/eigen/nsm2.example", false);
325  //test_eigen("../../examples/testdata/eigen/nsm3.example", false);
326  //test_eigen("../../examples/testdata/eigen/nsm4.example", false); //Note: This test suffers from round-off errors in single precision, hence disabled
327 
328  std::cout << std::endl;
329  std::cout << "------- Test completed --------" << std::endl;
330  std::cout << std::endl;
331 
332  return EXIT_SUCCESS;
333 }
void read_matrix_body(std::fstream &f, viennacl::matrix< NumericT, MatrixLayout > &A)
Definition: qr_method.cpp:59
Simple timer class based on gettimeofday (POSIX) or QueryPerformanceCounter (Windows).
Definition: timer.hpp:90
#define EPS
Definition: bisect.cpp:38
Helper class for checking whether a matrix has a row-major layout.
Definition: forwards.h:484
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
bool check_tridiag(viennacl::matrix< NumericT, MatrixLayout > &A_orig, NumericT EPS)
Definition: qr_method.cpp:93
A dense matrix class.
Definition: forwards.h:375
T max(const T &lhs, const T &rhs)
Maximum.
Definition: util.hpp:59
viennacl::ocl::device const & current_device()
Convenience function for returning the active device in the current context.
Definition: backend.hpp:351
float NumericT
Definition: bisect.cpp:40
NumericT vector_compare(std::vector< NumericT > &res, std::vector< NumericT > &ref)
Definition: qr_method.cpp:148
void read_matrix_size(std::fstream &f, std::size_t &sz)
Definition: qr_method.cpp:48
void qr_method_nsm(viennacl::matrix< SCALARTYPE > &A, viennacl::matrix< SCALARTYPE > &Q, std::vector< SCALARTYPE > &D, std::vector< SCALARTYPE > &E)
Definition: qr-method.hpp:796
bool check_hessenberg(viennacl::matrix< NumericT, MatrixLayout > &A_orig, NumericT EPS)
Definition: qr_method.cpp:111
bool double_support() const
ViennaCL convenience function: Returns true if the device supports double precision.
Definition: device.hpp:956
size_type size2() const
Returns the number of columns.
Definition: matrix_def.hpp:226
size_type size1() const
Returns the number of rows.
Definition: matrix_def.hpp:224
int main()
Definition: qr_method.cpp:295
A simple, yet (mostly) sufficiently accurate timer for benchmarking and profiling.
Implementation of the QR method for eigenvalue computations. Experimental.
void test_eigen(const std::string &fn, bool is_symm, NumericT EPS)
Definition: qr_method.cpp:176
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) ...
double get() const
Definition: timer.hpp:104
ScalarType diff(ScalarType &s1, viennacl::scalar< ScalarType > &s2)
Definition: blas3_solve.cpp:50
void read_vector_body(std::fstream &f, std::vector< NumericT > &v)
Definition: qr_method.cpp:80
void matrix_print(viennacl::matrix< NumericT, MatrixLayout > &A)
Definition: qr_method.cpp:166
NumericT matrix_compare(ublas::matrix< NumericT > &res, ublas::matrix< NumericT > &ref)
Definition: qr_method.cpp:129
void qr_method_sym(viennacl::matrix< SCALARTYPE > &A, viennacl::matrix< SCALARTYPE > &Q, std::vector< SCALARTYPE > &D)
Definition: qr-method.hpp:806