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 243 244 245 246 247 248 249 250 251 252 253 254
|
<!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: power-iter.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
 <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('power-iter_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"> </span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark"> </span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark"> </span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark"> </span>Files</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark"> </span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark"> </span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark"> </span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark"> </span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark"> </span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark"> </span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark"> </span>Macros</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(11)"><span class="SelectionMark"> </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">power-iter.cpp</div> </div>
</div><!--header-->
<div class="contents">
<p>This tutorial demonstrates the calculation of the eigenvalue with largest modulus using the power iteration method.</p>
<p>We start with including the necessary headers: </p>
<div class="fragment"><div class="line"><span class="comment">// Sparse matrices in uBLAS are *very* slow if debug mode is enabled. Disable it:</span></div>
<div class="line"><span class="preprocessor">#ifndef NDEBUG</span></div>
<div class="line"><span class="preprocessor"> #define BOOST_UBLAS_NDEBUG</span></div>
<div class="line"><span class="preprocessor">#endif</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 <iostream></span></div>
<div class="line"><span class="preprocessor">#include <fstream></span></div>
<div class="line"><span class="preprocessor">#include <limits></span></div>
<div class="line"><span class="preprocessor">#include <string></span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS</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 "<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="power__iter_8hpp.html">viennacl/linalg/power_iter.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>"</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 <iostream></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Boost includes:</span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix_proxy.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix_expression.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix_sparse.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/vector.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/operation.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/vector_expression.hpp></span></div>
</div><!-- fragment --><p> To run the power iteration, we set up a sparse matrix using Boost.uBLAS, transfer it over to a ViennaCL matrix, and then run the algorithm. </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">// This example relies on double precision to be available and will provide only poor results with single precision</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 uBLAS matrix and read the matrix from the matrix-market file: </p>
<div class="fragment"><div class="line">boost::numeric::ublas::compressed_matrix<ScalarType> ublas_A;</div>
<div class="line"></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>(ublas_A, <span class="stringliteral">"../examples/testdata/mat65k.mtx"</span>))</div>
<div class="line">{</div>
<div class="line"> std::cout << <span class="stringliteral">"Error reading Matrix file"</span> << std::endl;</div>
<div class="line"> <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">}</div>
</div><!-- fragment --><p> Transfer the data from the uBLAS matrix over to the ViennaCL sparse matrix: </p>
<div class="fragment"><div class="line"><a name="_a3"></a><a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<ScalarType></a> vcl_A(ublas_A.size1(), ublas_A.size2());</div>
<div class="line"><a name="a4"></a><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_A, vcl_A);</div>
</div><!-- fragment --><p> Run the power iteration up until the largest eigenvalue changes by less than the specified tolerance. Print the results of running with uBLAS as well as ViennaCL and exit. </p>
<div class="fragment"><div class="line"><a name="_a5"></a><a class="code" href="classviennacl_1_1linalg_1_1power__iter__tag.html">viennacl::linalg::power_iter_tag</a> ptag(1e-6);</div>
<div class="line"></div>
<div class="line">std::cout << <span class="stringliteral">"Starting computation of eigenvalue with largest modulus (might take about a minute)..."</span> << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">"Result of power iteration with ublas matrix (single-threaded): "</span> << <a name="a6"></a><a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(ublas_A, ptag) << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">"Result of power iteration with ViennaCL (OpenCL accelerated): "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(vcl_A, ptag) << std::endl;</div>
</div><!-- fragment --><p> You can also obtain the associated <em>approximated</em> eigenvector by passing it as a third argument to <a class="el" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba" title="Implementation of the calculation of eigenvalues using lanczos (with and without reorthogonalization)...">eig()</a> Tighten the tolerance passed to ptag above in order to obtain more accurate results. </p>
<div class="fragment"><div class="line"> <a name="_a7"></a><a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> eigenvector(vcl_A.size1());</div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(vcl_A, ptag, eigenvector);</div>
<div class="line"> std::cout << <span class="stringliteral">"First three entries in eigenvector: "</span> << eigenvector[0] << <span class="stringliteral">" "</span> << eigenvector[1] << <span class="stringliteral">" "</span> << eigenvector[2] << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> Ax = <a name="a8"></a><a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(vcl_A, eigenvector);</div>
<div class="line"> std::cout << <span class="stringliteral">"First three entries in A*eigenvector: "</span> << Ax[0] << <span class="stringliteral">" "</span> << Ax[1] << <span class="stringliteral">" "</span> << Ax[2] << std::endl;</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">// Sparse matrices in uBLAS are *very* slow if debug mode is enabled. Disable it:</span></div>
<div class="line"><span class="preprocessor">#ifndef NDEBUG</span></div>
<div class="line"><span class="preprocessor"> #define BOOST_UBLAS_NDEBUG</span></div>
<div class="line"><span class="preprocessor">#endif</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 <iostream></span></div>
<div class="line"><span class="preprocessor">#include <fstream></span></div>
<div class="line"><span class="preprocessor">#include <limits></span></div>
<div class="line"><span class="preprocessor">#include <string></span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS</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 "<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="power__iter_8hpp.html">viennacl/linalg/power_iter.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>"</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 <iostream></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Boost includes:</span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix_proxy.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix_expression.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/matrix_sparse.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/vector.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/operation.hpp></span></div>
<div class="line"><span class="preprocessor">#include <boost/numeric/ublas/vector_expression.hpp></span></div>
<div class="line"></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">// This example relies on double precision to be available and will provide only poor results with single precision</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"> boost::numeric::ublas::compressed_matrix<ScalarType> ublas_A;</div>
<div class="line"></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>(ublas_A, <span class="stringliteral">"../examples/testdata/mat65k.mtx"</span>))</div>
<div class="line"> {</div>
<div class="line"> std::cout << <span class="stringliteral">"Error reading Matrix file"</span> << 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<ScalarType></a> vcl_A(ublas_A.size1(), ublas_A.size2());</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_A, vcl_A);</div>
<div class="line"></div>
<div class="line"> <a class="code" href="classviennacl_1_1linalg_1_1power__iter__tag.html">viennacl::linalg::power_iter_tag</a> ptag(1e-6);</div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"Starting computation of eigenvalue with largest modulus (might take about a minute)..."</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"Result of power iteration with ublas matrix (single-threaded): "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(ublas_A, ptag) << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"Result of power iteration with ViennaCL (OpenCL accelerated): "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(vcl_A, ptag) << std::endl;</div>
<div class="line"></div>
<div class="line"> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> eigenvector(vcl_A.size1());</div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1linalg.html#af5002a88fa4cc3fbe65a904797a0ecba">viennacl::linalg::eig</a>(vcl_A, ptag, eigenvector);</div>
<div class="line"> std::cout << <span class="stringliteral">"First three entries in eigenvector: "</span> << eigenvector[0] << <span class="stringliteral">" "</span> << eigenvector[1] << <span class="stringliteral">" "</span> << eigenvector[2] << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> Ax = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(vcl_A, eigenvector);</div>
<div class="line"> std::cout << <span class="stringliteral">"First three entries in A*eigenvector: "</span> << Ax[0] << <span class="stringliteral">" "</span> << Ax[1] << <span class="stringliteral">" "</span> << Ax[2] << std::endl;</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>
|