File: lanczos_8cpp-example.html

package info (click to toggle)
viennacl 1.7.1%2Bdfsg1-6
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, sid, trixie
  • size: 114,428 kB
  • sloc: sh: 454,206; cpp: 109,088; ansic: 2,103; perl: 104; makefile: 22
file content (242 lines) | stat: -rw-r--r-- 22,999 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.8.6"/>
<title>ViennaCL - The Vienna Computing Library: lanczos.cpp</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<script type="text/javascript">
  $(document).ready(initResizable);
  $(window).load(resizeHeight);
</script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/search.js"></script>
<script type="text/javascript">
  $(document).ready(function() { searchBox.OnSelectItem(0); });
</script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
 <tbody>
 <tr style="height: 56px;">
  <td style="padding-left: 0.5em;">
   <div id="projectname">ViennaCL - The Vienna Computing Library
   &#160;<span id="projectnumber">1.7.1</span>
   </div>
   <div id="projectbrief">Free open-source GPU-accelerated linear algebra and solver library.</div>
  </td>
   <td>        <div id="MSearchBox" class="MSearchBoxInactive">
        <span class="left">
          <img id="MSearchSelect" src="search/mag_sel.png"
               onmouseover="return searchBox.OnSearchSelectShow()"
               onmouseout="return searchBox.OnSearchSelectHide()"
               alt=""/>
          <input type="text" id="MSearchField" value="Search" accesskey="S"
               onfocus="searchBox.OnSearchFieldFocus(true)" 
               onblur="searchBox.OnSearchFieldFocus(false)" 
               onkeyup="searchBox.OnSearchFieldChange(event)"/>
          </span><span class="right">
            <a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a>
          </span>
        </div>
</td>
 </tr>
 </tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.6 -->
<script type="text/javascript">
var searchBox = new SearchBox("searchBox", "search",false,'Search');
</script>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
  <div id="nav-tree">
    <div id="nav-tree-contents">
      <div id="nav-sync" class="sync"></div>
    </div>
  </div>
  <div id="splitbar" style="-moz-user-select:none;" 
       class="ui-resizable-handle">
  </div>
</div>
<script type="text/javascript">
$(document).ready(function(){initNavTree('lanczos_8cpp-example.html','');});
</script>
<div id="doc-content">
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
     onmouseover="return searchBox.OnSearchSelectShow()"
     onmouseout="return searchBox.OnSearchSelectHide()"
     onkeydown="return searchBox.OnSearchSelectKey(event)">
<a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark">&#160;</span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark">&#160;</span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark">&#160;</span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark">&#160;</span>Files</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark">&#160;</span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark">&#160;</span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark">&#160;</span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark">&#160;</span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark">&#160;</span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark">&#160;</span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark">&#160;</span>Macros</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(11)"><span class="SelectionMark">&#160;</span>Pages</a></div>

<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<iframe src="javascript:void(0)" frameborder="0" 
        name="MSearchResults" id="MSearchResults">
</iframe>
</div>

<div class="header">
  <div class="headertitle">
