File: manual-operations.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 (331 lines) | stat: -rw-r--r-- 30,438 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
<!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: Basic Operations</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('manual-operations.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">Basic Operations </div>  </div>
</div><!--header-->
<div class="contents">
<div class="textblock"><p>The basic types have been introduced in the previous chapter, so we move on with the description of the basic BLAS operations. Almost all operations supported by Boost.uBLAS are available, including element-wise operations on vectors. Thus, consider the <a href="http://www.boost.org/doc/libs/1_52_0/libs/numeric/ublas/doc/operations_overview.htm">Boost.uBLAS documentation</a> as a reference as well.</p>
<h1><a class="anchor" id="manual-operations-blas1"></a>
Elementary Vector and Matrix Operations (BLAS Level 1)</h1>
<p>ViennaCL provides all vector-vector operations defined at level 1 of BLAS as well as the most important element-wise operations:</p>
<center> <table class="doxtable">
<tr>
<th>Verbal </th><th>Mathematics </th><th>ViennaCL Code  </th></tr>
<tr>
<td>swap </td><td><img class="formulaInl" alt="$ x \leftrightarrow y $" src="form_40.png"/> </td><td><code>swap(x,y);</code>  </td></tr>
<tr>
<td>stretch </td><td><img class="formulaInl" alt="$ x \leftarrow \alpha x $" src="form_41.png"/> </td><td><code>x *= alpha;</code>  </td></tr>
<tr>
<td>assignment </td><td><img class="formulaInl" alt="$ y \leftarrow x $" src="form_42.png"/> </td><td><code>y = x;</code>  </td></tr>
<tr>
<td>multiply add </td><td><img class="formulaInl" alt="$ y \leftarrow \alpha x + y $" src="form_43.png"/> </td><td><code>y += alpha * x;</code>  </td></tr>
<tr>
<td>multiply subtract </td><td><img class="formulaInl" alt="$ y \leftarrow \alpha x - y $" src="form_44.png"/> </td><td><code>y -= alpha * x;</code>  </td></tr>
<tr>
<td>inner dot product </td><td><img class="formulaInl" alt="$ \alpha \leftarrow x^{\mathrm{T}} y $" src="form_45.png"/> </td><td><code>inner_prod(x,y);</code>  </td></tr>
<tr>
<td><img class="formulaInl" alt="$L^1$" src="form_46.png"/> norm </td><td><img class="formulaInl" alt="$ \alpha \leftarrow \Vert x \Vert_1 $" src="form_47.png"/> </td><td><code>alpha = norm_1(x);</code>  </td></tr>
<tr>
<td><img class="formulaInl" alt="$L^2$" src="form_48.png"/> norm </td><td><img class="formulaInl" alt="$ \alpha \leftarrow \Vert x \Vert_2 $" src="form_49.png"/> </td><td><code>alpha = norm_2(x);</code>  </td></tr>
<tr>
<td><img class="formulaInl" alt="$L^\infty$" src="form_50.png"/> norm </td><td><img class="formulaInl" alt="$ \alpha \leftarrow \Vert x \Vert_\infty $" src="form_51.png"/> </td><td><code>alpha = norm_inf(x);</code>  </td></tr>
<tr>
<td><img class="formulaInl" alt="$L^\infty$" src="form_50.png"/> norm index </td><td><img class="formulaInl" alt="$ i \leftarrow \max_i \vert x_i \vert $" src="form_52.png"/> </td><td><code>i = index_norm_inf(x);</code>  </td></tr>
<tr>
<td>plane rotation </td><td><img class="formulaInl" alt="$(x,y) \leftarrow (\alpha x + \beta y, -\beta x + \alpha y) $" src="form_53.png"/> </td><td><code>plane_rotation(a, b, x, y);</code>  </td></tr>
<tr>
<td></td><td></td><td></td></tr>
<tr>
<td>elementwise product </td><td><img class="formulaInl" alt="$ y_i \leftarrow x_i \cdot z_i $" src="form_54.png"/> </td><td><code>y = element_prod(x,z);</code>  </td></tr>
<tr>
<td>elementwise division </td><td><img class="formulaInl" alt="$ y_i \leftarrow x_i \cdot z_i $" src="form_54.png"/> </td><td><code>y = element_div(x,z);</code>  </td></tr>
<tr>
<td>elementwise power </td><td><img class="formulaInl" alt="$ y_i \leftarrow x_i^{z_i} $" src="form_55.png"/> </td><td><code>y = element_pow(x,z);</code>  </td></tr>
<tr>
<td></td><td></td><td></td></tr>
<tr>
<td>elementwise modulus (ints) </td><td><img class="formulaInl" alt="$ y_i \leftarrow `x_i` $" src="form_56.png"/> </td><td><code>y = element_abs(x);</code>  </td></tr>
<tr>
<td>elementwise modulus (floats) </td><td><img class="formulaInl" alt="$ y_i \leftarrow `x_i` $" src="form_56.png"/> </td><td><code>y = element_fabs(x);</code>  </td></tr>
<tr>
<td>elementwise acos </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{acos}(x_i) $" src="form_57.png"/> </td><td><code>y = element_acos(x);</code>  </td></tr>
<tr>
<td>elementwise asin </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{asin}(x_i) $" src="form_58.png"/> </td><td><code>y = element_asin(x);</code>  </td></tr>
<tr>
<td>elementwise atan </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{atan}(x_i) $" src="form_59.png"/> </td><td><code>y = element_atan(x);</code>  </td></tr>
<tr>
<td>elementwise ceil </td><td><img class="formulaInl" alt="$ y_i \leftarrow \lceil x_i \rceil $" src="form_60.png"/> </td><td><code>y = element_ceil(x);</code>  </td></tr>
<tr>
<td>elementwise cos </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{cos}(x_i) $" src="form_61.png"/> </td><td><code>y = element_cos(x);</code>  </td></tr>
<tr>
<td>elementwise cosh </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{cosh}(x_i) $" src="form_62.png"/> </td><td><code>y = element_cosh(x);</code>  </td></tr>
<tr>
<td>elementwise exp </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{exp}(x_i) $" src="form_63.png"/> </td><td><code>y = element_exp(x);</code>  </td></tr>
<tr>
<td>elementwise floor </td><td><img class="formulaInl" alt="$ y_i \leftarrow \lfloor x_i \rfloor $" src="form_64.png"/> </td><td><code>y = element_floor(x);</code>  </td></tr>
<tr>
<td>elementwise log (base e) </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{ln}(x_i) $" src="form_65.png"/> </td><td><code>y = element_log(x);</code>  </td></tr>
<tr>
<td>elementwise log (base 10) </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{log}_{10}(x_i)$" src="form_66.png"/> </td><td><code>y = element_log10(x);</code>  </td></tr>
<tr>
<td>elementwise sin </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{sin}(x_i) $" src="form_67.png"/> </td><td><code>y = element_sin(x);</code>  </td></tr>
<tr>
<td>elementwise sinh </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{sinh}(x_i) $" src="form_68.png"/> </td><td><code>y = element_sinh(x);</code>  </td></tr>
<tr>
<td>elementwise sqrt </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{sqrt}(x_i) $" src="form_69.png"/> </td><td><code>y = element_sqrt(x);</code>  </td></tr>
<tr>
<td>elementwise tan </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{tan}(x_i) $" src="form_70.png"/> </td><td><code>y = element_tan(x);</code>  </td></tr>
<tr>
<td>elementwise tanh </td><td><img class="formulaInl" alt="$ y_i \leftarrow \textrm{tanh}(x_i) $" src="form_71.png"/> </td><td><code>y = element_tanh(x);</code>  </td></tr>
</table>
<p><b>BLAS level 1 routines mapped to ViennaCL. Note that the free functions reside in namespace <code><a class="el" href="namespaceviennacl_1_1linalg.html" title="Provides all linear algebra operations which are not covered by operator overloads. ">viennacl::linalg</a></code></b> </center><p>The function interface is compatible with Boost.uBLAS, thus allowing quick code migration for Boost.uBLAS users. Element-wise operations and standard operator overloads are available for dense matrices as well. The only dense matrix norm provided is <code><a class="el" href="namespaceviennacl_1_1linalg.html#a67813d48ed274df9fc19e6aad99fc761">norm_frobenius()</a></code> for the Frobenius norm.</p>
<dl class="section note"><dt>Note</dt><dd>Mixing operations between objects of different scalar types is not supported. Convert the data manually on the host if needed.</dd></dl>
<dl class="section warning"><dt>Warning</dt><dd>The operator overloads make extensive use of expression templates. Do not use the C++11 keyword <code>auto</code> for the result type, as this might result in unexpected performance regressions or dangling references.</dd></dl>
<h1><a class="anchor" id="manual-operations-blas2"></a>
Matrix-Vector Operations (BLAS Level 2)</h1>
<p>The interface for level 2 BLAS functions in ViennaCL is similar to that of Boost.uBLAS:</p>
<dl class="section note"><dt>Note</dt><dd>Mixing operations between objects of different scalar types is not supported. Convert the data manually on the host if needed.</dd></dl>
<center> <table class="doxtable">
<tr>
<th>Verbal </th><th>Mathematics </th><th>ViennaCL Code  </th></tr>
<tr>
<td>matrix vector product </td><td><img class="formulaInl" alt="$ y \leftarrow A x $" src="form_72.png"/> </td><td><code>y = prod(A, x);</code>  </td></tr>
<tr>
<td>matrix vector product </td><td><img class="formulaInl" alt="$ y \leftarrow A^\mathrm{T} x $" src="form_73.png"/> </td><td><code>y = prod(trans(A), x);</code>  </td></tr>
<tr>
<td>inplace mv product </td><td><img class="formulaInl" alt="$ x \leftarrow A x $" src="form_74.png"/> </td><td><code>x = prod(A, x);</code>  </td></tr>
<tr>
<td>inplace mv product </td><td><img class="formulaInl" alt="$ x \leftarrow A^\mathrm{T} x $" src="form_75.png"/> </td><td><code>x = prod(trans(A), x);</code>  </td></tr>
<tr>
<td></td><td></td><td></td></tr>
<tr>
<td>scaled product add </td><td><img class="formulaInl" alt="$ y \leftarrow \alpha A x + \beta y $" src="form_76.png"/> </td><td><code>y = alpha * prod(A, x) + beta * y</code>  </td></tr>
<tr>
<td>scaled product add </td><td><img class="formulaInl" alt="$ y \leftarrow \alpha A^{\mathrm T} x + \beta y $" src="form_77.png"/> </td><td><code>y = alpha * prod(trans(A), x) + beta * y</code>  </td></tr>
<tr>
<td></td><td></td><td></td></tr>
<tr>
<td>tri. matrix solve </td><td><img class="formulaInl" alt="$ y \leftarrow A^{-1} x $" src="form_78.png"/> </td><td><code>y = solve(A, x, tag);</code>  </td></tr>
<tr>
<td>tri. matrix solve </td><td><img class="formulaInl" alt="$ y \leftarrow A^\mathrm{T^{-1}} x $" src="form_79.png"/> </td><td><code>y = solve(trans(A), x, tag);</code>  </td></tr>
<tr>
<td>inplace solve </td><td><img class="formulaInl" alt="$ x \leftarrow A^{-1} x $" src="form_80.png"/> </td><td><code>inplace_solve(A, x, tag);</code>  </td></tr>
<tr>
<td>inplace solve </td><td><img class="formulaInl" alt="$ x \leftarrow A^\mathrm{T^{-1}} x $" src="form_81.png"/> </td><td><code>inplace_solve(trans(A), x, tag);</code>  </td></tr>
<tr>
<td></td><td></td><td></td></tr>
<tr>
<td>rank 1 update </td><td><img class="formulaInl" alt="$ A \leftarrow \alpha x y^{\mathrm T} + A $" src="form_82.png"/> </td><td><code>A += alpha * outer_prod(x,y);</code>  </td></tr>
<tr>
<td>symm. rank 1 update </td><td><img class="formulaInl" alt="$ A \leftarrow \alpha x x^{\mathrm T} + A $" src="form_83.png"/> </td><td><code>A += alpha * outer_prod(x,x);</code>  </td></tr>
<tr>
<td>rank 2 update </td><td><img class="formulaInl" alt="$ A \leftarrow \alpha (x y^{\mathrm T} + y x^{\mathrm T}) + A $" src="form_84.png"/> </td><td><code>A += alpha * outer_prod(x,y);</code> <code>A += alpha * outer_prod(y,x);</code>  </td></tr>
</table>
<p><b>BLAS level 2 routines mapped to ViennaCL. Note that the free functions reside in namespace <code><a class="el" href="namespaceviennacl_1_1linalg.html" title="Provides all linear algebra operations which are not covered by operator overloads. ">viennacl::linalg</a></code>.<br/>
 <code>tag</code> is one out of <code>lower_tag</code>, <code>unit_lower_tag</code>, <code>upper_tag</code>, and <code>unit_upper_tag</code>.</b> </center><dl class="section warning"><dt>Warning</dt><dd>The operator overloads make extensive use of expression templates. Do not use the C++11 keyword <code>auto</code> for the result type, as this might result in unexpected performance regressions or dangling references.</dd></dl>
<h1><a class="anchor" id="manual-operations-blas3"></a>
Matrix-Matrix Operations (BLAS Level 3)</h1>
<p>While BLAS levels 1 and 2 are mostly memory-bandwidth-limited, BLAS level 3 is mostly limited by the available computational power of the respective device. Hence, matrix-matrix products regularly show impressive performance gains on mid- to high-end GPUs when compared to a single CPU core (which is apparently not a fair comparison...)</p>
<p>Again, the ViennaCL API is identical to that of Boost.uBLAS and comparisons can be carried out immediately, as is shown in the tutorial located in <code><a class="el" href="examples_2tutorial_2blas3_8cpp.html">examples/tutorial/blas3.cpp</a></code>.</p>
<p>As for performance, ViennaCL yields decent performance gains at BLAS level 3 on mid- to high-end GPUs. However, highest performance is usually obtained only with careful tuning to the respective target device. ViennaCL provides kernels that are device-specific and thus achieve high efficiency, without sacrificing the portability of OpenCL. Best performance is usually obtained with the OpenCL backend. The CUDA backend provides good performance for NVIDIA GPUs. The OpenMP backend, however, currently provides only moderate performance and is not recommended for BLAS level 3 operations (as well as algorithms depending on them).</p>
<dl class="section note"><dt>Note</dt><dd>Mixing operations between objects of different scalar types is not supported. Convert the data manually on the host if needed.</dd></dl>
<center> <table class="doxtable">
<tr>
<th>Verbal </th><th>Mathematics </th><th>ViennaCL  </th></tr>
<tr>
<td>matrix-matrix product </td><td><img class="formulaInl" alt="$ C \leftarrow A \times B $" src="form_85.png"/> </td><td><code>C = prod(A, B);</code>  </td></tr>
<tr>
<td>matrix-matrix product </td><td><img class="formulaInl" alt="$ C \leftarrow A \times B^\mathrm{T} $" src="form_86.png"/> </td><td><code>C = prod(A, trans(B));</code>  </td></tr>
<tr>
<td>matrix-matrix product </td><td><img class="formulaInl" alt="$ C \leftarrow A^\mathrm{T} \times B $" src="form_87.png"/> </td><td><code>C = prod(trans(A), B);</code>  </td></tr>
<tr>
<td>matrix-matrix product </td><td><img class="formulaInl" alt="$ C \leftarrow A^\mathrm{T} \times B^\mathrm{T} $" src="form_88.png"/> </td><td><code>C = prod(trans(A), trans(B));</code>  </td></tr>
<tr>
<td></td><td></td><td></td></tr>
<tr>
<td>tri. matrix solve </td><td><img class="formulaInl" alt="$ C \leftarrow A^{-1} B $" src="form_89.png"/> </td><td><code>C = solve(A, B, tag);</code>  </td></tr>
<tr>
<td>tri. matrix solve </td><td><img class="formulaInl" alt="$ C \leftarrow A^\mathrm{T^{-1}} B $" src="form_90.png"/> </td><td><code>C = solve(trans(A), B, tag);</code>  </td></tr>
<tr>
<td>tri. matrix solve </td><td><img class="formulaInl" alt="$ C \leftarrow A^{-1} B^\mathrm{T} $" src="form_91.png"/> </td><td><code>C = solve(A, trans(B), tag);</code>  </td></tr>
<tr>
<td>tri. matrix solve </td><td><img class="formulaInl" alt="$ C \leftarrow A^\mathrm{T^{-1}} B^\mathrm{T} $" src="form_92.png"/> </td><td><code>C = solve(trans(A), trans(B), tag);</code>  </td></tr>
<tr>
<td></td><td></td><td></td></tr>
<tr>
<td>inplace solve </td><td><img class="formulaInl" alt="$ B \leftarrow A^{-1} B $" src="form_93.png"/> </td><td><code>inplace_solve(A, trans(B), tag);</code>  </td></tr>
<tr>
<td>inplace solve </td><td><img class="formulaInl" alt="$ B \leftarrow A^\mathrm{T^{-1}} B $" src="form_94.png"/> </td><td><code>inplace_solve(trans(A), x, tag);</code>  </td></tr>
<tr>
<td>inplace solve </td><td><img class="formulaInl" alt="$ B \leftarrow A^{-1} B^\mathrm{T} $" src="form_95.png"/> </td><td><code>inplace_solve(A, trans(B), tag);</code>  </td></tr>
<tr>
<td>inplace solve </td><td><img class="formulaInl" alt="$ B \leftarrow A^\mathrm{T^{-1}} B^\mathrm{T} $" src="form_96.png"/> </td><td><code>inplace_solve(trans(A), x, tag);</code>  </td></tr>
</table>
<p><b>BLAS level 3 routines mapped to ViennaCL. Note that the free functions reside in namespace <code><a class="el" href="namespaceviennacl_1_1linalg.html" title="Provides all linear algebra operations which are not covered by operator overloads. ">viennacl::linalg</a></code></b> </center><dl class="section warning"><dt>Warning</dt><dd>The operator overloads make extensive use of expression templates. Do not use the C++11 keyword <code>auto</code> for the result type, as this might result in unexpected performance regressions or dangling references.</dd></dl>
<h1><a class="anchor" id="manual-operations-sparse"></a>
Sparse Matrix Operations</h1>
<h2><a class="anchor" id="manual-operations-sparse-spmv"></a>
Sparse Matrix-Vector Products</h2>
<p>The most important operation on sparse matrices is the sparse matrix-vector product, available for all sparse matrix types in ViennaCL. An example for <code>compressed_matrix&lt;T&gt;</code> is as follows: </p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> <a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a>() {</div>
<div class="line">  <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;T&gt;</a> A;</div>
<div class="line">  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;T&gt;</a> x;</div>
<div class="line"></div>
<div class="line">  <span class="comment">// fill A and x with data here</span></div>
<div class="line"></div>
<div class="line">  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;T&gt;</a> y = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(A, x);</div>
<div class="line">}</div>
</div><!-- fragment --><p>For best performance we recommend <code>compressed_matrix&lt;T&gt;</code>, <code>hyb_matrix&lt;T&gt;</code>, or <code>sliced_ell_matrix&lt;T&gt;</code>. Unfortunately it is not possible to predict the fastest matrix type in general, thus a certain amount of trial-and-error by the library user is required.</p>
<p>Have a look at the <a href="http://viennacl.sourceforge.net/viennacl-benchmarks.html">ViennaCL Benchmark Section</a> for a performance comparison of sparse matrix-vector products for different storage formats.</p>
<h2><a class="anchor" id="manual-operations-sparse-spgemm"></a>
Sparse Matrix-Matrix Products</h2>
<p>Several graph algorithms or multigrid methods require the multiplication of two sparse matrices. ViennaCL provides an implementation based on the row-merge algorithm developed by Gremse et al. for NVIDIA GPUs <a class="el" href="citelist.html#CITEREF_Gremse:SpGEMM">[15]</a>. We extended the row-merge algorithm and derived high-performance implemented for CPUs, OpenCL-enabled GPUs from AMD, and INTEL's Xeon Phi. The calling code, however, remains as simple as possible: </p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> <a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a>() {</div>
<div class="line">  <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;T&gt;</a> A;</div>
<div class="line">  <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;T&gt;</a> B;</div>
<div class="line"></div>
<div class="line">  <span class="comment">// fill A and B with data here</span></div>
<div class="line"></div>
<div class="line">  <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;T&gt;</a> C = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(A, B);</div>
<div class="line">}</div>
</div><!-- fragment --><p>Have a look at the <a href="http://viennacl.sourceforge.net/viennacl-benchmarks.html">ViennaCL Benchmark Section</a> for a performance comparison of sparse matrix-matrix products.</p>
<dl class="section note"><dt>Note</dt><dd>Sparse Matrix times Sparse Matrix products are only available for <code>compressed_matrix&lt;T&gt;</code>. Other sparse matrix formats do not allow for a high-performance implementation.</dd></dl>
<h1><a class="anchor" id="manual-operations-row-column-diagonal"></a>
Row, Column, and Diagonal Extraction</h1>
<p>For many algorithms it is of interest to extract a single row or column of a dense matrix, or to access the matrix diagonal. This is provided in the same way as for Boost.uBLAS through the free functions <code><a class="el" href="namespaceviennacl.html#a0a574e6cd04ca0e42298b4ab845700e4">row()</a></code>, <code><a class="el" href="namespaceviennacl.html#a7fca08f4a83edffe7f47666d298ca87d">column()</a></code>, and <code><a class="el" href="namespaceviennacl.html#a507d2ac469c79997f2bb6e82b37b7483">diag()</a></code>: </p>
<div class="fragment"><div class="line"><span class="comment">// A is a viennacl::matrix&lt;T&gt;</span></div>
<div class="line"><span class="comment">// Extract 5-th row of A, then overwrite with 6-th diagonal:</span></div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;T&gt;</a> r = <a class="code" href="namespaceviennacl.html#a0a574e6cd04ca0e42298b4ab845700e4">viennacl::row</a>(A, 4);</div>
<div class="line">r = <a class="code" href="namespaceviennacl.html#a0a574e6cd04ca0e42298b4ab845700e4">viennacl::row</a>(A, 5);</div>
<div class="line"></div>
<div class="line"><span class="comment">// Extract 4-th column of A, then overwrite with second column:</span></div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;T&gt;</a> c = <a class="code" href="namespaceviennacl.html#a7fca08f4a83edffe7f47666d298ca87d">viennacl::column</a>(A, 3);</div>
<div class="line">c = <a class="code" href="namespaceviennacl.html#a7fca08f4a83edffe7f47666d298ca87d">viennacl::column</a>(A, 1);</div>
<div class="line"></div>
<div class="line"><span class="comment">// Extract diagonal:</span></div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;T&gt;</a> d = <a class="code" href="namespaceviennacl.html#a507d2ac469c79997f2bb6e82b37b7483">viennacl::diag</a>(A);</div>
</div><!-- fragment --><p> The function <code>diag</code> can also be used to create a matrix which has the provided vector entries in the off-diagonal: </p>
<div class="fragment"><div class="line"><span class="comment">// Create the matrix</span></div>
<div class="line"><span class="comment">// 0 1 0 0</span></div>
<div class="line"><span class="comment">// 0 0 2 0</span></div>
<div class="line"><span class="comment">// 0 0 0 3</span></div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;float&gt;</a> v(3);</div>
<div class="line">v[0] = 1.0f; v[1] = 2.0f; v[2] = 3.0f;</div>
<div class="line"><a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix&lt;float&gt;</a> A = <a class="code" href="namespaceviennacl.html#a507d2ac469c79997f2bb6e82b37b7483">viennacl::diag</a>(v, 1);</div>
</div><!-- fragment --><p> This is similar to MATLAB's <code><a class="el" href="namespaceviennacl.html#a507d2ac469c79997f2bb6e82b37b7483">diag()</a></code> function.</p>
<h1><a class="anchor" id="manual-operations-conversion"></a>
Conversion Operations</h1>
<p>Operations such as the vector operation <code>x = y + z;</code> in ViennaCL require that the numeric type (i.e. the first template parameter for the type) is the same. Thus, it is not possible to e.g. mix single-precision floating point vectors with integer vectors within the same operation. However, it is possible to convert vectors and dense matrices by direct assignment. An example code snippet for vectors is as follows: </p>
<div class="fragment"><div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;long&gt;</a> x_long(10);</div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;int&gt;</a>  x_int = <a class="code" href="structviennacl_1_1scalar__vector.html">viennacl::scalar_vector&lt;int&gt;</a>(10, 1);</div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;float&gt;</a> x_float(10);</div>
<div class="line"></div>
<div class="line">x_long = x_int;   <span class="comment">// conversion from int to long</span></div>
<div class="line">x_float = x_long; <span class="comment">// conversion from long to float</span></div>
</div><!-- fragment --> </div></div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
  <ul>
    <li class="footer">Generated on Wed Jan 20 2016 22:32:44 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>