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 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522
|
<!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: viennacl/linalg/sparse_matrix_operations.hpp Source File</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('sparse__matrix__operations_8hpp_source.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">sparse_matrix_operations.hpp</div> </div>
</div><!--header-->
<div class="contents">
<a href="sparse__matrix__operations_8hpp.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="preprocessor">#ifndef VIENNACL_LINALG_SPARSE_MATRIX_OPERATIONS_HPP_</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="preprocessor">#define VIENNACL_LINALG_SPARSE_MATRIX_OPERATIONS_HPP_</span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span> </div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment">/* =========================================================================</span></div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment"> Copyright (c) 2010-2016, Institute for Microelectronics,</span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment"> Institute for Analysis and Scientific Computing,</span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment"> TU Wien.</span></div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment"> Portions of this software are copyright by UChicago Argonne, LLC.</span></div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="comment"></span></div>
<div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="comment"> -----------------</span></div>
<div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="comment"> ViennaCL - The Vienna Computing Library</span></div>
<div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="comment"> -----------------</span></div>
<div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="comment"></span></div>
<div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="comment"> Project Head: Karl Rupp rupp@iue.tuwien.ac.at</span></div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span> <span class="comment"></span></div>
<div class="line"><a name="l00016"></a><span class="lineno"> 16</span> <span class="comment"> (A list of authors and contributors can be found in the manual)</span></div>
<div class="line"><a name="l00017"></a><span class="lineno"> 17</span> <span class="comment"></span></div>
<div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="comment"> License: MIT (X11), see file LICENSE in the base directory</span></div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="comment">============================================================================= */</span></div>
<div class="line"><a name="l00020"></a><span class="lineno"> 20</span> </div>
<div class="line"><a name="l00025"></a><span class="lineno"> 25</span> <span class="preprocessor">#include "<a class="code" href="forwards_8h.html">viennacl/forwards.h</a>"</span></div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span> <span class="preprocessor">#include "<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>"</span></div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span> <span class="preprocessor">#include "<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>"</span></div>
<div class="line"><a name="l00028"></a><span class="lineno"> 28</span> <span class="preprocessor">#include "<a class="code" href="matrix_8hpp.html">viennacl/matrix.hpp</a>"</span></div>
<div class="line"><a name="l00029"></a><span class="lineno"> 29</span> <span class="preprocessor">#include "<a class="code" href="tools_8hpp.html">viennacl/tools/tools.hpp</a>"</span></div>
<div class="line"><a name="l00030"></a><span class="lineno"> 30</span> <span class="preprocessor">#include "<a class="code" href="host__based_2sparse__matrix__operations_8hpp.html">viennacl/linalg/host_based/sparse_matrix_operations.hpp</a>"</span></div>
<div class="line"><a name="l00031"></a><span class="lineno"> 31</span> </div>
<div class="line"><a name="l00032"></a><span class="lineno"> 32</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span> <span class="preprocessor"> #include "<a class="code" href="opencl_2sparse__matrix__operations_8hpp.html">viennacl/linalg/opencl/sparse_matrix_operations.hpp</a>"</span></div>
<div class="line"><a name="l00034"></a><span class="lineno"> 34</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00035"></a><span class="lineno"> 35</span> </div>
<div class="line"><a name="l00036"></a><span class="lineno"> 36</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span> <span class="preprocessor"> #include "<a class="code" href="cuda_2sparse__matrix__operations_8hpp.html">viennacl/linalg/cuda/sparse_matrix_operations.hpp</a>"</span></div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00039"></a><span class="lineno"> 39</span> </div>
<div class="line"><a name="l00040"></a><span class="lineno"> 40</span> <span class="keyword">namespace </span>viennacl</div>
<div class="line"><a name="l00041"></a><span class="lineno"> 41</span> {</div>
<div class="line"><a name="l00042"></a><span class="lineno"> 42</span>  <span class="keyword">namespace </span>linalg</div>
<div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  {</div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span> </div>
<div class="line"><a name="l00045"></a><span class="lineno"> 45</span>  <span class="keyword">namespace </span>detail</div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span>  {</div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span> </div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SparseMatrixType, <span class="keyword">typename</span> SCALARTYPE, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> VEC_ALIGNMENT></div>
<div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a> >::type</div>
<div class="line"><a name="l00050"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#a5ef8943054b39077573f94234980bd0f"> 50</a></span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a5ef8943054b39077573f94234980bd0f">row_info</a>(SparseMatrixType <span class="keyword">const</span> & mat,</div>
<div class="line"><a name="l00051"></a><span class="lineno"> 51</span>  <a class="code" href="classviennacl_1_1vector.html">vector<SCALARTYPE, VEC_ALIGNMENT></a> & vec,</div>
<div class="line"><a name="l00052"></a><span class="lineno"> 52</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#adc2e9ed3e4d9fb0258840837aaa722ac">row_info_types</a> info_selector)</div>
<div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  {</div>
<div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(mat).get_active_handle_id())</div>
<div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  {</div>
<div class="line"><a name="l00056"></a><span class="lineno"> 56</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based_1_1detail.html#a58dfcdc80252eb2bc06c50a1ef533fc9">viennacl::linalg::host_based::detail::row_info</a>(mat, vec, info_selector);</div>
<div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl_1_1detail.html#a256a599dedd093618c9936e85a75b3d4">viennacl::linalg::opencl::detail::row_info</a>(mat, vec, info_selector);</div>
<div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00063"></a><span class="lineno"> 63</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda_1_1detail.html#a45615bf3f3314b66269070ae06bc3386">viennacl::linalg::cuda::detail::row_info</a>(mat, vec, info_selector);</div>
<div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00068"></a><span class="lineno"> 68</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00071"></a><span class="lineno"> 71</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  }</div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  }</div>
<div class="line"><a name="l00075"></a><span class="lineno"> 75</span> </div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  }</div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span> </div>
<div class="line"><a name="l00078"></a><span class="lineno"> 78</span> </div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span> </div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  <span class="comment">// A * x</span></div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span> </div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SparseMatrixType, <span class="keyword">class</span> ScalarType></div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>>::type</div>
<div class="line"><a name="l00092"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#ac4ecc092e443847d1968d9d2cef07eda"> 92</a></span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">prod_impl</a>(<span class="keyword">const</span> SparseMatrixType & mat,</div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<ScalarType></a> & vec,</div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a> alpha,</div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<ScalarType></a> & result,</div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a> beta)</div>
<div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  {</div>
<div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  assert( (mat.size1() == result.<a class="code" href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">size</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for compressed matrix-vector product: size1(mat) != size(result)"</span>));</div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  assert( (mat.size2() == vec.<a class="code" href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">size</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for compressed matrix-vector product: size2(mat) != size(x)"</span>));</div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span> </div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(mat).<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#a5951aceb881068a77b4af8825b547ea3">get_active_handle_id</a>())</div>
<div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  {</div>
<div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based.html#a73ab507a304067a153e0ad741ed6a40d">viennacl::linalg::host_based::prod_impl</a>(mat, vec, alpha, result, beta);</div>
<div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00106"></a><span class="lineno"> 106</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl.html#a18a93d38ca1eff472149093974edb5cb">viennacl::linalg::opencl::prod_impl</a>(mat, vec, alpha, result, beta);</div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00110"></a><span class="lineno"> 110</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#afd7583ae441f66185764ae2fed763feb">viennacl::linalg::cuda::prod_impl</a>(mat, vec, alpha, result, beta);</div>
<div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  }</div>
<div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  }</div>
<div class="line"><a name="l00122"></a><span class="lineno"> 122</span> </div>
<div class="line"><a name="l00123"></a><span class="lineno"> 123</span> </div>
<div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  <span class="comment">// A * B</span></div>
<div class="line"><a name="l00133"></a><span class="lineno"> 133</span> <span class="comment"></span> <span class="keyword">template</span><<span class="keyword">typename</span> SparseMatrixType, <span class="keyword">class</span> ScalarType></div>
<div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>>::type</div>
<div class="line"><a name="l00135"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#ac6448c1d1c0cdb1010a146c0238c8d20"> 135</a></span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">prod_impl</a>(<span class="keyword">const</span> SparseMatrixType & sp_mat,</div>
<div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1matrix__base.html">viennacl::matrix_base<ScalarType></a> & d_mat,</div>
<div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  <a class="code" href="classviennacl_1_1matrix__base.html">viennacl::matrix_base<ScalarType></a> & result)</div>
<div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  {</div>
<div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  assert( (sp_mat.size1() == result.<a class="code" href="classviennacl_1_1matrix__base.html#a24e4f5fa27a1af5ad47e52a97c065d68">size1</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for compressed matrix - dense matrix product: size1(sp_mat) != size1(result)"</span>));</div>
<div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  assert( (sp_mat.size2() == d_mat.<a class="code" href="classviennacl_1_1matrix__base.html#a24e4f5fa27a1af5ad47e52a97c065d68">size1</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for compressed matrix - dense matrix product: size2(sp_mat) != size1(d_mat)"</span>));</div>
<div class="line"><a name="l00141"></a><span class="lineno"> 141</span> </div>
<div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(sp_mat).<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#a5951aceb881068a77b4af8825b547ea3">get_active_handle_id</a>())</div>
<div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  {</div>
<div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based.html#a73ab507a304067a153e0ad741ed6a40d">viennacl::linalg::host_based::prod_impl</a>(sp_mat, d_mat, result);</div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl.html#a18a93d38ca1eff472149093974edb5cb">viennacl::linalg::opencl::prod_impl</a>(sp_mat, d_mat, result);</div>
<div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#afd7583ae441f66185764ae2fed763feb">viennacl::linalg::cuda::prod_impl</a>(sp_mat, d_mat, result);</div>
<div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  }</div>
<div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  }</div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span> </div>
<div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  <span class="comment">// A * transpose(B)</span></div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span> <span class="comment"></span> <span class="keyword">template</span><<span class="keyword">typename</span> SparseMatrixType, <span class="keyword">class</span> ScalarType></div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>>::type</div>
<div class="line"><a name="l00175"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a07ddf85f2f61d198160ffd13bf7aebac"> 175</a></span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">prod_impl</a>(<span class="keyword">const</span> SparseMatrixType & sp_mat,</div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1matrix__expression.html">viennacl::matrix_expression</a><<span class="keyword">const</span> <a class="code" href="classviennacl_1_1matrix__base.html">viennacl::matrix_base<ScalarType></a>,</div>
<div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1matrix__base.html">viennacl::matrix_base<ScalarType></a>,</div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  <a class="code" href="structviennacl_1_1op__trans.html">viennacl::op_trans</a>>& d_mat,</div>
<div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  <a class="code" href="classviennacl_1_1matrix__base.html">viennacl::matrix_base<ScalarType></a> & result)</div>
<div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  {</div>
<div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  assert( (sp_mat.size1() == result.<a class="code" href="classviennacl_1_1matrix__base.html#a24e4f5fa27a1af5ad47e52a97c065d68">size1</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for compressed matrix - dense matrix product: size1(sp_mat) != size1(result)"</span>));</div>
<div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  assert( (sp_mat.size2() == d_mat.size1()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for compressed matrix - dense matrix product: size2(sp_mat) != size1(d_mat)"</span>));</div>
<div class="line"><a name="l00183"></a><span class="lineno"> 183</span> </div>
<div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(sp_mat).<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#a5951aceb881068a77b4af8825b547ea3">get_active_handle_id</a>())</div>
<div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  {</div>
<div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based.html#a73ab507a304067a153e0ad741ed6a40d">viennacl::linalg::host_based::prod_impl</a>(sp_mat, d_mat, result);</div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00189"></a><span class="lineno"> 189</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00191"></a><span class="lineno"> 191</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl.html#a18a93d38ca1eff472149093974edb5cb">viennacl::linalg::opencl::prod_impl</a>(sp_mat, d_mat, result);</div>
<div class="line"><a name="l00192"></a><span class="lineno"> 192</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00193"></a><span class="lineno"> 193</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00194"></a><span class="lineno"> 194</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#afd7583ae441f66185764ae2fed763feb">viennacl::linalg::cuda::prod_impl</a>(sp_mat, d_mat, result);</div>
<div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00198"></a><span class="lineno"> 198</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00199"></a><span class="lineno"> 199</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  }</div>
<div class="line"><a name="l00204"></a><span class="lineno"> 204</span>  }</div>
<div class="line"><a name="l00205"></a><span class="lineno"> 205</span> </div>
<div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  <span class="comment">// A * B with both A and B sparse</span></div>
<div class="line"><a name="l00207"></a><span class="lineno"> 207</span> </div>
<div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  <span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  <span class="keywordtype">void</span></div>
<div class="line"><a name="l00219"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#ac371515668e3ef88cd49a4adba597e94"> 219</a></span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">prod_impl</a>(<span class="keyword">const</span> <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<NumericT></a> & A,</div>
<div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<NumericT></a> & B,</div>
<div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<NumericT></a> & C)</div>
<div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  {</div>
<div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  assert( (A.<a class="code" href="classviennacl_1_1compressed__matrix.html#a4fc12fc4abfef4a1426575a2d73f18ab">size2</a>() == B.<a class="code" href="classviennacl_1_1compressed__matrix.html#a463cf1739f9cdd387aa185cb574db183">size1</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for sparse matrix-matrix product: size2(A) != size1(B)"</span>));</div>
<div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  assert( (C.<a class="code" href="classviennacl_1_1compressed__matrix.html#a463cf1739f9cdd387aa185cb574db183">size1</a>() == 0 || C.<a class="code" href="classviennacl_1_1compressed__matrix.html#a463cf1739f9cdd387aa185cb574db183">size1</a>() == A.<a class="code" href="classviennacl_1_1compressed__matrix.html#a463cf1739f9cdd387aa185cb574db183">size1</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for sparse matrix-matrix product: size1(A) != size1(C)"</span>));</div>
<div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  assert( (C.<a class="code" href="classviennacl_1_1compressed__matrix.html#a4fc12fc4abfef4a1426575a2d73f18ab">size2</a>() == 0 || C.<a class="code" href="classviennacl_1_1compressed__matrix.html#a4fc12fc4abfef4a1426575a2d73f18ab">size2</a>() == B.<a class="code" href="classviennacl_1_1compressed__matrix.html#a4fc12fc4abfef4a1426575a2d73f18ab">size2</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for sparse matrix-matrix product: size2(B) != size2(B)"</span>));</div>
<div class="line"><a name="l00226"></a><span class="lineno"> 226</span> </div>
<div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(A).<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#a5951aceb881068a77b4af8825b547ea3">get_active_handle_id</a>())</div>
<div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  {</div>
<div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based.html#a73ab507a304067a153e0ad741ed6a40d">viennacl::linalg::host_based::prod_impl</a>(A, B, C);</div>
<div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00232"></a><span class="lineno"> 232</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl.html#a18a93d38ca1eff472149093974edb5cb">viennacl::linalg::opencl::prod_impl</a>(A, B, C);</div>
<div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00236"></a><span class="lineno"> 236</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00237"></a><span class="lineno"> 237</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#afd7583ae441f66185764ae2fed763feb">viennacl::linalg::cuda::prod_impl</a>(A, B, C);</div>
<div class="line"><a name="l00240"></a><span class="lineno"> 240</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00241"></a><span class="lineno"> 241</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  }</div>
<div class="line"><a name="l00247"></a><span class="lineno"> 247</span>  }</div>
<div class="line"><a name="l00248"></a><span class="lineno"> 248</span> </div>
<div class="line"><a name="l00249"></a><span class="lineno"> 249</span> </div>
<div class="line"><a name="l00256"></a><span class="lineno"> 256</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SparseMatrixType, <span class="keyword">class</span> ScalarType, <span class="keyword">typename</span> SOLVERTAG></div>
<div class="line"><a name="l00257"></a><span class="lineno"> 257</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>>::type</div>
<div class="line"><a name="l00258"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a2cd904fcd01a99005ca9f83f3c696b41"> 258</a></span>  <a class="code" href="namespaceviennacl_1_1linalg.html#adf06f61be9418e5817eed66004f9dd2b">inplace_solve</a>(<span class="keyword">const</span> SparseMatrixType & mat,</div>
<div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<ScalarType></a> & vec,</div>
<div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  SOLVERTAG tag)</div>
<div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  {</div>
<div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  assert( (mat.size1() == mat.size2()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for triangular solve on compressed matrix: size1(mat) != size2(mat)"</span>));</div>
<div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  assert( (mat.size2() == vec.<a class="code" href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">size</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for compressed matrix-vector product: size2(mat) != size(x)"</span>));</div>
<div class="line"><a name="l00264"></a><span class="lineno"> 264</span> </div>
<div class="line"><a name="l00265"></a><span class="lineno"> 265</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(mat).<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#a5951aceb881068a77b4af8825b547ea3">get_active_handle_id</a>())</div>
<div class="line"><a name="l00266"></a><span class="lineno"> 266</span>  {</div>
<div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based.html#af83ff307c97a945a568c30253c1a7c61">viennacl::linalg::host_based::inplace_solve</a>(mat, vec, tag);</div>
<div class="line"><a name="l00269"></a><span class="lineno"> 269</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00270"></a><span class="lineno"> 270</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl.html#ac90a3dca0595b15ae04662bddb6cf57c">viennacl::linalg::opencl::inplace_solve</a>(mat, vec, tag);</div>
<div class="line"><a name="l00273"></a><span class="lineno"> 273</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00274"></a><span class="lineno"> 274</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00275"></a><span class="lineno"> 275</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00276"></a><span class="lineno"> 276</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00277"></a><span class="lineno"> 277</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#acb0b34e3ac22c72081eceaabeab80f22">viennacl::linalg::cuda::inplace_solve</a>(mat, vec, tag);</div>
<div class="line"><a name="l00278"></a><span class="lineno"> 278</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00279"></a><span class="lineno"> 279</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00280"></a><span class="lineno"> 280</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00282"></a><span class="lineno"> 282</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00283"></a><span class="lineno"> 283</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  }</div>
<div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  }</div>
<div class="line"><a name="l00286"></a><span class="lineno"> 286</span> </div>
<div class="line"><a name="l00287"></a><span class="lineno"> 287</span> </div>
<div class="line"><a name="l00294"></a><span class="lineno"> 294</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SparseMatrixType, <span class="keyword">class</span> ScalarType, <span class="keyword">typename</span> SOLVERTAG></div>
<div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>>::type</div>
<div class="line"><a name="l00296"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a02babae1af157516bd1e92d28c711957"> 296</a></span>  <a class="code" href="namespaceviennacl_1_1linalg.html#adf06f61be9418e5817eed66004f9dd2b">inplace_solve</a>(<span class="keyword">const</span> <a class="code" href="classviennacl_1_1matrix__expression.html">matrix_expression<const SparseMatrixType, const SparseMatrixType, op_trans></a> & mat,</div>
<div class="line"><a name="l00297"></a><span class="lineno"> 297</span>  <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<ScalarType></a> & vec,</div>
<div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  SOLVERTAG tag)</div>
<div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  {</div>
<div class="line"><a name="l00300"></a><span class="lineno"> 300</span>  assert( (mat.<a class="code" href="classviennacl_1_1matrix__expression.html#a3125cc5f209063de6a99dc67e19d42bb">size1</a>() == mat.<a class="code" href="classviennacl_1_1matrix__expression.html#ad8516f3a2ef4db4d5379fb64832fd39f">size2</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for triangular solve on transposed compressed matrix: size1(mat) != size2(mat)"</span>));</div>
<div class="line"><a name="l00301"></a><span class="lineno"> 301</span>  assert( (mat.<a class="code" href="classviennacl_1_1matrix__expression.html#a3125cc5f209063de6a99dc67e19d42bb">size1</a>() == vec.<a class="code" href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">size</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for transposed compressed matrix triangular solve: size1(mat) != size(x)"</span>));</div>
<div class="line"><a name="l00302"></a><span class="lineno"> 302</span> </div>
<div class="line"><a name="l00303"></a><span class="lineno"> 303</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(mat.<a class="code" href="classviennacl_1_1matrix__expression.html#aa2c582e8d2c1cafbf2ffdf05c9fdfaca">lhs</a>()).get_active_handle_id())</div>
<div class="line"><a name="l00304"></a><span class="lineno"> 304</span>  {</div>
<div class="line"><a name="l00305"></a><span class="lineno"> 305</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00306"></a><span class="lineno"> 306</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based.html#af83ff307c97a945a568c30253c1a7c61">viennacl::linalg::host_based::inplace_solve</a>(mat, vec, tag);</div>
<div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00308"></a><span class="lineno"> 308</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00309"></a><span class="lineno"> 309</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl.html#ac90a3dca0595b15ae04662bddb6cf57c">viennacl::linalg::opencl::inplace_solve</a>(mat, vec, tag);</div>
<div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00312"></a><span class="lineno"> 312</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00313"></a><span class="lineno"> 313</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00314"></a><span class="lineno"> 314</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#acb0b34e3ac22c72081eceaabeab80f22">viennacl::linalg::cuda::inplace_solve</a>(mat, vec, tag);</div>
<div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00317"></a><span class="lineno"> 317</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  }</div>
<div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  }</div>
<div class="line"><a name="l00324"></a><span class="lineno"> 324</span> </div>
<div class="line"><a name="l00325"></a><span class="lineno"> 325</span> </div>
<div class="line"><a name="l00326"></a><span class="lineno"> 326</span> </div>
<div class="line"><a name="l00327"></a><span class="lineno"> 327</span>  <span class="keyword">namespace </span>detail</div>
<div class="line"><a name="l00328"></a><span class="lineno"> 328</span>  {</div>
<div class="line"><a name="l00329"></a><span class="lineno"> 329</span> </div>
<div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SparseMatrixType, <span class="keyword">class</span> ScalarType, <span class="keyword">typename</span> SOLVERTAG></div>
<div class="line"><a name="l00331"></a><span class="lineno"> 331</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>>::type</div>
<div class="line"><a name="l00332"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#a926494a1bd9b3c66def0c860bf564b42"> 332</a></span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a926494a1bd9b3c66def0c860bf564b42">block_inplace_solve</a>(<span class="keyword">const</span> <a class="code" href="classviennacl_1_1matrix__expression.html">matrix_expression<const SparseMatrixType, const SparseMatrixType, op_trans></a> & mat,</div>
<div class="line"><a name="l00333"></a><span class="lineno"> 333</span>  <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">viennacl::backend::mem_handle</a> <span class="keyword">const</span> & block_index_array, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> num_blocks,</div>
<div class="line"><a name="l00334"></a><span class="lineno"> 334</span>  <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<ScalarType></a> <span class="keyword">const</span> & mat_diagonal,</div>
<div class="line"><a name="l00335"></a><span class="lineno"> 335</span>  <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<ScalarType></a> & vec,</div>
<div class="line"><a name="l00336"></a><span class="lineno"> 336</span>  SOLVERTAG tag)</div>
<div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  {</div>
<div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  assert( (mat.<a class="code" href="classviennacl_1_1matrix__expression.html#a3125cc5f209063de6a99dc67e19d42bb">size1</a>() == mat.<a class="code" href="classviennacl_1_1matrix__expression.html#ad8516f3a2ef4db4d5379fb64832fd39f">size2</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for triangular solve on transposed compressed matrix: size1(mat) != size2(mat)"</span>));</div>
<div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  assert( (mat.<a class="code" href="classviennacl_1_1matrix__expression.html#a3125cc5f209063de6a99dc67e19d42bb">size1</a>() == vec.<a class="code" href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">size</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size check failed for transposed compressed matrix triangular solve: size1(mat) != size(x)"</span>));</div>
<div class="line"><a name="l00340"></a><span class="lineno"> 340</span> </div>
<div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  <span class="keywordflow">switch</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(mat.<a class="code" href="classviennacl_1_1matrix__expression.html#aa2c582e8d2c1cafbf2ffdf05c9fdfaca">lhs</a>()).get_active_handle_id())</div>
<div class="line"><a name="l00342"></a><span class="lineno"> 342</span>  {</div>
<div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>:</div>
<div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1host__based_1_1detail.html#abc856cc45c42e8809e7e867904a4a7d7">viennacl::linalg::host_based::detail::block_inplace_solve</a>(mat, block_index_array, num_blocks, mat_diagonal, vec, tag);</div>
<div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00346"></a><span class="lineno"> 346</span> <span class="preprocessor"> #ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>:</div>
<div class="line"><a name="l00348"></a><span class="lineno"> 348</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1opencl_1_1detail.html#a144fd7dc25be2029356c24f07f281edc">viennacl::linalg::opencl::detail::block_inplace_solve</a>(mat, block_index_array, num_blocks, mat_diagonal, vec, tag);</div>
<div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00350"></a><span class="lineno"> 350</span> <span class="preprocessor"> #endif</span></div>
<div class="line"><a name="l00351"></a><span class="lineno"> 351</span> <span class="preprocessor"> #ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>:</div>
<div class="line"><a name="l00353"></a><span class="lineno"> 353</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda_1_1detail.html#ac588666828831dc470e4ed55da6d9a10">viennacl::linalg::cuda::detail::block_inplace_solve</a>(mat, block_index_array, num_blocks, mat_diagonal, vec, tag);</div>
<div class="line"><a name="l00354"></a><span class="lineno"> 354</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00355"></a><span class="lineno"> 355</span> <span class="preprocessor"> #endif</span></div>
<div class="line"><a name="l00356"></a><span class="lineno"> 356</span>  <span class="keywordflow">case</span> <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a>:</div>
<div class="line"><a name="l00357"></a><span class="lineno"> 357</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not initialised!"</span>);</div>
<div class="line"><a name="l00358"></a><span class="lineno"> 358</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00359"></a><span class="lineno"> 359</span>  <span class="keywordflow">throw</span> <a class="code" href="classviennacl_1_1memory__exception.html">memory_exception</a>(<span class="stringliteral">"not implemented"</span>);</div>
<div class="line"><a name="l00360"></a><span class="lineno"> 360</span>  }</div>
<div class="line"><a name="l00361"></a><span class="lineno"> 361</span>  }</div>
<div class="line"><a name="l00362"></a><span class="lineno"> 362</span> </div>
<div class="line"><a name="l00363"></a><span class="lineno"> 363</span> </div>
<div class="line"><a name="l00364"></a><span class="lineno"> 364</span>  }</div>
<div class="line"><a name="l00365"></a><span class="lineno"> 365</span> </div>
<div class="line"><a name="l00366"></a><span class="lineno"> 366</span> </div>
<div class="line"><a name="l00367"></a><span class="lineno"> 367</span> </div>
<div class="line"><a name="l00368"></a><span class="lineno"> 368</span>  } <span class="comment">//namespace linalg</span></div>
<div class="line"><a name="l00369"></a><span class="lineno"> 369</span> </div>
<div class="line"><a name="l00370"></a><span class="lineno"> 370</span> </div>
<div class="line"><a name="l00372"></a><span class="lineno"> 372</span>  <span class="keyword">template</span><<span class="keyword">typename</span> M1></div>
<div class="line"><a name="l00373"></a><span class="lineno"> 373</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if<viennacl::is_any_sparse_matrix<M1>::value</a>,</div>
<div class="line"><a name="l00374"></a><span class="lineno"> 374</span>  matrix_expression< const M1, const M1, op_trans></div>
<div class="line"><a name="l00375"></a><span class="lineno"> 375</span>  >::type</div>
<div class="line"><a name="l00376"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#a2ee5dd77d41040e0a937a60346475b84"> 376</a></span>  <a class="code" href="namespaceviennacl.html#a2ee5dd77d41040e0a937a60346475b84">trans</a>(<span class="keyword">const</span> M1 & mat)</div>
<div class="line"><a name="l00377"></a><span class="lineno"> 377</span>  {</div>
<div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  <span class="keywordflow">return</span> <a class="code" href="classviennacl_1_1matrix__expression.html">matrix_expression< const M1, const M1, op_trans></a>(mat, mat);</div>
<div class="line"><a name="l00379"></a><span class="lineno"> 379</span>  }</div>
<div class="line"><a name="l00380"></a><span class="lineno"> 380</span> </div>
<div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  <span class="comment">//free functions:</span></div>
<div class="line"><a name="l00387"></a><span class="lineno"> 387</span> <span class="comment"></span> <span class="keyword">template</span><<span class="keyword">typename</span> SCALARTYPE, <span class="keyword">typename</span> SparseMatrixType></div>
<div class="line"><a name="l00388"></a><span class="lineno"> 388</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>,</div>
<div class="line"><a name="l00389"></a><span class="lineno"> 389</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<SCALARTYPE></a> >::type</div>
<div class="line"><a name="l00390"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#ab85e3912a40d5d479835e7608808f07e"> 390</a></span>  <a class="code" href="namespaceviennacl.html#a4d3a8290e94b15b653b3a0f0a3f80496">operator+</a>(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<SCALARTYPE></a> & result,</div>
<div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1vector__expression.html">viennacl::vector_expression</a>< <span class="keyword">const</span> SparseMatrixType, <span class="keyword">const</span> <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<SCALARTYPE></a>, <a class="code" href="structviennacl_1_1op__prod.html">viennacl::op_prod</a>> & proxy)</div>
<div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  {</div>
<div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  assert(proxy.lhs().size1() == result.<a class="code" href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">size</a>() && bool(<span class="stringliteral">"Dimensions for addition of sparse matrix-vector product to vector don't match!"</span>));</div>
<div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  <a class="code" href="classviennacl_1_1vector.html">vector<SCALARTYPE></a> temp(proxy.lhs().size1());</div>
<div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(proxy.lhs(), proxy.rhs(), temp);</div>
<div class="line"><a name="l00396"></a><span class="lineno"> 396</span>  result += temp;</div>
<div class="line"><a name="l00397"></a><span class="lineno"> 397</span>  <span class="keywordflow">return</span> result;</div>
<div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  }</div>
<div class="line"><a name="l00399"></a><span class="lineno"> 399</span> </div>
<div class="line"><a name="l00405"></a><span class="lineno"> 405</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SCALARTYPE, <span class="keyword">typename</span> SparseMatrixType></div>
<div class="line"><a name="l00406"></a><span class="lineno"> 406</span>  <span class="keyword">typename</span> <a class="code" href="structviennacl_1_1enable__if.html">viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value</a>,</div>
<div class="line"><a name="l00407"></a><span class="lineno"> 407</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<SCALARTYPE></a> >::type</div>
<div class="line"><a name="l00408"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#a5309df539f189a889f45f5e25098930f"> 408</a></span>  <a class="code" href="namespaceviennacl.html#aeb3693cfc246adc889875c208ca23aff">operator-</a>(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<SCALARTYPE></a> & result,</div>
<div class="line"><a name="l00409"></a><span class="lineno"> 409</span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1vector__expression.html">viennacl::vector_expression</a>< <span class="keyword">const</span> SparseMatrixType, <span class="keyword">const</span> <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<SCALARTYPE></a>, <a class="code" href="structviennacl_1_1op__prod.html">viennacl::op_prod</a>> & proxy)</div>
<div class="line"><a name="l00410"></a><span class="lineno"> 410</span>  {</div>
<div class="line"><a name="l00411"></a><span class="lineno"> 411</span>  assert(proxy.lhs().size1() == result.<a class="code" href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">size</a>() && bool(<span class="stringliteral">"Dimensions for addition of sparse matrix-vector product to vector don't match!"</span>));</div>
<div class="line"><a name="l00412"></a><span class="lineno"> 412</span>  <a class="code" href="classviennacl_1_1vector.html">vector<SCALARTYPE></a> temp(proxy.lhs().size1());</div>
<div class="line"><a name="l00413"></a><span class="lineno"> 413</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(proxy.lhs(), proxy.rhs(), temp);</div>
<div class="line"><a name="l00414"></a><span class="lineno"> 414</span>  result += temp;</div>
<div class="line"><a name="l00415"></a><span class="lineno"> 415</span>  <span class="keywordflow">return</span> result;</div>
<div class="line"><a name="l00416"></a><span class="lineno"> 416</span>  }</div>
<div class="line"><a name="l00417"></a><span class="lineno"> 417</span> </div>
<div class="line"><a name="l00418"></a><span class="lineno"> 418</span> } <span class="comment">//namespace viennacl</span></div>
<div class="line"><a name="l00419"></a><span class="lineno"> 419</span> </div>
<div class="line"><a name="l00420"></a><span class="lineno"> 420</span> </div>
<div class="line"><a name="l00421"></a><span class="lineno"> 421</span> <span class="preprocessor">#endif</span></div>
<div class="ttc" id="classviennacl_1_1compressed__matrix_html_a4fc12fc4abfef4a1426575a2d73f18ab"><div class="ttname"><a href="classviennacl_1_1compressed__matrix.html#a4fc12fc4abfef4a1426575a2d73f18ab">viennacl::compressed_matrix::size2</a></div><div class="ttdeci">const vcl_size_t & size2() const </div><div class="ttdoc">Returns the number of columns. </div><div class="ttdef"><b>Definition:</b> <a href="compressed__matrix_8hpp_source.html#l00929">compressed_matrix.hpp:929</a></div></div>
<div class="ttc" id="structviennacl_1_1enable__if_html"><div class="ttname"><a href="structviennacl_1_1enable__if.html">viennacl::enable_if</a></div><div class="ttdoc">Simple enable-if variant that uses the SFINAE pattern. </div><div class="ttdef"><b>Definition:</b> <a href="enable__if_8hpp_source.html#l00030">enable_if.hpp:30</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1cuda_html_acb0b34e3ac22c72081eceaabeab80f22"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1cuda.html#acb0b34e3ac22c72081eceaabeab80f22">viennacl::linalg::cuda::inplace_solve</a></div><div class="ttdeci">void inplace_solve(matrix_base< NumericT > const &A, matrix_base< NumericT > &B, SolverTagT tag)</div><div class="ttdoc">Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...</div><div class="ttdef"><b>Definition:</b> <a href="cuda_2direct__solve_8hpp_source.html#l00253">direct_solve.hpp:253</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_adf06f61be9418e5817eed66004f9dd2b"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#adf06f61be9418e5817eed66004f9dd2b">viennacl::linalg::inplace_solve</a></div><div class="ttdeci">void inplace_solve(const matrix_base< NumericT > &A, matrix_base< NumericT > &B, SolverTagT)</div><div class="ttdoc">Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...</div><div class="ttdef"><b>Definition:</b> <a href="direct__solve_8hpp_source.html#l00217">direct_solve.hpp:217</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a2ee5dd77d41040e0a937a60346475b84"><div class="ttname"><a href="namespaceviennacl.html#a2ee5dd77d41040e0a937a60346475b84">viennacl::trans</a></div><div class="ttdeci">viennacl::enable_if< viennacl::is_any_sparse_matrix< M1 >::value, matrix_expression< const M1, const M1, op_trans > >::type trans(const M1 &mat)</div><div class="ttdoc">Returns an expression template class representing a transposed matrix. </div><div class="ttdef"><b>Definition:</b> <a href="sparse__matrix__operations_8hpp_source.html#l00376">sparse_matrix_operations.hpp:376</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_adc2e9ed3e4d9fb0258840837aaa722ac"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#adc2e9ed3e4d9fb0258840837aaa722ac">viennacl::linalg::detail::row_info_types</a></div><div class="ttdeci">row_info_types</div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00837">forwards.h:837</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1opencl_html_a18a93d38ca1eff472149093974edb5cb"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1opencl.html#a18a93d38ca1eff472149093974edb5cb">viennacl::linalg::opencl::prod_impl</a></div><div class="ttdeci">void prod_impl(const matrix_base< NumericT > &mat, bool trans_A, const vector_base< NumericT > &vec, vector_base< NumericT > &result)</div><div class="ttdoc">Carries out matrix-vector multiplication. </div><div class="ttdef"><b>Definition:</b> <a href="opencl_2matrix__operations_8hpp_source.html#l00620">matrix_operations.hpp:620</a></div></div>
<div class="ttc" id="classviennacl_1_1memory__exception_html"><div class="ttname"><a href="classviennacl_1_1memory__exception.html">viennacl::memory_exception</a></div><div class="ttdoc">Exception class in case of memory errors. </div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00572">forwards.h:572</a></div></div>
<div class="ttc" id="classviennacl_1_1compressed__matrix_html_a463cf1739f9cdd387aa185cb574db183"><div class="ttname"><a href="classviennacl_1_1compressed__matrix.html#a463cf1739f9cdd387aa185cb574db183">viennacl::compressed_matrix::size1</a></div><div class="ttdeci">const vcl_size_t & size1() const </div><div class="ttdoc">Returns the number of rows. </div><div class="ttdef"><b>Definition:</b> <a href="compressed__matrix_8hpp_source.html#l00927">compressed_matrix.hpp:927</a></div></div>
<div class="ttc" id="matrix_8hpp_html"><div class="ttname"><a href="matrix_8hpp.html">matrix.hpp</a></div><div class="ttdoc">Implementation of the dense matrix class. </div></div>
<div class="ttc" id="tools_8hpp_html"><div class="ttname"><a href="tools_8hpp.html">tools.hpp</a></div><div class="ttdoc">Various little tools used here and there in ViennaCL. </div></div>
<div class="ttc" id="host__based_2sparse__matrix__operations_8hpp_html"><div class="ttname"><a href="host__based_2sparse__matrix__operations_8hpp.html">sparse_matrix_operations.hpp</a></div><div class="ttdoc">Implementations of operations using sparse matrices on the CPU using a single thread or OpenMP...</div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1host__based_html_af83ff307c97a945a568c30253c1a7c61"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1host__based.html#af83ff307c97a945a568c30253c1a7c61">viennacl::linalg::host_based::inplace_solve</a></div><div class="ttdeci">void inplace_solve(matrix_base< NumericT > const &A, matrix_base< NumericT > &B, SolverTagT)</div><div class="ttdoc">Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...</div><div class="ttdef"><b>Definition:</b> <a href="host__based_2direct__solve_8hpp_source.html#l00130">direct_solve.hpp:130</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1"><div class="ttname"><a href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a></div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00349">forwards.h:349</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__base_html"><div class="ttname"><a href="classviennacl_1_1matrix__base.html">viennacl::matrix_base</a></div><div class="ttdef"><b>Definition:</b> <a href="matrix__def_8hpp_source.html#l00103">matrix_def.hpp:103</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__expression_html"><div class="ttname"><a href="classviennacl_1_1matrix__expression.html">viennacl::matrix_expression</a></div><div class="ttdoc">Expression template class for representing a tree of expressions which ultimately result in a matrix...</div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00341">forwards.h:341</a></div></div>
<div class="ttc" id="forwards_8h_html"><div class="ttname"><a href="forwards_8h.html">forwards.h</a></div><div class="ttdoc">This file provides the forward declarations for the main types used within ViennaCL. </div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1cuda_1_1detail_html_a45615bf3f3314b66269070ae06bc3386"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1cuda_1_1detail.html#a45615bf3f3314b66269070ae06bc3386">viennacl::linalg::cuda::detail::row_info</a></div><div class="ttdeci">void row_info(compressed_matrix< NumericT, AligmentV > const &mat, vector_base< NumericT > &vec, viennacl::linalg::detail::row_info_types info_selector)</div><div class="ttdef"><b>Definition:</b> <a href="cuda_2sparse__matrix__operations_8hpp_source.html#l00107">sparse_matrix_operations.hpp:107</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1cuda_html_afd7583ae441f66185764ae2fed763feb"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1cuda.html#afd7583ae441f66185764ae2fed763feb">viennacl::linalg::cuda::prod_impl</a></div><div class="ttdeci">void prod_impl(const matrix_base< NumericT > &mat, bool mat_transpose, const vector_base< NumericT > &vec, vector_base< NumericT > &result)</div><div class="ttdoc">Carries out matrix-vector multiplication. </div><div class="ttdef"><b>Definition:</b> <a href="cuda_2matrix__operations_8hpp_source.html#l01464">matrix_operations.hpp:1464</a></div></div>
<div class="ttc" id="classviennacl_1_1vector__expression_html"><div class="ttname"><a href="classviennacl_1_1vector__expression.html">viennacl::vector_expression</a></div><div class="ttdoc">An expression template class that represents a binary operation that yields a vector. </div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00239">forwards.h:239</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__expression_html_a3125cc5f209063de6a99dc67e19d42bb"><div class="ttname"><a href="classviennacl_1_1matrix__expression.html#a3125cc5f209063de6a99dc67e19d42bb">viennacl::matrix_expression::size1</a></div><div class="ttdeci">vcl_size_t size1() const </div><div class="ttdoc">Returns the size of the result vector. </div><div class="ttdef"><b>Definition:</b> <a href="matrix_8hpp_source.html#l00072">matrix.hpp:72</a></div></div>
<div class="ttc" id="namespaceviennacl_html_aeb3693cfc246adc889875c208ca23aff"><div class="ttname"><a href="namespaceviennacl.html#aeb3693cfc246adc889875c208ca23aff">viennacl::operator-</a></div><div class="ttdeci">viennacl::vector< NumericT > operator-(const vector_base< NumericT > &v1, const vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, op_prod > &proxy)</div><div class="ttdoc">Implementation of the operation 'result = v1 - A * v2', where A is a matrix. </div><div class="ttdef"><b>Definition:</b> <a href="matrix__operations_8hpp_source.html#l01200">matrix_operations.hpp:1200</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f"><div class="ttname"><a href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a020ede27e1975479bce748e0e4ea3c7f">viennacl::MEMORY_NOT_INITIALIZED</a></div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00347">forwards.h:347</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2"><div class="ttname"><a href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a></div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00350">forwards.h:350</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1host__based_html_a73ab507a304067a153e0ad741ed6a40d"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1host__based.html#a73ab507a304067a153e0ad741ed6a40d">viennacl::linalg::host_based::prod_impl</a></div><div class="ttdeci">void prod_impl(const matrix_base< NumericT > &mat, bool trans, const vector_base< NumericT > &vec, vector_base< NumericT > &result)</div><div class="ttdoc">Carries out matrix-vector multiplication. </div><div class="ttdef"><b>Definition:</b> <a href="host__based_2matrix__operations_8hpp_source.html#l00993">matrix_operations.hpp:993</a></div></div>
<div class="ttc" id="cuda_2sparse__matrix__operations_8hpp_html"><div class="ttname"><a href="cuda_2sparse__matrix__operations_8hpp.html">sparse_matrix_operations.hpp</a></div><div class="ttdoc">Implementations of operations using sparse matrices using CUDA. </div></div>
<div class="ttc" id="classviennacl_1_1matrix__expression_html_ad8516f3a2ef4db4d5379fb64832fd39f"><div class="ttname"><a href="classviennacl_1_1matrix__expression.html#ad8516f3a2ef4db4d5379fb64832fd39f">viennacl::matrix_expression::size2</a></div><div class="ttdeci">vcl_size_t size2() const </div><div class="ttdef"><b>Definition:</b> <a href="matrix_8hpp_source.html#l00073">matrix.hpp:73</a></div></div>
<div class="ttc" id="classviennacl_1_1vector__base_html"><div class="ttname"><a href="classviennacl_1_1vector__base.html">viennacl::vector_base</a></div><div class="ttdoc">Common base class for dense vectors, vector ranges, and vector slices. </div><div class="ttdef"><b>Definition:</b> <a href="vector__def_8hpp_source.html#l00104">vector_def.hpp:104</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a98a0afcc513111ffa0bd138f891930df"><div class="ttname"><a href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">viennacl::vcl_size_t</a></div><div class="ttdeci">std::size_t vcl_size_t</div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00075">forwards.h:75</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1host__based_1_1detail_html_abc856cc45c42e8809e7e867904a4a7d7"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1host__based_1_1detail.html#abc856cc45c42e8809e7e867904a4a7d7">viennacl::linalg::host_based::detail::block_inplace_solve</a></div><div class="ttdeci">void block_inplace_solve(const matrix_expression< const compressed_matrix< NumericT, AlignmentV >, const compressed_matrix< NumericT, AlignmentV >, op_trans > &L, viennacl::backend::mem_handle const &, vcl_size_t, vector_base< NumericT > const &, vector_base< NumericT > &vec, viennacl::linalg::unit_lower_tag)</div><div class="ttdef"><b>Definition:</b> <a href="host__based_2sparse__matrix__operations_8hpp_source.html#l00824">sparse_matrix_operations.hpp:824</a></div></div>
<div class="ttc" id="classviennacl_1_1vector_html"><div class="ttname"><a href="classviennacl_1_1vector.html">viennacl::vector</a></div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00266">forwards.h:266</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1cuda_1_1detail_html_ac588666828831dc470e4ed55da6d9a10"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1cuda_1_1detail.html#ac588666828831dc470e4ed55da6d9a10">viennacl::linalg::cuda::detail::block_inplace_solve</a></div><div class="ttdeci">void block_inplace_solve(const matrix_expression< const compressed_matrix< NumericT, AlignmentV >, const compressed_matrix< NumericT, AlignmentV >, op_trans > &L, viennacl::backend::mem_handle const &block_indices, vcl_size_t num_blocks, vector_base< NumericT > const &, vector_base< NumericT > &vec, viennacl::linalg::unit_lower_tag)</div><div class="ttdef"><b>Definition:</b> <a href="cuda_2sparse__matrix__operations_8hpp_source.html#l00970">sparse_matrix_operations.hpp:970</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__base_html_a24e4f5fa27a1af5ad47e52a97c065d68"><div class="ttname"><a href="classviennacl_1_1matrix__base.html#a24e4f5fa27a1af5ad47e52a97c065d68">viennacl::matrix_base::size1</a></div><div class="ttdeci">size_type size1() const </div><div class="ttdoc">Returns the number of rows. </div><div class="ttdef"><b>Definition:</b> <a href="matrix__def_8hpp_source.html#l00224">matrix_def.hpp:224</a></div></div>
<div class="ttc" id="opencl_2sparse__matrix__operations_8hpp_html"><div class="ttname"><a href="opencl_2sparse__matrix__operations_8hpp.html">sparse_matrix_operations.hpp</a></div><div class="ttdoc">Implementations of operations using sparse matrices and OpenCL. </div></div>
<div class="ttc" id="namespaceviennacl_html_a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4"><div class="ttname"><a href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a></div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00348">forwards.h:348</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1opencl_html_ac90a3dca0595b15ae04662bddb6cf57c"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1opencl.html#ac90a3dca0595b15ae04662bddb6cf57c">viennacl::linalg::opencl::inplace_solve</a></div><div class="ttdeci">void inplace_solve(matrix_base< NumericT > const &A, matrix_base< NumericT > &B, SolverTagT)</div><div class="ttdoc">Direct inplace solver for dense triangular systems. Matlab notation: A \ B. </div><div class="ttdef"><b>Definition:</b> <a href="opencl_2direct__solve_8hpp_source.html#l00077">direct_solve.hpp:77</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a4d3a8290e94b15b653b3a0f0a3f80496"><div class="ttname"><a href="namespaceviennacl.html#a4d3a8290e94b15b653b3a0f0a3f80496">viennacl::operator+</a></div><div class="ttdeci">viennacl::vector< NumericT > operator+(const vector_base< NumericT > &v1, const vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, op_prod > &proxy)</div><div class="ttdoc">Implementation of the operation 'result = v1 + A * v2', where A is a matrix. </div><div class="ttdef"><b>Definition:</b> <a href="matrix__operations_8hpp_source.html#l01182">matrix_operations.hpp:1182</a></div></div>
<div class="ttc" id="structviennacl_1_1op__prod_html"><div class="ttname"><a href="structviennacl_1_1op__prod.html">viennacl::op_prod</a></div><div class="ttdoc">A tag class representing matrix-vector products and element-wise multiplications. ...</div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00094">forwards.h:94</a></div></div>
<div class="ttc" id="vector_8hpp_html"><div class="ttname"><a href="vector_8hpp.html">vector.hpp</a></div><div class="ttdoc">The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...</div></div>
<div class="ttc" id="classviennacl_1_1vector__base_html_a15c47ae4326098aeaa4ed6b91fc6df9b"><div class="ttname"><a href="classviennacl_1_1vector__base.html#a15c47ae4326098aeaa4ed6b91fc6df9b">viennacl::vector_base::size</a></div><div class="ttdeci">size_type size() const </div><div class="ttdoc">Returns the length of the vector (cf. std::vector) </div><div class="ttdef"><b>Definition:</b> <a href="vector__def_8hpp_source.html#l00118">vector_def.hpp:118</a></div></div>
<div class="ttc" id="fft__1d_8cpp_html_ad5c19ca4f47d3f8ec21232a5af2624e5"><div class="ttname"><a href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a></div><div class="ttdeci">float ScalarType</div><div class="ttdef"><b>Definition:</b> <a href="fft__1d_8cpp_source.html#l00042">fft_1d.cpp:42</a></div></div>
<div class="ttc" id="classviennacl_1_1backend_1_1mem__handle_html"><div class="ttname"><a href="classviennacl_1_1backend_1_1mem__handle.html">viennacl::backend::mem_handle</a></div><div class="ttdoc">Main abstraction class for multiple memory domains. Represents a buffer in either main RAM...</div><div class="ttdef"><b>Definition:</b> <a href="mem__handle_8hpp_source.html#l00089">mem_handle.hpp:89</a></div></div>
<div class="ttc" id="structviennacl_1_1op__trans_html"><div class="ttname"><a href="structviennacl_1_1op__trans.html">viennacl::op_trans</a></div><div class="ttdoc">A tag class representing transposed matrices. </div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00220">forwards.h:220</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1host__based_1_1detail_html_a58dfcdc80252eb2bc06c50a1ef533fc9"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1host__based_1_1detail.html#a58dfcdc80252eb2bc06c50a1ef533fc9">viennacl::linalg::host_based::detail::row_info</a></div><div class="ttdeci">void row_info(compressed_matrix< NumericT, AlignmentV > const &mat, vector_base< NumericT > &vec, viennacl::linalg::detail::row_info_types info_selector)</div><div class="ttdef"><b>Definition:</b> <a href="host__based_2sparse__matrix__operations_8hpp_source.html#l00053">sparse_matrix_operations.hpp:53</a></div></div>
<div class="ttc" id="classviennacl_1_1compressed__matrix_html"><div class="ttname"><a href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix< NumericT ></a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_aafc6db8d806f67c24f93eaaded84b853"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a></div><div class="ttdeci">void prod_impl(const matrix_base< NumericT > &mat, const vector_base< NumericT > &vec, vector_base< NumericT > &result)</div><div class="ttdoc">Carries out matrix-vector multiplication. </div><div class="ttdef"><b>Definition:</b> <a href="matrix__operations_8hpp_source.html#l00438">matrix_operations.hpp:438</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1traits_html_ae39853b7f291a697e119a139439178fb"><div class="ttname"><a href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a></div><div class="ttdeci">viennacl::backend::mem_handle & handle(T &obj)</div><div class="ttdoc">Returns the generic memory handle of an object. Non-const version. </div><div class="ttdef"><b>Definition:</b> <a href="traits_2handle_8hpp_source.html#l00041">handle.hpp:41</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1opencl_1_1detail_html_a144fd7dc25be2029356c24f07f281edc"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1opencl_1_1detail.html#a144fd7dc25be2029356c24f07f281edc">viennacl::linalg::opencl::detail::block_inplace_solve</a></div><div class="ttdeci">void block_inplace_solve(const matrix_expression< const compressed_matrix< NumericT, AlignmentV >, const compressed_matrix< NumericT, AlignmentV >, op_trans > &L, viennacl::backend::mem_handle const &block_indices, vcl_size_t num_blocks, vector_base< NumericT > const &, vector_base< NumericT > &x, viennacl::linalg::unit_lower_tag)</div><div class="ttdef"><b>Definition:</b> <a href="opencl_2sparse__matrix__operations_8hpp_source.html#l00483">sparse_matrix_operations.hpp:483</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__expression_html_aa2c582e8d2c1cafbf2ffdf05c9fdfaca"><div class="ttname"><a href="classviennacl_1_1matrix__expression.html#aa2c582e8d2c1cafbf2ffdf05c9fdfaca">viennacl::matrix_expression::lhs</a></div><div class="ttdeci">LHS & lhs() const </div><div class="ttdoc">Get left hand side operand. </div><div class="ttdef"><b>Definition:</b> <a href="matrix_8hpp_source.html#l00066">matrix.hpp:66</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_a5ef8943054b39077573f94234980bd0f"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a5ef8943054b39077573f94234980bd0f">viennacl::linalg::detail::row_info</a></div><div class="ttdeci">viennacl::enable_if< viennacl::is_any_sparse_matrix< SparseMatrixType >::value >::type row_info(SparseMatrixType const &mat, vector< SCALARTYPE, VEC_ALIGNMENT > &vec, row_info_types info_selector)</div><div class="ttdef"><b>Definition:</b> <a href="sparse__matrix__operations_8hpp_source.html#l00050">sparse_matrix_operations.hpp:50</a></div></div>
<div class="ttc" id="scalar_8hpp_html"><div class="ttname"><a href="scalar_8hpp.html">scalar.hpp</a></div><div class="ttdoc">Implementation of the ViennaCL scalar class. </div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1opencl_1_1detail_html_a256a599dedd093618c9936e85a75b3d4"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1opencl_1_1detail.html#a256a599dedd093618c9936e85a75b3d4">viennacl::linalg::opencl::detail::row_info</a></div><div class="ttdeci">void row_info(compressed_matrix< NumericT, AlignmentV > const &A, vector_base< NumericT > &x, viennacl::linalg::detail::row_info_types info_selector)</div><div class="ttdef"><b>Definition:</b> <a href="opencl_2sparse__matrix__operations_8hpp_source.html#l00056">sparse_matrix_operations.hpp:56</a></div></div>
<div class="ttc" id="classviennacl_1_1backend_1_1mem__handle_html_a5951aceb881068a77b4af8825b547ea3"><div class="ttname"><a href="classviennacl_1_1backend_1_1mem__handle.html#a5951aceb881068a77b4af8825b547ea3">viennacl::backend::mem_handle::get_active_handle_id</a></div><div class="ttdeci">memory_types get_active_handle_id() const </div><div class="ttdoc">Returns an ID for the currently active memory buffer. Other memory buffers might contain old or no da...</div><div class="ttdef"><b>Definition:</b> <a href="mem__handle_8hpp_source.html#l00118">mem_handle.hpp:118</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_a926494a1bd9b3c66def0c860bf564b42"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a926494a1bd9b3c66def0c860bf564b42">viennacl::linalg::detail::block_inplace_solve</a></div><div class="ttdeci">viennacl::enable_if< viennacl::is_any_sparse_matrix< SparseMatrixType >::value >::type block_inplace_solve(const matrix_expression< const SparseMatrixType, const SparseMatrixType, op_trans > &mat, viennacl::backend::mem_handle const &block_index_array, vcl_size_t num_blocks, viennacl::vector_base< ScalarType > const &mat_diagonal, viennacl::vector_base< ScalarType > &vec, SOLVERTAG tag)</div><div class="ttdef"><b>Definition:</b> <a href="sparse__matrix__operations_8hpp_source.html#l00332">sparse_matrix_operations.hpp:332</a></div></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="navelem"><a class="el" href="dir_c82e3d11dd171600f4a6e0cab1ec1e0d.html">viennacl</a></li><li class="navelem"><a class="el" href="dir_63cde087767c4ed65c7901ffc6e293fe.html">linalg</a></li><li class="navelem"><a class="el" href="sparse__matrix__operations_8hpp.html">sparse_matrix_operations.hpp</a></li>
<li class="footer">Generated on Wed Jan 20 2016 22:32:41 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>
|