<div class="title">lanczos.cpp</div>  </div>
</div><!--header-->
<div class="contents">
<p>This tutorial shows how to calculate the largest eigenvalues of a matrix using Lanczos' method.</p>
<p>The Lanczos method is particularly attractive for use with large, sparse matrices, since the only requirement on the matrix is to provide a matrix-vector product. Although less common, the method is sometimes also used with dense matrices.</p>
<p>We start with including the necessary headers: </p>
<div class="fragment"><div class="line"><span class="comment">// include necessary system headers</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">//include basic scalar and vector types of ViennaCL</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix_8hpp.html">viennacl/matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__proxy_8hpp.html">viennacl/matrix_proxy.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="lanczos_8hpp.html">viennacl/linalg/lanczos.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Some helper functions for this tutorial:</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;string&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;iomanip&gt;</span></div>
</div><!-- fragment --><p> We read a sparse matrix from a matrix-market file, then run the Lanczos method. Finally, the computed eigenvalues are printed. </p>
<div class="fragment"><div class="line"><span class="keywordtype">int</span> <a name="a0"></a><a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a>()</div>
<div class="line">{</div>
<div class="line">  <span class="comment">// If you GPU does not support double precision, use `float` instead of `double`:</span></div>
<div class="line">  <span class="keyword">typedef</span> <span class="keywordtype">double</span>     <a name="a1"></a><a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>;</div>
</div><!-- fragment --><p> Create the sparse matrix and read data from a Matrix-Market file: </p>
<div class="fragment"><div class="line">std::vector&lt; std::map&lt;unsigned int, ScalarType&gt; &gt; host_A;</div>
<div class="line"><span class="keywordflow">if</span> (!<a name="a2"></a><a class="code" href="namespaceviennacl_1_1io.html#ad934125ed3dbe661e264bcd7d62b1048">viennacl::io::read_matrix_market_file</a>(host_A, <span class="stringliteral">&quot;../examples/testdata/mat65k.mtx&quot;</span>))</div>
<div class="line">{</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading Matrix file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">}</div>
<div class="line"></div>
<div class="line"><a name="_a3"></a><a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;ScalarType&gt;</a> A;</div>
<div class="line"><a name="a4"></a><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(host_A, A);</div>
</div><!-- fragment --><p> Create the configuration for the Lanczos method. All constructor arguments are optional, so feel free to use default settings. </p>
<div class="fragment"><div class="line"><a name="_a5"></a><a class="code" href="classviennacl_1_1linalg_1_1lanczos__tag.html">viennacl::linalg::lanczos_tag</a> ltag(0.75,    <span class="comment">// Select a power of 0.75 as the tolerance for the machine precision.</span></div>
<div class="line">                                   10,      <span class="comment">// Compute (approximations to) the 10 largest eigenvalues</span></div>
<div class="line">                                   <a name="a6"></a><a class="code" href="classviennacl_1_1linalg_1_1lanczos__tag.html#a8a743b9f474124a0f779ed35c6f6a684a937ee8052e51309672b2bcdba4bb015e">viennacl::linalg::lanczos_tag::partial_reorthogonalization</a>, <span class="comment">// use partial reorthogonalization</span></div>
<div class="line">                                   30);   <span class="comment">// Maximum size of the Krylov space</span></div>
</div><!-- fragment --><p> Run the Lanczos method for computing eigenvalues by passing the tag to the routine <a class="el" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba" title="Implementation of the calculation of eigenvalues using lanczos (with and without reorthogonalization)...">viennacl::linalg::eig()</a>. </p>
<div class="fragment"><div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;Running Lanczos algorithm (eigenvalues only). This might take a while...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">std::vector&lt;ScalarType&gt; lanczos_eigenvalues = <a name="a7"></a><a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(A, ltag);</div>
</div><!-- fragment --><p> Re-run the Lanczos method, this time also computing eigenvectors. To do so, we pass a dense matrix for which each column will ultimately hold the computed eigenvectors. </p>
<div class="fragment"><div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;Running Lanczos algorithm (with eigenvectors). This might take a while...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="_a8"></a><a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix&lt;ScalarType&gt;</a> approx_eigenvectors_A(A.<a name="a9"></a><a class="code" href="classviennacl_1_1compressed__matrix.html#a463cf1739f9cdd387aa185cb574db183">size1</a>(), ltag.num_eigenvalues());</div>
<div class="line">lanczos_eigenvalues = <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(A, approx_eigenvectors_A, ltag);</div>
</div><!-- fragment --><p> Print the computed eigenvalues and exit: </p>
<div class="fragment"><div class="line">  <span class="keywordflow">for</span> (std::size_t i = 0; i&lt; lanczos_eigenvalues.size(); i++)</div>
<div class="line">  {</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot;Approx. eigenvalue &quot;</span> &lt;&lt; std::setprecision(7) &lt;&lt; lanczos_eigenvalues[i];</div>
<div class="line"></div>
<div class="line">    <span class="comment">// test approximated eigenvector by computing A*v:</span></div>
<div class="line">    <a name="_a10"></a><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> approx_eigenvector = <a name="a11"></a><a class="code" href="namespaceviennacl.html#a7fca08f4a83edffe7f47666d298ca87d">viennacl::column</a>(approx_eigenvectors_A, static_cast&lt;unsigned int&gt;(i));</div>
<div class="line">    <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> Aq = <a name="a12"></a><a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(A, approx_eigenvector);</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot; (&quot;</span> &lt;&lt; <a name="a13"></a><a class="code" href="namespaceviennacl_1_1linalg.html#ab35950c4374eb3be08a03d852508c01a">viennacl::linalg::inner_prod</a>(Aq, approx_eigenvector) &lt;&lt; <span class="stringliteral">&quot; for &lt;Av,v&gt; with approx. eigenvector v)&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  }</div>
<div class="line"></div>
<div class="line">  <span class="keywordflow">return</span> EXIT_SUCCESS;</div>
<div class="line">}</div>
</div><!-- fragment --> <h2>Full Example Code</h2>
<div class="fragment"><div class="line"><span class="comment">/* =========================================================================</span></div>
<div class="line"><span class="comment">   Copyright (c) 2010-2016, Institute for Microelectronics,</span></div>
<div class="line"><span class="comment">                            Institute for Analysis and Scientific Computing,</span></div>
<div class="line"><span class="comment">                            TU Wien.</span></div>
<div class="line"><span class="comment">   Portions of this software are copyright by UChicago Argonne, LLC.</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment">                            -----------------</span></div>
<div class="line"><span class="comment">                  ViennaCL - The Vienna Computing Library</span></div>
<div class="line"><span class="comment">                            -----------------</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment">   Project Head:    Karl Rupp                   rupp@iue.tuwien.ac.at</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment">   (A list of authors and contributors can be found in the PDF manual)</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment">   License:         MIT (X11), see file LICENSE in the base directory</span></div>
<div class="line"><span class="comment">============================================================================= */</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// include necessary system headers</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">//include basic scalar and vector types of ViennaCL</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix_8hpp.html">viennacl/matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__proxy_8hpp.html">viennacl/matrix_proxy.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="lanczos_8hpp.html">viennacl/linalg/lanczos.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Some helper functions for this tutorial:</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;string&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;iomanip&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> <a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a>()</div>
<div class="line">{</div>
<div class="line">  <span class="comment">// If you GPU does not support double precision, use `float` instead of `double`:</span></div>
<div class="line">  <span class="keyword">typedef</span> <span class="keywordtype">double</span>     <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>;</div>
<div class="line"></div>
<div class="line">  std::vector&lt; std::map&lt;unsigned int, ScalarType&gt; &gt; host_A;</div>
<div class="line">  <span class="keywordflow">if</span> (!<a class="code" href="namespaceviennacl_1_1io.html#ad934125ed3dbe661e264bcd7d62b1048">viennacl::io::read_matrix_market_file</a>(host_A, <span class="stringliteral">&quot;../examples/testdata/mat65k.mtx&quot;</span>))</div>
<div class="line">  {</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading Matrix file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">    <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">  }</div>
<div class="line"></div>
<div class="line">  <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;ScalarType&gt;</a> A;</div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(host_A, A);</div>
<div class="line"></div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1lanczos__tag.html">viennacl::linalg::lanczos_tag</a> ltag(0.75,    <span class="comment">// Select a power of 0.75 as the tolerance for the machine precision.</span></div>
<div class="line">                                     10,      <span class="comment">// Compute (approximations to) the 10 largest eigenvalues</span></div>
<div class="line">                                     <a class="code" href="classviennacl_1_1linalg_1_1lanczos__tag.html#a8a743b9f474124a0f779ed35c6f6a684a937ee8052e51309672b2bcdba4bb015e">viennacl::linalg::lanczos_tag::partial_reorthogonalization</a>, <span class="comment">// use partial reorthogonalization</span></div>
<div class="line">                                     30);   <span class="comment">// Maximum size of the Krylov space</span></div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Running Lanczos algorithm (eigenvalues only). This might take a while...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  std::vector&lt;ScalarType&gt; lanczos_eigenvalues = <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(A, ltag);</div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Running Lanczos algorithm (with eigenvectors). This might take a while...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  <a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix&lt;ScalarType&gt;</a> approx_eigenvectors_A(A.<a class="code" href="classviennacl_1_1compressed__matrix.html#a463cf1739f9cdd387aa185cb574db183">size1</a>(), ltag.<a name="a14"></a><a class="code" href="classviennacl_1_1linalg_1_1lanczos__tag.html#ae744c43467774ee300f5fab0c607aed1">num_eigenvalues</a>());</div>
<div class="line">  lanczos_eigenvalues = <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(A, approx_eigenvectors_A, ltag);</div>
<div class="line"></div>
<div class="line">  <span class="keywordflow">for</span> (std::size_t i = 0; i&lt; lanczos_eigenvalues.size(); i++)</div>
<div class="line">  {</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot;Approx. eigenvalue &quot;</span> &lt;&lt; std::setprecision(7) &lt;&lt; lanczos_eigenvalues[i];</div>
<div class="line"></div>
<div class="line">    <span class="comment">// test approximated eigenvector by computing A*v:</span></div>
<div class="line">    <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> approx_eigenvector = <a class="code" href="namespaceviennacl.html#a7fca08f4a83edffe7f47666d298ca87d">viennacl::column</a>(approx_eigenvectors_A, static_cast&lt;unsigned int&gt;(i));</div>
<div class="line">    <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> Aq = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(A, approx_eigenvector);</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot; (&quot;</span> &lt;&lt; <a class="code" href="namespaceviennacl_1_1linalg.html#ab35950c4374eb3be08a03d852508c01a">viennacl::linalg::inner_prod</a>(Aq, approx_eigenvector) &lt;&lt; <span class="stringliteral">&quot; for &lt;Av,v&gt; with approx. eigenvector v)&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  }</div>
<div class="line"></div>
<div class="line">  <span class="keywordflow">return</span> EXIT_SUCCESS;</div>
<div class="line">}</div>
<div class="line"></div>
</div><!-- fragment --> </div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
  <ul>
    <li class="footer">Generated on Wed Jan 20 2016 22:32:38 for ViennaCL - The Vienna Computing Library by
    <a href="http://www.doxygen.org/index.html">
    <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.8.6 </li>
  </ul>
</div>
</body>
</html>