ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
qr_method.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2015, 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 
25 // System headers:
26 #include <iostream>
27 #include <fstream>
28 #include <stdexcept>
29 #include <vector>
30 
31 // ViennaCL headers:
32 #include "viennacl/matrix.hpp"
33 #include "viennacl/linalg/prod.hpp"
35 
36 
40 template <typename ScalarType>
41 void initialize(viennacl::matrix<ScalarType> & A, std::vector<ScalarType> & v)
42 {
43  // System matrix:
44  ScalarType M[9][9] = {{ 4, 1, -2, 2, -7, 3, 9, -6, -2},
45  { 1, -2, 0, 1, -1, 5, 4, 7, 3},
46  {-2, 0, 3, 2, 0, 3, 6, 1, -1},
47  { 2, 1, 2, 1, 4, 5, 6, 7, 8},
48  {-7, -1, 0, 4, 5, 4, 9, 1, -8},
49  { 3, 5, 3, 5, 4, 9, -3, 3, 3},
50  { 9, 4, 6, 6, 9, -3, 3, 6, -7},
51  {-6, 7, 1, 7, 1, 3, 6, 2, 6},
52  {-2, 3, -1, 8, -8, 3, -7, 6, 1}};
53 
54  for(std::size_t i = 0; i < 9; i++)
55  for(std::size_t j = 0; j < 9; j++)
56  A(i, j) = M[i][j];
57 
58  // Known eigenvalues:
59  ScalarType V[9] = {ScalarType(12.6005), ScalarType(19.5905), ScalarType(8.06067), ScalarType(2.95074), ScalarType(0.223506),
60  ScalarType(24.3642), ScalarType(-9.62084), ScalarType(-13.8374), ScalarType(-18.3319)};
61 
62  for(std::size_t i = 0; i < 9; i++)
63  v[i] = V[i];
64 }
65 
69 template <typename ScalarType>
70 void vector_print(std::vector<ScalarType>& v )
71 {
72  for (unsigned int i = 0; i < v.size(); i++)
73  std::cout << std::setprecision(6) << std::fixed << v[i] << "\t";
74  std::cout << std::endl;
75 }
76 
80 int main()
81 {
82  // Change to 'double' if supported by your hardware.
83  typedef float ScalarType;
84 
85  std::cout << "Testing matrix of size " << 9 << "-by-" << 9 << std::endl;
86 
87  viennacl::matrix<ScalarType> A_input(9,9);
89  std::vector<ScalarType> eigenvalues_ref(9);
90  std::vector<ScalarType> eigenvalues(9);
91 
92  viennacl::vector<ScalarType> vcl_eigenvalues(9);
93 
94  initialize(A_input, eigenvalues_ref);
95 
96  std::cout << std::endl <<"Input matrix: " << std::endl;
97  std::cout << A_input << std::endl;
98 
99  viennacl::matrix<ScalarType> A_input2(A_input); // duplicate to demonstrate the use with both std::vector and viennacl::vector
100 
101 
110  std::cout << "Calculation..." << std::endl;
111  viennacl::linalg::qr_method_sym(A_input, Q, eigenvalues);
112 
116  viennacl::linalg::qr_method_sym(A_input2, Q, vcl_eigenvalues);
117 
121  std::cout << std::endl << "Eigenvalues with std::vector<T>:" << std::endl;
122  vector_print(eigenvalues);
123  std::cout << "Eigenvalues with viennacl::vector<T>: " << std::endl << vcl_eigenvalues << std::endl;
124  std::cout << std::endl << "Reference eigenvalues:" << std::endl;
125  vector_print(eigenvalues_ref);
126  std::cout << std::endl;
127  std::cout << "Eigenvectors - each column is an eigenvector" << std::endl;
128  std::cout << Q << std::endl;
129 
133  std::cout << std::endl;
134  std::cout << "------- Tutorial completed --------" << std::endl;
135  std::cout << std::endl;
136 
137  return EXIT_SUCCESS;
138 }
139 
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
Implementation of the dense matrix class.
A dense matrix class.
Definition: forwards.h:375
int main()
Definition: qr_method.cpp:295
Implementation of the QR method for eigenvalue computations. Experimental.
float ScalarType
Definition: fft_1d.cpp:42
void vector_print(std::vector< ScalarType > &v)
void qr_method_sym(viennacl::matrix< SCALARTYPE > &A, viennacl::matrix< SCALARTYPE > &Q, std::vector< SCALARTYPE > &D)
Definition: qr-method.hpp:806