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
|
<!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 Types</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<script type="text/javascript">
$(document).ready(initResizable);
$(window).load(resizeHeight);
</script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/search.js"></script>
<script type="text/javascript">
$(document).ready(function() { searchBox.OnSelectItem(0); });
</script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr style="height: 56px;">
<td style="padding-left: 0.5em;">
<div id="projectname">ViennaCL - The Vienna Computing Library
 <span id="projectnumber">1.7.1</span>
</div>
<div id="projectbrief">Free open-source GPU-accelerated linear algebra and solver library.</div>
</td>
<td> <div id="MSearchBox" class="MSearchBoxInactive">
<span class="left">
<img id="MSearchSelect" src="search/mag_sel.png"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
alt=""/>
<input type="text" id="MSearchField" value="Search" accesskey="S"
onfocus="searchBox.OnSearchFieldFocus(true)"
onblur="searchBox.OnSearchFieldFocus(false)"
onkeyup="searchBox.OnSearchFieldChange(event)"/>
</span><span class="right">
<a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a>
</span>
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.6 -->
<script type="text/javascript">
var searchBox = new SearchBox("searchBox", "search",false,'Search');
</script>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
<div id="nav-tree">
<div id="nav-tree-contents">
<div id="nav-sync" class="sync"></div>
</div>
</div>
<div id="splitbar" style="-moz-user-select:none;"
class="ui-resizable-handle">
</div>
</div>
<script type="text/javascript">
$(document).ready(function(){initNavTree('manual-types.html','');});
</script>
<div id="doc-content">
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
<a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark"> </span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark"> </span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark"> </span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark"> </span>Files</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark"> </span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark"> </span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark"> </span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark"> </span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark"> </span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark"> </span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark"> </span>Macros</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(11)"><span class="SelectionMark"> </span>Pages</a></div>
<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<iframe src="javascript:void(0)" frameborder="0"
name="MSearchResults" id="MSearchResults">
</iframe>
</div>
<div class="header">
<div class="headertitle">
<div class="title">Basic Types </div> </div>
</div><!--header-->
<div class="contents">
<div class="textblock"><p>This chapter provides a brief overview of the basic interfaces and usage of the provided data types. Operations on the various types are explained in <a class="el" href="manual-operations.html">Basic Operations</a>.</p>
<h1><a class="anchor" id="manual-types-scalar"></a>
Scalar Type</h1>
<p>The scalar type <code>scalar<T></code> with template parameter <code>T</code> denoting the underlying CPU scalar type (<code>char</code>, <code>short</code>, <code>int</code>, <code>long</code>, <code>float</code> and <code>double</code>, if supported - see <a class="el" href="manual-introduction.html#manual-introduction-hardware-table">table of supported hardware</a>) and represents single scalar value on the computing device. <code>scalar<T></code> is designed to behave much like a scalar type on conventional host-based CPU processing, but library users have to keep in mind that every operation on <code>scalar<T>|</code> may require the launch of an appropriate compute kernel on the GPU, thus making the operation much slower then the conventional CPU equivalent. Even if the host-based computing backend of ViennaCL is used, some (small) overheads occur.</p>
<dl class="section note"><dt>Note</dt><dd>Be aware that operations between objects of type <code>scalar<T></code> (e.g.~additions, comparisons) have large overhead on GPU backends. A separate compute kernel launch is required for every operation in such case.</dd></dl>
<h2><a class="anchor" id="manual-types-scalar-usage"></a>
Example Usage</h2>
<p>The scalar type of ViennaCL can be used just like the built-in types, as the following snippet shows: </p>
<div class="fragment"><div class="line"><span class="keywordtype">float</span> cpu_float = 42.0f;</div>
<div class="line"><span class="keywordtype">double</span> cpu_double = 13.7603;</div>
<div class="line"><a class="code" href="classviennacl_1_1scalar.html">viennacl::scalar<float></a> gpu_float(3.1415f);</div>
<div class="line"><a class="code" href="classviennacl_1_1scalar.html">viennacl::scalar<double></a> gpu_double = 2.71828;</div>
<div class="line"></div>
<div class="line"><span class="comment">//conversions</span></div>
<div class="line">cpu_float = gpu_float;</div>
<div class="line">gpu_float = cpu_double; <span class="comment">//automatic transfer and conversion</span></div>
<div class="line"></div>
<div class="line">cpu_float = gpu_float * 2.0f;</div>
<div class="line">cpu_double = gpu_float - cpu_float;</div>
</div><!-- fragment --><p> Mixing built-in types with the ViennaCL scalar is usually not a problem. Nevertheless, since every operation requires OpenCL calls, such arithmetics should be used sparsingly.</p>
<dl class="section note"><dt>Note</dt><dd>It is not possible to assign a <code>scalar<float></code> to a <code>scalar<double></code> directly.</dd>
<dd>
Mixing operations between objects of different scalar types is not supported. Convert the data manually on the host if needed.</dd></dl>
<h2><a class="anchor" id="manual-types-scalar-members"></a>
Members</h2>
<p>Apart from suitably overloaded operators that mimic the behavior of the respective CPU counterparts, only a single public member function <code><a class="el" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb" title="Returns the generic memory handle of an object. Non-const version. ">handle()</a></code> is available:</p>
<center> <table class="doxtable">
<tr>
<th>Interface</th><th>Comment </th></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb" title="Returns the generic memory handle of an object. Non-const version. ">v.handle()</a></code></td><td>The memory handle (CPU, CUDA, or OpenCL) </td></tr>
</table>
<p><b>Interface of <code>scalar<T></code> in ViennaCL. Destructors and operator overloads for BLAS are not listed.</b> </center><h1><a class="anchor" id="manual-types-vector"></a>
Vector Type</h1>
<p>The main vector type in ViennaCL is <code>vector<T, alignment></code>, representing a chunk of memory on the compute device. <code>T</code> is the underlying scalar type (<code>char</code>, <code>short</code>, <code>int</code>, <code>long</code>, <code>float</code>, or <code>double</code> if supported, - see <a class="el" href="manual-introduction.html#manual-introduction-hardware-table">table of supported hardware</a>). Complex types are not supported in ViennaCL. The second template argument <code>alignment</code> is deprecated and should not be specified by the library user.</p>
<p>At construction, <code>vector<T, alignment></code> is initialized to have the supplied length and the memory is initialized to zero (similar to <code>std::vector<T></code>). A difference to CPU implementations is that accessing single vector elements is very costly, because every time an element is accessed, it has to be transferred from the CPU to the compute device or vice versa.</p>
<h2><a class="anchor" id="manual-types-vector-usage"></a>
Example Usage</h2>
<p>The following code snippet shows the typical use of the vector type provided by ViennaCL. The overloaded function <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code>, which is used similar to <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">std::copy()</a></code> from the C++ Standard Template Library (STL), should be used for writing vector entries: </p>
<div class="fragment"><div class="line">std::vector<ScalarType> stl_vec(10);</div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> vcl_vec(10);</div>
<div class="line"></div>
<div class="line"><span class="comment">//fill the STL vector:</span></div>
<div class="line"><span class="keywordflow">for</span> (<span class="keywordtype">size_t</span> i=0; i<stl_vec.size(); ++i)</div>
<div class="line"> stl_vec[i] = i;</div>
<div class="line"></div>
<div class="line"><span class="comment">//copy content to GPU vector (recommended initialization)</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(stl_vec.begin(), stl_vec.end(), vcl_vec.begin());</div>
<div class="line"></div>
<div class="line"><span class="comment">//manipulate GPU vector here</span></div>
<div class="line"></div>
<div class="line"><span class="comment">//copy content from GPU vector back to STL vector</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(vcl_vec.begin(), vcl_vec.end(), stl_vec.begin());</div>
</div><!-- fragment --><p>The function <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code> does not assume that the values of the supplied CPU object are located in a linear memory sequence. If this is the case, the function <code>fast_copy</code> provides better performance.</p>
<p>Once the vectors are set up on the GPU, they can be used like objects on the CPU (refer to <a class="el" href="manual-operations.html">Basic Operations</a> for more details): </p>
<div class="fragment"><div class="line"><span class="comment">// let vcl_vec1 and vcl_vec2 denote two vector on the GPU</span></div>
<div class="line">vcl_vec1 *= 2.0;</div>
<div class="line">vcl_vec2 += vcl_vec1;</div>
<div class="line">vcl_vec1 = vcl_vec1 - 3.0 * vcl_vec2;</div>
</div><!-- fragment --><h2><a class="anchor" id="manual-types-vector-members"></a>
Members</h2>
<p>At construction, <code>vector<T, alignment></code> is initialized to have the supplied length, but memory is not initialized. If initialization is desired, the memory can be initialized with zero values using the member function <a class="el" href="namespaceviennacl_1_1traits.html#a22ab64b1df12a9da0423e5cad52ea367" title="Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...">clear()</a>. Other member functions are as follows:</p>
<center> <table class="doxtable">
<tr>
<th>Interface </th><th>Comment </th></tr>
<tr>
<td><code>CTOR(n)</code> </td><td>Constructor with number of entries </td></tr>
<tr>
<td><code>v(i)</code> </td><td>Access to the <code>i</code>-th element of v (slow for GPUs!) </td></tr>
<tr>
<td><code>v[i]</code> </td><td>Access to the <code>i</code>-th element of v (slow for GPUs!) </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a22ab64b1df12a9da0423e5cad52ea367" title="Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...">v.clear()</a></code> </td><td>Initialize v with zeros </td></tr>
<tr>
<td><code>v.resize(n, bool preserve)</code> </td><td>Resize v to length n. Preserves old values if bool is true. </td></tr>
<tr>
<td><code>v.begin()</code> </td><td>Iterator to the begin of the matrix </td></tr>
<tr>
<td><code>v.end()</code> </td><td>Iterator to the end of the matrix </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0" title="Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.) ">v.size()</a></code> </td><td>Length of the vector </td></tr>
<tr>
<td><code>v.swap(v2)</code> </td><td>Swap the content of v with v2 </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a0cd530f8a76a5fb3d60a9230bea30bca" title="Helper routine for obtaining the buffer length of a ViennaCL vector. ">v.internal_size()</a></code> </td><td>Returns the number of entries allocated on the GPU (taking alignment into account) </td></tr>
<tr>
<td><code>v.empty()</code> </td><td>Shorthand notation for <code><a class="el" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0" title="Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.) ">v.size()</a> == 0</code> </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a22ab64b1df12a9da0423e5cad52ea367" title="Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...">v.clear()</a></code> </td><td>Sets all entries in v to zero </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb" title="Returns the generic memory handle of an object. Non-const version. ">v.handle()</a></code> </td><td>Returns the memory handle (needed for custom kernels, see <a class="el" href="manual-custom-kernels.html">Custom Compute Kernels</a>) </td></tr>
</table>
<p><b>Interface of <code>vector<T></code> in ViennaCL. Destructors and operator overloads for BLAS are not listed.</b> </center><dl class="section note"><dt>Note</dt><dd>Accessing single elements of a vector using operator() or operator[] is very slow for GPUs due to PCI-Express latency! Use with care!</dd></dl>
<p>One important difference to pure CPU implementations is that the bracket operator as well as the parenthesis operator are very slow, because for each access an OpenCL data transfer has to be initiated. The overhead of this transfer is orders of magnitude. For example: </p>
<div class="fragment"><div class="line"><span class="comment">// fill a vector on CPU</span></div>
<div class="line"><span class="keywordflow">for</span> (<span class="keywordtype">size_t</span> i=0; i<cpu_vector.size(); ++i)</div>
<div class="line"> cpu_vector(i) = 1e-3f;</div>
<div class="line"></div>
<div class="line"><span class="comment">// fill a ViennaCL vector - VERY SLOW with GPU backends!!</span></div>
<div class="line"><span class="keywordflow">for</span> (<span class="keywordtype">size_t</span> i=0; i<gpu_vector.size(); ++i)</div>
<div class="line"> vcl_vector(i) = 1e-3f;</div>
</div><!-- fragment --><p> The difference in execution speed is typically several orders of magnitude, therefore direct vector element access should be used only if a very small number of entries is accessed in this way. A much faster initialization is as follows: </p>
<div class="fragment"><div class="line"><span class="comment">// fill a vector on CPU</span></div>
<div class="line"><span class="keywordflow">for</span> (<span class="keywordtype">long</span> i=0; i<cpu_vector.size(); ++i)</div>
<div class="line"> cpu_vector(i) = 1e-3f;</div>
<div class="line"></div>
<div class="line"><span class="comment">// fill a vector on GPU with data from CPU - faster versions:</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(cpu_vector, vcl_vector); <span class="comment">//option 1</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(cpu_vector.begin(), cpu_vector.end(), vcl_vector.begin()); <span class="comment">//option 2</span></div>
</div><!-- fragment --><p> In this way, setup costs for the CPU vector and the ViennaCL vector are comparable.</p>
<h1><a class="anchor" id="manual-types-matrix"></a>
Dense Matrix Type</h1>
<p><code>matrix<T, F, alignment></code> represents a dense matrix. The second optional template argument <code>F</code> specifies the storage layout and defaults to <code>row_major</code>. As an alternative, a <code>column_major</code> memory layout can be used. The third template argument <code>alignment</code> denotes an alignment for the rows and columns for row-major and column-major memory layout and should no longer be specified by the user (cf. <code>alignment</code> for the <code>vector</code> type).</p>
<div class="image">
<img src="matrix-padding.svg" alt="matrix-padding.svg"/>
<div class="caption">
Memory layout of a row-major dense matrix. The rows and columns of the matrix may be padded, such that the internal buffer is of size `internal_size1()*internal_size2()` elements.</div></div>
<dl class="section note"><dt>Note</dt><dd>For efficiency purposes, the rows and columns of matrices in ViennaCL may be padded by additional zeros.</dd></dl>
<h2><a class="anchor" id="manual-types-matrix-usage"></a>
Example Usage</h2>
<p>The use of <code>matrix<T, F></code> is similar to that of the counterpart in Boost.uBLAS. The operators are overloaded similarly. </p>
<div class="fragment"><div class="line"><span class="comment">//set up a 3 by 5 matrix:</span></div>
<div class="line"><a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix<float></a> vcl_matrix(4, 5);</div>
<div class="line"></div>
<div class="line"><span class="comment">//fill it up:</span></div>
<div class="line">vcl_matrix(0,2) = 1.0;</div>
<div class="line">vcl_matrix(1,2) = -1.5;</div>
<div class="line">vcl_matrix(2,0) = 4.2;</div>
<div class="line">vcl_matrix(3,4) = 3.1415;</div>
</div><!-- fragment --><dl class="section note"><dt>Note</dt><dd>Accessing single elements of a matrix using <code>operator()</code> is very slow on GPU backends! Use with care!</dd></dl>
<p>A much better way is to initialize a dense matrix using the provided <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code> function: </p>
<div class="fragment"><div class="line"><span class="comment">//copy content from CPU matrix to GPU matrix</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(cpu_matrix, gpu_matrix);</div>
<div class="line"></div>
<div class="line"><span class="comment">//copy content from GPU matrix to CPU matrix</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(gpu_matrix, cpu_matrix);</div>
</div><!-- fragment --><p> The type requirement for a class instantiated in an object <code>cpu_matrix is that</code>operator() can be used for accessing entries, that a member function <code><a class="el" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2" title="Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.) ">size1()</a></code> returns the number of rows and that <code><a class="el" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98" title="Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...">size2()</a></code> returns the number of columns. Please refer to <a class="el" href="manual-interfacing.html">Interfacing Other Libraries</a> for an overview of other libraries for which an overload of <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code> is provided.</p>
<dl class="section note"><dt>Note</dt><dd>The internal memory buffer of a <code>matrix<></code> is by default padded with zeros so that the internal matrix size is a multiple of e.g. a power of two.</dd>
<dd>
When using <code><a class="el" href="namespaceviennacl.html#aaa5c8726b45bc89a523ca2fa8c42107a">fast_copy()</a></code> on a matrix, the padded zeros need to be taken into account correctly. Query <code><a class="el" href="namespaceviennacl_1_1traits.html#aedc426f055f1b4c5d00111cc8a46e50d" title="Helper routine for obtaining the internal number of entries per row of a ViennaCL matrix...">internal_size1()</a></code> and <code><a class="el" href="namespaceviennacl_1_1traits.html#a2cd7269b5d00b5ebadab1c4b5a94a7c1" title="Helper routine for obtaining the internal number of entries per column of a ViennaCL matrix...">internal_size2()</a></code> to do so.</dd></dl>
<h2><a class="anchor" id="manual-types-matrix-members"></a>
Members</h2>
<p>The members are as follows, with the usual operator overloads not listed explicitly:</p>
<center> <table class="doxtable">
<tr>
<th>Interface </th><th>Comment </th></tr>
<tr>
<td><code>CTOR(nrows, ncols)</code> </td><td>Constructor with number of rows and columns </td></tr>
<tr>
<td><code>mat(i,j)</code> </td><td>Access to the element in the <code>i</code>-th row and the <code>j</code>-th column of <code>mat</code> </td></tr>
<tr>
<td><code>mat.resize(m, n, bool preserve)</code> </td><td>Resize mat to m rows and n columns. Currently, the boolean flag is ignored and entries always discarded. </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2" title="Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.) ">mat.size1()</a></code> </td><td>Number of rows in mat </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#aedc426f055f1b4c5d00111cc8a46e50d" title="Helper routine for obtaining the internal number of entries per row of a ViennaCL matrix...">mat.internal_size1()</a></code> </td><td>Internal number of rows in mat </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98" title="Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...">mat.size2()</a></code> </td><td>Number of columns in mat </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a2cd7269b5d00b5ebadab1c4b5a94a7c1" title="Helper routine for obtaining the internal number of entries per column of a ViennaCL matrix...">mat.internal_size2()</a></code> </td><td>Internal number of columns in mat </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a22ab64b1df12a9da0423e5cad52ea367" title="Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...">mat.clear()</a></code> </td><td>Sets all entries in v to zero </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb" title="Returns the generic memory handle of an object. Non-const version. ">mat.handle()</a></code> </td><td>Returns the memory handle (needed for custom kernels, see <a class="el" href="manual-custom-kernels.html">Custom Compute Kernels</a>) </td></tr>
</table>
<p><b>Interface of the dense matrix type <code>matrix<T, F></code> in ViennaCL. Constructors, Destructors and operator overloads for BLAS are not listed.</b> </center><h1><a class="anchor" id="manual-operations-initializers"></a>
Initializer Types</h1>
<dl class="section note"><dt>Note</dt><dd>Initializer types in ViennaCL can currently only be used for initializing vectors and matrices, not for computations!</dd></dl>
<p>In order to initialize vectors, the following initializer types are provided, again similar to Boost.uBLAS: </p>
<center> <table class="doxtable">
<tr>
<td><code>unit_vector<T>(s, i)</code> </td><td>Unit vector of size <img class="formulaInl" alt="$ s $" src="form_101.png"/> with entry <img class="formulaInl" alt="$ 1 $" src="form_102.png"/> at index <img class="formulaInl" alt="$ i $" src="form_103.png"/>, zero elsewhere. </td></tr>
<tr>
<td><code>zero_vector<T>(s)</code> </td><td>Vector of size <img class="formulaInl" alt="$ s $" src="form_101.png"/> with all entries being zero. </td></tr>
<tr>
<td><code>scalar_vector<T>(s, v)</code> </td><td>Vector of size <img class="formulaInl" alt="$ s $" src="form_101.png"/> with all entries equal to <img class="formulaInl" alt="$ v $" src="form_104.png"/>. </td></tr>
<tr>
<td><code>random_vector<T>(s, d)</code> </td><td>Vector of size <img class="formulaInl" alt="$ s $" src="form_101.png"/> with all entries random according to the distribution specified by <img class="formulaInl" alt="$ d $" src="form_15.png"/>. </td></tr>
</table>
</center><p> For example, to initialize a vector <code>v1</code> with all <img class="formulaInl" alt="$ 42 $" src="form_105.png"/> entries being <img class="formulaInl" alt="$ 42.0 $" src="form_106.png"/>, use </p>
<pre class="fragment">viennacl::vector<float> v1 = viennacl::scalar_vector<float>(42, 42.0f);
</pre><p>Similarly the following initializer types are available for matrices: </p>
<center> <table class="doxtable">
<tr>
<td><code>identity_matrix<T>(s)</code> </td><td>Identity matrix of dimension <img class="formulaInl" alt="$ s \times s $" src="form_107.png"/>. </td></tr>
<tr>
<td><code>zero_matrix<T>(s1, s2)</code> </td><td>Matrix of size <img class="formulaInl" alt="$ s_1 \times s_2 $" src="form_108.png"/> with all entries being zero. </td></tr>
<tr>
<td><code>scalar_matrix<T>(s1, s2, v)</code> </td><td>Matrix of size <img class="formulaInl" alt="$ s_1 \times s_2 $" src="form_108.png"/> with all entries equal to <img class="formulaInl" alt="$ v $" src="form_104.png"/>. </td></tr>
<tr>
<td><code>random_matrix<T>(s1, s2, d)</code> </td><td>Vector of size <img class="formulaInl" alt="$ s $" src="form_101.png"/> with all entries random according to the distribution specified by <img class="formulaInl" alt="$ d $" src="form_15.png"/>. </td></tr>
</table>
</center><h1><a class="anchor" id="manual-types-sparse"></a>
Sparse Matrix Types</h1>
<p>Several different sparse matrix types are provided in ViennaCL, which will be discussed in the following.</p>
<h2><a class="anchor" id="manual-types-sparse-compressed"></a>
Compressed Matrix</h2>
<p><code>compressed_matrix<T, alignment></code> represents a sparse matrix using a compressed sparse row (CSR) scheme, for which a sparse matrix-vector multiplication kernel based on CSR-adaptive <a class="el" href="citelist.html#CITEREF_Greathouse-CSR-adaptive">[14]</a> is available. <code>T</code> is the floating point type. <code>alignment</code> is the alignment and defaults to <code>1</code> at present. In general, sparse matrices should be set up on the CPU and then be pushed to the compute device using <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code>, because dynamic memory management of sparse matrices is not provided on OpenCL compute devices such as GPUs.</p>
<p><a class="anchor" id="manual-types-sparse-compressed-table-members"></a></p>
<center> <table class="doxtable">
<tr>
<th>Interface</th><th>Comment </th></tr>
<tr>
<td><code>CTOR(nrows, ncols)</code></td><td>Constructor with number of rows and columns </td></tr>
<tr>
<td><code>mat.set()</code> </td><td>Initialize mat with the data provided as arguments </td></tr>
<tr>
<td><code>mat.reserve(num)</code> </td><td><p class="starttd">Reserve memory for up to <code>num</code> nonzero entries </p>
<p class="endtd"></p>
</td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2" title="Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.) ">mat.size1()</a></code> </td><td>Number of rows in <code>mat</code> </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98" title="Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...">mat.size2()</a></code> </td><td>Number of columns in <code>mat</code> </td></tr>
<tr>
<td><code>mat.nnz()</code> </td><td>Number of nonzeroes in <code>mat</code> </td></tr>
<tr>
<td><code>mat.resize(m, n, bool preserve)</code> </td><td>Resize mat to <code>m</code> rows and <code>n</code> columns. Currently, the boolean flag is ignored and entries always discarded. </td></tr>
<tr>
<td><code>mat.handle1()</code> </td><td>Returns the memory handle holding the row indices (needed for custom kernels, see <a class="el" href="manual-custom-kernels.html">Custom Compute Kernels</a>) </td></tr>
<tr>
<td><code>mat.handle2()</code> </td><td>Returns the memory handle holding the column indices (needed for custom kernels, see <a class="el" href="manual-custom-kernels.html">Custom Compute Kernels</a>) </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb" title="Returns the generic memory handle of an object. Non-const version. ">mat.handle()</a></code> </td><td>Returns the memory handle holding the entries (needed for custom kernels, see <a class="el" href="manual-custom-kernels.html">Custom Compute Kernels</a>) </td></tr>
</table>
<p><b>Interface of the sparse matrix type <code>compressed_matrix<T, F></code> in ViennaCL. Destructors and operator overloads for BLAS are not listed.</b> </center><h3><a class="anchor" id="manual-types-sparse-compressed-usage"></a>
Example Usage</h3>
<p>The use of <code>compressed_matrix<T, alignment></code> is similar to that of the counterpart in Boost.uBLAS. The operators are overloaded similarly. There is a direct interfacing with the standard implementation using a vector of maps from the STL: </p>
<div class="fragment"><div class="line"><span class="comment">//set up a sparse 3 by 5 matrix on the CPU:</span></div>
<div class="line">std::vector< std::map< unsigned int, float> > cpu_sparse_matrix(4);</div>
<div class="line"></div>
<div class="line"><span class="comment">//fill it up:</span></div>
<div class="line">cpu_sparse_matrix[0][2] = 1.0;</div>
<div class="line">cpu_sparse_matrix[1][2] = -1.5;</div>
<div class="line">cpu_sparse_matrix[3][0] = 4.2;</div>
<div class="line"></div>
<div class="line"><span class="comment">//set up a sparse ViennaCL matrix:</span></div>
<div class="line"><a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<float></a> vcl_sparse_matrix(4, 5);</div>
<div class="line"></div>
<div class="line"><span class="comment">//copy to OpenCL device:</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(cpu_sparse_matrix, vcl_sparse_matrix);</div>
<div class="line"></div>
<div class="line"><span class="comment">//copy back to CPU:</span></div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">copy</a>(vcl_sparse_matrix, cpu_sparse_matrix);</div>
</div><!-- fragment --><p> The <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code> functions can also be used with a generic sparse matrix data type fulfilling the following requirements:</p>
<ul>
<li>The <code>const_iterator1</code> type is provided for iteration along increasing row index</li>
<li>The <code>const_iterator2</code> type is provided for iteration along increasing column index</li>
<li><code>.begin1()</code> returns an iterator pointing to the element with indices <code>(0,0)</code>.</li>
<li><code>.end1()</code> returns an iterator pointing to the end of the first column</li>
<li>When copying to the cpu type: Write operation via <code>operator()</code></li>
<li>When copying to the cpu type: <code>resize(m,n,preserve)</code> member (cf. <a class="el" href="manual-types.html#manual-types-sparse-compressed-table-members">Table of members</a>)</li>
</ul>
<p>The iterator returned from the cpu sparse matrix type via <code>begin1()</code> has to fulfill the following requirements:</p>
<ul>
<li><code>.begin()</code> returns an column iterator pointing to the first nonzero element in the particular row.</li>
<li><code>.end()</code> returns an iterator pointing to the end of the row</li>
<li>Increment and dereference</li>
</ul>
<p>For the sparse matrix types in Boost.uBLAS, these requirements are all fulfilled. Please refer to <a class="el" href="manual-interfacing.html">Interfacing Other Libraries</a> for an overview of other libraries for which an overload of <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code> is provided.</p>
<h3><a class="anchor" id="manual-types-sparse-compressed-members"></a>
Members</h3>
<p>The interface is described in <a class="el" href="manual-types.html#manual-types-sparse-compressed-table-members">Table of members</a>.</p>
<h2><a class="anchor" id="manual-types-sparse-coordinate"></a>
Coordinate Matrix</h2>
<p>In the second sparse matrix type, <code>coordinate_matrix<T, alignment></code>, entries are stored as triplets <code>(i,j,val)</code>, where <code>i</code> is the row index, <code>j</code> is the column index, and <code>val</code> is the entry. <code>T</code> is the floating point type. The optional <code>alignment</code> defaults to <code>128</code> at present and should not be provided by the user. In general, sparse matrices should be set up on the CPU and then be pushed to the compute device using <code><a class="el" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453" title="Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...">copy()</a></code>, because dynamic memory management of sparse matrices is not provided on OpenCL compute devices such as GPUs.</p>
<p><a class="anchor" id="manual-types-sparse-coordinate-table-members"></a></p>
<center> <table class="doxtable">
<tr>
<th>Interface</th><th>Comment </th></tr>
<tr>
<td><code>CTOR(nrows, ncols)</code> </td><td>Constructor with number of rows and columns </td></tr>
<tr>
<td><code>mat.reserve(num)</code> </td><td>Reserve memory for <code>num</code> nonzero entries </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2" title="Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.) ">mat.size1()</a></code> </td><td>Number of rows in <code>mat</code> </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#a3658e7c29ac0f60a20cb5871f5b5fd98" title="Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...">mat.size2()</a></code> </td><td>Number of columns in <code>mat</code> </td></tr>
<tr>
<td><code>mat.nnz()</code> </td><td>Number of nonzeroes in <code>mat</code> </td></tr>
<tr>
<td><code>mat.resize(m, n, bool preserve)</code> </td><td>Resize <code>mat</code> to <code>m</code> rows and <code>n</code> columns. Currently, the boolean flag is ignored and entries always discarded. </td></tr>
<tr>
<td><code>mat.resize(m, n)</code> </td><td>Resize mat to <code>m</code> rows and <code>n</code> columns. Does not preserve old values. </td></tr>
<tr>
<td><code>mat.handle12()</code> </td><td>Returns the memory handle holding the row and column indices (needed for custom kernels, see <a class="el" href="manual-custom-kernels.html">Custom Compute Kernels</a>) </td></tr>
<tr>
<td><code><a class="el" href="namespaceviennacl_1_1traits.html#ae39853b7f291a697e119a139439178fb" title="Returns the generic memory handle of an object. Non-const version. ">mat.handle()</a></code> </td><td>Returns the memory handle holding the entries (needed for custom kernels, see <a class="el" href="manual-custom-kernels.html">Custom Compute Kernels</a>) </td></tr>
</table>
<p><b>Interface of the sparse matrix type <code>coordinate_matrix<T, A></code> in ViennaCL. Destructors and operator overloads for BLAS operations are not listed.</b> </center><h3><a class="anchor" id="manual-types-sparse-coordinate-usage"></a>
Example Usage</h3>
<p>The use of <code>coordinate_matrix<T, alignment></code> is similar to that of the first sparse matrix type <code>compressed_matrix<T, alignment></code>, thus we refer to <a class="el" href="manual-types.html#manual-types-sparse-coordinate-usage">the example usage</a> of <code>compressed_matrix<></code>.</p>
<h3><a class="anchor" id="manual-types-sparse-coordinate-members"></a>
Members</h3>
<p>The interface is described in <a class="el" href="manual-types.html#manual-types-sparse-coordinate-table-members">this table</a>.</p>
<dl class="section note"><dt>Note</dt><dd>Note that only a few preconditioners work with <code>coordinate_matrix</code> so far.</dd></dl>
<h2><a class="anchor" id="manual-types-sparse-ell"></a>
ELL Matrix</h2>
<p>A sparse matrix in ELL format of type <code>ell_matrix</code> is stored in a block of memory of size <img class="formulaInl" alt="$ N \times n_{\max} $" src="form_109.png"/>, where <code>N</code> is the number of rows of the matrix and <img class="formulaInl" alt="$ n_{\max} $" src="form_110.png"/> is the maximum number of nonzeros per row. Rows with less than <img class="formulaInl" alt="$ n_{\max} $" src="form_110.png"/> entries are padded with zeros. In a second memory block, the respective column indices are stored.</p>
<p>The ELL format is well suited for matrices where most rows have approximately the same number of nonzeros. This is often the case for matrices arising from the discretization of partial differential equations using e.g. the finite element method. On the other hand, the ELL format introduces substantial overhead if the number of nonzeros per row varies a lot.</p>
<p>For an example use of an <code>ell_matrix</code>, have a look at <a class="el" href="examples_2benchmarks_2sparse_8cpp.html">examples/benchmarks/sparse.cpp</a>.</p>
<dl class="section note"><dt>Note</dt><dd>Note that preconditioners do not work with <code>ell_matrix</code> yet.</dd></dl>
<h2><a class="anchor" id="manual-types-sparse-sliced-ell"></a>
Sliced ELL Matrix</h2>
<p>A variation of the ELL format was recently proposed by Kreutzer et al. for use on CPUs, GPUs, and Intel's MIC architecture. The implementation in ViennaCL does not reorder the rows of the matrix, but is otherwise as proposed in the paper.</p>
<p>For an example use of <code>sliced_ell_matrix</code>, have a look at <a class="el" href="examples_2benchmarks_2sparse_8cpp.html">examples/benchmarks/sparse.cpp</a>.</p>
<dl class="section note"><dt>Note</dt><dd>Note that preconditioners do not work with <code>sliced_ell_matrix</code> yet.</dd></dl>
<h2><a class="anchor" id="manual-types-sparse-hyb"></a>
Hybrid Matrix</h2>
<p>The higher performance of the ELL format for matrices with approximately the same number of entries per row and the higher flexibility of the CSR format is combined in the <code>hyb_matrix</code> type, where the main part of the system matrix is stored in ELL format and excess entries are stored in CSR format.</p>
<p>For an example use of an <code>hyb_matrix</code>, have a look at <a class="el" href="examples_2benchmarks_2sparse_8cpp.html">examples/benchmarks/sparse.cpp</a>.</p>
<dl class="section note"><dt>Note</dt><dd>Note that preconditioners do not work with <code>hyb_matrix</code> yet.</dd></dl>
<h2><a class="anchor" id="manual-types-sparse-compressed-compressed"></a>
Compressed Compressed Matrix</h2>
<p>If only a few rows of a sparse matrix are populated, then the previous sparse matrix formats are fairly expensive in terms of memory consumption. This is addressed by the <code>compressed_compressed_matrix<></code> format, which is similar to the standard CSR format, but only stores the rows containing nonzero elements. An additional array is used to store the global row index <code>r</code> in the sparse matrix <code>A</code> of the <code>i</code>-th nonzero row.</p>
<dl class="section note"><dt>Note</dt><dd>Note that preconditioners do not work with <code>compressed_compressed_matrix</code> yet.</dd></dl>
<h1><a class="anchor" id="manual-types-proxies"></a>
Proxies</h1>
<p>Similar to Boost.uBLAS, ViennaCL provides <code>range</code> and <code>slice</code> objects in order to conveniently manipulate dense submatrices and vectors. The functionality is provided in the headers <code><a class="el" href="vector__proxy_8hpp.html" title="Proxy classes for vectors. ">viennacl/vector_proxy.hpp</a></code> and <code><a class="el" href="matrix__proxy_8hpp.html" title="Proxy classes for matrices. ">viennacl/matrix_proxy.hpp</a></code> respectively. A range refers to a contiguous integer interval and is set up as </p>
<pre class="fragment">std::size_t lower_bound = 1;
std::size_t upper_bound = 7;
viennacl::range r(lower_bound, upper_bound);
</pre><p>A slice is similar to a range and allows in addition for arbitrary increments (<em>stride</em>). For example, to create a slice consisting of the indices <code>2, 5, 8, 11, 14</code>, use the code </p>
<pre class="fragment">std::size_t start = 2;
std::size_t stride = 3;
std::size_t size = 5
viennacl::slice s(start, stride, size);
</pre><p>In order to address a subvector of a vector <code>v</code> and a submatrix of a matrix <code>M</code>, the proxy objects <code>v_sub</code> and <code>M_sub</code> are created as follows: </p>
<div class="fragment"><div class="line"><span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> VectorType;</div>
<div class="line"><span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix<ScalarType, viennacl::row_major></a> MatrixType;</div>
<div class="line"></div>
<div class="line"><a class="code" href="classviennacl_1_1vector__range.html">viennacl::vector_range<VCLVectorType></a> v_sub(v, r);</div>
<div class="line"><a class="code" href="classviennacl_1_1matrix__range.html">viennacl::matrix_range<VCLMatrixType></a> M_sub(M, r, r);</div>
</div><!-- fragment --><p> As a shortcut, one may use the free function <code><a class="el" href="namespaceviennacl.html#adc45a895937fe299100e2b235a442748">project()</a></code> in order to avoid having to write the type explicitly: </p>
<pre class="fragment">project(v, r); //returns a vector_range as above
project(M, r, r); //returns a matrix_range as above
</pre><p>In the same way a <code>vector_slice</code> and a <code>matrix_slice</code> are set up.</p>
<p>The proxy objects can now be manipulated in the same way as vectors and dense matrices. In particular, operations such as vector proxy additions and matrix additions work as usual, e.g. </p>
<pre class="fragment">vcl_sub += vcl_sub; //or project(v, r) += project(v, r);
M_sub += M_sub; //or project(M, r, r) += project(M, r, r);
</pre><p>Submatrix-Submatrix products are computed in the same manner and are handy for many block-based linear algebra algorithms.</p>
<p>Example code can be found in <a class="el" href="vector-range_8cpp.html">examples/tutorial/vector-range.cpp</a> and <a class="el" href="matrix-range_8cpp.html">examples/tutorial/matrix-range.cpp</a> </p>
</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>
|