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 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442
|
<!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: eigen-with-viennacl.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('eigen-with-viennacl_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">eigen-with-viennacl.cpp</div> </div>
</div><!--header-->
<div class="contents">
<p>This tutorial shows how data can be directly transferred from the <a href="http://eigen.tuxfamily.org/">Eigen Library</a> to ViennaCL objects using the built-in convenience wrappers.</p>
<p>The first step is to include the necessary headers and activate the Eigen convenience functions in ViennaCL: </p>
<div class="fragment"><div class="line"><span class="comment">// System headers</span></div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Eigen headers</span></div>
<div class="line"><span class="preprocessor">#include <Eigen/Core></span></div>
<div class="line"><span class="preprocessor">#include <Eigen/Sparse></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// IMPORTANT: Must be set prior to any ViennaCL includes if you want to use ViennaCL algorithms on Eigen objects</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_EIGEN 1</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// ViennaCL includes</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="matrix_8hpp.html">viennacl/matrix.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"><span class="preprocessor">#include "<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// Helper functions for this tutorial:</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector-io_8hpp.html">vector-io.hpp</a>"</span></div>
</div><!-- fragment --><p> The following is a set of auxiliary dispatchers for obtaining the right Eigen types for a given floating point type. This is merely an implementation detail, so feel free to skip over it. </p>
<div class="fragment"><div class="line"><span class="comment">//dense matrix:</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> T></div>
<div class="line"><span class="keyword">struct </span>Eigen_dense_matrix</div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keyword">typename</span> T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_dense_matrix<float></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::MatrixXf type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_dense_matrix<double></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::MatrixXd type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">//sparse matrix</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> T></div>
<div class="line"><span class="keyword">struct </span>Eigen_vector</div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keyword">typename</span> T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_vector<float></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::VectorXf type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_vector<double></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::VectorXd type;</div>
<div class="line">};</div>
</div><!-- fragment --><p> The following function contains the main code for this tutorial. It consists of the following steps:</p>
<ul>
<li>Creates Eigen matrices and vectors</li>
<li>Initializes them with data</li>
<li>Create ViennaCL objects</li>
<li>Copy them over to the respective ViennaCL objects</li>
<li>Compute matrix-vector products in both Eigen and ViennaCL and compare results.</li>
</ul>
<div class="fragment"><div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> ScalarType></div>
<div class="line"><span class="keywordtype">void</span> run_tutorial()</div>
<div class="line">{</div>
</div><!-- fragment --><p> Get Eigen matrix and vector types for the provided ScalarType. Involves a little bit of template-metaprogramming. </p>
<div class="fragment"><div class="line"><span class="keyword">typedef</span> <span class="keyword">typename</span> Eigen_dense_matrix<ScalarType>::type EigenMatrix;</div>
<div class="line"><span class="keyword">typedef</span> <span class="keyword">typename</span> Eigen_vector<ScalarType>::type EigenVector;</div>
</div><!-- fragment --><p> Create and fill dense matrices from the Eigen library: </p>
<div class="fragment"><div class="line">EigenMatrix eigen_densemat(6, 5);</div>
<div class="line">EigenMatrix eigen_densemat2(6, 5);</div>
<div class="line">eigen_densemat(0,0) = 2.0; eigen_densemat(0,1) = -1.0;</div>
<div class="line">eigen_densemat(1,0) = -1.0; eigen_densemat(1,1) = 2.0; eigen_densemat(1,2) = -1.0;</div>
<div class="line">eigen_densemat(2,1) = -1.0; eigen_densemat(2,2) = -1.0; eigen_densemat(2,3) = -1.0;</div>
<div class="line">eigen_densemat(3,2) = -1.0; eigen_densemat(3,3) = 2.0; eigen_densemat(3,4) = -1.0;</div>
<div class="line"> eigen_densemat(5,4) = -1.0; eigen_densemat(4,4) = -1.0;</div>
<div class="line">Eigen::Map<EigenMatrix> eigen_densemat_map(eigen_densemat.data(), 6, 5); <span class="comment">// same as eigen_densemat, but emulating user-provided buffer</span></div>
</div><!-- fragment --><p> Create and fill sparse matrices from the Eigen library: </p>
<div class="fragment"><div class="line">Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat(6, 5);</div>
<div class="line">Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat2(6, 5);</div>
<div class="line">eigen_sparsemat.reserve(5*2);</div>
<div class="line">eigen_sparsemat.insert(0,0) = 2.0; eigen_sparsemat.insert(0,1) = -1.0;</div>
<div class="line">eigen_sparsemat.insert(1,1) = 2.0; eigen_sparsemat.insert(1,2) = -1.0;</div>
<div class="line">eigen_sparsemat.insert(2,2) = -1.0; eigen_sparsemat.insert(2,3) = -1.0;</div>
<div class="line">eigen_sparsemat.insert(3,3) = 2.0; eigen_sparsemat.insert(3,4) = -1.0;</div>
<div class="line">eigen_sparsemat.insert(5,4) = -1.0;</div>
<div class="line"><span class="comment">//eigen_sparsemat.endFill();</span></div>
</div><!-- fragment --><p> Create and fill a few vectors from the Eigen library: </p>
<div class="fragment"><div class="line">EigenVector eigen_rhs(5);</div>
<div class="line">Eigen::Map<EigenVector> eigen_rhs_map(eigen_rhs.data(), 5);</div>
<div class="line">EigenVector eigen_result(6);</div>
<div class="line">EigenVector eigen_temp(6);</div>
<div class="line"></div>
<div class="line">eigen_rhs(0) = 10.0;</div>
<div class="line">eigen_rhs(1) = 11.0;</div>
<div class="line">eigen_rhs(2) = 12.0;</div>
<div class="line">eigen_rhs(3) = 13.0;</div>
<div class="line">eigen_rhs(4) = 14.0;</div>
</div><!-- fragment --><p> Create the corresponding ViennaCL objects: </p>
<div class="fragment"><div class="line"><a name="_a0"></a><a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> vcl_rhs(5);</div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> vcl_result(6);</div>
<div class="line"><a name="_a1"></a><a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix<ScalarType></a> vcl_densemat(6, 5);</div>
<div class="line"><a name="_a2"></a><a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<ScalarType></a> vcl_sparsemat(6, 5);</div>
</div><!-- fragment --><p> Directly copy the Eigen objects to ViennaCL objects </p>
<div class="fragment"><div class="line"><a name="a3"></a><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(&(eigen_rhs[0]), &(eigen_rhs[0]) + 5, vcl_rhs.begin()); <span class="comment">// Method 1: via iterator interface (cf. std::copy())</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_rhs, vcl_rhs); <span class="comment">// Method 2: via built-in wrappers (convenience layer)</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_rhs_map, vcl_rhs); <span class="comment">// Same as method 2, but for a mapped vector</span></div>
<div class="line"></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_densemat, vcl_densemat);</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_densemat_map, vcl_densemat); <span class="comment">//same as above, using mapped matrix</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_sparsemat, vcl_sparsemat);</div>
<div class="line">std::cout << <span class="stringliteral">"VCL sparsematrix dimensions: "</span> << vcl_sparsemat.size1() << <span class="stringliteral">", "</span> << vcl_sparsemat.size2() << std::endl;</div>
<div class="line"></div>
<div class="line"><span class="comment">// For completeness: Copy matrices from ViennaCL back to Eigen:</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_densemat, eigen_densemat2);</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_sparsemat, eigen_sparsemat2);</div>
</div><!-- fragment --><p> Run dense matrix-vector products and compare results: </p>
<div class="fragment"><div class="line">eigen_result = eigen_densemat * eigen_rhs;</div>
<div class="line">vcl_result = <a name="a4"></a><a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(vcl_densemat, vcl_rhs);</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_result, eigen_temp);</div>
<div class="line">std::cout << <span class="stringliteral">"Difference for dense matrix-vector product: "</span> << (eigen_result - eigen_temp).norm() << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">"Difference for dense matrix-vector product (Eigen->ViennaCL->Eigen): "</span></div>
<div class="line"> << (eigen_densemat2 * eigen_rhs - eigen_temp).norm() << std::endl;</div>
</div><!-- fragment --><p> Run sparse matrix-vector products and compare results: </p>
<div class="fragment"><div class="line"> eigen_result = eigen_sparsemat * eigen_rhs;</div>
<div class="line"> vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(vcl_sparsemat, vcl_rhs);</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_result, eigen_temp);</div>
<div class="line"> std::cout << <span class="stringliteral">"Difference for sparse matrix-vector product: "</span> << (eigen_result - eigen_temp).norm() << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"Difference for sparse matrix-vector product (Eigen->ViennaCL->Eigen): "</span></div>
<div class="line"> << (eigen_sparsemat2 * eigen_rhs - eigen_temp).norm() << std::endl;</div>
<div class="line">}</div>
</div><!-- fragment --><p> In the <a class="el" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main()</a> routine we only call the worker function defined above with both single and double precision arithmetic. </p>
<div class="fragment"><div class="line"><span class="keywordtype">int</span> <a name="a5"></a><a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a>(<span class="keywordtype">int</span>, <span class="keywordtype">char</span> *[])</div>
<div class="line">{</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"## Single precision"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> run_tutorial<float>();</div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#ifdef VIENNACL_HAVE_OPENCL</span></div>
<div class="line"> <span class="keywordflow">if</span> ( <a name="a6"></a><a class="code" href="namespaceviennacl_1_1ocl.html#ac54d59a74deaccec81f64e738b856348">viennacl::ocl::current_device</a>().<a name="a7"></a><a class="code" href="classviennacl_1_1ocl_1_1device.html#a5182025f56c76de6770e26b14165b00a">double_support</a>() )</div>
<div class="line">#endif</div>
<div class="line"> {</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"## Double precision"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> run_tutorial<double>();</div>
<div class="line"> }</div>
</div><!-- fragment --><p> That's it. Print a success message and exit. </p>
<div class="fragment"><div class="line"> std::cout << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!"</span> << std::endl;</div>
<div class="line"> std::cout << std::endl;</div>
<div class="line"></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">// System headers</span></div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Eigen headers</span></div>
<div class="line"><span class="preprocessor">#include <Eigen/Core></span></div>
<div class="line"><span class="preprocessor">#include <Eigen/Sparse></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// IMPORTANT: Must be set prior to any ViennaCL includes if you want to use ViennaCL algorithms on Eigen objects</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_EIGEN 1</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// ViennaCL includes</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="matrix_8hpp.html">viennacl/matrix.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"><span class="preprocessor">#include "<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// Helper functions for this tutorial:</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector-io_8hpp.html">vector-io.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"><span class="comment">//dense matrix:</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> T></div>
<div class="line"><span class="keyword">struct </span>Eigen_dense_matrix</div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keyword">typename</span> T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_dense_matrix<float></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::MatrixXf type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_dense_matrix<double></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::MatrixXd type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">//sparse matrix</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> T></div>
<div class="line"><span class="keyword">struct </span>Eigen_vector</div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keyword">typename</span> T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_vector<float></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::VectorXf type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><></div>
<div class="line"><span class="keyword">struct </span>Eigen_vector<double></div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> Eigen::VectorXd type;</div>
<div class="line">};</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> ScalarType></div>
<div class="line"><span class="keywordtype">void</span> run_tutorial()</div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keyword">typename</span> Eigen_dense_matrix<ScalarType>::type EigenMatrix;</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keyword">typename</span> Eigen_vector<ScalarType>::type EigenVector;</div>
<div class="line"></div>
<div class="line"> EigenMatrix eigen_densemat(6, 5);</div>
<div class="line"> EigenMatrix eigen_densemat2(6, 5);</div>
<div class="line"> eigen_densemat(0,0) = 2.0; eigen_densemat(0,1) = -1.0;</div>
<div class="line"> eigen_densemat(1,0) = -1.0; eigen_densemat(1,1) = 2.0; eigen_densemat(1,2) = -1.0;</div>
<div class="line"> eigen_densemat(2,1) = -1.0; eigen_densemat(2,2) = -1.0; eigen_densemat(2,3) = -1.0;</div>
<div class="line"> eigen_densemat(3,2) = -1.0; eigen_densemat(3,3) = 2.0; eigen_densemat(3,4) = -1.0;</div>
<div class="line"> eigen_densemat(5,4) = -1.0; eigen_densemat(4,4) = -1.0;</div>
<div class="line"> Eigen::Map<EigenMatrix> eigen_densemat_map(eigen_densemat.data(), 6, 5); <span class="comment">// same as eigen_densemat, but emulating user-provided buffer</span></div>
<div class="line"></div>
<div class="line"> Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat(6, 5);</div>
<div class="line"> Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat2(6, 5);</div>
<div class="line"> eigen_sparsemat.reserve(5*2);</div>
<div class="line"> eigen_sparsemat.insert(0,0) = 2.0; eigen_sparsemat.insert(0,1) = -1.0;</div>
<div class="line"> eigen_sparsemat.insert(1,1) = 2.0; eigen_sparsemat.insert(1,2) = -1.0;</div>
<div class="line"> eigen_sparsemat.insert(2,2) = -1.0; eigen_sparsemat.insert(2,3) = -1.0;</div>
<div class="line"> eigen_sparsemat.insert(3,3) = 2.0; eigen_sparsemat.insert(3,4) = -1.0;</div>
<div class="line"> eigen_sparsemat.insert(5,4) = -1.0;</div>
<div class="line"> <span class="comment">//eigen_sparsemat.endFill();</span></div>
<div class="line"></div>
<div class="line"> EigenVector eigen_rhs(5);</div>
<div class="line"> Eigen::Map<EigenVector> eigen_rhs_map(eigen_rhs.data(), 5);</div>
<div class="line"> EigenVector eigen_result(6);</div>
<div class="line"> EigenVector eigen_temp(6);</div>
<div class="line"></div>
<div class="line"> eigen_rhs(0) = 10.0;</div>
<div class="line"> eigen_rhs(1) = 11.0;</div>
<div class="line"> eigen_rhs(2) = 12.0;</div>
<div class="line"> eigen_rhs(3) = 13.0;</div>
<div class="line"> eigen_rhs(4) = 14.0;</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> vcl_rhs(5);</div>
<div class="line"> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> vcl_result(6);</div>
<div class="line"> <a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix<ScalarType></a> vcl_densemat(6, 5);</div>
<div class="line"> <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<ScalarType></a> vcl_sparsemat(6, 5);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(&(eigen_rhs[0]), &(eigen_rhs[0]) + 5, vcl_rhs.begin()); <span class="comment">// Method 1: via iterator interface (cf. std::copy())</span></div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_rhs, vcl_rhs); <span class="comment">// Method 2: via built-in wrappers (convenience layer)</span></div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_rhs_map, vcl_rhs); <span class="comment">// Same as method 2, but for a mapped vector</span></div>
<div class="line"></div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_densemat, vcl_densemat);</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_densemat_map, vcl_densemat); <span class="comment">//same as above, using mapped matrix</span></div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(eigen_sparsemat, vcl_sparsemat);</div>
<div class="line"> std::cout << <span class="stringliteral">"VCL sparsematrix dimensions: "</span> << vcl_sparsemat.size1() << <span class="stringliteral">", "</span> << vcl_sparsemat.size2() << std::endl;</div>
<div class="line"></div>
<div class="line"> <span class="comment">// For completeness: Copy matrices from ViennaCL back to Eigen:</span></div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_densemat, eigen_densemat2);</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_sparsemat, eigen_sparsemat2);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"> eigen_result = eigen_densemat * eigen_rhs;</div>
<div class="line"> vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(vcl_densemat, vcl_rhs);</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_result, eigen_temp);</div>
<div class="line"> std::cout << <span class="stringliteral">"Difference for dense matrix-vector product: "</span> << (eigen_result - eigen_temp).norm() << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"Difference for dense matrix-vector product (Eigen->ViennaCL->Eigen): "</span></div>
<div class="line"> << (eigen_densemat2 * eigen_rhs - eigen_temp).norm() << std::endl;</div>
<div class="line"></div>
<div class="line"> eigen_result = eigen_sparsemat * eigen_rhs;</div>
<div class="line"> vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(vcl_sparsemat, vcl_rhs);</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_result, eigen_temp);</div>
<div class="line"> std::cout << <span class="stringliteral">"Difference for sparse matrix-vector product: "</span> << (eigen_result - eigen_temp).norm() << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"Difference for sparse matrix-vector product (Eigen->ViennaCL->Eigen): "</span></div>
<div class="line"> << (eigen_sparsemat2 * eigen_rhs - eigen_temp).norm() << std::endl;</div>
<div class="line">}</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>(<span class="keywordtype">int</span>, <span class="keywordtype">char</span> *[])</div>
<div class="line">{</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"## Single precision"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> run_tutorial<float>();</div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#ifdef VIENNACL_HAVE_OPENCL</span></div>
<div class="line"> <span class="keywordflow">if</span> ( <a class="code" href="namespaceviennacl_1_1ocl.html#ac54d59a74deaccec81f64e738b856348">viennacl::ocl::current_device</a>().<a class="code" href="classviennacl_1_1ocl_1_1device.html#a5182025f56c76de6770e26b14165b00a">double_support</a>() )</div>
<div class="line">#endif</div>
<div class="line"> {</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"## Double precision"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"----------------------------------------------"</span> << std::endl;</div>
<div class="line"> run_tutorial<double>();</div>
<div class="line"> }</div>
<div class="line"></div>
<div class="line"> std::cout << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!"</span> << std::endl;</div>
<div class="line"> std::cout << std::endl;</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>
|