File: sliced__ell__matrix_8hpp_source.html

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

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

<div class="header">
  <div class="headertitle">
<div class="title">sliced_ell_matrix.hpp</div>  </div>
</div><!--header-->
<div class="contents">
<a href="sliced__ell__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>&#160;<span class="preprocessor">#ifndef VIENNACL_SLICED_ELL_MATRIX_HPP_</span></div>
<div class="line"><a name="l00002"></a><span class="lineno">    2</span>&#160;<span class="preprocessor">#define VIENNACL_SLICED_ELL_MATRIX_HPP_</span></div>
<div class="line"><a name="l00003"></a><span class="lineno">    3</span>&#160;</div>
<div class="line"><a name="l00004"></a><span class="lineno">    4</span>&#160;<span class="comment">/* =========================================================================</span></div>
<div class="line"><a name="l00005"></a><span class="lineno">    5</span>&#160;<span class="comment">   Copyright (c) 2010-2016, Institute for Microelectronics,</span></div>
<div class="line"><a name="l00006"></a><span class="lineno">    6</span>&#160;<span class="comment">                            Institute for Analysis and Scientific Computing,</span></div>
<div class="line"><a name="l00007"></a><span class="lineno">    7</span>&#160;<span class="comment">                            TU Wien.</span></div>
<div class="line"><a name="l00008"></a><span class="lineno">    8</span>&#160;<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>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00010"></a><span class="lineno">   10</span>&#160;<span class="comment">                            -----------------</span></div>
<div class="line"><a name="l00011"></a><span class="lineno">   11</span>&#160;<span class="comment">                  ViennaCL - The Vienna Computing Library</span></div>
<div class="line"><a name="l00012"></a><span class="lineno">   12</span>&#160;<span class="comment">                            -----------------</span></div>
<div class="line"><a name="l00013"></a><span class="lineno">   13</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00014"></a><span class="lineno">   14</span>&#160;<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>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00016"></a><span class="lineno">   16</span>&#160;<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>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00018"></a><span class="lineno">   18</span>&#160;<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>&#160;<span class="comment">============================================================================= */</span></div>
<div class="line"><a name="l00020"></a><span class="lineno">   20</span>&#160;</div>
<div class="line"><a name="l00028"></a><span class="lineno">   28</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="forwards_8h.html">viennacl/forwards.h</a>&quot;</span></div>
<div class="line"><a name="l00029"></a><span class="lineno">   29</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>&quot;</span></div>
<div class="line"><a name="l00030"></a><span class="lineno">   30</span>&#160;</div>
<div class="line"><a name="l00031"></a><span class="lineno">   31</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="tools_8hpp.html">viennacl/tools/tools.hpp</a>&quot;</span></div>
<div class="line"><a name="l00032"></a><span class="lineno">   32</span>&#160;</div>
<div class="line"><a name="l00033"></a><span class="lineno">   33</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="sparse__matrix__operations_8hpp.html">viennacl/linalg/sparse_matrix_operations.hpp</a>&quot;</span></div>
<div class="line"><a name="l00034"></a><span class="lineno">   34</span>&#160;</div>
<div class="line"><a name="l00035"></a><span class="lineno">   35</span>&#160;<span class="keyword">namespace </span>viennacl</div>
<div class="line"><a name="l00036"></a><span class="lineno">   36</span>&#160;{</div>
<div class="line"><a name="l00045"></a><span class="lineno">   45</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT <span class="comment">/* see forwards.h = unsigned int */</span>&gt;</div>
<div class="line"><a name="l00046"></a><span class="lineno">   46</span>&#160;<span class="keyword">class </span>sliced_ell_matrix</div>
<div class="line"><a name="l00047"></a><span class="lineno">   47</span>&#160;{</div>
<div class="line"><a name="l00048"></a><span class="lineno">   48</span>&#160;<span class="keyword">public</span>:</div>
<div class="line"><a name="l00049"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a2cdf3aeccfe9a10b89c48b3430cb4c93">   49</a></span>&#160;  <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_1sliced__ell__matrix.html#a2cdf3aeccfe9a10b89c48b3430cb4c93">handle_type</a>;</div>
<div class="line"><a name="l00050"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#af0672b686ad27ddcf948116bae059329">   50</a></span>&#160;  <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1scalar.html">scalar&lt;typename viennacl::tools::CHECK_SCALAR_TEMPLATE_ARGUMENT&lt;ScalarT&gt;::ResultType</a>&gt;   <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#af0672b686ad27ddcf948116bae059329">value_type</a>;</div>
<div class="line"><a name="l00051"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#adf0fb80bb38ab0256a622960c0c18699">   51</a></span>&#160;  <span class="keyword">typedef</span> <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a>                                                                              <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#adf0fb80bb38ab0256a622960c0c18699">size_type</a>;</div>
<div class="line"><a name="l00052"></a><span class="lineno">   52</span>&#160;</div>
<div class="line"><a name="l00053"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a362dc0c5f0d644cbcf6514f3b3ee6529">   53</a></span>&#160;  <span class="keyword">explicit</span> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a362dc0c5f0d644cbcf6514f3b3ee6529">sliced_ell_matrix</a>() : rows_(0), cols_(0), rows_per_block_(0) {}</div>
<div class="line"><a name="l00054"></a><span class="lineno">   54</span>&#160;</div>
<div class="line"><a name="l00059"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a9f244082eff24e93495d3b4abec735b0">   59</a></span>&#160;  <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a9f244082eff24e93495d3b4abec735b0">sliced_ell_matrix</a>(<a class="code" href="classviennacl_1_1sliced__ell__matrix.html#adf0fb80bb38ab0256a622960c0c18699">size_type</a> num_rows,</div>
<div class="line"><a name="l00060"></a><span class="lineno">   60</span>&#160;                    <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#adf0fb80bb38ab0256a622960c0c18699">size_type</a> num_cols,</div>
<div class="line"><a name="l00061"></a><span class="lineno">   61</span>&#160;                    <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#adf0fb80bb38ab0256a622960c0c18699">size_type</a> num_rows_per_block_ = 0)</div>
<div class="line"><a name="l00062"></a><span class="lineno">   62</span>&#160;    : rows_(num_rows),</div>
<div class="line"><a name="l00063"></a><span class="lineno">   63</span>&#160;      cols_(num_cols),</div>
<div class="line"><a name="l00064"></a><span class="lineno">   64</span>&#160;      rows_per_block_(num_rows_per_block_) {}</div>
<div class="line"><a name="l00065"></a><span class="lineno">   65</span>&#160;</div>
<div class="line"><a name="l00066"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#aa725bf2d5b213f119bd4b9fe9258c201">   66</a></span>&#160;  <span class="keyword">explicit</span> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#aa725bf2d5b213f119bd4b9fe9258c201">sliced_ell_matrix</a>(<a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx) : rows_(0), cols_(0), rows_per_block_(0)</div>
<div class="line"><a name="l00067"></a><span class="lineno">   67</span>&#160;  {</div>
<div class="line"><a name="l00068"></a><span class="lineno">   68</span>&#160;    columns_per_block_.<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="l00069"></a><span class="lineno">   69</span>&#160;    column_indices_.<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="l00070"></a><span class="lineno">   70</span>&#160;    block_start_.<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="l00071"></a><span class="lineno">   71</span>&#160;    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="l00072"></a><span class="lineno">   72</span>&#160;</div>
<div class="line"><a name="l00073"></a><span class="lineno">   73</span>&#160;<span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00074"></a><span class="lineno">   74</span>&#160;    <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="l00075"></a><span class="lineno">   75</span>&#160;    {</div>
<div class="line"><a name="l00076"></a><span class="lineno">   76</span>&#160;      columns_per_block_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00077"></a><span class="lineno">   77</span>&#160;      column_indices_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00078"></a><span class="lineno">   78</span>&#160;      block_start_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00079"></a><span class="lineno">   79</span>&#160;      elements_.opencl_handle().context(ctx.opencl_context());</div>
<div class="line"><a name="l00080"></a><span class="lineno">   80</span>&#160;    }</div>
<div class="line"><a name="l00081"></a><span class="lineno">   81</span>&#160;<span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00082"></a><span class="lineno">   82</span>&#160;  }</div>
<div class="line"><a name="l00083"></a><span class="lineno">   83</span>&#160;</div>
<div class="line"><a name="l00085"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#ac681e3c8d41b686282711644cbdcadfb">   85</a></span>&#160;  <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#ac681e3c8d41b686282711644cbdcadfb">clear</a>()</div>
<div class="line"><a name="l00086"></a><span class="lineno">   86</span>&#160;  {</div>
<div class="line"><a name="l00087"></a><span class="lineno">   87</span>&#160;    <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array&lt;IndexT&gt;</a> host_columns_per_block_buffer(columns_per_block_, rows_ / rows_per_block_ + 1);</div>
<div class="line"><a name="l00088"></a><span class="lineno">   88</span>&#160;    <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array&lt;IndexT&gt;</a> host_column_buffer(column_indices_, <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#ab8a2233384ac8dcdb85b5385706548a1">internal_size1</a>());</div>
<div class="line"><a name="l00089"></a><span class="lineno">   89</span>&#160;    <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array&lt;IndexT&gt;</a> host_block_start_buffer(block_start_, (rows_ - 1) / rows_per_block_ + 1);</div>
<div class="line"><a name="l00090"></a><span class="lineno">   90</span>&#160;    std::vector&lt;ScalarT&gt; host_elements(1);</div>
<div class="line"><a name="l00091"></a><span class="lineno">   91</span>&#160;</div>
<div class="line"><a name="l00092"></a><span class="lineno">   92</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(columns_per_block_, host_columns_per_block_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac8646f4841c47048a1d2e004b9536af8">element_size</a>() * (rows_ / rows_per_block_ + 1), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(columns_per_block_), host_columns_per_block_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#a1727eeb845f50cf5bd51ef54c5938b8c">get</a>());</div>
<div class="line"><a name="l00093"></a><span class="lineno">   93</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(column_indices_,    host_column_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac8646f4841c47048a1d2e004b9536af8">element_size</a>() * <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#ab8a2233384ac8dcdb85b5385706548a1">internal_size1</a>(),                         <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(column_indices_),    host_column_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#a1727eeb845f50cf5bd51ef54c5938b8c">get</a>());</div>
<div class="line"><a name="l00094"></a><span class="lineno">   94</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(block_start_,       host_block_start_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac8646f4841c47048a1d2e004b9536af8">element_size</a>() * ((rows_ - 1) / rows_per_block_ + 1), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(block_start_),       host_block_start_buffer.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#a1727eeb845f50cf5bd51ef54c5938b8c">get</a>());</div>
<div class="line"><a name="l00095"></a><span class="lineno">   95</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(elements_,          <span class="keyword">sizeof</span>(ScalarT) * 1,                                                          <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(elements_),          &amp;(host_elements[0]));</div>
<div class="line"><a name="l00096"></a><span class="lineno">   96</span>&#160;  }</div>
<div class="line"><a name="l00097"></a><span class="lineno">   97</span>&#160;</div>
<div class="line"><a name="l00098"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#ab8a2233384ac8dcdb85b5385706548a1">   98</a></span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#ab8a2233384ac8dcdb85b5385706548a1">internal_size1</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> viennacl::tools::align_to_multiple&lt;vcl_size_t&gt;(rows_, rows_per_block_); }</div>
<div class="line"><a name="l00099"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#ad4fd3359611ca475b85d2041bdec1377">   99</a></span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#ad4fd3359611ca475b85d2041bdec1377">internal_size2</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> cols_; }</div>
<div class="line"><a name="l00100"></a><span class="lineno">  100</span>&#160;</div>
<div class="line"><a name="l00101"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a16cd8194ca5b263de32f131362498d77">  101</a></span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a16cd8194ca5b263de32f131362498d77">size1</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> rows_; }</div>
<div class="line"><a name="l00102"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#ae0afefc19f6dc7692a87cd8a16b51d30">  102</a></span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#ae0afefc19f6dc7692a87cd8a16b51d30">size2</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> cols_; }</div>
<div class="line"><a name="l00103"></a><span class="lineno">  103</span>&#160;</div>
<div class="line"><a name="l00104"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#aa9a7f4003687cae9a2af97f9fa328db5">  104</a></span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#aa9a7f4003687cae9a2af97f9fa328db5">rows_per_block</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> rows_per_block_; }</div>
<div class="line"><a name="l00105"></a><span class="lineno">  105</span>&#160;</div>
<div class="line"><a name="l00106"></a><span class="lineno">  106</span>&#160;  <span class="comment">//vcl_size_t nnz() const { return rows_ * maxnnz_; }</span></div>
<div class="line"><a name="l00107"></a><span class="lineno">  107</span>&#160;  <span class="comment">//vcl_size_t internal_nnz() const { return internal_size1() * internal_maxnnz(); }</span></div>
<div class="line"><a name="l00108"></a><span class="lineno">  108</span>&#160;</div>
<div class="line"><a name="l00109"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a94043424c01974b2db8062f74f625825">  109</a></span>&#160;  <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a94043424c01974b2db8062f74f625825">handle1</a>()       { <span class="keywordflow">return</span> columns_per_block_; }</div>
<div class="line"><a name="l00110"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a727fd3a25232545c59204885ff3e32d4">  110</a></span>&#160;  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a727fd3a25232545c59204885ff3e32d4">handle1</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> columns_per_block_; }</div>
<div class="line"><a name="l00111"></a><span class="lineno">  111</span>&#160;</div>
<div class="line"><a name="l00112"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a31d5fa9ec25dcfeb260f7e11c4db6571">  112</a></span>&#160;  <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a31d5fa9ec25dcfeb260f7e11c4db6571">handle2</a>()       { <span class="keywordflow">return</span> column_indices_; }</div>
<div class="line"><a name="l00113"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a45249ad5b44d0f506adc27e1e18aecdf">  113</a></span>&#160;  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a45249ad5b44d0f506adc27e1e18aecdf">handle2</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> column_indices_; }</div>
<div class="line"><a name="l00114"></a><span class="lineno">  114</span>&#160;</div>
<div class="line"><a name="l00115"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#aa59ae796d592a679a4c6036ba84f61b0">  115</a></span>&#160;  <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#aa59ae796d592a679a4c6036ba84f61b0">handle3</a>()       { <span class="keywordflow">return</span> block_start_; }</div>
<div class="line"><a name="l00116"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#aec497ed037ef674d36729bda4db7e93f">  116</a></span>&#160;  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#aec497ed037ef674d36729bda4db7e93f">handle3</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> block_start_; }</div>
<div class="line"><a name="l00117"></a><span class="lineno">  117</span>&#160;</div>
<div class="line"><a name="l00118"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#abbc671358c44037b7c2b22a01f17c3a7">  118</a></span>&#160;  <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#abbc671358c44037b7c2b22a01f17c3a7">handle</a>()       { <span class="keywordflow">return</span> elements_; }</div>
<div class="line"><a name="l00119"></a><span class="lineno"><a class="line" href="classviennacl_1_1sliced__ell__matrix.html#a1b9a6e9229a8694123562166eb2d4e32">  119</a></span>&#160;  <span class="keyword">const</span> <a class="code" href="classviennacl_1_1backend_1_1mem__handle.html">handle_type</a> &amp; <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a1b9a6e9229a8694123562166eb2d4e32">handle</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> elements_; }</div>
<div class="line"><a name="l00120"></a><span class="lineno">  120</span>&#160;</div>
<div class="line"><a name="l00121"></a><span class="lineno">  121</span>&#160;<span class="preprocessor">#if defined(_MSC_VER) &amp;&amp; _MSC_VER &lt; 1500          //Visual Studio 2005 needs special treatment</span></div>
<div class="line"><a name="l00122"></a><span class="lineno">  122</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> CPUMatrixT&gt;</div>
<div class="line"><a name="l00123"></a><span class="lineno">  123</span>&#160;  <span class="keyword">friend</span> <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a40c479e880ef7e7f50ed67cda91bb8e7">copy</a>(CPUMatrixT <span class="keyword">const</span> &amp; cpu_matrix, <a class="code" href="classviennacl_1_1sliced__ell__matrix.html">sliced_ell_matrix</a> &amp; gpu_matrix );</div>
<div class="line"><a name="l00124"></a><span class="lineno">  124</span>&#160;<span class="preprocessor">#else</span></div>
<div class="line"><a name="l00125"></a><span class="lineno">  125</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> CPUMatrixT, <span class="keyword">typename</span> ScalarT2, <span class="keyword">typename</span> IndexT2&gt;</div>
<div class="line"><a name="l00126"></a><span class="lineno">  126</span>&#160;  <span class="keyword">friend</span> <span class="keywordtype">void</span> <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a40c479e880ef7e7f50ed67cda91bb8e7">copy</a>(CPUMatrixT <span class="keyword">const</span> &amp; cpu_matrix, <a class="code" href="classviennacl_1_1sliced__ell__matrix.html">sliced_ell_matrix&lt;ScalarT2, IndexT2&gt;</a> &amp; gpu_matrix );</div>
<div class="line"><a name="l00127"></a><span class="lineno">  127</span>&#160;<span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00128"></a><span class="lineno">  128</span>&#160;</div>
<div class="line"><a name="l00129"></a><span class="lineno">  129</span>&#160;<span class="keyword">private</span>:</div>
<div class="line"><a name="l00130"></a><span class="lineno">  130</span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> rows_;</div>
<div class="line"><a name="l00131"></a><span class="lineno">  131</span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> cols_;</div>
<div class="line"><a name="l00132"></a><span class="lineno">  132</span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> rows_per_block_; <span class="comment">//parameter C in the paper by Kreutzer et al.</span></div>
<div class="line"><a name="l00133"></a><span class="lineno">  133</span>&#160;</div>
<div class="line"><a name="l00134"></a><span class="lineno">  134</span>&#160;  <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a2cdf3aeccfe9a10b89c48b3430cb4c93">handle_type</a> columns_per_block_;</div>
<div class="line"><a name="l00135"></a><span class="lineno">  135</span>&#160;  <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a2cdf3aeccfe9a10b89c48b3430cb4c93">handle_type</a> column_indices_;</div>
<div class="line"><a name="l00136"></a><span class="lineno">  136</span>&#160;  <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a2cdf3aeccfe9a10b89c48b3430cb4c93">handle_type</a> block_start_;</div>
<div class="line"><a name="l00137"></a><span class="lineno">  137</span>&#160;  <a class="code" href="classviennacl_1_1sliced__ell__matrix.html#a2cdf3aeccfe9a10b89c48b3430cb4c93">handle_type</a> elements_;</div>
<div class="line"><a name="l00138"></a><span class="lineno">  138</span>&#160;};</div>
<div class="line"><a name="l00139"></a><span class="lineno">  139</span>&#160;</div>
<div class="line"><a name="l00140"></a><span class="lineno">  140</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> CPUMatrixT, <span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT&gt;</div>
<div class="line"><a name="l00141"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#a40fb810a020b6bd780ee1e12e6d02763">  141</a></span>&#160;<span class="keywordtype">void</span> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(CPUMatrixT <span class="keyword">const</span> &amp; cpu_matrix, <a class="code" href="classviennacl_1_1sliced__ell__matrix.html">sliced_ell_matrix&lt;ScalarT, IndexT&gt;</a> &amp; gpu_matrix )</div>
<div class="line"><a name="l00142"></a><span class="lineno">  142</span>&#160;{</div>
<div class="line"><a name="l00143"></a><span class="lineno">  143</span>&#160;  assert( (gpu_matrix.size1() == 0 || <a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) == gpu_matrix.size1()) &amp;&amp; <span class="keywordtype">bool</span>(<span class="stringliteral">&quot;Size mismatch&quot;</span>) );</div>
<div class="line"><a name="l00144"></a><span class="lineno">  144</span>&#160;  assert( (gpu_matrix.size2() == 0 || <a class="code" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98">viennacl::traits::size2</a>(cpu_matrix) == gpu_matrix.size2()) &amp;&amp; <span class="keywordtype">bool</span>(<span class="stringliteral">&quot;Size mismatch&quot;</span>) );</div>
<div class="line"><a name="l00145"></a><span class="lineno">  145</span>&#160;</div>
<div class="line"><a name="l00146"></a><span class="lineno">  146</span>&#160;  <span class="keywordflow">if</span> (gpu_matrix.rows_per_block() == 0) <span class="comment">// not yet initialized by user. Set default: 32 is perfect for NVIDIA GPUs and older AMD GPUs. Still okay for newer AMD GPUs.</span></div>
<div class="line"><a name="l00147"></a><span class="lineno">  147</span>&#160;    gpu_matrix.rows_per_block_ = 32;</div>
<div class="line"><a name="l00148"></a><span class="lineno">  148</span>&#160;</div>
<div class="line"><a name="l00149"></a><span class="lineno">  149</span>&#160;  if (<a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) &gt; 0 &amp;&amp; <a class="code" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98">viennacl::traits::size2</a>(cpu_matrix) &gt; 0)</div>
<div class="line"><a name="l00150"></a><span class="lineno">  150</span>&#160;  {</div>
<div class="line"><a name="l00151"></a><span class="lineno">  151</span>&#160;    <span class="comment">//determine max capacity for row</span></div>
<div class="line"><a name="l00152"></a><span class="lineno">  152</span>&#160;    IndexT columns_in_current_block = 0;</div>
<div class="line"><a name="l00153"></a><span class="lineno">  153</span>&#160;    <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> total_element_buffer_size = 0;</div>
<div class="line"><a name="l00154"></a><span class="lineno">  154</span>&#160;    <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array&lt;IndexT&gt;</a> columns_in_block_buffer(gpu_matrix.handle1(), (<a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) - 1) / gpu_matrix.rows_per_block() + 1);</div>
<div class="line"><a name="l00155"></a><span class="lineno">  155</span>&#160;    <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="l00156"></a><span class="lineno">  156</span>&#160;    {</div>
<div class="line"><a name="l00157"></a><span class="lineno">  157</span>&#160;      <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> entries_in_row = 0;</div>
<div class="line"><a name="l00158"></a><span class="lineno">  158</span>&#160;      <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="l00159"></a><span class="lineno">  159</span>&#160;        ++entries_in_row;</div>
<div class="line"><a name="l00160"></a><span class="lineno">  160</span>&#160;</div>
<div class="line"><a name="l00161"></a><span class="lineno">  161</span>&#160;      columns_in_current_block = <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a5d46fe9558b0e462f10fd44942ad4fc6">std::max</a>(columns_in_current_block, static_cast&lt;IndexT&gt;(entries_in_row));</div>
<div class="line"><a name="l00162"></a><span class="lineno">  162</span>&#160;</div>
<div class="line"><a name="l00163"></a><span class="lineno">  163</span>&#160;      <span class="comment">// check for end of block</span></div>
<div class="line"><a name="l00164"></a><span class="lineno">  164</span>&#160;      <span class="keywordflow">if</span> ( (row_it.index1() % gpu_matrix.rows_per_block() == gpu_matrix.rows_per_block() - 1)</div>
<div class="line"><a name="l00165"></a><span class="lineno">  165</span>&#160;           || row_it.index1() == <a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) - 1)</div>
<div class="line"><a name="l00166"></a><span class="lineno">  166</span>&#160;      {</div>
<div class="line"><a name="l00167"></a><span class="lineno">  167</span>&#160;        total_element_buffer_size += columns_in_current_block * gpu_matrix.rows_per_block();</div>
<div class="line"><a name="l00168"></a><span class="lineno">  168</span>&#160;        columns_in_block_buffer.set(row_it.index1() / gpu_matrix.rows_per_block(), columns_in_current_block);</div>
<div class="line"><a name="l00169"></a><span class="lineno">  169</span>&#160;        columns_in_current_block = 0;</div>
<div class="line"><a name="l00170"></a><span class="lineno">  170</span>&#160;      }</div>
<div class="line"><a name="l00171"></a><span class="lineno">  171</span>&#160;    }</div>
<div class="line"><a name="l00172"></a><span class="lineno">  172</span>&#160;</div>
<div class="line"><a name="l00173"></a><span class="lineno">  173</span>&#160;    <span class="comment">//setup GPU matrix</span></div>
<div class="line"><a name="l00174"></a><span class="lineno">  174</span>&#160;    gpu_matrix.rows_ = cpu_matrix.size1();</div>
<div class="line"><a name="l00175"></a><span class="lineno">  175</span>&#160;    gpu_matrix.cols_ = cpu_matrix.size2();</div>
<div class="line"><a name="l00176"></a><span class="lineno">  176</span>&#160;</div>
<div class="line"><a name="l00177"></a><span class="lineno">  177</span>&#160;    <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array&lt;IndexT&gt;</a> coords(gpu_matrix.handle2(), total_element_buffer_size);</div>
<div class="line"><a name="l00178"></a><span class="lineno">  178</span>&#160;    <a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html">viennacl::backend::typesafe_host_array&lt;IndexT&gt;</a> block_start(gpu_matrix.handle3(), (<a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) - 1) / gpu_matrix.rows_per_block() + 1);</div>
<div class="line"><a name="l00179"></a><span class="lineno">  179</span>&#160;    std::vector&lt;ScalarT&gt; elements(total_element_buffer_size, 0);</div>
<div class="line"><a name="l00180"></a><span class="lineno">  180</span>&#160;</div>
<div class="line"><a name="l00181"></a><span class="lineno">  181</span>&#160;    <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> block_offset = 0;</div>
<div class="line"><a name="l00182"></a><span class="lineno">  182</span>&#160;    <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> block_index  = 0;</div>
<div class="line"><a name="l00183"></a><span class="lineno">  183</span>&#160;    <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> row_in_block = 0;</div>
<div class="line"><a name="l00184"></a><span class="lineno">  184</span>&#160;    <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="l00185"></a><span class="lineno">  185</span>&#160;    {</div>
<div class="line"><a name="l00186"></a><span class="lineno">  186</span>&#160;      <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> entry_in_row = 0;</div>
<div class="line"><a name="l00187"></a><span class="lineno">  187</span>&#160;</div>
<div class="line"><a name="l00188"></a><span class="lineno">  188</span>&#160;      <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="l00189"></a><span class="lineno">  189</span>&#160;      {</div>
<div class="line"><a name="l00190"></a><span class="lineno">  190</span>&#160;        <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> buffer_index = block_offset + entry_in_row * gpu_matrix.rows_per_block() + row_in_block;</div>
<div class="line"><a name="l00191"></a><span class="lineno">  191</span>&#160;        coords.set(buffer_index, col_it.index2());</div>
<div class="line"><a name="l00192"></a><span class="lineno">  192</span>&#160;        elements[buffer_index] = *col_it;</div>
<div class="line"><a name="l00193"></a><span class="lineno">  193</span>&#160;        entry_in_row++;</div>
<div class="line"><a name="l00194"></a><span class="lineno">  194</span>&#160;      }</div>
<div class="line"><a name="l00195"></a><span class="lineno">  195</span>&#160;</div>
<div class="line"><a name="l00196"></a><span class="lineno">  196</span>&#160;      ++row_in_block;</div>
<div class="line"><a name="l00197"></a><span class="lineno">  197</span>&#160;</div>
<div class="line"><a name="l00198"></a><span class="lineno">  198</span>&#160;      <span class="comment">// check for end of block:</span></div>
<div class="line"><a name="l00199"></a><span class="lineno">  199</span>&#160;      <span class="keywordflow">if</span> ( (row_it.index1() % gpu_matrix.rows_per_block() == gpu_matrix.rows_per_block() - 1)</div>
<div class="line"><a name="l00200"></a><span class="lineno">  200</span>&#160;           || row_it.index1() == <a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">viennacl::traits::size1</a>(cpu_matrix) - 1)</div>
<div class="line"><a name="l00201"></a><span class="lineno">  201</span>&#160;      {</div>
<div class="line"><a name="l00202"></a><span class="lineno">  202</span>&#160;        block_start.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#ac31bfde381eb40e0349678252c444af9">set</a>(block_index, static_cast&lt;IndexT&gt;(block_offset));</div>
<div class="line"><a name="l00203"></a><span class="lineno">  203</span>&#160;        block_offset += columns_in_block_buffer[block_index] * gpu_matrix.rows_per_block();</div>
<div class="line"><a name="l00204"></a><span class="lineno">  204</span>&#160;        ++block_index;</div>
<div class="line"><a name="l00205"></a><span class="lineno">  205</span>&#160;        row_in_block = 0;</div>
<div class="line"><a name="l00206"></a><span class="lineno">  206</span>&#160;      }</div>
<div class="line"><a name="l00207"></a><span class="lineno">  207</span>&#160;    }</div>
<div class="line"><a name="l00208"></a><span class="lineno">  208</span>&#160;</div>
<div class="line"><a name="l00209"></a><span class="lineno">  209</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(gpu_matrix.handle1(), columns_in_block_buffer.raw_size(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">traits::context</a>(gpu_matrix.handle1()), columns_in_block_buffer.get());</div>
<div class="line"><a name="l00210"></a><span class="lineno">  210</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(gpu_matrix.handle2(), coords.raw_size(),                  <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">traits::context</a>(gpu_matrix.handle2()), coords.get());</div>
<div class="line"><a name="l00211"></a><span class="lineno">  211</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(gpu_matrix.handle3(), block_start.raw_size(),             <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">traits::context</a>(gpu_matrix.handle3()), block_start.get());</div>
<div class="line"><a name="l00212"></a><span class="lineno">  212</span>&#160;    <a class="code" href="namespaceviennacl_1_1backend.html#a1499f19634964e2c7c8aeeefc6206126">viennacl::backend::memory_create</a>(gpu_matrix.handle(),  <span class="keyword">sizeof</span>(ScalarT) * elements.<a class="code" href="classviennacl_1_1backend_1_1typesafe__host__array.html#a3943ea09faa2df9b7bee130d3a2d5f20">size</a>(),  <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">traits::context</a>(gpu_matrix.handle()), &amp;(elements[0]));</div>
<div class="line"><a name="l00213"></a><span class="lineno">  213</span>&#160;  }</div>
<div class="line"><a name="l00214"></a><span class="lineno">  214</span>&#160;}</div>
<div class="line"><a name="l00215"></a><span class="lineno">  215</span>&#160;</div>
<div class="line"><a name="l00216"></a><span class="lineno">  216</span>&#160;</div>
<div class="line"><a name="l00217"></a><span class="lineno">  217</span>&#160;</div>
<div class="line"><a name="l00223"></a><span class="lineno">  223</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> IndexT, <span class="keyword">typename</span> NumericT, <span class="keyword">typename</span> IndexT2&gt;</div>
<div class="line"><a name="l00224"></a><span class="lineno"><a class="line" href="namespaceviennacl.html#ab1bf9bf2d8e0714db63fc999000c1510">  224</a></span>&#160;<span class="keywordtype">void</span> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(std::vector&lt; std::map&lt;IndexT, NumericT&gt; &gt; <span class="keyword">const</span> &amp; cpu_matrix,</div>
<div class="line"><a name="l00225"></a><span class="lineno">  225</span>&#160;          <a class="code" href="classviennacl_1_1sliced__ell__matrix.html">sliced_ell_matrix&lt;NumericT, IndexT2&gt;</a> &amp; gpu_matrix)</div>
<div class="line"><a name="l00226"></a><span class="lineno">  226</span>&#160;{</div>
<div class="line"><a name="l00227"></a><span class="lineno">  227</span>&#160;  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> max_col = 0;</div>
<div class="line"><a name="l00228"></a><span class="lineno">  228</span>&#160;  <span class="keywordflow">for</span> (<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> i=0; i&lt;cpu_matrix.size(); ++i)</div>
<div class="line"><a name="l00229"></a><span class="lineno">  229</span>&#160;  {</div>
<div class="line"><a name="l00230"></a><span class="lineno">  230</span>&#160;    <span class="keywordflow">if</span> (cpu_matrix[i].<a class="code" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">size</a>() &gt; 0)</div>
<div class="line"><a name="l00231"></a><span class="lineno">  231</span>&#160;      max_col = std::max&lt;vcl_size_t&gt;(max_col, (cpu_matrix[i].rbegin())-&gt;first);</div>
<div class="line"><a name="l00232"></a><span class="lineno">  232</span>&#160;  }</div>
<div class="line"><a name="l00233"></a><span class="lineno">  233</span>&#160;</div>
<div class="line"><a name="l00234"></a><span class="lineno">  234</span>&#160;  <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&lt;NumericT, IndexT&gt;</a>(cpu_matrix, cpu_matrix.size(), max_col + 1), gpu_matrix);</div>
<div class="line"><a name="l00235"></a><span class="lineno">  235</span>&#160;}</div>
<div class="line"><a name="l00236"></a><span class="lineno">  236</span>&#160;</div>
<div class="line"><a name="l00237"></a><span class="lineno">  237</span>&#160;</div>
<div class="line"><a name="l00238"></a><span class="lineno">  238</span>&#160;<span class="comment">/*</span></div>
<div class="line"><a name="l00239"></a><span class="lineno">  239</span>&#160;<span class="comment">template&lt;typename CPUMatrixT, typename ScalarT, typename IndexT&gt;</span></div>
<div class="line"><a name="l00240"></a><span class="lineno">  240</span>&#160;<span class="comment">void copy(sliced_ell_matrix&lt;ScalarT, IndexT&gt; const &amp; gpu_matrix, CPUMatrixT &amp; cpu_matrix )</span></div>
<div class="line"><a name="l00241"></a><span class="lineno">  241</span>&#160;<span class="comment">{</span></div>
<div class="line"><a name="l00242"></a><span class="lineno">  242</span>&#160;<span class="comment">  assert( (viennacl::traits::size1(cpu_matrix) == gpu_matrix.size1()) &amp;&amp; bool(&quot;Size mismatch&quot;) );</span></div>
<div class="line"><a name="l00243"></a><span class="lineno">  243</span>&#160;<span class="comment">  assert( (viennacl::traits::size2(cpu_matrix) == gpu_matrix.size2()) &amp;&amp; bool(&quot;Size mismatch&quot;) );</span></div>
<div class="line"><a name="l00244"></a><span class="lineno">  244</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00245"></a><span class="lineno">  245</span>&#160;<span class="comment">  if (gpu_matrix.size1() &gt; 0 &amp;&amp; gpu_matrix.size2() &gt; 0)</span></div>
<div class="line"><a name="l00246"></a><span class="lineno">  246</span>&#160;<span class="comment">  {</span></div>
<div class="line"><a name="l00247"></a><span class="lineno">  247</span>&#160;<span class="comment">    std::vector&lt;NumericT&gt; elements(gpu_matrix.internal_nnz());</span></div>
<div class="line"><a name="l00248"></a><span class="lineno">  248</span>&#160;<span class="comment">    viennacl::backend::typesafe_host_array&lt;unsigned int&gt; coords(gpu_matrix.handle2(), gpu_matrix.internal_nnz());</span></div>
<div class="line"><a name="l00249"></a><span class="lineno">  249</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00250"></a><span class="lineno">  250</span>&#160;<span class="comment">    viennacl::backend::memory_read(gpu_matrix.handle(), 0, sizeof(NumericT) * elements.size(), &amp;(elements[0]));</span></div>
<div class="line"><a name="l00251"></a><span class="lineno">  251</span>&#160;<span class="comment">    viennacl::backend::memory_read(gpu_matrix.handle2(), 0, coords.raw_size(), coords.get());</span></div>
<div class="line"><a name="l00252"></a><span class="lineno">  252</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00253"></a><span class="lineno">  253</span>&#160;<span class="comment">    for (vcl_size_t row = 0; row &lt; gpu_matrix.size1(); row++)</span></div>
<div class="line"><a name="l00254"></a><span class="lineno">  254</span>&#160;<span class="comment">    {</span></div>
<div class="line"><a name="l00255"></a><span class="lineno">  255</span>&#160;<span class="comment">      for (vcl_size_t ind = 0; ind &lt; gpu_matrix.internal_maxnnz(); ind++)</span></div>
<div class="line"><a name="l00256"></a><span class="lineno">  256</span>&#160;<span class="comment">      {</span></div>
<div class="line"><a name="l00257"></a><span class="lineno">  257</span>&#160;<span class="comment">        vcl_size_t offset = gpu_matrix.internal_size1() * ind + row;</span></div>
<div class="line"><a name="l00258"></a><span class="lineno">  258</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00259"></a><span class="lineno">  259</span>&#160;<span class="comment">        if (elements[offset] == static_cast&lt;NumericT&gt;(0.0))</span></div>
<div class="line"><a name="l00260"></a><span class="lineno">  260</span>&#160;<span class="comment">            continue;</span></div>
<div class="line"><a name="l00261"></a><span class="lineno">  261</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00262"></a><span class="lineno">  262</span>&#160;<span class="comment">        if (coords[offset] &gt;= gpu_matrix.size2())</span></div>
<div class="line"><a name="l00263"></a><span class="lineno">  263</span>&#160;<span class="comment">        {</span></div>
<div class="line"><a name="l00264"></a><span class="lineno">  264</span>&#160;<span class="comment">            std::cerr &lt;&lt; &quot;ViennaCL encountered invalid data &quot; &lt;&lt; offset &lt;&lt; &quot; &quot; &lt;&lt; ind &lt;&lt; &quot; &quot; &lt;&lt; row &lt;&lt; &quot; &quot; &lt;&lt; coords[offset] &lt;&lt; &quot; &quot; &lt;&lt; gpu_matrix.size2() &lt;&lt; std::endl;</span></div>
<div class="line"><a name="l00265"></a><span class="lineno">  265</span>&#160;<span class="comment">            return;</span></div>
<div class="line"><a name="l00266"></a><span class="lineno">  266</span>&#160;<span class="comment">        }</span></div>
<div class="line"><a name="l00267"></a><span class="lineno">  267</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00268"></a><span class="lineno">  268</span>&#160;<span class="comment">        cpu_matrix(row, coords[offset]) = elements[offset];</span></div>
<div class="line"><a name="l00269"></a><span class="lineno">  269</span>&#160;<span class="comment">      }</span></div>
<div class="line"><a name="l00270"></a><span class="lineno">  270</span>&#160;<span class="comment">    }</span></div>
<div class="line"><a name="l00271"></a><span class="lineno">  271</span>&#160;<span class="comment">  }</span></div>
<div class="line"><a name="l00272"></a><span class="lineno">  272</span>&#160;<span class="comment">} */</span></div>
<div class="line"><a name="l00273"></a><span class="lineno">  273</span>&#160;</div>
<div class="line"><a name="l00274"></a><span class="lineno">  274</span>&#160;</div>
<div class="line"><a name="l00275"></a><span class="lineno">  275</span>&#160;<span class="comment">//</span></div>
<div class="line"><a name="l00276"></a><span class="lineno">  276</span>&#160;<span class="comment">// Specify available operations:</span></div>
<div class="line"><a name="l00277"></a><span class="lineno">  277</span>&#160;<span class="comment">//</span></div>
<div class="line"><a name="l00278"></a><span class="lineno">  278</span>&#160;</div>
<div class="line"><a name="l00281"></a><span class="lineno">  281</span>&#160;<span class="keyword">namespace </span>linalg</div>
<div class="line"><a name="l00282"></a><span class="lineno">  282</span>&#160;{</div>
<div class="line"><a name="l00283"></a><span class="lineno">  283</span>&#160;<span class="keyword">namespace </span>detail</div>
<div class="line"><a name="l00284"></a><span class="lineno">  284</span>&#160;{</div>
<div class="line"><a name="l00285"></a><span class="lineno">  285</span>&#160;  <span class="comment">// x = A * y</span></div>
<div class="line"><a name="l00286"></a><span class="lineno">  286</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT&gt;</div>
<div class="line"><a name="l00287"></a><span class="lineno">  287</span>&#160;  <span class="keyword">struct </span>op_executor&lt;vector_base&lt;ScalarT&gt;, <a class="code" href="structop__assign.html">op_assign</a>, vector_expression&lt;const sliced_ell_matrix&lt;ScalarT, IndexT&gt;, const vector_base&lt;ScalarT&gt;, op_prod&gt; &gt;</div>
<div class="line"><a name="l00288"></a><span class="lineno">  288</span>&#160;  {</div>
<div class="line"><a name="l00289"></a><span class="lineno">  289</span>&#160;    <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base&lt;ScalarT&gt; &amp; lhs, vector_expression&lt;<span class="keyword">const</span> sliced_ell_matrix&lt;ScalarT, IndexT&gt;, <span class="keyword">const</span> vector_base&lt;ScalarT&gt;, op_prod&gt; <span class="keyword">const</span> &amp; rhs)</div>
<div class="line"><a name="l00290"></a><span class="lineno">  290</span>&#160;    {</div>
<div class="line"><a name="l00291"></a><span class="lineno">  291</span>&#160;      <span class="comment">// check for the special case x = A * x</span></div>
<div class="line"><a name="l00292"></a><span class="lineno">  292</span>&#160;      <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="l00293"></a><span class="lineno">  293</span>&#160;      {</div>
<div class="line"><a name="l00294"></a><span class="lineno">  294</span>&#160;        <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp(lhs);</div>
<div class="line"><a name="l00295"></a><span class="lineno">  295</span>&#160;        <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), ScalarT(1), temp, ScalarT(0));</div>
<div class="line"><a name="l00296"></a><span class="lineno">  296</span>&#160;        lhs = temp;</div>
<div class="line"><a name="l00297"></a><span class="lineno">  297</span>&#160;      }</div>
<div class="line"><a name="l00298"></a><span class="lineno">  298</span>&#160;      <span class="keywordflow">else</span></div>
<div class="line"><a name="l00299"></a><span class="lineno">  299</span>&#160;        <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), ScalarT(1), lhs, ScalarT(0));</div>
<div class="line"><a name="l00300"></a><span class="lineno">  300</span>&#160;    }</div>
<div class="line"><a name="l00301"></a><span class="lineno">  301</span>&#160;  };</div>
<div class="line"><a name="l00302"></a><span class="lineno">  302</span>&#160;</div>
<div class="line"><a name="l00303"></a><span class="lineno">  303</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT&gt;</div>
<div class="line"><a name="l00304"></a><span class="lineno">  304</span>&#160;  <span class="keyword">struct </span>op_executor&lt;vector_base&lt;ScalarT&gt;, op_inplace_add, vector_expression&lt;const sliced_ell_matrix&lt;ScalarT, IndexT&gt;, const vector_base&lt;ScalarT&gt;, op_prod&gt; &gt;</div>
<div class="line"><a name="l00305"></a><span class="lineno">  305</span>&#160;  {</div>
<div class="line"><a name="l00306"></a><span class="lineno">  306</span>&#160;    <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base&lt;ScalarT&gt; &amp; lhs, vector_expression&lt;<span class="keyword">const</span> sliced_ell_matrix&lt;ScalarT, IndexT&gt;, <span class="keyword">const</span> vector_base&lt;ScalarT&gt;, op_prod&gt; <span class="keyword">const</span> &amp; rhs)</div>
<div class="line"><a name="l00307"></a><span class="lineno">  307</span>&#160;    {</div>
<div class="line"><a name="l00308"></a><span class="lineno">  308</span>&#160;      <span class="comment">// check for the special case x += A * x</span></div>
<div class="line"><a name="l00309"></a><span class="lineno">  309</span>&#160;      <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="l00310"></a><span class="lineno">  310</span>&#160;      {</div>
<div class="line"><a name="l00311"></a><span class="lineno">  311</span>&#160;        <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp(lhs);</div>
<div class="line"><a name="l00312"></a><span class="lineno">  312</span>&#160;        <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), ScalarT(1), temp, ScalarT(0));</div>
<div class="line"><a name="l00313"></a><span class="lineno">  313</span>&#160;        lhs += temp;</div>
<div class="line"><a name="l00314"></a><span class="lineno">  314</span>&#160;      }</div>
<div class="line"><a name="l00315"></a><span class="lineno">  315</span>&#160;      <span class="keywordflow">else</span></div>
<div class="line"><a name="l00316"></a><span class="lineno">  316</span>&#160;        <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), ScalarT(1), lhs, ScalarT(1));</div>
<div class="line"><a name="l00317"></a><span class="lineno">  317</span>&#160;    }</div>
<div class="line"><a name="l00318"></a><span class="lineno">  318</span>&#160;  };</div>
<div class="line"><a name="l00319"></a><span class="lineno">  319</span>&#160;</div>
<div class="line"><a name="l00320"></a><span class="lineno">  320</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT&gt;</div>
<div class="line"><a name="l00321"></a><span class="lineno">  321</span>&#160;  <span class="keyword">struct </span>op_executor&lt;vector_base&lt;ScalarT&gt;, op_inplace_sub, vector_expression&lt;const sliced_ell_matrix&lt;ScalarT, IndexT&gt;, const vector_base&lt;ScalarT&gt;, op_prod&gt; &gt;</div>
<div class="line"><a name="l00322"></a><span class="lineno">  322</span>&#160;  {</div>
<div class="line"><a name="l00323"></a><span class="lineno">  323</span>&#160;    <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base&lt;ScalarT&gt; &amp; lhs, vector_expression&lt;<span class="keyword">const</span> sliced_ell_matrix&lt;ScalarT, IndexT&gt;, <span class="keyword">const</span> vector_base&lt;ScalarT&gt;, op_prod&gt; <span class="keyword">const</span> &amp; rhs)</div>
<div class="line"><a name="l00324"></a><span class="lineno">  324</span>&#160;    {</div>
<div class="line"><a name="l00325"></a><span class="lineno">  325</span>&#160;      <span class="comment">// check for the special case x -= A * x</span></div>
<div class="line"><a name="l00326"></a><span class="lineno">  326</span>&#160;      <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="l00327"></a><span class="lineno">  327</span>&#160;      {</div>
<div class="line"><a name="l00328"></a><span class="lineno">  328</span>&#160;        <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp(lhs);</div>
<div class="line"><a name="l00329"></a><span class="lineno">  329</span>&#160;        <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), ScalarT(1), temp, ScalarT(0));</div>
<div class="line"><a name="l00330"></a><span class="lineno">  330</span>&#160;        lhs -= temp;</div>
<div class="line"><a name="l00331"></a><span class="lineno">  331</span>&#160;      }</div>
<div class="line"><a name="l00332"></a><span class="lineno">  332</span>&#160;      <span class="keywordflow">else</span></div>
<div class="line"><a name="l00333"></a><span class="lineno">  333</span>&#160;        <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), rhs.rhs(), ScalarT(-1), lhs, ScalarT(1));</div>
<div class="line"><a name="l00334"></a><span class="lineno">  334</span>&#160;    }</div>
<div class="line"><a name="l00335"></a><span class="lineno">  335</span>&#160;  };</div>
<div class="line"><a name="l00336"></a><span class="lineno">  336</span>&#160;</div>
<div class="line"><a name="l00337"></a><span class="lineno">  337</span>&#160;</div>
<div class="line"><a name="l00338"></a><span class="lineno">  338</span>&#160;  <span class="comment">// x = A * vec_op</span></div>
<div class="line"><a name="l00339"></a><span class="lineno">  339</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT, <span class="keyword">typename</span> LHS, <span class="keyword">typename</span> RHS, <span class="keyword">typename</span> OP&gt;</div>
<div class="line"><a name="l00340"></a><span class="lineno">  340</span>&#160;  <span class="keyword">struct </span>op_executor&lt;vector_base&lt;ScalarT&gt;, <a class="code" href="structop__assign.html">op_assign</a>, vector_expression&lt;const sliced_ell_matrix&lt;ScalarT, IndexT&gt;, const vector_expression&lt;const LHS, const RHS, OP&gt;, op_prod&gt; &gt;</div>
<div class="line"><a name="l00341"></a><span class="lineno">  341</span>&#160;  {</div>
<div class="line"><a name="l00342"></a><span class="lineno">  342</span>&#160;    <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base&lt;ScalarT&gt; &amp; lhs, vector_expression&lt;<span class="keyword">const</span> sliced_ell_matrix&lt;ScalarT, IndexT&gt;, <span class="keyword">const</span> vector_expression&lt;const LHS, const RHS, OP&gt;, op_prod&gt; <span class="keyword">const</span> &amp; rhs)</div>
<div class="line"><a name="l00343"></a><span class="lineno">  343</span>&#160;    {</div>
<div class="line"><a name="l00344"></a><span class="lineno">  344</span>&#160;      <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp(rhs.rhs(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(rhs));</div>
<div class="line"><a name="l00345"></a><span class="lineno">  345</span>&#160;      <a class="code" href="namespaceviennacl_1_1linalg.html#aafc6db8d806f67c24f93eaaded84b853">viennacl::linalg::prod_impl</a>(rhs.lhs(), temp, lhs);</div>
<div class="line"><a name="l00346"></a><span class="lineno">  346</span>&#160;    }</div>
<div class="line"><a name="l00347"></a><span class="lineno">  347</span>&#160;  };</div>
<div class="line"><a name="l00348"></a><span class="lineno">  348</span>&#160;</div>
<div class="line"><a name="l00349"></a><span class="lineno">  349</span>&#160;  <span class="comment">// x = A * vec_op</span></div>
<div class="line"><a name="l00350"></a><span class="lineno">  350</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT, <span class="keyword">typename</span> LHS, <span class="keyword">typename</span> RHS, <span class="keyword">typename</span> OP&gt;</div>
<div class="line"><a name="l00351"></a><span class="lineno">  351</span>&#160;  <span class="keyword">struct </span>op_executor&lt;vector_base&lt;ScalarT&gt;, op_inplace_add, vector_expression&lt;const sliced_ell_matrix&lt;ScalarT, IndexT&gt;, const vector_expression&lt;const LHS, const RHS, OP&gt;, op_prod&gt; &gt;</div>
<div class="line"><a name="l00352"></a><span class="lineno">  352</span>&#160;  {</div>
<div class="line"><a name="l00353"></a><span class="lineno">  353</span>&#160;    <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base&lt;ScalarT&gt; &amp; lhs, vector_expression&lt;<span class="keyword">const</span> sliced_ell_matrix&lt;ScalarT, IndexT&gt;, <span class="keyword">const</span> vector_expression&lt;const LHS, const RHS, OP&gt;, op_prod&gt; <span class="keyword">const</span> &amp; rhs)</div>
<div class="line"><a name="l00354"></a><span class="lineno">  354</span>&#160;    {</div>
<div class="line"><a name="l00355"></a><span class="lineno">  355</span>&#160;      <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp(rhs.rhs(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(rhs));</div>
<div class="line"><a name="l00356"></a><span class="lineno">  356</span>&#160;      <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp_result(lhs);</div>
<div class="line"><a name="l00357"></a><span class="lineno">  357</span>&#160;      <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="l00358"></a><span class="lineno">  358</span>&#160;      lhs += temp_result;</div>
<div class="line"><a name="l00359"></a><span class="lineno">  359</span>&#160;    }</div>
<div class="line"><a name="l00360"></a><span class="lineno">  360</span>&#160;  };</div>
<div class="line"><a name="l00361"></a><span class="lineno">  361</span>&#160;</div>
<div class="line"><a name="l00362"></a><span class="lineno">  362</span>&#160;  <span class="comment">// x = A * vec_op</span></div>
<div class="line"><a name="l00363"></a><span class="lineno">  363</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> ScalarT, <span class="keyword">typename</span> IndexT, <span class="keyword">typename</span> LHS, <span class="keyword">typename</span> RHS, <span class="keyword">typename</span> OP&gt;</div>
<div class="line"><a name="l00364"></a><span class="lineno">  364</span>&#160;  <span class="keyword">struct </span>op_executor&lt;vector_base&lt;ScalarT&gt;, op_inplace_sub, vector_expression&lt;const sliced_ell_matrix&lt;ScalarT, IndexT&gt;, const vector_expression&lt;const LHS, const RHS, OP&gt;, op_prod&gt; &gt;</div>
<div class="line"><a name="l00365"></a><span class="lineno">  365</span>&#160;  {</div>
<div class="line"><a name="l00366"></a><span class="lineno">  366</span>&#160;    <span class="keyword">static</span> <span class="keywordtype">void</span> apply(vector_base&lt;ScalarT&gt; &amp; lhs, vector_expression&lt;<span class="keyword">const</span> sliced_ell_matrix&lt;ScalarT, IndexT&gt;, <span class="keyword">const</span> vector_expression&lt;const LHS, const RHS, OP&gt;, op_prod&gt; <span class="keyword">const</span> &amp; rhs)</div>
<div class="line"><a name="l00367"></a><span class="lineno">  367</span>&#160;    {</div>
<div class="line"><a name="l00368"></a><span class="lineno">  368</span>&#160;      <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp(rhs.rhs(), <a class="code" href="namespaceviennacl_1_1traits.html#a6707f5dab8f482170d2046a605f46ef8">viennacl::traits::context</a>(rhs));</div>
<div class="line"><a name="l00369"></a><span class="lineno">  369</span>&#160;      <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarT&gt;</a> temp_result(lhs);</div>
<div class="line"><a name="l00370"></a><span class="lineno">  370</span>&#160;      <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="l00371"></a><span class="lineno">  371</span>&#160;      lhs -= temp_result;</div>
<div class="line"><a name="l00372"></a><span class="lineno">  372</span>&#160;    }</div>
<div class="line"><a name="l00373"></a><span class="lineno">  373</span>&#160;  };</div>
<div class="line"><a name="l00374"></a><span class="lineno">  374</span>&#160;</div>
<div class="line"><a name="l00375"></a><span class="lineno">  375</span>&#160;} <span class="comment">// namespace detail</span></div>
<div class="line"><a name="l00376"></a><span class="lineno">  376</span>&#160;} <span class="comment">// namespace linalg</span></div>
<div class="line"><a name="l00377"></a><span class="lineno">  377</span>&#160;</div>
<div class="line"><a name="l00379"></a><span class="lineno">  379</span>&#160;}</div>
<div class="line"><a name="l00380"></a><span class="lineno">  380</span>&#160;</div>
<div class="line"><a name="l00381"></a><span class="lineno">  381</span>&#160;<span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00382"></a><span class="lineno">  382</span>&#160;</div>
<div class="line"><a name="l00383"></a><span class="lineno">  383</span>&#160;</div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_a31d5fa9ec25dcfeb260f7e11c4db6571"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a31d5fa9ec25dcfeb260f7e11c4db6571">viennacl::sliced_ell_matrix::handle2</a></div><div class="ttdeci">handle_type &amp; handle2()</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00112">sliced_ell_matrix.hpp:112</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_ac681e3c8d41b686282711644cbdcadfb"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#ac681e3c8d41b686282711644cbdcadfb">viennacl::sliced_ell_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="sliced__ell__matrix_8hpp_source.html#l00085">sliced_ell_matrix.hpp:85</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_aec497ed037ef674d36729bda4db7e93f"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#aec497ed037ef674d36729bda4db7e93f">viennacl::sliced_ell_matrix::handle3</a></div><div class="ttdeci">const handle_type &amp; handle3() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00116">sliced_ell_matrix.hpp:116</a></div></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_1sliced__ell__matrix_html_a362dc0c5f0d644cbcf6514f3b3ee6529"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a362dc0c5f0d644cbcf6514f3b3ee6529">viennacl::sliced_ell_matrix::sliced_ell_matrix</a></div><div class="ttdeci">sliced_ell_matrix()</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00053">sliced_ell_matrix.hpp:53</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_1sliced__ell__matrix_html_adf0fb80bb38ab0256a622960c0c18699"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#adf0fb80bb38ab0256a622960c0c18699">viennacl::sliced_ell_matrix::size_type</a></div><div class="ttdeci">vcl_size_t size_type</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00051">sliced_ell_matrix.hpp:51</a></div></div>
<div class="ttc" id="classviennacl_1_1backend_1_1typesafe__host__array_html_a3943ea09faa2df9b7bee130d3a2d5f20"><div class="ttname"><a href="classviennacl_1_1backend_1_1typesafe__host__array.html#a3943ea09faa2df9b7bee130d3a2d5f20">viennacl::backend::typesafe_host_array::size</a></div><div class="ttdeci">vcl_size_t size() const </div><div class="ttdef"><b>Definition:</b> <a href="backend_2util_8hpp_source.html#l00113">util.hpp:113</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_a727fd3a25232545c59204885ff3e32d4"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a727fd3a25232545c59204885ff3e32d4">viennacl::sliced_ell_matrix::handle1</a></div><div class="ttdeci">const handle_type &amp; handle1() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00110">sliced_ell_matrix.hpp:110</a></div></div>
<div class="ttc" id="tools_8hpp_html"><div class="ttname"><a href="tools_8hpp.html">tools.hpp</a></div><div class="ttdoc">Various little tools used here and there in ViennaCL. </div></div>
<div class="ttc" id="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 &amp;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="namespaceviennacl_html_a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1"><div class="ttname"><a href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a></div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00349">forwards.h:349</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_a45249ad5b44d0f506adc27e1e18aecdf"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a45249ad5b44d0f506adc27e1e18aecdf">viennacl::sliced_ell_matrix::handle2</a></div><div class="ttdeci">const handle_type &amp; handle2() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00113">sliced_ell_matrix.hpp:113</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_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_1linalg_1_1detail_html_a5d46fe9558b0e462f10fd44942ad4fc6"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a5d46fe9558b0e462f10fd44942ad4fc6">viennacl::linalg::detail::max</a></div><div class="ttdeci">T max(const T &amp;lhs, const T &amp;rhs)</div><div class="ttdoc">Maximum. </div><div class="ttdef"><b>Definition:</b> <a href="linalg_2detail_2bisect_2util_8hpp_source.html#l00059">util.hpp:59</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_aa9a7f4003687cae9a2af97f9fa328db5"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#aa9a7f4003687cae9a2af97f9fa328db5">viennacl::sliced_ell_matrix::rows_per_block</a></div><div class="ttdeci">vcl_size_t rows_per_block() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00104">sliced_ell_matrix.hpp:104</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_a40c479e880ef7e7f50ed67cda91bb8e7"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a40c479e880ef7e7f50ed67cda91bb8e7">viennacl::sliced_ell_matrix::copy</a></div><div class="ttdeci">friend void copy(CPUMatrixT const &amp;cpu_matrix, sliced_ell_matrix&lt; ScalarT2, IndexT2 &gt; &amp;gpu_matrix)</div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_af0672b686ad27ddcf948116bae059329"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#af0672b686ad27ddcf948116bae059329">viennacl::sliced_ell_matrix::value_type</a></div><div class="ttdeci">scalar&lt; typename viennacl::tools::CHECK_SCALAR_TEMPLATE_ARGUMENT&lt; ScalarT &gt;::ResultType &gt; value_type</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00050">sliced_ell_matrix.hpp:50</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&lt; MatrixType &gt;::type size2(MatrixType const &amp;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="classviennacl_1_1sliced__ell__matrix_html_abbc671358c44037b7c2b22a01f17c3a7"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#abbc671358c44037b7c2b22a01f17c3a7">viennacl::sliced_ell_matrix::handle</a></div><div class="ttdeci">handle_type &amp; handle()</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00118">sliced_ell_matrix.hpp:118</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 &#39;context&#39; 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_1sliced__ell__matrix_html_aa725bf2d5b213f119bd4b9fe9258c201"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#aa725bf2d5b213f119bd4b9fe9258c201">viennacl::sliced_ell_matrix::sliced_ell_matrix</a></div><div class="ttdeci">sliced_ell_matrix(viennacl::context ctx)</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00066">sliced_ell_matrix.hpp:66</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 &amp;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_1sliced__ell__matrix_html"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html">viennacl::sliced_ell_matrix</a></div><div class="ttdoc">Sparse matrix class using the sliced ELLPACK with parameters C, . </div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00403">forwards.h:403</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_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="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_1sliced__ell__matrix_html_a16cd8194ca5b263de32f131362498d77"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a16cd8194ca5b263de32f131362498d77">viennacl::sliced_ell_matrix::size1</a></div><div class="ttdeci">vcl_size_t size1() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00101">sliced_ell_matrix.hpp:101</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_a9f244082eff24e93495d3b4abec735b0"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a9f244082eff24e93495d3b4abec735b0">viennacl::sliced_ell_matrix::sliced_ell_matrix</a></div><div class="ttdeci">sliced_ell_matrix(size_type num_rows, size_type num_cols, size_type num_rows_per_block_=0)</div><div class="ttdoc">Standard constructor for setting the row and column sizes as well as the block size. </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00059">sliced_ell_matrix.hpp:59</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_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_1sliced__ell__matrix_html_a1b9a6e9229a8694123562166eb2d4e32"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a1b9a6e9229a8694123562166eb2d4e32">viennacl::sliced_ell_matrix::handle</a></div><div class="ttdeci">const handle_type &amp; handle() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00119">sliced_ell_matrix.hpp:119</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_ad4fd3359611ca475b85d2041bdec1377"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#ad4fd3359611ca475b85d2041bdec1377">viennacl::sliced_ell_matrix::internal_size2</a></div><div class="ttdeci">vcl_size_t internal_size2() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00099">sliced_ell_matrix.hpp:99</a></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_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 &amp;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="vector_8hpp_html"><div class="ttname"><a href="vector_8hpp.html">vector.hpp</a></div><div class="ttdoc">The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...</div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_aa59ae796d592a679a4c6036ba84f61b0"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#aa59ae796d592a679a4c6036ba84f61b0">viennacl::sliced_ell_matrix::handle3</a></div><div class="ttdeci">handle_type &amp; handle3()</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00115">sliced_ell_matrix.hpp:115</a></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&lt; NumericT &gt; &amp;cpu_vec, circulant_matrix&lt; NumericT, AlignmentV &gt; &amp;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_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="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 &amp;handle, vcl_size_t size_in_bytes, viennacl::context const &amp;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&lt; NumericT &gt; &amp;mat, const vector_base&lt; NumericT &gt; &amp;vec, vector_base&lt; NumericT &gt; &amp;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 &amp; handle(T &amp;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_1sliced__ell__matrix_html_ae0afefc19f6dc7692a87cd8a16b51d30"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#ae0afefc19f6dc7692a87cd8a16b51d30">viennacl::sliced_ell_matrix::size2</a></div><div class="ttdeci">vcl_size_t size2() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00102">sliced_ell_matrix.hpp:102</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_ab8a2233384ac8dcdb85b5385706548a1"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#ab8a2233384ac8dcdb85b5385706548a1">viennacl::sliced_ell_matrix::internal_size1</a></div><div class="ttdeci">vcl_size_t internal_size1() const </div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00098">sliced_ell_matrix.hpp:98</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_a94043424c01974b2db8062f74f625825"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a94043424c01974b2db8062f74f625825">viennacl::sliced_ell_matrix::handle1</a></div><div class="ttdeci">handle_type &amp; handle1()</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00109">sliced_ell_matrix.hpp:109</a></div></div>
<div class="ttc" id="classviennacl_1_1sliced__ell__matrix_html_a2cdf3aeccfe9a10b89c48b3430cb4c93"><div class="ttname"><a href="classviennacl_1_1sliced__ell__matrix.html#a2cdf3aeccfe9a10b89c48b3430cb4c93">viennacl::sliced_ell_matrix::handle_type</a></div><div class="ttdeci">viennacl::backend::mem_handle handle_type</div><div class="ttdef"><b>Definition:</b> <a href="sliced__ell__matrix_8hpp_source.html#l00049">sliced_ell_matrix.hpp:49</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><!-- 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="sliced__ell__matrix_8hpp.html">sliced_ell_matrix.hpp</a></li>
    <li class="footer">Generated on Wed Jan 20 2016 22:32:42 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>