File: iterative_8cpp-example.html

package info (click to toggle)
viennacl 1.7.1%2Bdfsg1-6
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, sid, trixie
  • size: 114,428 kB
  • sloc: sh: 454,206; cpp: 109,088; ansic: 2,103; perl: 104; makefile: 22
file content (462 lines) | stat: -rw-r--r-- 58,277 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
<!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: iterative.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
   &#160;<span id="projectnumber">1.7.1</span>
   </div>
   <div id="projectbrief">Free open-source GPU-accelerated linear algebra and solver library.</div>
  </td>
   <td>        <div id="MSearchBox" class="MSearchBoxInactive">
        <span class="left">
          <img id="MSearchSelect" src="search/mag_sel.png"
               onmouseover="return searchBox.OnSearchSelectShow()"
               onmouseout="return searchBox.OnSearchSelectHide()"
               alt=""/>
          <input type="text" id="MSearchField" value="Search" accesskey="S"
               onfocus="searchBox.OnSearchFieldFocus(true)" 
               onblur="searchBox.OnSearchFieldFocus(false)" 
               onkeyup="searchBox.OnSearchFieldChange(event)"/>
          </span><span class="right">
            <a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a>
          </span>
        </div>
</td>
 </tr>
 </tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.6 -->
<script type="text/javascript">
var searchBox = new SearchBox("searchBox", "search",false,'Search');
</script>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
  <div id="nav-tree">
    <div id="nav-tree-contents">
      <div id="nav-sync" class="sync"></div>
    </div>
  </div>
  <div id="splitbar" style="-moz-user-select:none;" 
       class="ui-resizable-handle">
  </div>
