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
|
<!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: spai.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('spai_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">spai.cpp</div> </div>
</div><!--header-->
<div class="contents">
<p>This tutorial shows how to use the sparse approximate inverse (SPAI) preconditioner.</p>
<dl class="section warning"><dt>Warning</dt><dd>SPAI is currently only available with the OpenCL backend and is experimental. API-changes may happen any time in the future.</dd></dl>
<p>We start with including the necessary headers: </p>
<div class="fragment"><div class="line"><span class="comment">// enable Boost.uBLAS support</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#ifndef NDEBUG</span></div>
<div class="line"><span class="preprocessor"> #define BOOST_UBLAS_NDEBUG</span></div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// System headers:</span></div>
<div class="line"><span class="preprocessor">#include <utility></span></div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"><span class="preprocessor">#include <fstream></span></div>
<div class="line"><span class="preprocessor">#include <string></span></div>
<div class="line"><span class="preprocessor">#include <cmath></span></div>
<div class="line"><span class="preprocessor">#include <algorithm></span></div>
<div class="line"><span class="preprocessor">#include <stdio.h></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// ViennaCL headers:</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="matrix_8hpp.html">viennacl/matrix.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.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="prod_8hpp.html">viennacl/linalg/prod.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="inner__prod_8hpp.html">viennacl/linalg/inner_prod.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="ilu_8hpp.html">viennacl/linalg/ilu.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="norm__2_8hpp.html">viennacl/linalg/norm_2.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="spai_8hpp.html">viennacl/linalg/spai.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Boost headers:</span></div>
<div class="line"><span class="preprocessor">#include "boost/numeric/ublas/vector.hpp"</span></div>
<div class="line"><span class="preprocessor">#include "boost/numeric/ublas/matrix.hpp"</span></div>
<div class="line"><span class="preprocessor">#include "boost/numeric/ublas/io.hpp"</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// auxiliary functionality:</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector-io_8hpp.html">vector-io.hpp</a>"</span></div>
</div><!-- fragment --><p> The following helper routine is used to run a solver with the provided preconditioner and to print the resulting residual norm. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> MatrixT, <span class="keyword">typename</span> VectorT, <span class="keyword">typename</span> SolverTagT, <span class="keyword">typename</span> PreconditionerT></div>
<div class="line"><span class="keywordtype">void</span> <a name="a0"></a><a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(MatrixT <span class="keyword">const</span> & A, VectorT <span class="keyword">const</span> & b, SolverTagT <span class="keyword">const</span> & solver_tag, PreconditionerT <span class="keyword">const</span> & precond)</div>
<div class="line">{</div>
<div class="line"> VectorT result = <a name="a1"></a><a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(A, b, solver_tag, precond);</div>
<div class="line"> std::cout << <span class="stringliteral">" * Solver iterations: "</span> << solver_tag.iters() << std::endl;</div>
<div class="line"> VectorT residual = <a name="a2"></a><a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(A, result);</div>
<div class="line"> residual -= b;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Rel. Residual: "</span> << <a name="a3"></a><a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(residual) / <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(b) << std::endl;</div>
<div class="line">}</div>
</div><!-- fragment --><p> The main steps in this tutorial are the following:</p>
<ul>
<li>Setup the systems</li>
<li>Run solvers without preconditioner and with ILUT preconditioner for comparison</li>
<li>Run solver with SPAI preconditioner on CPU</li>
<li>Run solver with SPAI preconditioner on GPU</li>
<li>Run solver with factored SPAI preconditioner on CPU</li>
<li>Run solver with factored SPAI preconditioner on GPU</li>
</ul>
<div class="fragment"><div class="line"><span class="keywordtype">int</span> <a name="a4"></a><a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a> (<span class="keywordtype">int</span>, <span class="keyword">const</span> <span class="keywordtype">char</span> **)</div>
<div class="line">{</div>
<div class="line"> <span class="keyword">typedef</span> <span class="keywordtype">float</span> <a name="a5"></a><a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>;</div>
<div class="line"> <span class="keyword">typedef</span> boost::numeric::ublas::compressed_matrix<ScalarType> MatrixType;</div>
<div class="line"> <span class="keyword">typedef</span> boost::numeric::ublas::vector<ScalarType> VectorType;</div>
<div class="line"> <span class="keyword">typedef</span> <a name="_a6"></a><a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<ScalarType></a> GPUMatrixType;</div>
<div class="line"> <span class="keyword">typedef</span> <a name="_a7"></a><a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> GPUVectorType;</div>
</div><!-- fragment --><p> If you have multiple OpenCL-capable devices in your system, we pick the second device for this tutorial. </p>
<div class="fragment"><div class="line"><span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"> <span class="comment">// Optional: Customize OpenCL backend</span></div>
<div class="line"> <a name="_a8"></a><a class="code" href="classviennacl_1_1ocl_1_1platform.html">viennacl::ocl::platform</a> pf = <a name="a9"></a><a class="code" href="namespaceviennacl_1_1ocl.html#ae1e931f6efd155240b33ca440dfcfef5">viennacl::ocl::get_platforms</a>()[0];</div>
<div class="line"> std::vector<viennacl::ocl::device> <span class="keyword">const</span> & devices = pf.<a name="a10"></a><a class="code" href="classviennacl_1_1ocl_1_1platform.html#afaa563522ebe9ce7b80ef02c40e7fe31">devices</a>();</div>
<div class="line"></div>
<div class="line"> <span class="comment">// Optional: Set first device to first context:</span></div>
<div class="line"> <a name="a11"></a><a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(0, devices[0]);</div>
<div class="line"></div>
<div class="line"> <span class="comment">// Optional: Set second device for second context (use the same device for the second context if only one device available):</span></div>
<div class="line"> <span class="keywordflow">if</span> (devices.size() > 1)</div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(1, devices[1]);</div>
<div class="line"> <span class="keywordflow">else</span></div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(1, devices[0]);</div>
<div class="line"></div>
<div class="line"> std::cout << <a name="a12"></a><a class="code" href="namespaceviennacl_1_1ocl.html#ac54d59a74deaccec81f64e738b856348">viennacl::ocl::current_device</a>().<a name="a13"></a><a class="code" href="classviennacl_1_1ocl_1_1device.html#a26a424782f982725453808a94763bd19">info</a>() << std::endl;</div>
<div class="line"> <a name="_a14"></a><a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx(<a name="a15"></a><a class="code" href="namespaceviennacl_1_1ocl.html#a82c1aba632a7ee0991eee480d5340966">viennacl::ocl::get_context</a>(1));</div>
<div class="line"><span class="preprocessor">#else</span></div>
<div class="line"> <a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
</div><!-- fragment --><p> Create uBLAS-based sparse matrix and read system matrix from file </p>
<div class="fragment"><div class="line">MatrixType M;</div>
<div class="line"></div>
<div class="line"><span class="keywordflow">if</span> (!<a name="a16"></a><a class="code" href="namespaceviennacl_1_1io.html#ad934125ed3dbe661e264bcd7d62b1048">viennacl::io::read_matrix_market_file</a>(M, <span class="stringliteral">"../examples/testdata/mat65k.mtx"</span>))</div>
<div class="line">{</div>
<div class="line"> std::cerr<<<span class="stringliteral">"ERROR: Could not read matrix file "</span> << std::endl;</div>
<div class="line"> exit(EXIT_FAILURE);</div>
<div class="line">}</div>
<div class="line"></div>
<div class="line">std::cout << <span class="stringliteral">"Size of matrix: "</span> << M.size1() << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">"Avg. Entries per row: "</span> << double(M.nnz()) / static_cast<double>(M.size1()) << std::endl;</div>
</div><!-- fragment --><p> Use a constant load vector for simplicity </p>
<div class="fragment"><div class="line">VectorType rhs(M.size2());</div>
<div class="line"><span class="keywordflow">for</span> (std::size_t i=0; i<rhs.size(); ++i)</div>
<div class="line"> rhs(i) = <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>(1);</div>
</div><!-- fragment --><p> Create the ViennaCL matrix and vector and initialize with uBLAS data: </p>
<div class="fragment"><div class="line">GPUMatrixType gpu_M(M.size1(), M.size2(), ctx);</div>
<div class="line">GPUVectorType gpu_rhs(M.size1(), ctx);</div>
<div class="line"><a name="a17"></a><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(M, gpu_M);</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(rhs, gpu_rhs);</div>
</div><!-- fragment --> <h2>Solver Runs</h2>
<p>We use a relative tolerance of <img class="formulaInl" alt="$ 10^{-10} $" src="form_23.png"/> with a maximum of 50 iterations for each use case. Usually more than 50 solver iterations are required for convergence, but this choice ensures shorter execution times and suffices for this tutorial. </p>
<div class="fragment"><div class="line"><a name="_a18"></a><a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a> solver_tag(1e-10, 50); <span class="comment">//for simplicity and reasonably short execution times we use only 50 iterations here</span></div>
</div><!-- fragment --><p> The first reference is to use no preconditioner (CPU and GPU): </p>
<div class="fragment"><div class="line">std::cout << <span class="stringliteral">"--- Reference 1: Pure BiCGStab on CPU ---"</span> << std::endl;</div>
<div class="line">VectorType result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(M, rhs, solver_tag);</div>
<div class="line">std::cout << <span class="stringliteral">" * Solver iterations: "</span> << solver_tag.iters() << std::endl;</div>
<div class="line">VectorType residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(M, result) - rhs;</div>
<div class="line">std::cout << <span class="stringliteral">" * Rel. Residual: "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(residual) / <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(rhs) << std::endl;</div>
<div class="line"></div>
<div class="line">std::cout << <span class="stringliteral">"--- Reference 2: Pure BiCGStab on GPU ---"</span> << std::endl;</div>
<div class="line">GPUVectorType gpu_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(gpu_M, gpu_rhs, solver_tag);</div>
<div class="line">std::cout << <span class="stringliteral">" * Solver iterations: "</span> << solver_tag.iters() << std::endl;</div>
<div class="line">GPUVectorType gpu_residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(gpu_M, gpu_result);</div>
<div class="line">gpu_residual -= gpu_rhs;</div>
<div class="line">std::cout << <span class="stringliteral">" * Rel. Residual: "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(gpu_residual) / <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(gpu_rhs) << std::endl;</div>
</div><!-- fragment --><p> The second reference is a standard ILUT preconditioner (only CPU): </p>
<div class="fragment"><div class="line">std::cout << <span class="stringliteral">"--- Reference 2: BiCGStab with ILUT on CPU ---"</span> << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"><a name="_a19"></a><a class="code" href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond<MatrixType></a> ilut(M, <a name="_a20"></a><a class="code" href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a>());</div>
<div class="line">std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"><a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, ilut);</div>
</div><!-- fragment --> <h2>Step 1: SPAI with CPU</h2>
<div class="fragment"><div class="line">std::cout << <span class="stringliteral">"--- Test 1: CPU-based SPAI ---"</span> << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"><a name="_a21"></a><a class="code" href="classviennacl_1_1linalg_1_1spai__precond.html">viennacl::linalg::spai_precond<MatrixType></a> spai_cpu(M, <a name="_a22"></a><a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1spai__tag.html">viennacl::linalg::spai_tag</a>(1e-3, 3, 5e-2));</div>
<div class="line">std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"><a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, spai_cpu);</div>
</div><!-- fragment --> <h2>Step 2: FSPAI with CPU</h2>
<div class="fragment"><div class="line">std::cout << <span class="stringliteral">"--- Test 2: CPU-based FSPAI ---"</span> << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"><a name="_a23"></a><a class="code" href="classviennacl_1_1linalg_1_1fspai__precond.html">viennacl::linalg::fspai_precond<MatrixType></a> fspai_cpu(M, <a name="_a24"></a><a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1fspai__tag.html">viennacl::linalg::fspai_tag</a>());</div>
<div class="line">std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"><a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, fspai_cpu);</div>
</div><!-- fragment --> <h2>Step 3: SPAI with GPU</h2>
<div class="fragment"><div class="line">std::cout << <span class="stringliteral">"--- Test 3: GPU-based SPAI ---"</span> << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"><a class="code" href="classviennacl_1_1linalg_1_1spai__precond.html">viennacl::linalg::spai_precond<GPUMatrixType></a> spai_gpu(gpu_M, <a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1spai__tag.html">viennacl::linalg::spai_tag</a>(1e-3, 3, 5e-2));</div>
<div class="line">std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"><a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(gpu_M, gpu_rhs, solver_tag, spai_gpu);</div>
</div><!-- fragment --> <h2>Step 4: FSPAI with GPU</h2>
<div class="fragment"><div class="line">std::cout << <span class="stringliteral">"--- Test 4: GPU-based FSPAI ---"</span> << std::endl;</div>
<div class="line">std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"><a class="code" href="classviennacl_1_1linalg_1_1fspai__precond.html">viennacl::linalg::fspai_precond<GPUMatrixType></a> fspai_gpu(gpu_M, <a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1fspai__tag.html">viennacl::linalg::fspai_tag</a>());</div>
<div class="line">std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"><a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(gpu_M, gpu_rhs, solver_tag, fspai_gpu);</div>
</div><!-- fragment --><p> That's it! Print success message and exit. </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>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">// enable Boost.uBLAS support</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS</span></div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#ifndef NDEBUG</span></div>
<div class="line"><span class="preprocessor"> #define BOOST_UBLAS_NDEBUG</span></div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// System headers:</span></div>
<div class="line"><span class="preprocessor">#include <utility></span></div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"><span class="preprocessor">#include <fstream></span></div>
<div class="line"><span class="preprocessor">#include <string></span></div>
<div class="line"><span class="preprocessor">#include <cmath></span></div>
<div class="line"><span class="preprocessor">#include <algorithm></span></div>
<div class="line"><span class="preprocessor">#include <stdio.h></span></div>
<div class="line"></div>
<div class="line"><span class="comment">// ViennaCL headers:</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="matrix_8hpp.html">viennacl/matrix.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.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="prod_8hpp.html">viennacl/linalg/prod.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="inner__prod_8hpp.html">viennacl/linalg/inner_prod.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="ilu_8hpp.html">viennacl/linalg/ilu.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="norm__2_8hpp.html">viennacl/linalg/norm_2.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>"</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="spai_8hpp.html">viennacl/linalg/spai.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Boost headers:</span></div>
<div class="line"><span class="preprocessor">#include "boost/numeric/ublas/vector.hpp"</span></div>
<div class="line"><span class="preprocessor">#include "boost/numeric/ublas/matrix.hpp"</span></div>
<div class="line"><span class="preprocessor">#include "boost/numeric/ublas/io.hpp"</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// auxiliary functionality:</span></div>
<div class="line"><span class="preprocessor">#include "<a class="code" href="vector-io_8hpp.html">vector-io.hpp</a>"</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span><<span class="keyword">typename</span> MatrixT, <span class="keyword">typename</span> VectorT, <span class="keyword">typename</span> SolverTagT, <span class="keyword">typename</span> PreconditionerT></div>
<div class="line"><span class="keywordtype">void</span> <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(MatrixT <span class="keyword">const</span> & A, VectorT <span class="keyword">const</span> & b, SolverTagT <span class="keyword">const</span> & solver_tag, PreconditionerT <span class="keyword">const</span> & precond)</div>
<div class="line">{</div>
<div class="line"> VectorT result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(A, b, solver_tag, precond);</div>
<div class="line"> std::cout << <span class="stringliteral">" * Solver iterations: "</span> << solver_tag.iters() << std::endl;</div>
<div class="line"> VectorT residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(A, result);</div>
<div class="line"> residual -= b;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Rel. Residual: "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(residual) / <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(b) << std::endl;</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> (<span class="keywordtype">int</span>, <span class="keyword">const</span> <span class="keywordtype">char</span> **)</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>;</div>
<div class="line"> <span class="keyword">typedef</span> boost::numeric::ublas::compressed_matrix<ScalarType> MatrixType;</div>
<div class="line"> <span class="keyword">typedef</span> boost::numeric::ublas::vector<ScalarType> VectorType;</div>
<div class="line"> <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix<ScalarType></a> GPUMatrixType;</div>
<div class="line"> <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector<ScalarType></a> GPUVectorType;</div>
<div class="line"></div>
<div class="line"><span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"> <span class="comment">// Optional: Customize OpenCL backend</span></div>
<div class="line"> <a class="code" href="classviennacl_1_1ocl_1_1platform.html">viennacl::ocl::platform</a> pf = <a class="code" href="namespaceviennacl_1_1ocl.html#ae1e931f6efd155240b33ca440dfcfef5">viennacl::ocl::get_platforms</a>()[0];</div>
<div class="line"> std::vector<viennacl::ocl::device> <span class="keyword">const</span> & devices = pf.<a class="code" href="classviennacl_1_1ocl_1_1platform.html#afaa563522ebe9ce7b80ef02c40e7fe31">devices</a>();</div>
<div class="line"></div>
<div class="line"> <span class="comment">// Optional: Set first device to first context:</span></div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(0, devices[0]);</div>
<div class="line"></div>
<div class="line"> <span class="comment">// Optional: Set second device for second context (use the same device for the second context if only one device available):</span></div>
<div class="line"> <span class="keywordflow">if</span> (devices.size() > 1)</div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(1, devices[1]);</div>
<div class="line"> <span class="keywordflow">else</span></div>
<div class="line"> <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(1, devices[0]);</div>
<div class="line"></div>
<div class="line"> std::cout << <a class="code" href="namespaceviennacl_1_1ocl.html#ac54d59a74deaccec81f64e738b856348">viennacl::ocl::current_device</a>().<a class="code" href="classviennacl_1_1ocl_1_1device.html#a26a424782f982725453808a94763bd19">info</a>() << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx(<a class="code" href="namespaceviennacl_1_1ocl.html#a82c1aba632a7ee0991eee480d5340966">viennacl::ocl::get_context</a>(1));</div>
<div class="line"><span class="preprocessor">#else</span></div>
<div class="line"> <a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx;</div>
<div class="line"><span class="preprocessor">#endif</span></div>
<div class="line"></div>
<div class="line"> MatrixType M;</div>
<div class="line"></div>
<div class="line"> <span class="keywordflow">if</span> (!<a class="code" href="namespaceviennacl_1_1io.html#ad934125ed3dbe661e264bcd7d62b1048">viennacl::io::read_matrix_market_file</a>(M, <span class="stringliteral">"../examples/testdata/mat65k.mtx"</span>))</div>
<div class="line"> {</div>
<div class="line"> std::cerr<<<span class="stringliteral">"ERROR: Could not read matrix file "</span> << std::endl;</div>
<div class="line"> exit(EXIT_FAILURE);</div>
<div class="line"> }</div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"Size of matrix: "</span> << M.size1() << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">"Avg. Entries per row: "</span> << double(M.nnz()) / static_cast<double>(M.size1()) << std::endl;</div>
<div class="line"></div>
<div class="line"> VectorType rhs(M.size2());</div>
<div class="line"> <span class="keywordflow">for</span> (std::size_t i=0; i<rhs.size(); ++i)</div>
<div class="line"> rhs(i) = <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>(1);</div>
<div class="line"></div>
<div class="line"> GPUMatrixType gpu_M(M.size1(), M.size2(), ctx);</div>
<div class="line"> GPUVectorType gpu_rhs(M.size1(), ctx);</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(M, gpu_M);</div>
<div class="line"> <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(rhs, gpu_rhs);</div>
<div class="line"></div>
<div class="line"> <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a> solver_tag(1e-10, 50); <span class="comment">//for simplicity and reasonably short execution times we use only 50 iterations here</span></div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"--- Reference 1: Pure BiCGStab on CPU ---"</span> << std::endl;</div>
<div class="line"> VectorType result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(M, rhs, solver_tag);</div>
<div class="line"> std::cout << <span class="stringliteral">" * Solver iterations: "</span> << solver_tag.iters() << std::endl;</div>
<div class="line"> VectorType residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(M, result) - rhs;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Rel. Residual: "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(residual) / <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(rhs) << std::endl;</div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"--- Reference 2: Pure BiCGStab on GPU ---"</span> << std::endl;</div>
<div class="line"> GPUVectorType gpu_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(gpu_M, gpu_rhs, solver_tag);</div>
<div class="line"> std::cout << <span class="stringliteral">" * Solver iterations: "</span> << solver_tag.iters() << std::endl;</div>
<div class="line"> GPUVectorType gpu_residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(gpu_M, gpu_result);</div>
<div class="line"> gpu_residual -= gpu_rhs;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Rel. Residual: "</span> << <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(gpu_residual) / <a class="code" href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a>(gpu_rhs) << std::endl;</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"--- Reference 2: BiCGStab with ILUT on CPU ---"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond<MatrixType></a> ilut(M, <a class="code" href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a>());</div>
<div class="line"> std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, ilut);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"--- Test 1: CPU-based SPAI ---"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1linalg_1_1spai__precond.html">viennacl::linalg::spai_precond<MatrixType></a> spai_cpu(M, <a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1spai__tag.html">viennacl::linalg::spai_tag</a>(1e-3, 3, 5e-2));</div>
<div class="line"> std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, spai_cpu);</div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"--- Test 2: CPU-based FSPAI ---"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1linalg_1_1fspai__precond.html">viennacl::linalg::fspai_precond<MatrixType></a> fspai_cpu(M, <a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1fspai__tag.html">viennacl::linalg::fspai_tag</a>());</div>
<div class="line"> std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, fspai_cpu);</div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"--- Test 3: GPU-based SPAI ---"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1linalg_1_1spai__precond.html">viennacl::linalg::spai_precond<GPUMatrixType></a> spai_gpu(gpu_M, <a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1spai__tag.html">viennacl::linalg::spai_tag</a>(1e-3, 3, 5e-2));</div>
<div class="line"> std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(gpu_M, gpu_rhs, solver_tag, spai_gpu);</div>
<div class="line"></div>
<div class="line"> std::cout << <span class="stringliteral">"--- Test 4: GPU-based FSPAI ---"</span> << std::endl;</div>
<div class="line"> std::cout << <span class="stringliteral">" * Preconditioner setup..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="classviennacl_1_1linalg_1_1fspai__precond.html">viennacl::linalg::fspai_precond<GPUMatrixType></a> fspai_gpu(gpu_M, <a class="code" href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1fspai__tag.html">viennacl::linalg::fspai_tag</a>());</div>
<div class="line"> std::cout << <span class="stringliteral">" * Iterative solver run..."</span> << std::endl;</div>
<div class="line"> <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(gpu_M, gpu_rhs, solver_tag, fspai_gpu);</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><!-- 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>
|