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 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594
|
<!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/coordinate_matrix.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('coordinate__matrix_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">coordinate_matrix.hpp</div> </div>
</div><!--header-->
<div class="contents">
<a href="coordinate__matrix_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_COORDINATE_MATRIX_HPP_</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="preprocessor">#define VIENNACL_COORDINATE_MATRIX_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 <map></span></div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span> <span class="preprocessor">#include <vector></span></div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span> <span class="preprocessor">#include <list></span></div>
<div class="line"><a name="l00028"></a><span class="lineno"> 28</span> </div>
<div class="line"><a name="l00029"></a><span class="lineno"> 29</span> <span class="preprocessor">#include "<a class="code" href="forwards_8h.html">viennacl/forwards.h</a>"</span></div>
<div class="line"><a name="l00030"></a><span class="lineno"> 30</span> <span class="preprocessor">#include "<a class="code" href="vector_8hpp.html">viennacl/vector.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">#include "<a class="code" href="sparse__matrix__operations_8hpp.html">viennacl/linalg/sparse_matrix_operations.hpp</a>"</span></div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span> </div>
<div class="line"><a name="l00034"></a><span class="lineno"> 34</span> <span class="keyword">namespace </span>viennacl</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> </div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span> </div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span> <span class="comment">//provide copy-operation:</span></div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span> <span class="comment"></span><span class="keyword">template</span><<span class="keyword">typename</span> CPUMatrixT, <span class="keyword">typename</span> NumericT, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> AlignmentV></div>
<div class="line"><a name="l00047"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#a3947c4a75b60a14b7051746805e1bc47"> 47</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(<span class="keyword">const</span> CPUMatrixT & cpu_matrix,</div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix<NumericT, AlignmentV></a> & gpu_matrix )</div>
<div class="line"><a name="l00049"></a><span class="lineno"> 49</span> {</div>
<div class="line"><a name="l00050"></a><span class="lineno"> 50</span>  assert( (gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>() == 0 || <a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) == gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size mismatch"</span>) );</div>
<div class="line"><a name="l00051"></a><span class="lineno"> 51</span>  assert( (gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99">size2</a>() == 0 || <a class="code" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98">viennacl::traits::size2</a>(cpu_matrix) == gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99">size2</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size mismatch"</span>) );</div>
<div class="line"><a name="l00052"></a><span class="lineno"> 52</span> </div>
<div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> group_num = 64;</div>
<div class="line"><a name="l00054"></a><span class="lineno"> 54</span> </div>
<div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  <span class="comment">// Step 1: Determine nonzeros:</span></div>
<div class="line"><a name="l00056"></a><span class="lineno"> 56</span>  <span class="keywordflow">if</span> ( cpu_matrix.size1() > 0 && cpu_matrix.size2() > 0 )</div>
<div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  {</div>
<div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> num_entries = 0;</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  <span class="keywordflow">for</span> (<span class="keyword">typename</span> CPUMatrixT::const_iterator1 row_it = cpu_matrix.begin1(); row_it != cpu_matrix.end1(); ++row_it)</div>
<div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  <span class="keywordflow">for</span> (<span class="keyword">typename</span> CPUMatrixT::const_iterator2 col_it = row_it.begin(); col_it != row_it.end(); ++col_it)</div>
<div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  ++num_entries;</div>
<div class="line"><a name="l00062"></a><span class="lineno"> 62</span> </div>
<div class="line"><a name="l00063"></a><span class="lineno"> 63</span>  <span class="comment">// Step 2: Set up matrix data:</span></div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  gpu_matrix.nonzeros_ = num_entries;</div>
<div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  gpu_matrix.rows_ = cpu_matrix.size1();</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  gpu_matrix.cols_ = cpu_matrix.size2();</div>
<div class="line"><a name="l00067"></a><span class="lineno"> 67</span> </div>
<div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a> group_boundaries(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a66b7f5c5487b719c91bc76984f89f97f">handle3</a>(), group_num + 1);</div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a> coord_buffer(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a875b6ccacbf0df8ef48ab7f7c46244f1">handle12</a>(), 2*gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a2035ac2e7acc39c21cf13bb6f9c2dc08">internal_nnz</a>());</div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  std::vector<NumericT> elements(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a2035ac2e7acc39c21cf13bb6f9c2dc08">internal_nnz</a>());</div>
<div class="line"><a name="l00071"></a><span class="lineno"> 71</span> </div>
<div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> data_index = 0;</div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> current_fraction = 0;</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>  group_boundaries.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac31bfde381eb40e0349678252c444af9">set</a>(0, 0);</div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  <span class="keywordflow">for</span> (<span class="keyword">typename</span> CPUMatrixT::const_iterator1 row_it = cpu_matrix.begin1(); row_it != cpu_matrix.end1(); ++row_it)</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>  <span class="keywordflow">for</span> (<span class="keyword">typename</span> CPUMatrixT::const_iterator2 col_it = row_it.begin(); col_it != row_it.end(); ++col_it)</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>  coord_buffer.set(2*data_index, col_it.index1());</div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  coord_buffer.set(2*data_index + 1, col_it.index2());</div>
<div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  elements[data_index] = *col_it;</div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  ++data_index;</div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  }</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span> </div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  <span class="keywordflow">while</span> (data_index > <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a>(static_cast<double>(current_fraction + 1) / static_cast<double>(group_num)) * num_entries) <span class="comment">//split data equally over 64 groups</span></div>
<div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  group_boundaries.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac31bfde381eb40e0349678252c444af9">set</a>(++current_fraction, data_index);</div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  }</div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span> </div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="comment">//write end of last group:</span></div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  group_boundaries.set(group_num, data_index);</div>
<div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  <span class="comment">//group_boundaries[1] = data_index; //for one compute unit</span></div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span> </div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  <span class="comment">//std::cout << "Group boundaries: " << std::endl;</span></div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  <span class="comment">//for (vcl_size_t i=0; i<group_boundaries.size(); ++i)</span></div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <span class="comment">// std::cout << group_boundaries[i] << std::endl;</span></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>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(gpu_matrix.group_boundaries_, group_boundaries.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac8373f0d899b89c843e14de4cb7a1c4a">raw_size</a>(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">traits::context</a>(gpu_matrix.group_boundaries_), group_boundaries.get());</div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(gpu_matrix.coord_buffer_, coord_buffer.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac8373f0d899b89c843e14de4cb7a1c4a">raw_size</a>(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">traits::context</a>(gpu_matrix.coord_buffer_), coord_buffer.get());</div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(gpu_matrix.elements_, <span class="keyword">sizeof</span>(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>)*elements.size(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">traits::context</a>(gpu_matrix.elements_), &(elements[0]));</div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  }</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> </div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span> <span class="keyword">template</span><<span class="keyword">typename</span> NumericT, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> AlignmentV></div>
<div class="line"><a name="l00110"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#a9285f6f62eebc19c0756c5a31f3806f1"> 110</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(<span class="keyword">const</span> std::vector< std::map<unsigned int, NumericT> > & cpu_matrix,</div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix<NumericT, AlignmentV></a> & gpu_matrix )</div>
<div class="line"><a name="l00112"></a><span class="lineno"> 112</span> {</div>
<div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> max_col = 0;</div>
<div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <span class="keywordflow">for</span> (<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> i=0; i<cpu_matrix.size(); ++i)</div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  {</div>
<div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  <span class="keywordflow">if</span> (cpu_matrix[i].<a class="code" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">size</a>() > 0)</div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  max_col = std::max<vcl_size_t>(max_col, (cpu_matrix[i].rbegin())->first);</div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  }</div>
<div class="line"><a name="l00119"></a><span class="lineno"> 119</span> </div>
<div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(<a class="code" href="classviennacl_1_1tools_1_1const__sparse__matrix__adapter.html">tools::const_sparse_matrix_adapter<NumericT></a>(cpu_matrix, cpu_matrix.size(), max_col + 1), gpu_matrix);</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> <span class="comment">//gpu to cpu:</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> CPUMatrixT, <span class="keyword">typename</span> NumericT, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> AlignmentV></div>
<div class="line"><a name="l00134"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#a6f355f732e924a4233950b0e77a14db5"> 134</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(<span class="keyword">const</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix<NumericT, AlignmentV></a> & gpu_matrix,</div>
<div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  CPUMatrixT & cpu_matrix )</div>
<div class="line"><a name="l00136"></a><span class="lineno"> 136</span> {</div>
<div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  assert( (<a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) == gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size mismatch"</span>) );</div>
<div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  assert( (<a class="code" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98">viennacl::traits::size2</a>(cpu_matrix) == gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99">size2</a>()) && <span class="keywordtype">bool</span>(<span class="stringliteral">"Size mismatch"</span>) );</div>
<div class="line"><a name="l00139"></a><span class="lineno"> 139</span> </div>
<div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  <span class="keywordflow">if</span> ( gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>() > 0 && gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99">size2</a>() > 0 )</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="comment">//get raw data from memory:</span></div>
<div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a> coord_buffer(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a875b6ccacbf0df8ef48ab7f7c46244f1">handle12</a>(), 2*gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a04c51b826e4bae9c97a3749ee2161b87">nnz</a>());</div>
<div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  std::vector<NumericT> elements(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a04c51b826e4bae9c97a3749ee2161b87">nnz</a>());</div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span> </div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  <span class="comment">//std::cout << "GPU nonzeros: " << gpu_matrix.nnz() << std::endl;</span></div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span> </div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a62854cfd6f04404b274f8ede36f63e2d">viennacl::backend::memory_read</a>(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a875b6ccacbf0df8ef48ab7f7c46244f1">handle12</a>(), 0, coord_buffer.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac8373f0d899b89c843e14de4cb7a1c4a">raw_size</a>(), coord_buffer.get());</div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a62854cfd6f04404b274f8ede36f63e2d">viennacl::backend::memory_read</a>(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#aeb74dd2ab673d7912a17f60777e9bf6c">handle</a>(), 0, <span class="keyword">sizeof</span>(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>) * elements.size(), &(elements[0]));</div>
<div class="line"><a name="l00150"></a><span class="lineno"> 150</span> </div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  <span class="comment">//fill the cpu_matrix:</span></div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  <span class="keywordflow">for</span> (<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> index = 0; index < gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a04c51b826e4bae9c97a3749ee2161b87">nnz</a>(); ++index)</div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  cpu_matrix(coord_buffer[2*index], coord_buffer[2*index+1]) = elements[index];</div>
<div class="line"><a name="l00154"></a><span class="lineno"> 154</span> </div>
<div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  }</div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span> }</div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</span> </div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span> <span class="keyword">template</span><<span class="keyword">typename</span> NumericT, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> AlignmentV></div>
<div class="line"><a name="l00164"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#a09970b6d02bfa2dd17a7d63d1fcbe18f"> 164</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(<span class="keyword">const</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix<NumericT, AlignmentV></a> & gpu_matrix,</div>
<div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  std::vector< std::map<unsigned int, NumericT> > & cpu_matrix)</div>
<div class="line"><a name="l00166"></a><span class="lineno"> 166</span> {</div>
<div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  <span class="keywordflow">if</span> (cpu_matrix.size() == 0)</div>
<div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  cpu_matrix.resize(gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>());</div>
<div class="line"><a name="l00169"></a><span class="lineno"> 169</span> </div>
<div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  assert(cpu_matrix.size() == gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>() && bool(<span class="stringliteral">"Matrix dimension mismatch!"</span>));</div>
<div class="line"><a name="l00171"></a><span class="lineno"> 171</span> </div>
<div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  <a class="code" href="classviennacl_1_1tools_1_1sparse__matrix__adapter.html">tools::sparse_matrix_adapter<NumericT></a> temp(cpu_matrix, gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>(), gpu_matrix.<a class="code" href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99">size2</a>());</div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(gpu_matrix, temp);</div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span> }</div>
<div class="line"><a name="l00175"></a><span class="lineno"> 175</span> </div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span> </div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span> </div>
<div class="line"><a name="l00185"></a><span class="lineno"> 185</span> <span class="keyword">template</span><<span class="keyword">class </span><a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> AlignmentV <span class="comment">/* see forwards.h */</span> ></div>
<div class="line"><a name="l00186"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html"> 186</a></span> <span class="keyword">class </span><a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix</a></div>
<div class="line"><a name="l00187"></a><span class="lineno"> 187</span> {</div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span> <span class="keyword">public</span>:</div>
<div class="line"><a name="l00189"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a5497d1e7fdce125b7a5bd12530002b4f"> 189</a></span>  <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">viennacl::backend::mem_handle</a> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a5497d1e7fdce125b7a5bd12530002b4f">handle_type</a>;</div>
<div class="line"><a name="l00190"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a9fd9168600c330faa69cae2626509602"> 190</a></span>  <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1scalar.html">scalar<typename viennacl::tools::CHECK_SCALAR_TEMPLATE_ARGUMENT<NumericT>::ResultType</a>> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a9fd9168600c330faa69cae2626509602">value_type</a>;</div>
<div class="line"><a name="l00191"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#ad4218882f3bd97eeff6d81b575ea1e3c"> 191</a></span>  <span class="keyword">typedef</span> <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1coordinate__matrix.html#ad4218882f3bd97eeff6d81b575ea1e3c">size_type</a>;</div>
<div class="line"><a name="l00192"></a><span class="lineno"> 192</span> </div>
<div class="line"><a name="l00194"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a626ac84303ae01ac4f7e964f2c2972b6"> 194</a></span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html#a626ac84303ae01ac4f7e964f2c2972b6">coordinate_matrix</a>() : rows_(0), cols_(0), nonzeros_(0), group_num_(64) {}</div>
<div class="line"><a name="l00195"></a><span class="lineno"> 195</span> </div>
<div class="line"><a name="l00196"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a9cde8a8c5cb1516c5503d011388daed4"> 196</a></span>  <span class="keyword">explicit</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a9cde8a8c5cb1516c5503d011388daed4">coordinate_matrix</a>(<a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx) : rows_(0), cols_(0), nonzeros_(0), group_num_(64)</div>
<div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  {</div>
<div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  group_boundaries_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>());</div>
<div class="line"><a name="l00199"></a><span class="lineno"> 199</span>  coord_buffer_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>());</div>
<div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  elements_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>());</div>
<div class="line"><a name="l00201"></a><span class="lineno"> 201</span> </div>
<div class="line"><a name="l00202"></a><span class="lineno"> 202</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  <span class="keywordflow">if</span> (ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>() == <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">OPENCL_MEMORY</a>)</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>  group_boundaries_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  coord_buffer_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  elements_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  }</div>
<div class="line"><a name="l00209"></a><span class="lineno"> 209</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  }</div>
<div class="line"><a name="l00211"></a><span class="lineno"> 211</span> </div>
<div class="line"><a name="l00219"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#ad5b57c8ab8bcc8d839aba53be7cb1463"> 219</a></span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html#ad5b57c8ab8bcc8d839aba53be7cb1463">coordinate_matrix</a>(<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> rows, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> cols, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> nonzeros = 0, <a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx = <a class="code" href="classviennacl_1_1context.html">viennacl::context</a>()) :</div>
<div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  rows_(rows), cols_(cols), nonzeros_(nonzeros)</div>
<div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  {</div>
<div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  <span class="keywordflow">if</span> (nonzeros > 0)</div>
<div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  {</div>
<div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(group_boundaries_, <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a>().<a class="code" href="namespaceviennacl_1_1backend_1_1detail.html#af9f9cc6c50debc5b1d93e458d1fd392a">element_size</a>() * (group_num_ + 1), ctx);</div>
<div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(coord_buffer_, <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a>().<a class="code" href="namespaceviennacl_1_1backend_1_1detail.html#af9f9cc6c50debc5b1d93e458d1fd392a">element_size</a>() * 2 * <a class="code" href="classviennacl_1_1coordinate__matrix.html#a2035ac2e7acc39c21cf13bb6f9c2dc08">internal_nnz</a>(), ctx);</div>
<div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(elements_, <span class="keyword">sizeof</span>(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>) * <a class="code" href="classviennacl_1_1coordinate__matrix.html#a2035ac2e7acc39c21cf13bb6f9c2dc08">internal_nnz</a>(), ctx);</div>
<div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  }</div>
<div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  {</div>
<div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  group_boundaries_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.memory_type());</div>
<div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  coord_buffer_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.memory_type());</div>
<div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  elements_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.memory_type());</div>
<div class="line"><a name="l00233"></a><span class="lineno"> 233</span> </div>
<div class="line"><a name="l00234"></a><span class="lineno"> 234</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  <span class="keywordflow">if</span> (ctx.memory_type() == <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">OPENCL_MEMORY</a>)</div>
<div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  {</div>
<div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  group_boundaries_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  coord_buffer_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  elements_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00240"></a><span class="lineno"> 240</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>  }</div>
<div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  }</div>
<div class="line"><a name="l00244"></a><span class="lineno"> 244</span> </div>
<div class="line"><a name="l00251"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a472604b5e8e7a818fcee171898df1c92"> 251</a></span>  <span class="keyword">explicit</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a472604b5e8e7a818fcee171898df1c92">coordinate_matrix</a>(<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> rows, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> cols, <a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx)</div>
<div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  : rows_(rows), cols_(cols), nonzeros_(0)</div>
<div class="line"><a name="l00253"></a><span class="lineno"> 253</span>  {</div>
<div class="line"><a name="l00254"></a><span class="lineno"> 254</span>  group_boundaries_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>());</div>
<div class="line"><a name="l00255"></a><span class="lineno"> 255</span>  coord_buffer_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>());</div>
<div class="line"><a name="l00256"></a><span class="lineno"> 256</span>  elements_.<a class="code" href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">switch_active_handle_id</a>(ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>());</div>
<div class="line"><a name="l00257"></a><span class="lineno"> 257</span> </div>
<div class="line"><a name="l00258"></a><span class="lineno"> 258</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  <span class="keywordflow">if</span> (ctx.<a class="code" href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">memory_type</a>() == <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">OPENCL_MEMORY</a>)</div>
<div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  {</div>
<div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  group_boundaries_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  coord_buffer_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  elements_.opencl_handle().context(ctx.opencl_context());</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="preprocessor">#endif</span></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> </div>
<div class="line"><a name="l00268"></a><span class="lineno"> 268</span> </div>
<div class="line"><a name="l00270"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#ad091412ca52bb131234464ed3f5b103c"> 270</a></span>  <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html#ad091412ca52bb131234464ed3f5b103c">reserve</a>(<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> new_nonzeros)</div>
<div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  {</div>
<div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  <span class="keywordflow">if</span> (new_nonzeros > nonzeros_) <span class="comment">//TODO: Do we need to initialize new memory with zero?</span></div>
<div class="line"><a name="l00273"></a><span class="lineno"> 273</span>  {</div>
<div class="line"><a name="l00274"></a><span class="lineno"> 274</span>  <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> coord_buffer_old;</div>
<div class="line"><a name="l00275"></a><span class="lineno"> 275</span>  <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> elements_old;</div>
<div class="line"><a name="l00276"></a><span class="lineno"> 276</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a8548a37535eac92e886cb183182b618f">viennacl::backend::memory_shallow_copy</a>(coord_buffer_, coord_buffer_old);</div>
<div class="line"><a name="l00277"></a><span class="lineno"> 277</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a8548a37535eac92e886cb183182b618f">viennacl::backend::memory_shallow_copy</a>(elements_, elements_old);</div>
<div class="line"><a name="l00278"></a><span class="lineno"> 278</span> </div>
<div class="line"><a name="l00279"></a><span class="lineno"> 279</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> internal_new_nnz = viennacl::tools::align_to_multiple<vcl_size_t>(new_nonzeros, AlignmentV);</div>
<div class="line"><a name="l00280"></a><span class="lineno"> 280</span>  <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a> size_deducer(coord_buffer_);</div>
<div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(coord_buffer_, size_deducer.element_size() * 2 * internal_new_nnz, <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(coord_buffer_));</div>
<div class="line"><a name="l00282"></a><span class="lineno"> 282</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(elements_, <span class="keyword">sizeof</span>(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>) * internal_new_nnz, <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(elements_));</div>
<div class="line"><a name="l00283"></a><span class="lineno"> 283</span> </div>
<div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a39a91fd7aa68ebf64bd7bfa409732ec0">viennacl::backend::memory_copy</a>(coord_buffer_old, coord_buffer_, 0, 0, size_deducer.element_size() * 2 * nonzeros_);</div>
<div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a39a91fd7aa68ebf64bd7bfa409732ec0">viennacl::backend::memory_copy</a>(elements_old, elements_, 0, 0, <span class="keyword">sizeof</span>(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>) * nonzeros_);</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>  nonzeros_ = new_nonzeros;</div>
<div class="line"><a name="l00288"></a><span class="lineno"> 288</span>  }</div>
<div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  }</div>
<div class="line"><a name="l00290"></a><span class="lineno"> 290</span> </div>
<div class="line"><a name="l00297"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a5363ed7b82aaad8549c85d880755f846"> 297</a></span>  <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a5363ed7b82aaad8549c85d880755f846">resize</a>(<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> new_size1, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> new_size2, <span class="keywordtype">bool</span> preserve = <span class="keyword">true</span>)</div>
<div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  {</div>
<div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  assert (new_size1 > 0 && new_size2 > 0);</div>
<div class="line"><a name="l00300"></a><span class="lineno"> 300</span> </div>
<div class="line"><a name="l00301"></a><span class="lineno"> 301</span>  <span class="keywordflow">if</span> (new_size1 < rows_ || new_size2 < cols_) <span class="comment">//enlarge buffer</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>  std::vector<std::map<unsigned int, NumericT> > stl_sparse_matrix;</div>
<div class="line"><a name="l00304"></a><span class="lineno"> 304</span>  <span class="keywordflow">if</span> (rows_ > 0)</div>
<div class="line"><a name="l00305"></a><span class="lineno"> 305</span>  stl_sparse_matrix.resize(rows_);</div>
<div class="line"><a name="l00306"></a><span class="lineno"> 306</span> </div>
<div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  <span class="keywordflow">if</span> (preserve && rows_ > 0)</div>
<div class="line"><a name="l00308"></a><span class="lineno"> 308</span>  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(*<span class="keyword">this</span>, stl_sparse_matrix);</div>
<div class="line"><a name="l00309"></a><span class="lineno"> 309</span> </div>
<div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  stl_sparse_matrix.resize(new_size1);</div>
<div class="line"><a name="l00311"></a><span class="lineno"> 311</span> </div>
<div class="line"><a name="l00312"></a><span class="lineno"> 312</span>  <span class="comment">//std::cout << "Cropping STL matrix of size " << stl_sparse_matrix.size() << std::endl;</span></div>
<div class="line"><a name="l00313"></a><span class="lineno"> 313</span>  <span class="keywordflow">if</span> (new_size2 < cols_ && rows_ > 0)</div>
<div class="line"><a name="l00314"></a><span class="lineno"> 314</span>  {</div>
<div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  <span class="keywordflow">for</span> (<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> i=0; i<stl_sparse_matrix.size(); ++i)</div>
<div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  {</div>
<div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  std::list<unsigned int> to_delete;</div>
<div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  <span class="keywordflow">for</span> (<span class="keyword">typename</span> std::map<unsigned int, NumericT>::iterator it = stl_sparse_matrix[i].begin();</div>
<div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  it != stl_sparse_matrix[i].end();</div>
<div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  ++it)</div>
<div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  {</div>
<div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  <span class="keywordflow">if</span> (it->first >= new_size2)</div>
<div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  to_delete.push_back(it->first);</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>  <span class="keywordflow">for</span> (std::list<unsigned int>::iterator it = to_delete.begin(); it != to_delete.end(); ++it)</div>
<div class="line"><a name="l00327"></a><span class="lineno"> 327</span>  stl_sparse_matrix[i].erase(*it);</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>  <span class="comment">//std::cout << "Cropping done..." << std::endl;</span></div>
<div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  }</div>
<div class="line"><a name="l00331"></a><span class="lineno"> 331</span> </div>
<div class="line"><a name="l00332"></a><span class="lineno"> 332</span>  rows_ = new_size1;</div>
<div class="line"><a name="l00333"></a><span class="lineno"> 333</span>  cols_ = new_size2;</div>
<div class="line"><a name="l00334"></a><span class="lineno"> 334</span>  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(stl_sparse_matrix, *<span class="keyword">this</span>);</div>
<div class="line"><a name="l00335"></a><span class="lineno"> 335</span>  }</div>
<div class="line"><a name="l00336"></a><span class="lineno"> 336</span> </div>
<div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  rows_ = new_size1;</div>
<div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  cols_ = new_size2;</div>
<div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  }</div>
<div class="line"><a name="l00340"></a><span class="lineno"> 340</span> </div>
<div class="line"><a name="l00342"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a4fccbe4157aa8833d23b0541803e06af"> 342</a></span>  <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a4fccbe4157aa8833d23b0541803e06af">clear</a>()</div>
<div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  {</div>
<div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a> host_group_buffer(group_boundaries_, 65);</div>
<div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array<unsigned int></a> host_coord_buffer(coord_buffer_, 2);</div>
<div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  std::vector<NumericT> host_elements(1);</div>
<div class="line"><a name="l00347"></a><span class="lineno"> 347</span> </div>
<div class="line"><a name="l00348"></a><span class="lineno"> 348</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(group_boundaries_, host_group_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac8646f4841c47048a1d2e004b9536af8">element_size</a>() * 65, <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(group_boundaries_), host_group_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#a1727eeb845f50cf5bd51ef54c5938b8c">get</a>());</div>
<div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(coord_buffer_, host_coord_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac8646f4841c47048a1d2e004b9536af8">element_size</a>() * 2, <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(coord_buffer_), host_coord_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#a1727eeb845f50cf5bd51ef54c5938b8c">get</a>());</div>
<div class="line"><a name="l00350"></a><span class="lineno"> 350</span>  <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(elements_, <span class="keyword">sizeof</span>(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>) * 1, <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(elements_), &(host_elements[0]));</div>
<div class="line"><a name="l00351"></a><span class="lineno"> 351</span> </div>
<div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  nonzeros_ = 0;</div>
<div class="line"><a name="l00353"></a><span class="lineno"> 353</span>  group_num_ = 64;</div>
<div class="line"><a name="l00354"></a><span class="lineno"> 354</span>  }</div>
<div class="line"><a name="l00355"></a><span class="lineno"> 355</span> </div>
<div class="line"><a name="l00357"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377"> 357</a></span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">size1</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> rows_; }</div>
<div class="line"><a name="l00359"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99"> 359</a></span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99">size2</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> cols_; }</div>
<div class="line"><a name="l00361"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a04c51b826e4bae9c97a3749ee2161b87"> 361</a></span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a04c51b826e4bae9c97a3749ee2161b87">nnz</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> nonzeros_; }</div>
<div class="line"><a name="l00363"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a2035ac2e7acc39c21cf13bb6f9c2dc08"> 363</a></span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a2035ac2e7acc39c21cf13bb6f9c2dc08">internal_nnz</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> viennacl::tools::align_to_multiple<vcl_size_t>(nonzeros_, AlignmentV); }</div>
<div class="line"><a name="l00364"></a><span class="lineno"> 364</span> </div>
<div class="line"><a name="l00366"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a875b6ccacbf0df8ef48ab7f7c46244f1"> 366</a></span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> & <a class="code" href="classviennacl_1_1coordinate__matrix.html#a875b6ccacbf0df8ef48ab7f7c46244f1">handle12</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> coord_buffer_; }</div>
<div class="line"><a name="l00368"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#aeb74dd2ab673d7912a17f60777e9bf6c"> 368</a></span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> & <a class="code" href="classviennacl_1_1coordinate__matrix.html#aeb74dd2ab673d7912a17f60777e9bf6c">handle</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> elements_; }</div>
<div class="line"><a name="l00370"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a66b7f5c5487b719c91bc76984f89f97f"> 370</a></span>  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> & <a class="code" href="classviennacl_1_1coordinate__matrix.html#a66b7f5c5487b719c91bc76984f89f97f">handle3</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> group_boundaries_; }</div>
<div class="line"><a name="l00371"></a><span class="lineno"> 371</span> </div>
<div class="line"><a name="l00372"></a><span class="lineno"><a class="line" href="classviennacl_1_1coordinate__matrix.html#a7f83ad82b3bab755d76e7f9d665367f7"> 372</a></span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a7f83ad82b3bab755d76e7f9d665367f7">groups</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> group_num_; }</div>
<div class="line"><a name="l00373"></a><span class="lineno"> 373</span> </div>
<div class="line"><a name="l00374"></a><span class="lineno"> 374</span> <span class="preprocessor">#if defined(_MSC_VER) && _MSC_VER < 1500 //Visual Studio 2005 needs special treatment</span></div>
<div class="line"><a name="l00375"></a><span class="lineno"> 375</span>  <span class="keyword">template</span><<span class="keyword">typename</span> CPUMatrixT></div>
<div class="line"><a name="l00376"></a><span class="lineno"> 376</span>  <span class="keyword">friend</span> <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a3cd8e9989ee1b547dbcba2fe87f7b29d">copy</a>(<span class="keyword">const</span> CPUMatrixT & cpu_matrix, <a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix</a> & gpu_matrix );</div>
<div class="line"><a name="l00377"></a><span class="lineno"> 377</span> <span class="preprocessor">#else</span></div>
<div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  <span class="keyword">template</span><<span class="keyword">typename</span> CPUMatrixT, <span class="keyword">typename</span> NumericT2, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> AlignmentV2></div>
<div class="line"><a name="l00379"></a><span class="lineno"> 379</span>  <span class="keyword">friend</span> <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1coordinate__matrix.html#a3cd8e9989ee1b547dbcba2fe87f7b29d">copy</a>(<span class="keyword">const</span> CPUMatrixT & cpu_matrix, <a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix<NumericT2, AlignmentV2></a> & gpu_matrix );</div>
<div class="line"><a name="l00380"></a><span class="lineno"> 380</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00381"></a><span class="lineno"> 381</span> </div>
<div class="line"><a name="l00382"></a><span class="lineno"> 382</span> <span class="keyword">private</span>:</div>
<div class="line"><a name="l00384"></a><span class="lineno"> 384</span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html#a626ac84303ae01ac4f7e964f2c2972b6">coordinate_matrix</a>(<a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix</a> <span class="keyword">const</span> &);</div>
<div class="line"><a name="l00385"></a><span class="lineno"> 385</span> </div>
<div class="line"><a name="l00387"></a><span class="lineno"> 387</span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix</a> & operator=(<a class="code" href="classviennacl_1_1coordinate__matrix.html">coordinate_matrix</a> <span class="keyword">const</span> &);</div>
<div class="line"><a name="l00388"></a><span class="lineno"> 388</span> </div>
<div class="line"><a name="l00389"></a><span class="lineno"> 389</span> </div>
<div class="line"><a name="l00390"></a><span class="lineno"> 390</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> rows_;</div>
<div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> cols_;</div>
<div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> nonzeros_;</div>
<div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> group_num_;</div>
<div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html#a5497d1e7fdce125b7a5bd12530002b4f">handle_type</a> coord_buffer_;</div>
<div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html#a5497d1e7fdce125b7a5bd12530002b4f">handle_type</a> elements_;</div>
<div class="line"><a name="l00396"></a><span class="lineno"> 396</span>  <a class="code" href="classviennacl_1_1coordinate__matrix.html#a5497d1e7fdce125b7a5bd12530002b4f">handle_type</a> group_boundaries_;</div>
<div class="line"><a name="l00397"></a><span class="lineno"> 397</span> };</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="l00400"></a><span class="lineno"> 400</span> <span class="comment">//</span></div>
<div class="line"><a name="l00401"></a><span class="lineno"> 401</span> <span class="comment">// Specify available operations:</span></div>
<div class="line"><a name="l00402"></a><span class="lineno"> 402</span> <span class="comment">//</span></div>
<div class="line"><a name="l00403"></a><span class="lineno"> 403</span> </div>
<div class="line"><a name="l00406"></a><span class="lineno"> 406</span> <span class="keyword">namespace </span>linalg</div>
<div class="line"><a name="l00407"></a><span class="lineno"> 407</span> {</div>
<div class="line"><a name="l00408"></a><span class="lineno"> 408</span> <span class="keyword">namespace </span>detail</div>
<div class="line"><a name="l00409"></a><span class="lineno"> 409</span> {</div>
<div class="line"><a name="l00410"></a><span class="lineno"> 410</span>  <span class="comment">// x = A * y</span></div>
<div class="line"><a name="l00411"></a><span class="lineno"> 411</span>  <span class="keyword">template</span><<span class="keyword">typename</span> T, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> A></div>
<div class="line"><a name="l00412"></a><span class="lineno"> 412</span>  <span class="keyword">struct </span>op_executor<vector_base<T>, <a class="code" href="structop__assign.html">op_assign</a>, vector_expression<const coordinate_matrix<T, A>, const vector_base<T>, op_prod> ></div>
<div class="line"><a name="l00413"></a><span class="lineno"> 413</span>  {</div>
<div class="line"><a name="l00414"></a><span class="lineno"> 414</span>  <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base<T> & lhs, vector_expression<<span class="keyword">const</span> coordinate_matrix<T, A>, <span class="keyword">const</span> vector_base<T>, op_prod> <span class="keyword">const</span> & rhs)</div>
<div class="line"><a name="l00415"></a><span class="lineno"> 415</span>  {</div>
<div class="line"><a name="l00416"></a><span class="lineno"> 416</span>  <span class="comment">// check for the special case x = A * x</span></div>
<div class="line"><a name="l00417"></a><span class="lineno"> 417</span>  <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(lhs) == <a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(rhs.rhs()))</div>
<div class="line"><a name="l00418"></a><span class="lineno"> 418</span>  {</div>
<div class="line"><a name="l00419"></a><span class="lineno"> 419</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp(lhs);</div>
<div class="line"><a name="l00420"></a><span class="lineno"> 420</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), T(1), temp, T(0));</div>
<div class="line"><a name="l00421"></a><span class="lineno"> 421</span>  lhs = temp;</div>
<div class="line"><a name="l00422"></a><span class="lineno"> 422</span>  }</div>
<div class="line"><a name="l00423"></a><span class="lineno"> 423</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00424"></a><span class="lineno"> 424</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), T(1), lhs, T(0));</div>
<div class="line"><a name="l00425"></a><span class="lineno"> 425</span>  }</div>
<div class="line"><a name="l00426"></a><span class="lineno"> 426</span>  };</div>
<div class="line"><a name="l00427"></a><span class="lineno"> 427</span> </div>
<div class="line"><a name="l00428"></a><span class="lineno"> 428</span>  <span class="keyword">template</span><<span class="keyword">typename</span> T, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> A></div>
<div class="line"><a name="l00429"></a><span class="lineno"> 429</span>  <span class="keyword">struct </span>op_executor<vector_base<T>, op_inplace_add, vector_expression<const coordinate_matrix<T, A>, const vector_base<T>, op_prod> ></div>
<div class="line"><a name="l00430"></a><span class="lineno"> 430</span>  {</div>
<div class="line"><a name="l00431"></a><span class="lineno"> 431</span>  <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base<T> & lhs, vector_expression<<span class="keyword">const</span> coordinate_matrix<T, A>, <span class="keyword">const</span> vector_base<T>, op_prod> <span class="keyword">const</span> & rhs)</div>
<div class="line"><a name="l00432"></a><span class="lineno"> 432</span>  {</div>
<div class="line"><a name="l00433"></a><span class="lineno"> 433</span>  <span class="comment">// check for the special case x += A * x</span></div>
<div class="line"><a name="l00434"></a><span class="lineno"> 434</span>  <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(lhs) == <a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(rhs.rhs()))</div>
<div class="line"><a name="l00435"></a><span class="lineno"> 435</span>  {</div>
<div class="line"><a name="l00436"></a><span class="lineno"> 436</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp(lhs);</div>
<div class="line"><a name="l00437"></a><span class="lineno"> 437</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), T(1), temp, T(0));</div>
<div class="line"><a name="l00438"></a><span class="lineno"> 438</span>  lhs += temp;</div>
<div class="line"><a name="l00439"></a><span class="lineno"> 439</span>  }</div>
<div class="line"><a name="l00440"></a><span class="lineno"> 440</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00441"></a><span class="lineno"> 441</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), T(1), lhs, T(1));</div>
<div class="line"><a name="l00442"></a><span class="lineno"> 442</span>  }</div>
<div class="line"><a name="l00443"></a><span class="lineno"> 443</span>  };</div>
<div class="line"><a name="l00444"></a><span class="lineno"> 444</span> </div>
<div class="line"><a name="l00445"></a><span class="lineno"> 445</span>  <span class="keyword">template</span><<span class="keyword">typename</span> T, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> A></div>
<div class="line"><a name="l00446"></a><span class="lineno"> 446</span>  <span class="keyword">struct </span>op_executor<vector_base<T>, op_inplace_sub, vector_expression<const coordinate_matrix<T, A>, const vector_base<T>, op_prod> ></div>
<div class="line"><a name="l00447"></a><span class="lineno"> 447</span>  {</div>
<div class="line"><a name="l00448"></a><span class="lineno"> 448</span>  <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base<T> & lhs, vector_expression<<span class="keyword">const</span> coordinate_matrix<T, A>, <span class="keyword">const</span> vector_base<T>, op_prod> <span class="keyword">const</span> & rhs)</div>
<div class="line"><a name="l00449"></a><span class="lineno"> 449</span>  {</div>
<div class="line"><a name="l00450"></a><span class="lineno"> 450</span>  <span class="comment">// check for the special case x -= A * x</span></div>
<div class="line"><a name="l00451"></a><span class="lineno"> 451</span>  <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(lhs) == <a class="code" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb">viennacl::traits::handle</a>(rhs.rhs()))</div>
<div class="line"><a name="l00452"></a><span class="lineno"> 452</span>  {</div>
<div class="line"><a name="l00453"></a><span class="lineno"> 453</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp(lhs);</div>
<div class="line"><a name="l00454"></a><span class="lineno"> 454</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), T(1), temp, T(0));</div>
<div class="line"><a name="l00455"></a><span class="lineno"> 455</span>  lhs -= temp;</div>
<div class="line"><a name="l00456"></a><span class="lineno"> 456</span>  }</div>
<div class="line"><a name="l00457"></a><span class="lineno"> 457</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00458"></a><span class="lineno"> 458</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), T(-1), lhs, T(1));</div>
<div class="line"><a name="l00459"></a><span class="lineno"> 459</span>  }</div>
<div class="line"><a name="l00460"></a><span class="lineno"> 460</span>  };</div>
<div class="line"><a name="l00461"></a><span class="lineno"> 461</span> </div>
<div class="line"><a name="l00462"></a><span class="lineno"> 462</span> </div>
<div class="line"><a name="l00463"></a><span class="lineno"> 463</span>  <span class="comment">// x = A * vec_op</span></div>
<div class="line"><a name="l00464"></a><span class="lineno"> 464</span>  <span class="keyword">template</span><<span class="keyword">typename</span> T, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> A, <span class="keyword">typename</span> LHS, <span class="keyword">typename</span> RHS, <span class="keyword">typename</span> OP></div>
<div class="line"><a name="l00465"></a><span class="lineno"> 465</span>  <span class="keyword">struct </span>op_executor<vector_base<T>, <a class="code" href="structop__assign.html">op_assign</a>, vector_expression<const coordinate_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> ></div>
<div class="line"><a name="l00466"></a><span class="lineno"> 466</span>  {</div>
<div class="line"><a name="l00467"></a><span class="lineno"> 467</span>  <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base<T> & lhs, vector_expression<<span class="keyword">const</span> coordinate_matrix<T, A>, <span class="keyword">const</span> vector_expression<const LHS, const RHS, OP>, op_prod> <span class="keyword">const</span> & rhs)</div>
<div class="line"><a name="l00468"></a><span class="lineno"> 468</span>  {</div>
<div class="line"><a name="l00469"></a><span class="lineno"> 469</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp(rhs.rhs(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(rhs));</div>
<div class="line"><a name="l00470"></a><span class="lineno"> 470</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), temp, lhs);</div>
<div class="line"><a name="l00471"></a><span class="lineno"> 471</span>  }</div>
<div class="line"><a name="l00472"></a><span class="lineno"> 472</span>  };</div>
<div class="line"><a name="l00473"></a><span class="lineno"> 473</span> </div>
<div class="line"><a name="l00474"></a><span class="lineno"> 474</span>  <span class="comment">// x += A * vec_op</span></div>
<div class="line"><a name="l00475"></a><span class="lineno"> 475</span>  <span class="keyword">template</span><<span class="keyword">typename</span> T, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> A, <span class="keyword">typename</span> LHS, <span class="keyword">typename</span> RHS, <span class="keyword">typename</span> OP></div>
<div class="line"><a name="l00476"></a><span class="lineno"> 476</span>  <span class="keyword">struct </span>op_executor<vector_base<T>, op_inplace_add, vector_expression<const coordinate_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> ></div>
<div class="line"><a name="l00477"></a><span class="lineno"> 477</span>  {</div>
<div class="line"><a name="l00478"></a><span class="lineno"> 478</span>  <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base<T> & lhs, vector_expression<<span class="keyword">const</span> coordinate_matrix<T, A>, <span class="keyword">const</span> vector_expression<const LHS, const RHS, OP>, op_prod> <span class="keyword">const</span> & rhs)</div>
<div class="line"><a name="l00479"></a><span class="lineno"> 479</span>  {</div>
<div class="line"><a name="l00480"></a><span class="lineno"> 480</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp(rhs.rhs(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(rhs));</div>
<div class="line"><a name="l00481"></a><span class="lineno"> 481</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp_result(lhs);</div>
<div class="line"><a name="l00482"></a><span class="lineno"> 482</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), temp, temp_result);</div>
<div class="line"><a name="l00483"></a><span class="lineno"> 483</span>  lhs += temp_result;</div>
<div class="line"><a name="l00484"></a><span class="lineno"> 484</span>  }</div>
<div class="line"><a name="l00485"></a><span class="lineno"> 485</span>  };</div>
<div class="line"><a name="l00486"></a><span class="lineno"> 486</span> </div>
<div class="line"><a name="l00487"></a><span class="lineno"> 487</span>  <span class="comment">// x -= A * vec_op</span></div>
<div class="line"><a name="l00488"></a><span class="lineno"> 488</span>  <span class="keyword">template</span><<span class="keyword">typename</span> T, <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> A, <span class="keyword">typename</span> LHS, <span class="keyword">typename</span> RHS, <span class="keyword">typename</span> OP></div>
<div class="line"><a name="l00489"></a><span class="lineno"> 489</span>  <span class="keyword">struct </span>op_executor<vector_base<T>, op_inplace_sub, vector_expression<const coordinate_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> ></div>
<div class="line"><a name="l00490"></a><span class="lineno"> 490</span>  {</div>
<div class="line"><a name="l00491"></a><span class="lineno"> 491</span>  <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base<T> & lhs, vector_expression<<span class="keyword">const</span> coordinate_matrix<T, A>, <span class="keyword">const</span> vector_expression<const LHS, const RHS, OP>, op_prod> <span class="keyword">const</span> & rhs)</div>
<div class="line"><a name="l00492"></a><span class="lineno"> 492</span>  {</div>
<div class="line"><a name="l00493"></a><span class="lineno"> 493</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp(rhs.rhs(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(rhs));</div>
<div class="line"><a name="l00494"></a><span class="lineno"> 494</span>  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<T></a> temp_result(lhs);</div>
<div class="line"><a name="l00495"></a><span class="lineno"> 495</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), temp, temp_result);</div>
<div class="line"><a name="l00496"></a><span class="lineno"> 496</span>  lhs -= temp_result;</div>
<div class="line"><a name="l00497"></a><span class="lineno"> 497</span>  }</div>
<div class="line"><a name="l00498"></a><span class="lineno"> 498</span>  };</div>
<div class="line"><a name="l00499"></a><span class="lineno"> 499</span> </div>
<div class="line"><a name="l00500"></a><span class="lineno"> 500</span> } <span class="comment">// namespace detail</span></div>
<div class="line"><a name="l00501"></a><span class="lineno"> 501</span> } <span class="comment">// namespace linalg</span></div>
<div class="line"><a name="l00502"></a><span class="lineno"> 502</span> </div>
<div class="line"><a name="l00504"></a><span class="lineno"> 504</span> }</div>
<div class="line"><a name="l00505"></a><span class="lineno"> 505</span> </div>
<div class="line"><a name="l00506"></a><span class="lineno"> 506</span> <span class="preprocessor">#endif</span></div>
<div class="ttc" id="classviennacl_1_1backend_1_1typesafe__host__array_html"><div class="ttname"><a href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array</a></div><div class="ttdoc">Helper class implementing an array on the host. Default case: No conversion necessary. </div><div class="ttdef"><b>Definition:</b> <a href="backend_2util_8hpp_source.html#l00092">util.hpp:92</a></div></div>
<div class="ttc" id="classviennacl_1_1backend_1_1typesafe__host__array_html_ac8646f4841c47048a1d2e004b9536af8"><div class="ttname"><a href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac8646f4841c47048a1d2e004b9536af8">viennacl::backend::typesafe_host_array::element_size</a></div><div class="ttdeci">vcl_size_t element_size() const </div><div class="ttdef"><b>Definition:</b> <a href="backend_2util_8hpp_source.html#l00112">util.hpp:112</a></div></div>
<div class="ttc" id="classviennacl_1_1scalar_html"><div class="ttname"><a href="classviennacl_1_1scalar.html">viennacl::scalar</a></div><div class="ttdoc">This class represents a single scalar value on the GPU and behaves mostly like a built-in scalar type...</div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00227">forwards.h:227</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a875b6ccacbf0df8ef48ab7f7c46244f1"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a875b6ccacbf0df8ef48ab7f7c46244f1">viennacl::coordinate_matrix::handle12</a></div><div class="ttdeci">const handle_type & handle12() const </div><div class="ttdoc">Returns the OpenCL handle to the (row, column) index array. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00366">coordinate_matrix.hpp:366</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1traits_html_aa756f5d6820722094cae0d8b9bb6d5e2"><div class="ttname"><a href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a></div><div class="ttdeci">vcl_size_t size1(MatrixType const &mat)</div><div class="ttdoc">Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.) </div><div class="ttdef"><b>Definition:</b> <a href="size_8hpp_source.html#l00163">size.hpp:163</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a9cde8a8c5cb1516c5503d011388daed4"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a9cde8a8c5cb1516c5503d011388daed4">viennacl::coordinate_matrix::coordinate_matrix</a></div><div class="ttdeci">coordinate_matrix(viennacl::context ctx)</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00196">coordinate_matrix.hpp:196</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a626ac84303ae01ac4f7e964f2c2972b6"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a626ac84303ae01ac4f7e964f2c2972b6">viennacl::coordinate_matrix::coordinate_matrix</a></div><div class="ttdeci">coordinate_matrix()</div><div class="ttdoc">Default construction of a coordinate matrix. No memory is allocated. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00194">coordinate_matrix.hpp:194</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="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="classviennacl_1_1coordinate__matrix_html_a04c51b826e4bae9c97a3749ee2161b87"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a04c51b826e4bae9c97a3749ee2161b87">viennacl::coordinate_matrix::nnz</a></div><div class="ttdeci">vcl_size_t nnz() const </div><div class="ttdoc">Returns the number of nonzero entries. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00361">coordinate_matrix.hpp:361</a></div></div>
<div class="ttc" id="classviennacl_1_1backend_1_1typesafe__host__array_html_a1727eeb845f50cf5bd51ef54c5938b8c"><div class="ttname"><a href="classviennacl_1_1backend_1_1typesafe__host__array.html#a1727eeb845f50cf5bd51ef54c5938b8c">viennacl::backend::typesafe_host_array::get</a></div><div class="ttdeci">void * get()</div><div class="ttdef"><b>Definition:</b> <a href="backend_2util_8hpp_source.html#l00110">util.hpp:110</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1backend_html_a62854cfd6f04404b274f8ede36f63e2d"><div class="ttname"><a href="namespaceviennacl_1_1backend.html#a62854cfd6f04404b274f8ede36f63e2d">viennacl::backend::memory_read</a></div><div class="ttdeci">void memory_read(mem_handle const &src_buffer, vcl_size_t src_offset, vcl_size_t bytes_to_read, void *ptr, bool async=false)</div><div class="ttdoc">Reads data from a buffer back to main RAM. </div><div class="ttdef"><b>Definition:</b> <a href="memory_8hpp_source.html#l00261">memory.hpp:261</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a5363ed7b82aaad8549c85d880755f846"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a5363ed7b82aaad8549c85d880755f846">viennacl::coordinate_matrix::resize</a></div><div class="ttdeci">void resize(vcl_size_t new_size1, vcl_size_t new_size2, bool preserve=true)</div><div class="ttdoc">Resize the matrix. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00297">coordinate_matrix.hpp:297</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1traits_html_a3658e7c29ac0f60a20cb5871f5b5fd98"><div class="ttname"><a href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98">viennacl::traits::size2</a></div><div class="ttdeci">result_of::size_type< MatrixType >::type size2(MatrixType const &mat)</div><div class="ttdoc">Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...</div><div class="ttdef"><b>Definition:</b> <a href="size_8hpp_source.html#l00201">size.hpp:201</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1backend_1_1detail_html_af9f9cc6c50debc5b1d93e458d1fd392a"><div class="ttname"><a href="namespaceviennacl_1_1backend_1_1detail.html#af9f9cc6c50debc5b1d93e458d1fd392a">viennacl::backend::detail::element_size</a></div><div class="ttdeci">vcl_size_t element_size(memory_types)</div><div class="ttdef"><b>Definition:</b> <a href="memory_8hpp_source.html#l00299">memory.hpp:299</a></div></div>
<div class="ttc" id="tests_2src_2bisect_8cpp_html_a52b5d30a2d7b064678644a3bf49b7f6c"><div class="ttname"><a href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a></div><div class="ttdeci">float NumericT</div><div class="ttdef"><b>Definition:</b> <a href="tests_2src_2bisect_8cpp_source.html#l00040">bisect.cpp:40</a></div></div>
<div class="ttc" id="classviennacl_1_1context_html"><div class="ttname"><a href="classviennacl_1_1context.html">viennacl::context</a></div><div class="ttdoc">Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...</div><div class="ttdef"><b>Definition:</b> <a href="context_8hpp_source.html#l00039">context.hpp:39</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a2035ac2e7acc39c21cf13bb6f9c2dc08"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a2035ac2e7acc39c21cf13bb6f9c2dc08">viennacl::coordinate_matrix::internal_nnz</a></div><div class="ttdeci">vcl_size_t internal_nnz() const </div><div class="ttdoc">Returns the number of internal nonzero entries. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00363">coordinate_matrix.hpp:363</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_ad5b57c8ab8bcc8d839aba53be7cb1463"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#ad5b57c8ab8bcc8d839aba53be7cb1463">viennacl::coordinate_matrix::coordinate_matrix</a></div><div class="ttdeci">coordinate_matrix(vcl_size_t rows, vcl_size_t cols, vcl_size_t nonzeros=0, viennacl::context ctx=viennacl::context())</div><div class="ttdoc">Construction of a coordinate matrix with the supplied number of rows and columns. If the number of no...</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00219">coordinate_matrix.hpp:219</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1traits_html_aa2344ea20469f55fbc15a8e9526494d0"><div class="ttname"><a href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">viennacl::traits::size</a></div><div class="ttdeci">vcl_size_t size(VectorType const &vec)</div><div class="ttdoc">Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.) </div><div class="ttdef"><b>Definition:</b> <a href="size_8hpp_source.html#l00239">size.hpp:239</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a4a8183e8f8e3a91e379d12ada4943b99"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a4a8183e8f8e3a91e379d12ada4943b99">viennacl::coordinate_matrix::size2</a></div><div class="ttdeci">vcl_size_t size2() const </div><div class="ttdoc">Returns the number of columns. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00359">coordinate_matrix.hpp:359</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a9fd9168600c330faa69cae2626509602"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a9fd9168600c330faa69cae2626509602">viennacl::coordinate_matrix::value_type</a></div><div class="ttdeci">scalar< typename viennacl::tools::CHECK_SCALAR_TEMPLATE_ARGUMENT< NumericT >::ResultType > value_type</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00190">coordinate_matrix.hpp:190</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_ad4218882f3bd97eeff6d81b575ea1e3c"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#ad4218882f3bd97eeff6d81b575ea1e3c">viennacl::coordinate_matrix::size_type</a></div><div class="ttdeci">vcl_size_t size_type</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00191">coordinate_matrix.hpp:191</a></div></div>
<div class="ttc" id="sparse__matrix__operations_8hpp_html"><div class="ttname"><a href="sparse__matrix__operations_8hpp.html">sparse_matrix_operations.hpp</a></div><div class="ttdoc">Implementations of operations using sparse matrices. </div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a5497d1e7fdce125b7a5bd12530002b4f"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a5497d1e7fdce125b7a5bd12530002b4f">viennacl::coordinate_matrix::handle_type</a></div><div class="ttdeci">viennacl::backend::mem_handle handle_type</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00189">coordinate_matrix.hpp:189</a></div></div>
<div class="ttc" id="classviennacl_1_1tools_1_1const__sparse__matrix__adapter_html"><div class="ttname"><a href="classviennacl_1_1tools_1_1const__sparse__matrix__adapter.html">viennacl::tools::const_sparse_matrix_adapter</a></div><div class="ttdoc">Adapts a constant sparse matrix type made up from std::vector<std::map<SizeT, NumericT> > to basic ub...</div><div class="ttdef"><b>Definition:</b> <a href="adapter_8hpp_source.html#l00183">adapter.hpp:183</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a4fccbe4157aa8833d23b0541803e06af"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a4fccbe4157aa8833d23b0541803e06af">viennacl::coordinate_matrix::clear</a></div><div class="ttdeci">void clear()</div><div class="ttdoc">Resets all entries in the matrix back to zero without changing the matrix size. Resets the sparsity p...</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00342">coordinate_matrix.hpp:342</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="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="classviennacl_1_1coordinate__matrix_html_a66b7f5c5487b719c91bc76984f89f97f"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a66b7f5c5487b719c91bc76984f89f97f">viennacl::coordinate_matrix::handle3</a></div><div class="ttdeci">const handle_type & handle3() const </div><div class="ttdoc">Returns the OpenCL handle to the group start index array. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00370">coordinate_matrix.hpp:370</a></div></div>
<div class="ttc" id="classviennacl_1_1context_html_aadf575d49ed144e8fa4707cdad69c349"><div class="ttname"><a href="classviennacl_1_1context.html#aadf575d49ed144e8fa4707cdad69c349">viennacl::context::memory_type</a></div><div class="ttdeci">viennacl::memory_types memory_type() const </div><div class="ttdef"><b>Definition:</b> <a href="context_8hpp_source.html#l00076">context.hpp:76</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a3cd8e9989ee1b547dbcba2fe87f7b29d"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a3cd8e9989ee1b547dbcba2fe87f7b29d">viennacl::coordinate_matrix::copy</a></div><div class="ttdeci">friend void copy(const CPUMatrixT &cpu_matrix, coordinate_matrix< NumericT2, AlignmentV2 > &gpu_matrix)</div></div>
<div class="ttc" id="classviennacl_1_1backend_1_1mem__handle_html_ac13e99f562746d6682116dec1e13a7da"><div class="ttname"><a href="classviennacl_1_1backend_1_1mem__handle.html#ac13e99f562746d6682116dec1e13a7da">viennacl::backend::mem_handle::switch_active_handle_id</a></div><div class="ttdeci">void switch_active_handle_id(memory_types new_id)</div><div class="ttdoc">Switches the currently active handle. If no support for that backend is provided, an exception is thr...</div><div class="ttdef"><b>Definition:</b> <a href="mem__handle_8hpp_source.html#l00121">mem_handle.hpp:121</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1backend_html_a39a91fd7aa68ebf64bd7bfa409732ec0"><div class="ttname"><a href="namespaceviennacl_1_1backend.html#a39a91fd7aa68ebf64bd7bfa409732ec0">viennacl::backend::memory_copy</a></div><div class="ttdeci">void memory_copy(mem_handle const &src_buffer, mem_handle &dst_buffer, vcl_size_t src_offset, vcl_size_t dst_offset, vcl_size_t bytes_to_copy)</div><div class="ttdoc">Copies 'bytes_to_copy' bytes from address 'src_buffer + src_offset' to memory starting at address 'ds...</div><div class="ttdef"><b>Definition:</b> <a href="memory_8hpp_source.html#l00140">memory.hpp:140</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_ad091412ca52bb131234464ed3f5b103c"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#ad091412ca52bb131234464ed3f5b103c">viennacl::coordinate_matrix::reserve</a></div><div class="ttdeci">void reserve(vcl_size_t new_nonzeros)</div><div class="ttdoc">Allocate memory for the supplied number of nonzeros in the matrix. Old values are preserved...</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00270">coordinate_matrix.hpp:270</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1traits_html_a6707f5dab8f482170d2046a605f46ef8"><div class="ttname"><a href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a></div><div class="ttdeci">viennacl::context context(T const &t)</div><div class="ttdoc">Returns an ID for the currently active memory domain of an object. </div><div class="ttdef"><b>Definition:</b> <a href="traits_2context_8hpp_source.html#l00040">context.hpp:40</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_aeb74dd2ab673d7912a17f60777e9bf6c"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#aeb74dd2ab673d7912a17f60777e9bf6c">viennacl::coordinate_matrix::handle</a></div><div class="ttdeci">const handle_type & handle() const </div><div class="ttdoc">Returns the OpenCL handle to the matrix entry array. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00368">coordinate_matrix.hpp:368</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="namespaceviennacl_html_a10b7f8cf6b8864a7aa196d670481a453"><div class="ttname"><a href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a></div><div class="ttdeci">void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)</div><div class="ttdoc">Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...</div><div class="ttdef"><b>Definition:</b> <a href="circulant__matrix_8hpp_source.html#l00150">circulant_matrix.hpp:150</a></div></div>
<div class="ttc" id="classviennacl_1_1backend_1_1typesafe__host__array_html_ac31bfde381eb40e0349678252c444af9"><div class="ttname"><a href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac31bfde381eb40e0349678252c444af9">viennacl::backend::typesafe_host_array::set</a></div><div class="ttdeci">void set(vcl_size_t index, U value)</div><div class="ttdef"><b>Definition:</b> <a href="backend_2util_8hpp_source.html#l00115">util.hpp:115</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a7f83ad82b3bab755d76e7f9d665367f7"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a7f83ad82b3bab755d76e7f9d665367f7">viennacl::coordinate_matrix::groups</a></div><div class="ttdeci">vcl_size_t groups() const </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00372">coordinate_matrix.hpp:372</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="classviennacl_1_1backend_1_1mem__handle_html_ac8373f0d899b89c843e14de4cb7a1c4a"><div class="ttname"><a href="classviennacl_1_1backend_1_1mem__handle.html#ac8373f0d899b89c843e14de4cb7a1c4a">viennacl::backend::mem_handle::raw_size</a></div><div class="ttdeci">vcl_size_t raw_size() const </div><div class="ttdoc">Returns the number of bytes of the currently active buffer. </div><div class="ttdef"><b>Definition:</b> <a href="mem__handle_8hpp_source.html#l00230">mem_handle.hpp:230</a></div></div>
<div class="ttc" id="classviennacl_1_1tools_1_1sparse__matrix__adapter_html"><div class="ttname"><a href="classviennacl_1_1tools_1_1sparse__matrix__adapter.html">viennacl::tools::sparse_matrix_adapter</a></div><div class="ttdoc">Adapts a non-const sparse matrix type made up from std::vector<std::map<SizeT, NumericT> > to basic u...</div><div class="ttdef"><b>Definition:</b> <a href="adapter_8hpp_source.html#l00357">adapter.hpp:357</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html_a472604b5e8e7a818fcee171898df1c92"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a472604b5e8e7a818fcee171898df1c92">viennacl::coordinate_matrix::coordinate_matrix</a></div><div class="ttdeci">coordinate_matrix(vcl_size_t rows, vcl_size_t cols, viennacl::context ctx)</div><div class="ttdoc">Construction of a coordinate matrix with the supplied number of rows and columns in the supplied cont...</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00251">coordinate_matrix.hpp:251</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1backend_html_a1499f19634964e2c7c8aeeefc6206126"><div class="ttname"><a href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a></div><div class="ttdeci">void memory_create(mem_handle &handle, vcl_size_t size_in_bytes, viennacl::context const &ctx, const void *host_ptr=NULL)</div><div class="ttdoc">Creates an array of the specified size. If the second argument is provided, the buffer is initialized...</div><div class="ttdef"><b>Definition:</b> <a href="memory_8hpp_source.html#l00087">memory.hpp:87</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="classviennacl_1_1coordinate__matrix_html_a734ed32ff1858ef5453c170326fda377"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html#a734ed32ff1858ef5453c170326fda377">viennacl::coordinate_matrix::size1</a></div><div class="ttdeci">vcl_size_t size1() const </div><div class="ttdoc">Returns the number of rows. </div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00357">coordinate_matrix.hpp:357</a></div></div>
<div class="ttc" id="structop__assign_html"><div class="ttname"><a href="structop__assign.html">op_assign</a></div><div class="ttdef"><b>Definition:</b> <a href="self__assign_8cpp_source.html#l00132">self_assign.cpp:132</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1backend_html_a8548a37535eac92e886cb183182b618f"><div class="ttname"><a href="namespaceviennacl_1_1backend.html#a8548a37535eac92e886cb183182b618f">viennacl::backend::memory_shallow_copy</a></div><div class="ttdeci">void memory_shallow_copy(mem_handle const &src_buffer, mem_handle &dst_buffer)</div><div class="ttdoc">A 'shallow' copy operation from an initialized buffer to an uninitialized buffer. The uninitialized b...</div><div class="ttdef"><b>Definition:</b> <a href="memory_8hpp_source.html#l00177">memory.hpp:177</a></div></div>
<div class="ttc" id="classviennacl_1_1coordinate__matrix_html"><div class="ttname"><a href="classviennacl_1_1coordinate__matrix.html">viennacl::coordinate_matrix</a></div><div class="ttdoc">A sparse square matrix, where entries are stored as triplets (i,j, val), where i and j are the row an...</div><div class="ttdef"><b>Definition:</b> <a href="coordinate__matrix_8hpp_source.html#l00186">coordinate_matrix.hpp:186</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="coordinate__matrix_8hpp.html">coordinate_matrix.hpp</a></li>
<li class="footer">Generated on Wed Jan 20 2016 22:32:39 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>
|