</div>
<script type="text/javascript">
$(document).ready(function(){initNavTree('iterative_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">&#160;</span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark">&#160;</span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark">&#160;</span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark">&#160;</span>Files</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark">&#160;</span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark">&#160;</span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark">&#160;</span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark">&#160;</span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark">&#160;</span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark">&#160;</span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark">&#160;</span>Macros</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(11)"><span class="SelectionMark">&#160;</span>Pages</a></div>

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

<div class="header">
  <div class="headertitle">
<div class="title">iterative.cpp</div>  </div>
</div><!--header-->
<div class="contents">
<p>This tutorial explains the use of iterative solvers in ViennaCL.</p>
<dl class="section note"><dt>Note</dt><dd><a class="el" href="iterative_8cpp.html">iterative.cpp</a> and iterative.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 &lt;iostream&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="comment">// Necessary to obtain a suitable performance in ublas</span></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">//</span></div>
<div class="line"><span class="comment">// ublas includes</span></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/io.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/triangular.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/matrix_sparse.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/matrix.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/matrix_proxy.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/operation.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/operation_sparse.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/io.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/lu.hpp&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Must be set if you want to use ViennaCL algorithms on ublas objects</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS 1</span></div>
<div class="line"></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 &quot;<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="coordinate__matrix_8hpp.html">viennacl/coordinate_matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="ilu_8hpp.html">viennacl/linalg/ilu.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="jacobi__precond_8hpp.html">viennacl/linalg/jacobi_precond.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="cg_8hpp.html">viennacl/linalg/cg.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="bicgstab_8hpp.html">viennacl/linalg/bicgstab.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="gmres_8hpp.html">viennacl/linalg/gmres.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// Some helper functions for this tutorial:</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="vector-io_8hpp.html">vector-io.hpp</a>&quot;</span></div>
</div><!-- fragment --><p> A shortcut for convience (use <code>ublas::</code> instead of <code>boost::numeric::ublas::</code>) </p>
<div class="fragment"><div class="line"><span class="keyword">using namespace </span>boost::numeric;</div>
</div><!-- fragment --><p> In the <a class="el" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main()</a> routine we create matrices and vectors using Boost.uBLAS types and fill them with data. Then, ViennaCL types are created and the data is copied from the respective Boost.uBLAS objects. Various preconditioners are set up and finally the iterative solvers get called. </p>
<div class="fragment"><div class="line"><span class="keywordtype">int</span> <a name="a0"></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="a1"></a><a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>;</div>
</div><!-- fragment --><p> Set up some ublas objects </p>
<div class="fragment"><div class="line">ublas::vector&lt;ScalarType&gt; rhs;</div>
<div class="line">ublas::vector&lt;ScalarType&gt; ref_result;</div>
<div class="line">ublas::vector&lt;ScalarType&gt; result;</div>
<div class="line">ublas::compressed_matrix&lt;ScalarType&gt; ublas_matrix;</div>
</div><!-- fragment --><p> Read system matrix from a matrix-market file </p>
<div class="fragment"><div class="line"><span class="keywordflow">if</span> (!<a name="a2"></a><a class="code" href="namespaceviennacl_1_1io.html#ad934125ed3dbe661e264bcd7d62b1048">viennacl::io::read_matrix_market_file</a>(ublas_matrix, <span class="stringliteral">&quot;../examples/testdata/mat65k.mtx&quot;</span>))</div>
<div class="line">{</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading Matrix file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">}</div>
<div class="line"></div>
<div class="line"><span class="keywordflow">if</span> (!<a name="a3"></a><a class="code" href="vector-io_8hpp.html#a303954a45833b84192a6dcf3065f85e3">readVectorFromFile</a>(<span class="stringliteral">&quot;../examples/testdata/rhs65025.txt&quot;</span>, rhs))</div>
<div class="line">{</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading RHS file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">}</div>
<div class="line"></div>
<div class="line"><span class="keywordflow">if</span> (!<a class="code" href="vector-io_8hpp.html#a303954a45833b84192a6dcf3065f85e3">readVectorFromFile</a>(<span class="stringliteral">&quot;../examples/testdata/result65025.txt&quot;</span>, ref_result))</div>
<div class="line">{</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading Result file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">}</div>
</div><!-- fragment --><p> Set up some ViennaCL objects </p>
<div class="fragment"><div class="line">std::size_t vcl_size = rhs.size();</div>
<div class="line"><a name="_a4"></a><a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;ScalarType&gt;</a> vcl_compressed_matrix;</div>
<div class="line"><a name="_a5"></a><a class="code" href="classviennacl_1_1coordinate__matrix.html">viennacl::coordinate_matrix&lt;ScalarType&gt;</a> vcl_coordinate_matrix;</div>
<div class="line"><a name="_a6"></a><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> vcl_rhs(vcl_size);</div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> vcl_result(vcl_size);</div>
<div class="line"><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> vcl_ref_result(vcl_size);</div>
<div class="line"></div>
<div class="line"><a name="a7"></a><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(rhs.begin(), rhs.end(), vcl_rhs.begin());</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ref_result.begin(), ref_result.end(), vcl_ref_result.begin());</div>
</div><!-- fragment --><p> Transfer ublas-matrix to GPU: </p>
<div class="fragment"><div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_matrix, vcl_compressed_matrix);</div>
</div><!-- fragment --><p> An alternative way is to use the C++ STL. In such case, the sparse matrix is represented as a std::vector&lt; std::map&lt; unsigned int, ScalarType&gt; &gt; After the instantiation, the data from the uBLAS matrix is copied over to the STL-matrix, which is then copied over to a ViennaCL matrix. This step is just for demonstration purposes. </p>
<div class="fragment"><div class="line">std::vector&lt; std::map&lt; unsigned int, ScalarType&gt; &gt; stl_matrix(rhs.size());</div>
<div class="line"><span class="keywordflow">for</span> (ublas::compressed_matrix&lt;ScalarType&gt;::iterator1 iter1 = ublas_matrix.begin1();</div>
<div class="line">     iter1 != ublas_matrix.end1();</div>
<div class="line">     ++iter1)</div>
<div class="line">{</div>
<div class="line">  <span class="keywordflow">for</span> (ublas::compressed_matrix&lt;ScalarType&gt;::iterator2 iter2 = iter1.begin();</div>
<div class="line">       iter2 != iter1.end();</div>
<div class="line">       ++iter2)</div>
<div class="line">       stl_matrix[iter2.index1()][<span class="keyword">static_cast&lt;</span><span class="keywordtype">unsigned</span> <span class="keywordtype">int</span><span class="keyword">&gt;</span>(iter2.index2())] = *iter2;</div>
<div class="line">}</div>
<div class="line">std::vector&lt;ScalarType&gt; stl_rhs(rhs.size()), stl_result(result.size());</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">std::copy</a>(rhs.begin(), rhs.end(), stl_rhs.begin());</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(stl_matrix, vcl_coordinate_matrix);</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_coordinate_matrix, stl_matrix);</div>
</div><!-- fragment --><p> Set up ILUT preconditioners for ublas and ViennaCL objects: </p>
<div class="fragment"><div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;Setting up preconditioners for uBLAS-matrix...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="_a8"></a><a class="code" href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a> &gt;    ublas_ilut(ublas_matrix, <a name="_a9"></a><a class="code" href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a>());</div>
<div class="line"><a name="_a10"></a><a class="code" href="classviennacl_1_1linalg_1_1ilu0__precond.html">viennacl::linalg::ilu0_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a> &gt;    ublas_ilu0(ublas_matrix, <a name="_a11"></a><a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
<div class="line"><a name="_a12"></a><a class="code" href="classviennacl_1_1linalg_1_1block__ilu__precond.html">viennacl::linalg::block_ilu_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a>,</div>
<div class="line">                                     <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>&gt;          ublas_block_ilu0(ublas_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
<div class="line"></div>
<div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;Setting up preconditioners for ViennaCL-matrix...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a class="code" href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a> &gt; vcl_ilut(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a>());</div>
<div class="line"><a class="code" href="classviennacl_1_1linalg_1_1ilu0__precond.html">viennacl::linalg::ilu0_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a> &gt; vcl_ilu0(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
<div class="line"><a class="code" href="classviennacl_1_1linalg_1_1block__ilu__precond.html">viennacl::linalg::block_ilu_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a>,</div>
<div class="line">                                     <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>&gt;          vcl_block_ilu0(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
</div><!-- fragment --><p> set up Jacobi preconditioners for ViennaCL and ublas objects: </p>
<div class="fragment"><div class="line"><a name="_a13"></a><a class="code" href="classviennacl_1_1linalg_1_1jacobi__precond.html">viennacl::linalg::jacobi_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a> &gt;    ublas_jacobi(ublas_matrix, <a name="_a14"></a><a class="code" href="classviennacl_1_1linalg_1_1jacobi__tag.html">viennacl::linalg::jacobi_tag</a>());</div>
<div class="line"><a class="code" href="classviennacl_1_1linalg_1_1jacobi__precond.html">viennacl::linalg::jacobi_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a> &gt; vcl_jacobi(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1jacobi__tag.html">viennacl::linalg::jacobi_tag</a>());</div>
</div><!-- fragment --> <h2>Conjugate Gradient Solver</h2>
<div class="fragment"><div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;----- CG Method -----&quot;</span> &lt;&lt; std::endl;</div>
</div><!-- fragment --><p> Run the CG method with uBLAS objects. Use a relative tolerance of <img class="formulaInl" alt="$ 10^{-6} $" src="form_111.png"/> and a maximum of 20 iterations when using ILUT or Jacobi preconditioners. You might need higher iteration counts for poorly conditioned systems. </p>
<div class="fragment"><div class="line">result = <a name="a15"></a><a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a name="_a16"></a><a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
<div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), ublas_ilut);</div>
<div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), ublas_jacobi);</div>
</div><!-- fragment --><p> Run the CG method for ViennaCL objects. The results here should be the same as with uBLAS objects (at least up to round-off error). </p>
<div class="fragment"><div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
<div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_jacobi);</div>
</div><!-- fragment --><p> Convenience option: Run the CG method by passing STL types. This will use appropriate ViennaCL objects internally. You need to include <a class="el" href="compressed__matrix_8hpp.html" title="Implementation of the compressed_matrix class. ">viennacl/compressed_matrix.hpp</a> and <a class="el" href="vector_8hpp.html" title="The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...">viennacl/vector.hpp</a> before <a class="el" href="cg_8hpp.html" title="The conjugate gradient method is implemented here. ">viennacl/linalg/cg.hpp</a> for this to work! </p>
<div class="fragment"><div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
<div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_jacobi);</div>
</div><!-- fragment --> <h2>Stabilized BiConjugate Gradient Solver</h2>
<div class="fragment"><div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;----- BiCGStab Method -----&quot;</span> &lt;&lt; std::endl;</div>
</div><!-- fragment --><p> Run BiCGStab for Boost.uBLAS objects. Parameters are the same as for the CG method. </p>
<div class="fragment"><div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a name="_a17"></a><a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>());          <span class="comment">//without preconditioner</span></div>
<div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), ublas_ilut); <span class="comment">//with preconditioner</span></div>
<div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), ublas_jacobi); <span class="comment">//with preconditioner</span></div>
</div><!-- fragment --><p> Run the BiCGStab method for ViennaCL objects. The results here should be the same as with uBLAS objects (at least up to round-off error). </p>
<div class="fragment"><div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>());   <span class="comment">//without preconditioner</span></div>
<div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_ilut); <span class="comment">//with preconditioner</span></div>
<div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_jacobi); <span class="comment">//with preconditioner</span></div>
</div><!-- fragment --><p> Convenience option: Run the BiCGStab method by passing STL types. This will use appropriate ViennaCL objects internally. You need to include <a class="el" href="compressed__matrix_8hpp.html" title="Implementation of the compressed_matrix class. ">viennacl/compressed_matrix.hpp</a> and <a class="el" href="vector_8hpp.html" title="The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...">viennacl/vector.hpp</a> before <a class="el" href="bicgstab_8hpp.html" title="The stabilized bi-conjugate gradient method is implemented here. ">viennacl/linalg/bicgstab.hpp</a> for this to work! </p>
<div class="fragment"><div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>());</div>
<div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_jacobi);</div>
</div><!-- fragment --> <h2>GMRES Solver</h2>
<div class="fragment"><div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;----- GMRES Method -----&quot;</span> &lt;&lt; std::endl;</div>
</div><!-- fragment --><p> Run GMRES for Boost.uBLAS objects. Parameters are the same as for the CG method. </p>
<div class="fragment"><div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a name="_a18"></a><a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>());   <span class="comment">//without preconditioner</span></div>
<div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), ublas_ilut);<span class="comment">//with preconditioner</span></div>
<div class="line">result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), ublas_jacobi);<span class="comment">//with preconditioner</span></div>
</div><!-- fragment --><p> Run the GMRES method for ViennaCL objects. The results here should be the same as with uBLAS objects (at least up to round-off error). </p>
<div class="fragment"><div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>());   <span class="comment">//without preconditioner</span></div>
<div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_ilut);<span class="comment">//with preconditioner</span></div>
<div class="line">vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_jacobi);<span class="comment">//with preconditioner</span></div>
</div><!-- fragment --><p> Convenience option: Run the GMRES method by passing STL types. This will use appropriate ViennaCL objects internally. You need to include <a class="el" href="compressed__matrix_8hpp.html" title="Implementation of the compressed_matrix class. ">viennacl/compressed_matrix.hpp</a> and <a class="el" href="vector_8hpp.html" title="The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...">viennacl/vector.hpp</a> before <a class="el" href="gmres_8hpp.html" title="Implementations of the generalized minimum residual method are in this file. ">viennacl/linalg/gmres.hpp</a> for this to work! </p>
<div class="fragment"><div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>());</div>
<div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_jacobi);</div>
</div><!-- fragment --><p> That's it, the tutorial is completed. </p>
<div class="fragment"><div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!&quot;</span> &lt;&lt; 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">//</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 &lt;iostream&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="comment">// Necessary to obtain a suitable performance in ublas</span></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">//</span></div>
<div class="line"><span class="comment">// ublas includes</span></div>
<div class="line"><span class="comment">//</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/io.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/triangular.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/matrix_sparse.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/matrix.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/matrix_proxy.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/operation.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/operation_sparse.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/io.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/lu.hpp&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="comment">// Must be set if you want to use ViennaCL algorithms on ublas objects</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS 1</span></div>
<div class="line"></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 &quot;<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="coordinate__matrix_8hpp.html">viennacl/coordinate_matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="ilu_8hpp.html">viennacl/linalg/ilu.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="jacobi__precond_8hpp.html">viennacl/linalg/jacobi_precond.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="cg_8hpp.html">viennacl/linalg/cg.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="bicgstab_8hpp.html">viennacl/linalg/bicgstab.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="gmres_8hpp.html">viennacl/linalg/gmres.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__market_8hpp.html">viennacl/io/matrix_market.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// Some helper functions for this tutorial:</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="vector-io_8hpp.html">vector-io.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>boost::numeric;</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>;</div>
<div class="line"></div>
<div class="line">  ublas::vector&lt;ScalarType&gt; rhs;</div>
<div class="line">  ublas::vector&lt;ScalarType&gt; ref_result;</div>
<div class="line">  ublas::vector&lt;ScalarType&gt; result;</div>
<div class="line">  ublas::compressed_matrix&lt;ScalarType&gt; ublas_matrix;</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>(ublas_matrix, <span class="stringliteral">&quot;../examples/testdata/mat65k.mtx&quot;</span>))</div>
<div class="line">  {</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading Matrix file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">    <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">  }</div>
<div class="line"></div>
<div class="line">  <span class="keywordflow">if</span> (!<a class="code" href="vector-io_8hpp.html#a303954a45833b84192a6dcf3065f85e3">readVectorFromFile</a>(<span class="stringliteral">&quot;../examples/testdata/rhs65025.txt&quot;</span>, rhs))</div>
<div class="line">  {</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading RHS file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">    <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">  }</div>
<div class="line"></div>
<div class="line">  <span class="keywordflow">if</span> (!<a class="code" href="vector-io_8hpp.html#a303954a45833b84192a6dcf3065f85e3">readVectorFromFile</a>(<span class="stringliteral">&quot;../examples/testdata/result65025.txt&quot;</span>, ref_result))</div>
<div class="line">  {</div>
<div class="line">    std::cout &lt;&lt; <span class="stringliteral">&quot;Error reading Result file&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">    <span class="keywordflow">return</span> EXIT_FAILURE;</div>
<div class="line">  }</div>
<div class="line"></div>
<div class="line">  std::size_t vcl_size = rhs.size();</div>
<div class="line">  <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;ScalarType&gt;</a> vcl_compressed_matrix;</div>
<div class="line">  <a class="code" href="classviennacl_1_1coordinate__matrix.html">viennacl::coordinate_matrix&lt;ScalarType&gt;</a> vcl_coordinate_matrix;</div>
<div class="line">  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> vcl_rhs(vcl_size);</div>
<div class="line">  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> vcl_result(vcl_size);</div>
<div class="line">  <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a> vcl_ref_result(vcl_size);</div>
<div class="line"></div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(rhs.begin(), rhs.end(), vcl_rhs.<a name="a19"></a><a class="code" href="classviennacl_1_1vector__base.html#a762d98e2f912fc534951f25555b6077f">begin</a>());</div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ref_result.begin(), ref_result.end(), vcl_ref_result.<a class="code" href="classviennacl_1_1vector__base.html#a762d98e2f912fc534951f25555b6077f">begin</a>());</div>
<div class="line"></div>
<div class="line"></div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_matrix, vcl_compressed_matrix);</div>
<div class="line"></div>
<div class="line">  std::vector&lt; std::map&lt; unsigned int, ScalarType&gt; &gt; stl_matrix(rhs.size());</div>
<div class="line">  <span class="keywordflow">for</span> (ublas::compressed_matrix&lt;ScalarType&gt;::iterator1 iter1 = ublas_matrix.begin1();</div>
<div class="line">       iter1 != ublas_matrix.end1();</div>
<div class="line">       ++iter1)</div>
<div class="line">  {</div>
<div class="line">    <span class="keywordflow">for</span> (ublas::compressed_matrix&lt;ScalarType&gt;::iterator2 iter2 = iter1.begin();</div>
<div class="line">         iter2 != iter1.end();</div>
<div class="line">         ++iter2)</div>
<div class="line">         stl_matrix[iter2.index1()][<span class="keyword">static_cast&lt;</span><span class="keywordtype">unsigned</span> <span class="keywordtype">int</span><span class="keyword">&gt;</span>(iter2.index2())] = *iter2;</div>
<div class="line">  }</div>
<div class="line">  std::vector&lt;ScalarType&gt; stl_rhs(rhs.size()), stl_result(result.size());</div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">std::copy</a>(rhs.begin(), rhs.end(), stl_rhs.begin());</div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(stl_matrix, vcl_coordinate_matrix);</div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(vcl_coordinate_matrix, stl_matrix);</div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Setting up preconditioners for uBLAS-matrix...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a> &gt;    ublas_ilut(ublas_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a>());</div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1ilu0__precond.html">viennacl::linalg::ilu0_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a> &gt;    ublas_ilu0(ublas_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1block__ilu__precond.html">viennacl::linalg::block_ilu_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a>,</div>
<div class="line">                                       <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>&gt;          ublas_block_ilu0(ublas_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Setting up preconditioners for ViennaCL-matrix...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a> &gt; vcl_ilut(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a>());</div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1ilu0__precond.html">viennacl::linalg::ilu0_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a> &gt; vcl_ilu0(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1block__ilu__precond.html">viennacl::linalg::block_ilu_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a>,</div>
<div class="line">                                       <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>&gt;          vcl_block_ilu0(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1ilu0__tag.html">viennacl::linalg::ilu0_tag</a>());</div>
<div class="line"></div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1jacobi__precond.html">viennacl::linalg::jacobi_precond&lt; ublas::compressed_matrix&lt;ScalarType&gt;</a> &gt;    ublas_jacobi(ublas_matrix, <a class="code" href="classviennacl_1_1linalg_1_1jacobi__tag.html">viennacl::linalg::jacobi_tag</a>());</div>
<div class="line">  <a class="code" href="classviennacl_1_1linalg_1_1jacobi__precond.html">viennacl::linalg::jacobi_precond&lt; viennacl::compressed_matrix&lt;ScalarType&gt;</a> &gt; vcl_jacobi(vcl_compressed_matrix, <a class="code" href="classviennacl_1_1linalg_1_1jacobi__tag.html">viennacl::linalg::jacobi_tag</a>());</div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;----- CG Method -----&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"></div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), ublas_ilut);</div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), ublas_jacobi);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_jacobi);</div>
<div class="line"></div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>());</div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1cg__tag.html">viennacl::linalg::cg_tag</a>(1e-6, 20), vcl_jacobi);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;----- BiCGStab Method -----&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"></div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>());          <span class="comment">//without preconditioner</span></div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), ublas_ilut); <span class="comment">//with preconditioner</span></div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), ublas_jacobi); <span class="comment">//with preconditioner</span></div>
<div class="line"></div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>());   <span class="comment">//without preconditioner</span></div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_ilut); <span class="comment">//with preconditioner</span></div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_jacobi); <span class="comment">//with preconditioner</span></div>
<div class="line"></div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>());</div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a>(1e-6, 20), vcl_jacobi);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;----- GMRES Method -----&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"></div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>());   <span class="comment">//without preconditioner</span></div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), ublas_ilut);<span class="comment">//with preconditioner</span></div>
<div class="line">  result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(ublas_matrix, rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), ublas_jacobi);<span class="comment">//with preconditioner</span></div>
<div class="line"></div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>());   <span class="comment">//without preconditioner</span></div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_ilut);<span class="comment">//with preconditioner</span></div>
<div class="line">  vcl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(vcl_compressed_matrix, vcl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_jacobi);<span class="comment">//with preconditioner</span></div>
<div class="line"></div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>());</div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_ilut);</div>
<div class="line">  stl_result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(stl_matrix, stl_rhs, <a class="code" href="classviennacl_1_1linalg_1_1gmres__tag.html">viennacl::linalg::gmres_tag</a>(1e-6, 20), vcl_jacobi);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!&quot;</span> &lt;&lt; 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>