1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516
|
<!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: matrix-free.cpp</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('matrix-free_8cpp-example.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">matrix-free.cpp</div> </div>
</div><!--header-->
<div class="contents">
<p>This tutorial explains how to use the iterative solvers in ViennaCL in a matrix-free manner, i.e. without explicitly assembling a matrix.</p>
<p>We consider the solution of the Poisson equation <img class="formulaInl" alt="$ \Delta \varphi = -1 $" src="form_138.png"/> on the unit square <img class="formulaInl" alt="$ [0,1] \times [0,1] $" src="form_139.png"/> with homogeneous Dirichlet boundary conditions using a finite-difference discretization. A <img class="formulaInl" alt="$ N \times N $" src="form_140.png"/> grid is used, where the first and the last points per dimensions represent the boundary. For simplicity we only consider the host-backend here. Have a look at custom-kernels.hpp and custom-cuda.cu on how to use custom kernels in such a matrix-free setting.</p>
<dl class="section note"><dt>Note</dt><dd><a class="el" href="matrix-free_8cpp.html">matrix-free.cpp</a> and matrix-free.cu are identical, the latter being required for compilation using CUDA nvcc</dd></dl>
<p>We start with including the necessary system headers: </p>
<div class="fragment"><div class="line"><span class="comment">//</span></div>
<div class="line"><span class="comment">// include necessary system headers</span></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="comment">// ViennaCL includes</span></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="cg_8hpp.html">viennacl/linalg/cg.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="bicgstab_8hpp.html">viennacl/linalg/bicgstab.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="gmres_8hpp.html">viennacl/linalg/gmres.hpp</a>"</span></div>
</div><!-- fragment --><p> ViennaCL imposes two type requirements on a user-provided operator to compute <code>y = prod(A, x)</code> for the iterative solvers:</p>
<ul>
<li>A member function <code>apply()</code>, taking two ViennaCL base vectors <code>x</code> and <code>y</code> as arguments. This member function carries out the action of the matrix to the vector.</li>
<li>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> returning the length of the result vectors. Keep in mind that you can always wrap your existing classes accordingly to fit ViennaCL's interface requirements.</li>
</ul>
<p>We define a simple class for dealing with the <img class="formulaInl" alt="$ N \times N $" src="form_140.png"/> grid for solving Poisson's equation. It only holds the number of grid points per coordinate direction and implements the <code>apply()</code> and <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> routines. Depending on whether the host, OpenCL, or CUDA is used for the computation, the respective implementation is called. We skip the details for now and discuss (and implement) them at the end of this tutorial. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keyword">class </span>MyOperator</div>
<div class="line">{</div>
<div class="line"><span class="keyword">public</span>:</div>
<div class="line"> MyOperator(std::size_t N) : N_(N) {}</div>
<div class="line"></div>
<div class="line"> <span class="comment">// Dispatcher for y = Ax</span></div>
<div class="line"> <span class="keywordtype">void</span> apply(<a name="_a0"></a><a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"> </span>{</div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"> <span class="keywordflow">if</span> (<a name="a1"></a><a class="code" href="namespaceviennacl_1_1traits.html#aad443a796634a9648e19b4156fcd880d">viennacl::traits::active_handle_id</a>(x) == <a name="a2"></a><a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>)</div>
<div class="line"> apply_cuda(x, y);</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"> <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#aad443a796634a9648e19b4156fcd880d">viennacl::traits::active_handle_id</a>(x) == <a name="a3"></a><a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>)</div>
<div class="line"> apply_opencl(x, y);</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"> <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#aad443a796634a9648e19b4156fcd880d">viennacl::traits::active_handle_id</a>(x) == <a name="a4"></a><a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>)</div>
<div class="line"> apply_host(x, y);</div>
<div class="line"> }</div>
<div class="line"></div>
<div class="line"> std::size_t <a name="a5"></a><a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">size1</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> N_ * N_; }</div>
<div class="line"></div>
<div class="line"><span class="keyword">private</span>:</div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"> <span class="keywordtype">void</span> apply_cuda(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y) <span class="keyword">const</span>;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"> <span class="keywordtype">void</span> apply_opencl(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y) <span class="keyword">const</span>;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"> <span class="keywordtype">void</span> apply_host(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y) <span class="keyword">const</span>;</div>
<div class="line"></div>
<div class="line"> std::size_t N_;</div>
<div class="line">};</div>
</div><!-- fragment --> <h2>Main Program</h2>
<p>In the <code><a class="el" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main()</a></code> routine we create the right hand side vector, instantiate the operator, and then call the solver. </p>
<div class="fragment"><div class="line"><span class="keywordtype">int</span> <a name="a6"></a><a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a>()</div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keywordtype">float</span> <a name="a7"></a><a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>; <span class="comment">// feel free to change to double (and change OpenCL kernel argument types accordingly)</span></div>
<div class="line"></div>
<div class="line"> std::size_t N = 10;</div>
<div class="line"> <a name="_a8"></a><a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> rhs = <a name="_a9"></a><a class="code" href="structviennacl_1_1scalar__vector.html">viennacl::scalar_vector<ScalarType></a>(N*N, <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>(-1));</div>
<div class="line"> MyOperator<ScalarType> op(N);</div>
</div><!-- fragment --><p> Run the CG method with our on-the-fly operator. Use <code><a class="el" href="classviennacl_1_1linalg_1_1bicgstab__tag.html" title="A tag for the stabilized Bi-conjugate gradient solver. Used for supplying solver parameters and for d...">viennacl::linalg::bicgstab_tag()</a></code> or <code><a class="el" href="classviennacl_1_1linalg_1_1gmres__tag.html" title="A tag for the solver GMRES. Used for supplying solver parameters and for dispatching the solve() func...">viennacl::linalg::gmres_tag()</a></code> instead of <code><a class="el" href="classviennacl_1_1linalg_1_1cg__tag.html" title="A tag for the conjugate gradient Used for supplying solver parameters and for dispatching the solve()...">viennacl::linalg::cg_tag()</a></code> to solve using BiCGStab or GMRES, respectively. </p>
<div class="fragment"><div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> result = <a name="a10"></a><a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(op, rhs, <a name="_a11"></a><a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
</div><!-- fragment --><p> Pretty-Print solution vector to verify solution. (We use a slow direct element-access via <code>operator[]</code> here for convenience.) </p>
<div class="fragment"><div class="line">std::cout.precision(3);</div>
<div class="line">std::cout << std::fixed;</div>
<div class="line">std::cout << <span class="stringliteral">"Result value map: "</span> << std::endl;</div>
<div class="line">std::cout << std::endl << <span class="stringliteral">"^ y "</span> << std::endl;</div>
<div class="line"><span class="keywordflow">for</span> (std::size_t i=0; i<N; ++i)</div>
<div class="line">{</div>
<div class="line"> std::cout << <span class="stringliteral">"| "</span>;</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t j=0; j<N; ++j)</div>
<div class="line"> std::cout << result[i * N + j] << <span class="stringliteral">" "</span>;</div>
<div class="line"> std::cout << std::endl;</div>
<div class="line">}</div>
<div class="line">std::cout << <span class="stringliteral">"*---------------------------------------------> x"</span> << std::endl;</div>
</div><!-- fragment --><p> That's it, print a completion message. Read on for details on how to implement the actual compute kernels. </p>
<div class="fragment"><div class="line"> std::cout << <span class="stringliteral">"!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!"</span> << std::endl;</div>
<div class="line"></div>
<div class="line"> <span class="keywordflow">return</span> EXIT_SUCCESS;</div>
<div class="line">}</div>
</div><!-- fragment --> <h2>Implementation Details </h2>
<p>So far we have only looked at the code for the control logic. In the following we define the actual 'worker' code for the matrix-free implementations.</p>
<h3>Execution on Host </h3>
<p>Since the execution on the host has the smallest amount of boilerplate code surrounding it, we use this case as a starting point. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keywordtype">void</span> MyOperator<NumericT>::apply_host(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"></span>{</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> <span class="keyword">const</span> * values_x = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(x.<a name="a12"></a><a class="code" href="classviennacl_1_1vector__base.html#a64e98aea5aa298ad63e3832a04c19648">handle</a>());</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> * values_y = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(y.<a class="code" href="classviennacl_1_1vector__base.html#a64e98aea5aa298ad63e3832a04c19648">handle</a>());</div>
<div class="line"></div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dx = <a name="a13"></a><a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(N_ + 1);</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dy = <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(N_ + 1);</div>
</div><!-- fragment --><p> In the following we iterate over all <img class="formulaInl" alt="$ N \times N $" src="form_140.png"/> points and apply the five-point stencil directly. This is done in a straightforward manner for illustration purposes. Multi-threaded execution via OpenMP can be obtained by uncommenting the pragma below.</p>
<p>Feel free to apply additional optimizations with respect to data access patterns and the like. </p>
<div class="fragment"><div class="line"> <span class="comment">// feel free to use</span></div>
<div class="line"> <span class="comment">// #pragma omp parallel for</span></div>
<div class="line"> <span class="comment">// here</span></div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t i=0; i<N_; ++i)</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t j=0; j<N_; ++j)</div>
<div class="line"> {</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_right = (j < N_ - 1) ? values_x[ i *N_ + j + 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_left = (j > 0 ) ? values_x[ i *N_ + j - 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_top = (i < N_ - 1) ? values_x[(i+1)*N_ + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_bottom = (i > 0 ) ? values_x[(i-1)*N_ + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_center = values_x[i*N_ + j];</div>
<div class="line"></div>
<div class="line"> values_y[i*N_ + j] = ((value_right - value_center) / dx - (value_center - value_left) / dx) / dx</div>
<div class="line"> + ((value_top - value_center) / dy - (value_center - value_bottom) / dy) / dy;</div>
<div class="line"> }</div>
<div class="line">}</div>
</div><!-- fragment --> <h3>Execution via CUDA </h3>
<p>The host-based kernel code serves as a basis for the CUDA kernel. The only thing we have to adjust are the array bounds: We assign one CUDA threadblock per index <code>i</code>. For a fixed <code>i</code>, parallelization across all threads in the block is obtained with respect to <code>j</code>.</p>
<p>Again, feel free to apply additional optimizations with respect to data access patterns and the like... </p>
<div class="fragment"><div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line">__global__ <span class="keywordtype">void</span> apply_cuda_kernel(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> <span class="keyword">const</span> * values_x,</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> * values_y,</div>
<div class="line"> std::size_t N)</div>
<div class="line">{</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dx = <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / (N + 1);</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dy = <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / (N + 1);</div>
<div class="line"></div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t i = blockIdx.x; i < N; i += gridDim.x)</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t j = threadIdx.x; j < N; j += blockDim.x)</div>
<div class="line"> {</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_right = (j < N - 1) ? values_x[ i *N + j + 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_left = (j > 0 ) ? values_x[ i *N + j - 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_top = (i < N - 1) ? values_x[(i+1)*N + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_bottom = (i > 0 ) ? values_x[(i-1)*N + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_center = values_x[i*N + j];</div>
<div class="line"></div>
<div class="line"> values_y[i*N + j] = ((value_right - value_center) / dx - (value_center - value_left) / dx) / dx</div>
<div class="line"> + ((value_top - value_center) / dy - (value_center - value_bottom) / dy) / dy;</div>
<div class="line"> }</div>
<div class="line">}</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keywordtype">void</span> MyOperator<NumericT>::apply_cuda(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"></span>{</div>
<div class="line"> apply_cuda_kernel<<<128, 128>>>(<a name="a14"></a><a class="code" href="namespaceviennacl.html#ae7d5db0c2c91be75218db5b52c4d13da">viennacl::cuda_arg</a>(x), <a class="code" href="namespaceviennacl.html#ae7d5db0c2c91be75218db5b52c4d13da">viennacl::cuda_arg</a>(y), N_);</div>
<div class="line">}</div>
<div class="line"><span class="preprocessor">#endif</span></div>
</div><!-- fragment --> <h3>Execution via OpenCL </h3>
<p>The OpenCL kernel is almost identical to the CUDA kernel: Only a couple of keywords need to be replaced. Also, the sources need to be packed into a string: </p>
<div class="fragment"><div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"><span class="keyword">static</span> <span class="keyword">const</span> <span class="keywordtype">char</span> * my_compute_program =</div>
<div class="line"><span class="stringliteral">"typedef float NumericT; \n"</span></div>
<div class="line"><span class="stringliteral">"__kernel void apply_opencl_kernel(__global NumericT const * values_x, \n"</span></div>
<div class="line"><span class="stringliteral">" __global NumericT * values_y, \n"</span></div>
<div class="line"><span class="stringliteral">" unsigned int N) {\n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" NumericT dx = (NumericT)1 / (N + 1); \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT dy = (NumericT)1 / (N + 1); \n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" for (unsigned int i = get_group_id(0); i < N; i += get_num_groups(0)) \n"</span></div>
<div class="line"><span class="stringliteral">" for (unsigned int j = get_local_id(0); j < N; j += get_local_size(0)) { \n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" NumericT value_right = (j < N - 1) ? values_x[ i *N + j + 1] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_left = (j > 0 ) ? values_x[ i *N + j - 1] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_top = (i < N - 1) ? values_x[(i+1)*N + j ] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_bottom = (i > 0 ) ? values_x[(i-1)*N + j ] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_center = values_x[i*N + j]; \n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" values_y[i*N + j] = ((value_right - value_center) / dx - (value_center - value_left) / dx) / dx \n"</span></div>
<div class="line"><span class="stringliteral">" + ((value_top - value_center) / dy - (value_center - value_bottom) / dy) / dy; \n"</span></div>
<div class="line"><span class="stringliteral">" } \n"</span></div>
<div class="line"><span class="stringliteral">" } \n"</span>;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
</div><!-- fragment --><p> Before the kernel is called for the first time, the OpenCL program containing the kernel needs to be compiled. We use a simple singleton using a static variable to achieve that.</p>
<p>Except for the kernel compilation at the first invocation, the OpenCL kernel launch is just one line of code just like for CUDA. Refer to <a class="el" href="custom-kernels_8cpp.html">custom-kernels.cpp</a> for some more details. </p>
<div class="fragment"><div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keywordtype">void</span> MyOperator<NumericT>::apply_opencl(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"> </span>{</div>
<div class="line"> <a name="_a15"></a><a class="code" href="classviennacl_1_1ocl_1_1context.html">viennacl::ocl::context</a> & ctx = <span class="keyword">const_cast<</span><a class="code" href="classviennacl_1_1ocl_1_1context.html">viennacl::ocl::context</a> &<span class="keyword">></span>(viennacl::traits::opencl_handle(x).context());</div>
<div class="line"> <span class="keyword">static</span> <span class="keywordtype">bool</span> first_run = <span class="keyword">true</span>;</div>
<div class="line"> <span class="keywordflow">if</span> (first_run) {</div>
<div class="line"> ctx.<a name="a16"></a><a class="code" href="classviennacl_1_1ocl_1_1context.html#a9094ac71f0cdf80df698c0c84ebb483d">add_program</a>(my_compute_program, <span class="stringliteral">"my_compute_program"</span>);</div>
<div class="line"> first_run = <span class="keyword">false</span>;</div>
<div class="line"> }</div>
<div class="line"> <a name="_a17"></a><a class="code" href="classviennacl_1_1ocl_1_1kernel.html">viennacl::ocl::kernel</a> & my_kernel = ctx.<a name="a18"></a><a class="code" href="classviennacl_1_1ocl_1_1context.html#afe68ad8c2df449f40f45db6605038946">get_kernel</a>(<span class="stringliteral">"my_compute_program"</span>, <span class="stringliteral">"apply_opencl_kernel"</span>);</div>
<div class="line"></div>
<div class="line"> <a name="a19"></a><a class="code" href="namespaceviennacl_1_1ocl.html#a5f2022f653ea1cf364d20e3ff84dcada">viennacl::ocl::enqueue</a>(my_kernel(x, y, static_cast<cl_uint>(N_)));</div>
<div class="line"> }</div>
<div class="line"><span class="preprocessor">#endif</span></div>
</div><!-- fragment --> <h2>Full Example Code</h2>
<div class="fragment"><div class="line"><span class="comment">/* =========================================================================</span></div>
<div class="line"><span class="comment"> Copyright (c) 2010-2016, Institute for Microelectronics,</span></div>
<div class="line"><span class="comment"> Institute for Analysis and Scientific Computing,</span></div>
<div class="line"><span class="comment"> TU Wien.</span></div>
<div class="line"><span class="comment"> Portions of this software are copyright by UChicago Argonne, LLC.</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment"> -----------------</span></div>
<div class="line"><span class="comment"> ViennaCL - The Vienna Computing Library</span></div>
<div class="line"><span class="comment"> -----------------</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment"> Project Head: Karl Rupp rupp@iue.tuwien.ac.at</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment"> (A list of authors and contributors can be found in the PDF manual)</span></div>
<div class="line"><span class="comment"></span></div>
<div class="line"><span class="comment"> License: MIT (X11), see file LICENSE in the base directory</span></div>
<div class="line"><span class="comment">============================================================================= */</span></div>
<div class="line"></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="comment">// include necessary system headers</span></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="comment">// ViennaCL includes</span></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="cg_8hpp.html">viennacl/linalg/cg.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="bicgstab_8hpp.html">viennacl/linalg/bicgstab.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="gmres_8hpp.html">viennacl/linalg/gmres.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keyword">class </span>MyOperator</div>
<div class="line">{</div>
<div class="line"><span class="keyword">public</span>:</div>
<div class="line"> MyOperator(std::size_t N) : N_(N) {}</div>
<div class="line"></div>
<div class="line"> <span class="comment">// Dispatcher for y = Ax</span></div>
<div class="line"> <span class="keywordtype">void</span> apply(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"> </span>{</div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"> <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#aad443a796634a9648e19b4156fcd880d">viennacl::traits::active_handle_id</a>(x) == <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8ab58facda25e2c7e20d9fe1b5e62f46d2">viennacl::CUDA_MEMORY</a>)</div>
<div class="line"> apply_cuda(x, y);</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"> <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#aad443a796634a9648e19b4156fcd880d">viennacl::traits::active_handle_id</a>(x) == <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8adb37af613f34867568e4f6cf720c68b1">viennacl::OPENCL_MEMORY</a>)</div>
<div class="line"> apply_opencl(x, y);</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"> <span class="keywordflow">if</span> (<a class="code" href="namespaceviennacl_1_1traits.html#aad443a796634a9648e19b4156fcd880d">viennacl::traits::active_handle_id</a>(x) == <a class="code" href="namespaceviennacl.html#a00b40450b6b2fd87aad3527d9b2084b8a427356f0fb1b8d32b28f37e36b272df4">viennacl::MAIN_MEMORY</a>)</div>
<div class="line"> apply_host(x, y);</div>
<div class="line"> }</div>
<div class="line"></div>
<div class="line"> std::size_t <a class="code" href="namespaceviennacl_1_1traits.html#aa756f5d6820722094cae0d8b9bb6d5e2">size1</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> N_ * N_; }</div>
<div class="line"></div>
<div class="line"><span class="keyword">private</span>:</div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"> <span class="keywordtype">void</span> apply_cuda(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y) <span class="keyword">const</span>;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"> <span class="keywordtype">void</span> apply_opencl(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y) <span class="keyword">const</span>;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"> <span class="keywordtype">void</span> apply_host(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y) <span class="keyword">const</span>;</div>
<div class="line"></div>
<div class="line"> std::size_t N_;</div>
<div class="line">};</div>
<div class="line"></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">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keywordtype">float</span> <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>; <span class="comment">// feel free to change to double (and change OpenCL kernel argument types accordingly)</span></div>
<div class="line"></div>
<div class="line"> std::size_t N = 10;</div>
<div class="line"> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> rhs = <a class="code" href="structviennacl_1_1scalar__vector.html">viennacl::scalar_vector<ScalarType></a>(N*N, <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>(-1));</div>
<div class="line"> MyOperator<ScalarType> op(N);</div>
<div class="line"></div>
<div class="line"> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(op, rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
<div class="line"></div>
<div class="line"> std::cout.precision(3);</div>
<div class="line"> std::cout << std::fixed;</div>
<div class="line"> std::cout << <span class="stringliteral">"Result value map: "</span> << std::endl;</div>
<div class="line"> std::cout << std::endl << <span class="stringliteral">"^ y "</span> << std::endl;</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t i=0; i<N; ++i)</div>
<div class="line"> {</div>
<div class="line"> std::cout << <span class="stringliteral">"| "</span>;</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t j=0; j<N; ++j)</div>
<div class="line"> std::cout << result[i * N + j] << <span class="stringliteral">" "</span>;</div>
<div class="line"> std::cout << std::endl;</div>
<div class="line"> }</div>
<div class="line"> std::cout << <span class="stringliteral">"*---------------------------------------------> x"</span> << std::endl;</div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!"</span> << std::endl;</div>
<div class="line"></div>
<div class="line"> <span class="keywordflow">return</span> EXIT_SUCCESS;</div>
<div class="line">}</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keywordtype">void</span> MyOperator<NumericT>::apply_host(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"></span>{</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> <span class="keyword">const</span> * values_x = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(x.<a class="code" href="classviennacl_1_1vector__base.html#a64e98aea5aa298ad63e3832a04c19648">handle</a>());</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> * values_y = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(y.<a class="code" href="classviennacl_1_1vector__base.html#a64e98aea5aa298ad63e3832a04c19648">handle</a>());</div>
<div class="line"></div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dx = <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(N_ + 1);</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dy = <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(N_ + 1);</div>
<div class="line"></div>
<div class="line"> <span class="comment">// feel free to use</span></div>
<div class="line"> <span class="comment">// #pragma omp parallel for</span></div>
<div class="line"> <span class="comment">// here</span></div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t i=0; i<N_; ++i)</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t j=0; j<N_; ++j)</div>
<div class="line"> {</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_right = (j < N_ - 1) ? values_x[ i *N_ + j + 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_left = (j > 0 ) ? values_x[ i *N_ + j - 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_top = (i < N_ - 1) ? values_x[(i+1)*N_ + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_bottom = (i > 0 ) ? values_x[(i-1)*N_ + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_center = values_x[i*N_ + j];</div>
<div class="line"></div>
<div class="line"> values_y[i*N_ + j] = ((value_right - value_center) / dx - (value_center - value_left) / dx) / dx</div>
<div class="line"> + ((value_top - value_center) / dy - (value_center - value_bottom) / dy) / dy;</div>
<div class="line"> }</div>
<div class="line">}</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line">__global__ <span class="keywordtype">void</span> apply_cuda_kernel(<a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> <span class="keyword">const</span> * values_x,</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> * values_y,</div>
<div class="line"> std::size_t N)</div>
<div class="line">{</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dx = <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / (N + 1);</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> dy = <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a>(1) / (N + 1);</div>
<div class="line"></div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t i = blockIdx.x; i < N; i += gridDim.x)</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t j = threadIdx.x; j < N; j += blockDim.x)</div>
<div class="line"> {</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_right = (j < N - 1) ? values_x[ i *N + j + 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_left = (j > 0 ) ? values_x[ i *N + j - 1] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_top = (i < N - 1) ? values_x[(i+1)*N + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_bottom = (i > 0 ) ? values_x[(i-1)*N + j ] : 0;</div>
<div class="line"> <a class="code" href="tests_2src_2bisect_8cpp.html#a52b5d30a2d7b064678644a3bf49b7f6c">NumericT</a> value_center = values_x[i*N + j];</div>
<div class="line"></div>
<div class="line"> values_y[i*N + j] = ((value_right - value_center) / dx - (value_center - value_left) / dx) / dx</div>
<div class="line"> + ((value_top - value_center) / dy - (value_center - value_bottom) / dy) / dy;</div>
<div class="line"> }</div>
<div class="line">}</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_CUDA)</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keywordtype">void</span> MyOperator<NumericT>::apply_cuda(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"></span>{</div>
<div class="line"> apply_cuda_kernel<<<128, 128>>>(<a class="code" href="namespaceviennacl.html#ae7d5db0c2c91be75218db5b52c4d13da">viennacl::cuda_arg</a>(x), <a class="code" href="namespaceviennacl.html#ae7d5db0c2c91be75218db5b52c4d13da">viennacl::cuda_arg</a>(y), N_);</div>
<div class="line">}</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"><span class="keyword">static</span> <span class="keyword">const</span> <span class="keywordtype">char</span> * my_compute_program =</div>
<div class="line"><span class="stringliteral">"typedef float NumericT; \n"</span></div>
<div class="line"><span class="stringliteral">"__kernel void apply_opencl_kernel(__global NumericT const * values_x, \n"</span></div>
<div class="line"><span class="stringliteral">" __global NumericT * values_y, \n"</span></div>
<div class="line"><span class="stringliteral">" unsigned int N) {\n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" NumericT dx = (NumericT)1 / (N + 1); \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT dy = (NumericT)1 / (N + 1); \n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" for (unsigned int i = get_group_id(0); i < N; i += get_num_groups(0)) \n"</span></div>
<div class="line"><span class="stringliteral">" for (unsigned int j = get_local_id(0); j < N; j += get_local_size(0)) { \n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" NumericT value_right = (j < N - 1) ? values_x[ i *N + j + 1] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_left = (j > 0 ) ? values_x[ i *N + j - 1] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_top = (i < N - 1) ? values_x[(i+1)*N + j ] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_bottom = (i > 0 ) ? values_x[(i-1)*N + j ] : 0; \n"</span></div>
<div class="line"><span class="stringliteral">" NumericT value_center = values_x[i*N + j]; \n"</span></div>
<div class="line"></div>
<div class="line"><span class="stringliteral">" values_y[i*N + j] = ((value_right - value_center) / dx - (value_center - value_left) / dx) / dx \n"</span></div>
<div class="line"><span class="stringliteral">" + ((value_top - value_center) / dy - (value_center - value_bottom) / dy) / dy; \n"</span></div>
<div class="line"><span class="stringliteral">" } \n"</span></div>
<div class="line"><span class="stringliteral">" } \n"</span>;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#if defined(VIENNACL_WITH_OPENCL)</span></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> NumericT></div>
<div class="line"><span class="keywordtype">void</span> MyOperator<NumericT>::apply_opencl(<a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> <span class="keyword">const</span> & x, <a class="code" href="classviennacl_1_1vector__base.html">viennacl::vector_base<NumericT></a> & y)<span class="keyword"> const</span></div>
<div class="line"><span class="keyword"> </span>{</div>
<div class="line"> <a class="code" href="classviennacl_1_1ocl_1_1context.html">viennacl::ocl::context</a> & ctx = <span class="keyword">const_cast<</span><a class="code" href="classviennacl_1_1ocl_1_1context.html">viennacl::ocl::context</a> &<span class="keyword">></span>(viennacl::traits::opencl_handle(x).context());</div>
<div class="line"> <span class="keyword">static</span> <span class="keywordtype">bool</span> first_run = <span class="keyword">true</span>;</div>
<div class="line"> <span class="keywordflow">if</span> (first_run) {</div>
<div class="line"> ctx.<a class="code" href="classviennacl_1_1ocl_1_1context.html#a9094ac71f0cdf80df698c0c84ebb483d">add_program</a>(my_compute_program, <span class="stringliteral">"my_compute_program"</span>);</div>
<div class="line"> first_run = <span class="keyword">false</span>;</div>
<div class="line"> }</div>
<div class="line"> <a class="code" href="classviennacl_1_1ocl_1_1kernel.html">viennacl::ocl::kernel</a> & my_kernel = ctx.<a class="code" href="classviennacl_1_1ocl_1_1context.html#afe68ad8c2df449f40f45db6605038946">get_kernel</a>(<span class="stringliteral">"my_compute_program"</span>, <span class="stringliteral">"apply_opencl_kernel"</span>);</div>
<div class="line"></div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1ocl.html#a5f2022f653ea1cf364d20e3ff84dcada">viennacl::ocl::enqueue</a>(my_kernel(x, y, static_cast<cl_uint>(N_)));</div>
<div class="line"> }</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"></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="footer">Generated on Wed Jan 20 2016 22:32:38 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>
|