File: spai_8cpp_source.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 (308 lines) | stat: -rw-r--r-- 48,838 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
<!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: examples/tutorial/spai.cpp Source File</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('spai_8cpp_source.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">spai.cpp</div>  </div>
</div><!--header-->
<div class="contents">
<a href="spai_8cpp.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno">    1</span>&#160;<span class="comment">/* =========================================================================</span></div>
<div class="line"><a name="l00002"></a><span class="lineno">    2</span>&#160;<span class="comment">   Copyright (c) 2010-2016, Institute for Microelectronics,</span></div>
<div class="line"><a name="l00003"></a><span class="lineno">    3</span>&#160;<span class="comment">                            Institute for Analysis and Scientific Computing,</span></div>
<div class="line"><a name="l00004"></a><span class="lineno">    4</span>&#160;<span class="comment">                            TU Wien.</span></div>
<div class="line"><a name="l00005"></a><span class="lineno">    5</span>&#160;<span class="comment">   Portions of this software are copyright by UChicago Argonne, LLC.</span></div>
<div class="line"><a name="l00006"></a><span class="lineno">    6</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00007"></a><span class="lineno">    7</span>&#160;<span class="comment">                            -----------------</span></div>
<div class="line"><a name="l00008"></a><span class="lineno">    8</span>&#160;<span class="comment">                  ViennaCL - The Vienna Computing Library</span></div>
<div class="line"><a name="l00009"></a><span class="lineno">    9</span>&#160;<span class="comment">                            -----------------</span></div>
<div class="line"><a name="l00010"></a><span class="lineno">   10</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00011"></a><span class="lineno">   11</span>&#160;<span class="comment">   Project Head:    Karl Rupp                   rupp@iue.tuwien.ac.at</span></div>
<div class="line"><a name="l00012"></a><span class="lineno">   12</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00013"></a><span class="lineno">   13</span>&#160;<span class="comment">   (A list of authors and contributors can be found in the PDF manual)</span></div>
<div class="line"><a name="l00014"></a><span class="lineno">   14</span>&#160;<span class="comment"></span></div>
<div class="line"><a name="l00015"></a><span class="lineno">   15</span>&#160;<span class="comment">   License:         MIT (X11), see file LICENSE in the base directory</span></div>
<div class="line"><a name="l00016"></a><span class="lineno">   16</span>&#160;<span class="comment">============================================================================= */</span></div>
<div class="line"><a name="l00017"></a><span class="lineno">   17</span>&#160;</div>
<div class="line"><a name="l00027"></a><span class="lineno">   27</span>&#160;<span class="comment">// enable Boost.uBLAS support</span></div>
<div class="line"><a name="l00028"></a><span class="lineno">   28</span>&#160;<span class="preprocessor">#define VIENNACL_WITH_UBLAS</span></div>
<div class="line"><a name="l00029"></a><span class="lineno">   29</span>&#160;</div>
<div class="line"><a name="l00030"></a><span class="lineno">   30</span>&#160;<span class="preprocessor">#ifndef NDEBUG</span></div>
<div class="line"><a name="l00031"></a><span class="lineno">   31</span>&#160;<span class="preprocessor"> #define BOOST_UBLAS_NDEBUG</span></div>
<div class="line"><a name="l00032"></a><span class="lineno">   32</span>&#160;<span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00033"></a><span class="lineno">   33</span>&#160;</div>
<div class="line"><a name="l00034"></a><span class="lineno">   34</span>&#160;<span class="comment">// System headers:</span></div>
<div class="line"><a name="l00035"></a><span class="lineno">   35</span>&#160;<span class="preprocessor">#include &lt;utility&gt;</span></div>
<div class="line"><a name="l00036"></a><span class="lineno">   36</span>&#160;<span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><a name="l00037"></a><span class="lineno">   37</span>&#160;<span class="preprocessor">#include &lt;fstream&gt;</span></div>
<div class="line"><a name="l00038"></a><span class="lineno">   38</span>&#160;<span class="preprocessor">#include &lt;string&gt;</span></div>
<div class="line"><a name="l00039"></a><span class="lineno">   39</span>&#160;<span class="preprocessor">#include &lt;cmath&gt;</span></div>
<div class="line"><a name="l00040"></a><span class="lineno">   40</span>&#160;<span class="preprocessor">#include &lt;algorithm&gt;</span></div>
<div class="line"><a name="l00041"></a><span class="lineno">   41</span>&#160;<span class="preprocessor">#include &lt;stdio.h&gt;</span></div>
<div class="line"><a name="l00042"></a><span class="lineno">   42</span>&#160;</div>
<div class="line"><a name="l00043"></a><span class="lineno">   43</span>&#160;<span class="comment">// ViennaCL headers:</span></div>
<div class="line"><a name="l00044"></a><span class="lineno">   44</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="scalar_8hpp.html">viennacl/scalar.hpp</a>&quot;</span></div>
<div class="line"><a name="l00045"></a><span class="lineno">   45</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="matrix_8hpp.html">viennacl/matrix.hpp</a>&quot;</span></div>
<div class="line"><a name="l00046"></a><span class="lineno">   46</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="compressed__matrix_8hpp.html">viennacl/compressed_matrix.hpp</a>&quot;</span></div>
<div class="line"><a name="l00047"></a><span class="lineno">   47</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="cg_8hpp.html">viennacl/linalg/cg.hpp</a>&quot;</span></div>
<div class="line"><a name="l00048"></a><span class="lineno">   48</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="bicgstab_8hpp.html">viennacl/linalg/bicgstab.hpp</a>&quot;</span></div>
<div class="line"><a name="l00049"></a><span class="lineno">   49</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="prod_8hpp.html">viennacl/linalg/prod.hpp</a>&quot;</span></div>
<div class="line"><a name="l00050"></a><span class="lineno">   50</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="inner__prod_8hpp.html">viennacl/linalg/inner_prod.hpp</a>&quot;</span></div>
<div class="line"><a name="l00051"></a><span class="lineno">   51</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="ilu_8hpp.html">viennacl/linalg/ilu.hpp</a>&quot;</span></div>
<div class="line"><a name="l00052"></a><span class="lineno">   52</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="norm__2_8hpp.html">viennacl/linalg/norm_2.hpp</a>&quot;</span></div>
<div class="line"><a name="l00053"></a><span class="lineno">   53</span>&#160;<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"><a name="l00054"></a><span class="lineno">   54</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="spai_8hpp.html">viennacl/linalg/spai.hpp</a>&quot;</span></div>
<div class="line"><a name="l00055"></a><span class="lineno">   55</span>&#160;</div>
<div class="line"><a name="l00056"></a><span class="lineno">   56</span>&#160;<span class="comment">// Boost headers:</span></div>
<div class="line"><a name="l00057"></a><span class="lineno">   57</span>&#160;<span class="preprocessor">#include &quot;boost/numeric/ublas/vector.hpp&quot;</span></div>
<div class="line"><a name="l00058"></a><span class="lineno">   58</span>&#160;<span class="preprocessor">#include &quot;boost/numeric/ublas/matrix.hpp&quot;</span></div>
<div class="line"><a name="l00059"></a><span class="lineno">   59</span>&#160;<span class="preprocessor">#include &quot;boost/numeric/ublas/io.hpp&quot;</span></div>
<div class="line"><a name="l00060"></a><span class="lineno">   60</span>&#160;</div>
<div class="line"><a name="l00061"></a><span class="lineno">   61</span>&#160;<span class="comment">// auxiliary functionality:</span></div>
<div class="line"><a name="l00062"></a><span class="lineno">   62</span>&#160;<span class="preprocessor">#include &quot;<a class="code" href="vector-io_8hpp.html">vector-io.hpp</a>&quot;</span></div>
<div class="line"><a name="l00063"></a><span class="lineno">   63</span>&#160;</div>
<div class="line"><a name="l00067"></a><span class="lineno">   67</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixT, <span class="keyword">typename</span> VectorT, <span class="keyword">typename</span> SolverTagT, <span class="keyword">typename</span> PreconditionerT&gt;</div>
<div class="line"><a name="l00068"></a><span class="lineno">   68</span>&#160;<span class="keywordtype">void</span> <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(MatrixT <span class="keyword">const</span> &amp; A, VectorT <span class="keyword">const</span> &amp; b, SolverTagT <span class="keyword">const</span> &amp; solver_tag, PreconditionerT <span class="keyword">const</span> &amp; precond)</div>
<div class="line"><a name="l00069"></a><span class="lineno">   69</span>&#160;{</div>
<div class="line"><a name="l00070"></a><span class="lineno">   70</span>&#160;  VectorT result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(A, b, solver_tag, precond);</div>
<div class="line"><a name="l00071"></a><span class="lineno">   71</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Solver iterations: &quot;</span> &lt;&lt; solver_tag.iters() &lt;&lt; std::endl;</div>
<div class="line"><a name="l00072"></a><span class="lineno">   72</span>&#160;  VectorT residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(A, result);</div>
<div class="line"><a name="l00073"></a><span class="lineno">   73</span>&#160;  residual -= b;</div>
<div class="line"><a name="l00074"></a><span class="lineno">   74</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Rel. Residual: &quot;</span> &lt;&lt; <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) &lt;&lt; std::endl;</div>
<div class="line"><a name="l00075"></a><span class="lineno">   75</span>&#160;}</div>
<div class="line"><a name="l00076"></a><span class="lineno">   76</span>&#160;</div>
<div class="line"><a name="l00087"></a><span class="lineno">   87</span>&#160;<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"><a name="l00088"></a><span class="lineno">   88</span>&#160;{</div>
<div class="line"><a name="l00089"></a><span class="lineno">   89</span>&#160;  <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"><a name="l00090"></a><span class="lineno">   90</span>&#160;  <span class="keyword">typedef</span> boost::numeric::ublas::compressed_matrix&lt;ScalarType&gt;        MatrixType;</div>
<div class="line"><a name="l00091"></a><span class="lineno">   91</span>&#160;  <span class="keyword">typedef</span> boost::numeric::ublas::vector&lt;ScalarType&gt;                   VectorType;</div>
<div class="line"><a name="l00092"></a><span class="lineno">   92</span>&#160;  <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix&lt;ScalarType&gt;</a>                     GPUMatrixType;</div>
<div class="line"><a name="l00093"></a><span class="lineno">   93</span>&#160;  <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a>                                GPUVectorType;</div>
<div class="line"><a name="l00094"></a><span class="lineno">   94</span>&#160;</div>
<div class="line"><a name="l00098"></a><span class="lineno">   98</span>&#160;<span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00099"></a><span class="lineno">   99</span>&#160;  <span class="comment">// Optional: Customize OpenCL backend</span></div>
<div class="line"><a name="l00100"></a><span class="lineno">  100</span>&#160;  <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"><a name="l00101"></a><span class="lineno">  101</span>&#160;  std::vector&lt;viennacl::ocl::device&gt; <span class="keyword">const</span> &amp; devices = pf.<a class="code" href="classviennacl_1_1ocl_1_1platform.html#afaa563522ebe9ce7b80ef02c40e7fe31">devices</a>();</div>
<div class="line"><a name="l00102"></a><span class="lineno">  102</span>&#160;</div>
<div class="line"><a name="l00103"></a><span class="lineno">  103</span>&#160;  <span class="comment">// Optional: Set first device to first context:</span></div>
<div class="line"><a name="l00104"></a><span class="lineno">  104</span>&#160;  <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(0, devices[0]);</div>
<div class="line"><a name="l00105"></a><span class="lineno">  105</span>&#160;</div>
<div class="line"><a name="l00106"></a><span class="lineno">  106</span>&#160;  <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"><a name="l00107"></a><span class="lineno">  107</span>&#160;  <span class="keywordflow">if</span> (devices.size() &gt; 1)</div>
<div class="line"><a name="l00108"></a><span class="lineno">  108</span>&#160;    <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(1, devices[1]);</div>
<div class="line"><a name="l00109"></a><span class="lineno">  109</span>&#160;  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00110"></a><span class="lineno">  110</span>&#160;    <a class="code" href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a>(1, devices[0]);</div>
<div class="line"><a name="l00111"></a><span class="lineno">  111</span>&#160;</div>
<div class="line"><a name="l00112"></a><span class="lineno">  112</span>&#160;  std::cout &lt;&lt; <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>() &lt;&lt; std::endl;</div>
<div class="line"><a name="l00113"></a><span class="lineno">  113</span>&#160;  <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"><a name="l00114"></a><span class="lineno">  114</span>&#160;<span class="preprocessor">#else</span></div>
<div class="line"><a name="l00115"></a><span class="lineno">  115</span>&#160;  <a class="code" href="classviennacl_1_1context.html">viennacl::context</a> ctx;</div>
<div class="line"><a name="l00116"></a><span class="lineno">  116</span>&#160;<span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00117"></a><span class="lineno">  117</span>&#160;</div>
<div class="line"><a name="l00121"></a><span class="lineno">  121</span>&#160;  MatrixType M;</div>
<div class="line"><a name="l00122"></a><span class="lineno">  122</span>&#160;</div>
<div class="line"><a name="l00123"></a><span class="lineno">  123</span>&#160;  <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">&quot;../examples/testdata/mat65k.mtx&quot;</span>))</div>
<div class="line"><a name="l00124"></a><span class="lineno">  124</span>&#160;  {</div>
<div class="line"><a name="l00125"></a><span class="lineno">  125</span>&#160;    std::cerr&lt;&lt;<span class="stringliteral">&quot;ERROR: Could not read matrix file &quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00126"></a><span class="lineno">  126</span>&#160;    exit(EXIT_FAILURE);</div>
<div class="line"><a name="l00127"></a><span class="lineno">  127</span>&#160;  }</div>
<div class="line"><a name="l00128"></a><span class="lineno">  128</span>&#160;</div>
<div class="line"><a name="l00129"></a><span class="lineno">  129</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;Size of matrix: &quot;</span> &lt;&lt; M.size1() &lt;&lt; std::endl;</div>
<div class="line"><a name="l00130"></a><span class="lineno">  130</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;Avg. Entries per row: &quot;</span> &lt;&lt; double(M.nnz()) / static_cast&lt;double&gt;(M.size1()) &lt;&lt; std::endl;</div>
<div class="line"><a name="l00131"></a><span class="lineno">  131</span>&#160;</div>
<div class="line"><a name="l00135"></a><span class="lineno">  135</span>&#160;  VectorType rhs(M.size2());</div>
<div class="line"><a name="l00136"></a><span class="lineno">  136</span>&#160;  <span class="keywordflow">for</span> (std::size_t i=0; i&lt;rhs.size(); ++i)</div>
<div class="line"><a name="l00137"></a><span class="lineno">  137</span>&#160;    rhs(i) = <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>(1);</div>
<div class="line"><a name="l00138"></a><span class="lineno">  138</span>&#160;</div>
<div class="line"><a name="l00142"></a><span class="lineno">  142</span>&#160;  GPUMatrixType  gpu_M(M.size1(), M.size2(), ctx);</div>
<div class="line"><a name="l00143"></a><span class="lineno">  143</span>&#160;  GPUVectorType  gpu_rhs(M.size1(), ctx);</div>
<div class="line"><a name="l00144"></a><span class="lineno">  144</span>&#160;  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(M, gpu_M);</div>
<div class="line"><a name="l00145"></a><span class="lineno">  145</span>&#160;  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(rhs, gpu_rhs);</div>
<div class="line"><a name="l00146"></a><span class="lineno">  146</span>&#160;</div>
<div class="line"><a name="l00153"></a><span class="lineno">  153</span>&#160;  <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"><a name="l00154"></a><span class="lineno">  154</span>&#160;</div>
<div class="line"><a name="l00158"></a><span class="lineno">  158</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Reference 1: Pure BiCGStab on CPU ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00159"></a><span class="lineno">  159</span>&#160;  VectorType result = <a class="code" href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a>(M, rhs, solver_tag);</div>
<div class="line"><a name="l00160"></a><span class="lineno">  160</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Solver iterations: &quot;</span> &lt;&lt; solver_tag.iters() &lt;&lt; std::endl;</div>
<div class="line"><a name="l00161"></a><span class="lineno">  161</span>&#160;  VectorType residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(M, result) - rhs;</div>
<div class="line"><a name="l00162"></a><span class="lineno">  162</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Rel. Residual: &quot;</span> &lt;&lt; <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) &lt;&lt; std::endl;</div>
<div class="line"><a name="l00163"></a><span class="lineno">  163</span>&#160;</div>
<div class="line"><a name="l00164"></a><span class="lineno">  164</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Reference 2: Pure BiCGStab on GPU ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00165"></a><span class="lineno">  165</span>&#160;  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"><a name="l00166"></a><span class="lineno">  166</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Solver iterations: &quot;</span> &lt;&lt; solver_tag.iters() &lt;&lt; std::endl;</div>
<div class="line"><a name="l00167"></a><span class="lineno">  167</span>&#160;  GPUVectorType gpu_residual = <a class="code" href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a>(gpu_M, gpu_result);</div>
<div class="line"><a name="l00168"></a><span class="lineno">  168</span>&#160;  gpu_residual -= gpu_rhs;</div>
<div class="line"><a name="l00169"></a><span class="lineno">  169</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Rel. Residual: &quot;</span> &lt;&lt; <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) &lt;&lt; std::endl;</div>
<div class="line"><a name="l00170"></a><span class="lineno">  170</span>&#160;</div>
<div class="line"><a name="l00171"></a><span class="lineno">  171</span>&#160;</div>
<div class="line"><a name="l00175"></a><span class="lineno">  175</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Reference 2: BiCGStab with ILUT on CPU ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00176"></a><span class="lineno">  176</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Preconditioner setup...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00177"></a><span class="lineno">  177</span>&#160;  <a class="code" href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond&lt;MatrixType&gt;</a> ilut(M, <a class="code" href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a>());</div>
<div class="line"><a name="l00178"></a><span class="lineno">  178</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Iterative solver run...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00179"></a><span class="lineno">  179</span>&#160;  <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, ilut);</div>
<div class="line"><a name="l00180"></a><span class="lineno">  180</span>&#160;</div>
<div class="line"><a name="l00181"></a><span class="lineno">  181</span>&#160;</div>
<div class="line"><a name="l00185"></a><span class="lineno">  185</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Test 1: CPU-based SPAI ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00186"></a><span class="lineno">  186</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Preconditioner setup...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00187"></a><span class="lineno">  187</span>&#160;  <a class="code" href="classviennacl_1_1linalg_1_1spai__precond.html">viennacl::linalg::spai_precond&lt;MatrixType&gt;</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"><a name="l00188"></a><span class="lineno">  188</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Iterative solver run...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00189"></a><span class="lineno">  189</span>&#160;  <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, spai_cpu);</div>
<div class="line"><a name="l00190"></a><span class="lineno">  190</span>&#160;</div>
<div class="line"><a name="l00194"></a><span class="lineno">  194</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Test 2: CPU-based FSPAI ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00195"></a><span class="lineno">  195</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Preconditioner setup...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00196"></a><span class="lineno">  196</span>&#160;  <a class="code" href="classviennacl_1_1linalg_1_1fspai__precond.html">viennacl::linalg::fspai_precond&lt;MatrixType&gt;</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"><a name="l00197"></a><span class="lineno">  197</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Iterative solver run...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00198"></a><span class="lineno">  198</span>&#160;  <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(M, rhs, solver_tag, fspai_cpu);</div>
<div class="line"><a name="l00199"></a><span class="lineno">  199</span>&#160;</div>
<div class="line"><a name="l00203"></a><span class="lineno">  203</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Test 3: GPU-based SPAI ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00204"></a><span class="lineno">  204</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Preconditioner setup...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00205"></a><span class="lineno">  205</span>&#160;  <a class="code" href="classviennacl_1_1linalg_1_1spai__precond.html">viennacl::linalg::spai_precond&lt;GPUMatrixType&gt;</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"><a name="l00206"></a><span class="lineno">  206</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Iterative solver run...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00207"></a><span class="lineno">  207</span>&#160;  <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(gpu_M, gpu_rhs, solver_tag, spai_gpu);</div>
<div class="line"><a name="l00208"></a><span class="lineno">  208</span>&#160;</div>
<div class="line"><a name="l00212"></a><span class="lineno">  212</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Test 4: GPU-based FSPAI ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00213"></a><span class="lineno">  213</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Preconditioner setup...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00214"></a><span class="lineno">  214</span>&#160;  <a class="code" href="classviennacl_1_1linalg_1_1fspai__precond.html">viennacl::linalg::fspai_precond&lt;GPUMatrixType&gt;</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"><a name="l00215"></a><span class="lineno">  215</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot; * Iterative solver run...&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00216"></a><span class="lineno">  216</span>&#160;  <a class="code" href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a>(gpu_M, gpu_rhs, solver_tag, fspai_gpu);</div>
<div class="line"><a name="l00217"></a><span class="lineno">  217</span>&#160;</div>
<div class="line"><a name="l00221"></a><span class="lineno">  221</span>&#160;  std::cout &lt;&lt; <span class="stringliteral">&quot;!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"><a name="l00222"></a><span class="lineno">  222</span>&#160;</div>
<div class="line"><a name="l00223"></a><span class="lineno">  223</span>&#160;  <span class="keywordflow">return</span> EXIT_SUCCESS;</div>
<div class="line"><a name="l00224"></a><span class="lineno">  224</span>&#160;}</div>
<div class="line"><a name="l00225"></a><span class="lineno">  225</span>&#160;</div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_ae46f15d01c01f92a153b3f555a15096b"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#ae46f15d01c01f92a153b3f555a15096b">viennacl::linalg::norm_2</a></div><div class="ttdeci">T norm_2(std::vector&lt; T, A &gt; const &amp;v1)</div><div class="ttdef"><b>Definition:</b> <a href="norm__2_8hpp_source.html#l00096">norm_2.hpp:96</a></div></div>
<div class="ttc" id="matrix__market_8hpp_html"><div class="ttname"><a href="matrix__market_8hpp.html">matrix_market.hpp</a></div><div class="ttdoc">A reader and writer for the matrix market format is implemented here. </div></div>
<div class="ttc" id="namespaceviennacl_1_1ocl_html_ae1e931f6efd155240b33ca440dfcfef5"><div class="ttname"><a href="namespaceviennacl_1_1ocl.html#ae1e931f6efd155240b33ca440dfcfef5">viennacl::ocl::get_platforms</a></div><div class="ttdeci">std::vector&lt; platform &gt; get_platforms()</div><div class="ttdef"><b>Definition:</b> <a href="platform_8hpp_source.html#l00124">platform.hpp:124</a></div></div>
<div class="ttc" id="norm__2_8hpp_html"><div class="ttname"><a href="norm__2_8hpp.html">norm_2.hpp</a></div><div class="ttdoc">Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...</div></div>
<div class="ttc" id="classviennacl_1_1ocl_1_1platform_html"><div class="ttname"><a href="classviennacl_1_1ocl_1_1platform.html">viennacl::ocl::platform</a></div><div class="ttdoc">Wrapper class for an OpenCL platform. </div><div class="ttdef"><b>Definition:</b> <a href="platform_8hpp_source.html#l00045">platform.hpp:45</a></div></div>
<div class="ttc" id="prod_8hpp_html"><div class="ttname"><a href="prod_8hpp.html">prod.hpp</a></div><div class="ttdoc">Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...</div></div>
<div class="ttc" id="matrix_8hpp_html"><div class="ttname"><a href="matrix_8hpp.html">matrix.hpp</a></div><div class="ttdoc">Implementation of the dense matrix class. </div></div>
<div class="ttc" id="tests_2src_2bisect_8cpp_html_ae66f6b31b5ad750f1fe042a706a4e3d4"><div class="ttname"><a href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a></div><div class="ttdeci">int main()</div><div class="ttdef"><b>Definition:</b> <a href="tests_2src_2bisect_8cpp_source.html#l00091">bisect.cpp:91</a></div></div>
<div class="ttc" id="classviennacl_1_1ocl_1_1platform_html_afaa563522ebe9ce7b80ef02c40e7fe31"><div class="ttname"><a href="classviennacl_1_1ocl_1_1platform.html#afaa563522ebe9ce7b80ef02c40e7fe31">viennacl::ocl::platform::devices</a></div><div class="ttdeci">std::vector&lt; device &gt; devices(cl_device_type dtype=CL_DEVICE_TYPE_DEFAULT)</div><div class="ttdoc">Returns the available devices of the supplied device type. </div><div class="ttdef"><b>Definition:</b> <a href="platform_8hpp_source.html#l00091">platform.hpp:91</a></div></div>
<div class="ttc" id="bicgstab_8hpp_html"><div class="ttname"><a href="bicgstab_8hpp.html">bicgstab.hpp</a></div><div class="ttdoc">The stabilized bi-conjugate gradient method is implemented here. </div></div>
<div class="ttc" id="vector-io_8hpp_html"><div class="ttname"><a href="vector-io_8hpp.html">vector-io.hpp</a></div></div>
<div class="ttc" id="classviennacl_1_1linalg_1_1detail_1_1spai_1_1spai__tag_html"><div class="ttname"><a href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1spai__tag.html">viennacl::linalg::detail::spai::spai_tag</a></div><div class="ttdoc">A tag for SPAI. </div><div class="ttdef"><b>Definition:</b> <a href="spai__tag_8hpp_source.html#l00064">spai_tag.hpp:64</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1ocl_html_ac54d59a74deaccec81f64e738b856348"><div class="ttname"><a href="namespaceviennacl_1_1ocl.html#ac54d59a74deaccec81f64e738b856348">viennacl::ocl::current_device</a></div><div class="ttdeci">viennacl::ocl::device const &amp; current_device()</div><div class="ttdoc">Convenience function for returning the active device in the current context. </div><div class="ttdef"><b>Definition:</b> <a href="backend_8hpp_source.html#l00351">backend.hpp:351</a></div></div>
<div class="ttc" id="inner__prod_8hpp_html"><div class="ttname"><a href="inner__prod_8hpp.html">inner_prod.hpp</a></div><div class="ttdoc">Generic interface for the computation of inner products. See viennacl/linalg/vector_operations.hpp for implementations. </div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a6e9b329b64ac782e6a5687ad2fc47a2a"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a6e9b329b64ac782e6a5687ad2fc47a2a">viennacl::linalg::solve</a></div><div class="ttdeci">VectorT solve(MatrixT const &amp;matrix, VectorT const &amp;rhs, bicgstab_tag const &amp;tag, PreconditionerT const &amp;precond)</div><div class="ttdef"><b>Definition:</b> <a href="bicgstab_8hpp_source.html#l00496">bicgstab.hpp:496</a></div></div>
<div class="ttc" id="classviennacl_1_1ocl_1_1device_html_a26a424782f982725453808a94763bd19"><div class="ttname"><a href="classviennacl_1_1ocl_1_1device.html#a26a424782f982725453808a94763bd19">viennacl::ocl::device::info</a></div><div class="ttdeci">std::string info(vcl_size_t indent=0, char indent_char= ' ') const </div><div class="ttdoc">Returns an info string with a few properties of the device. Use full_info() to get all details...</div><div class="ttdef"><b>Definition:</b> <a href="device_8hpp_source.html#l00995">device.hpp:995</a></div></div>
<div class="ttc" id="classviennacl_1_1context_html"><div class="ttname"><a href="classviennacl_1_1context.html">viennacl::context</a></div><div class="ttdoc">Represents a generic &#39;context&#39; similar to an OpenCL context, but is backend-agnostic and thus also su...</div><div class="ttdef"><b>Definition:</b> <a href="context_8hpp_source.html#l00039">context.hpp:39</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_aa18d10f8a90e38bd9ff43c650fc670ef"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#aa18d10f8a90e38bd9ff43c650fc670ef">viennacl::linalg::prod</a></div><div class="ttdeci">VectorT prod(std::vector&lt; std::vector&lt; T, A1 &gt;, A2 &gt; const &amp;matrix, VectorT const &amp;vector)</div><div class="ttdef"><b>Definition:</b> <a href="prod_8hpp_source.html#l00102">prod.hpp:102</a></div></div>
<div class="ttc" id="classviennacl_1_1linalg_1_1detail_1_1spai_1_1fspai__tag_html"><div class="ttname"><a href="classviennacl_1_1linalg_1_1detail_1_1spai_1_1fspai__tag.html">viennacl::linalg::detail::spai::fspai_tag</a></div><div class="ttdoc">A tag for FSPAI. Experimental. </div><div class="ttdef"><b>Definition:</b> <a href="fspai_8hpp_source.html#l00071">fspai.hpp:71</a></div></div>
<div class="ttc" id="ilu_8hpp_html"><div class="ttname"><a href="ilu_8hpp.html">ilu.hpp</a></div><div class="ttdoc">Implementations of incomplete factorization preconditioners. Convenience header file. </div></div>
<div class="ttc" id="classviennacl_1_1linalg_1_1fspai__precond_html"><div class="ttname"><a href="classviennacl_1_1linalg_1_1fspai__precond.html">viennacl::linalg::fspai_precond</a></div><div class="ttdoc">Implementation of the Factored SParse Approximate Inverse Algorithm for a generic, uBLAS-compatible matrix type. </div><div class="ttdef"><b>Definition:</b> <a href="spai_8hpp_source.html#l00189">spai.hpp:189</a></div></div>
<div class="ttc" id="classviennacl_1_1linalg_1_1ilut__tag_html"><div class="ttname"><a href="classviennacl_1_1linalg_1_1ilut__tag.html">viennacl::linalg::ilut_tag</a></div><div class="ttdoc">A tag for incomplete LU factorization with threshold (ILUT) </div><div class="ttdef"><b>Definition:</b> <a href="ilut_8hpp_source.html#l00045">ilut.hpp:45</a></div></div>
<div class="ttc" id="compressed__matrix_8hpp_html"><div class="ttname"><a href="compressed__matrix_8hpp.html">compressed_matrix.hpp</a></div><div class="ttdoc">Implementation of the compressed_matrix class. </div></div>
<div class="ttc" id="classviennacl_1_1linalg_1_1spai__precond_html"><div class="ttname"><a href="classviennacl_1_1linalg_1_1spai__precond.html">viennacl::linalg::spai_precond</a></div><div class="ttdoc">Implementation of the SParse Approximate Inverse Algorithm for a generic, uBLAS-compatible matrix typ...</div><div class="ttdef"><b>Definition:</b> <a href="spai_8hpp_source.html#l00075">spai.hpp:75</a></div></div>
<div class="ttc" id="solver_8cpp_html_afadbea4317c8f42f54d4622cf0a1829f"><div class="ttname"><a href="solver_8cpp.html#afadbea4317c8f42f54d4622cf0a1829f">run_solver</a></div><div class="ttdeci">void run_solver(MatrixType const &amp;matrix, VectorType const &amp;rhs, VectorType const &amp;ref_result, SolverTag const &amp;solver, PrecondTag const &amp;precond, long ops)</div><div class="ttdef"><b>Definition:</b> <a href="solver_8cpp_source.html#l00101">solver.cpp:101</a></div></div>
<div class="ttc" id="classviennacl_1_1linalg_1_1ilut__precond_html"><div class="ttname"><a href="classviennacl_1_1linalg_1_1ilut__precond.html">viennacl::linalg::ilut_precond</a></div><div class="ttdoc">ILUT preconditioner class, can be supplied to solve()-routines. </div><div class="ttdef"><b>Definition:</b> <a href="ilut_8hpp_source.html#l00352">ilut.hpp:352</a></div></div>
<div class="ttc" id="cg_8hpp_html"><div class="ttname"><a href="cg_8hpp.html">cg.hpp</a></div><div class="ttdoc">The conjugate gradient method is implemented here. </div></div>
<div class="ttc" id="classviennacl_1_1vector_html"><div class="ttname"><a href="classviennacl_1_1vector.html">viennacl::vector&lt; ScalarType &gt;</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a10b7f8cf6b8864a7aa196d670481a453"><div class="ttname"><a href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a></div><div class="ttdeci">void copy(std::vector&lt; NumericT &gt; &amp;cpu_vec, circulant_matrix&lt; NumericT, AlignmentV &gt; &amp;gpu_mat)</div><div class="ttdoc">Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...</div><div class="ttdef"><b>Definition:</b> <a href="circulant__matrix_8hpp_source.html#l00150">circulant_matrix.hpp:150</a></div></div>
<div class="ttc" id="fft__1d_8cpp_html_ad5c19ca4f47d3f8ec21232a5af2624e5"><div class="ttname"><a href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a></div><div class="ttdeci">float ScalarType</div><div class="ttdef"><b>Definition:</b> <a href="fft__1d_8cpp_source.html#l00042">fft_1d.cpp:42</a></div></div>
<div class="ttc" id="spai_8hpp_html"><div class="ttname"><a href="spai_8hpp.html">spai.hpp</a></div><div class="ttdoc">Main include file for the sparse approximate inverse preconditioner family (SPAI and FSPAI)...</div></div>
<div class="ttc" id="classviennacl_1_1compressed__matrix_html"><div class="ttname"><a href="classviennacl_1_1compressed__matrix.html">viennacl::compressed_matrix</a></div><div class="ttdoc">A sparse square matrix in compressed sparse rows format. </div><div class="ttdef"><b>Definition:</b> <a href="compressed__matrix_8hpp_source.html#l00559">compressed_matrix.hpp:559</a></div></div>
<div class="ttc" id="classviennacl_1_1linalg_1_1bicgstab__tag_html"><div class="ttname"><a href="classviennacl_1_1linalg_1_1bicgstab__tag.html">viennacl::linalg::bicgstab_tag</a></div><div class="ttdoc">A tag for the stabilized Bi-conjugate gradient solver. Used for supplying solver parameters and for d...</div><div class="ttdef"><b>Definition:</b> <a href="bicgstab_8hpp_source.html#l00047">bicgstab.hpp:47</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1io_html_ad934125ed3dbe661e264bcd7d62b1048"><div class="ttname"><a href="namespaceviennacl_1_1io.html#ad934125ed3dbe661e264bcd7d62b1048">viennacl::io::read_matrix_market_file</a></div><div class="ttdeci">long read_matrix_market_file(MatrixT &amp;mat, const char *file, long index_base=1)</div><div class="ttdoc">Reads a sparse matrix from a file (MatrixMarket format) </div><div class="ttdef"><b>Definition:</b> <a href="matrix__market_8hpp_source.html#l00339">matrix_market.hpp:339</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1ocl_html_a82c1aba632a7ee0991eee480d5340966"><div class="ttname"><a href="namespaceviennacl_1_1ocl.html#a82c1aba632a7ee0991eee480d5340966">viennacl::ocl::get_context</a></div><div class="ttdeci">viennacl::ocl::context &amp; get_context(long i)</div><div class="ttdoc">Convenience function for returning the current context. </div><div class="ttdef"><b>Definition:</b> <a href="backend_8hpp_source.html#l00225">backend.hpp:225</a></div></div>
<div class="ttc" id="scalar_8hpp_html"><div class="ttname"><a href="scalar_8hpp.html">scalar.hpp</a></div><div class="ttdoc">Implementation of the ViennaCL scalar class. </div></div>
<div class="ttc" id="namespaceviennacl_1_1ocl_html_add1725d48cfd159ce187e287369d1cdb"><div class="ttname"><a href="namespaceviennacl_1_1ocl.html#add1725d48cfd159ce187e287369d1cdb">viennacl::ocl::setup_context</a></div><div class="ttdeci">void setup_context(long i, std::vector&lt; cl_device_id &gt; const &amp;devices)</div><div class="ttdoc">Convenience function for setting devices for a context. </div><div class="ttdef"><b>Definition:</b> <a href="backend_8hpp_source.html#l00231">backend.hpp:231</a></div></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="navelem"><a class="el" href="dir_e6cbf65cd48c4de0cbf6207c238fcdc2.html">examples</a></li><li class="navelem"><a class="el" href="dir_a5d754652166e9a81bc9a52ffaf77760.html">tutorial</a></li><li class="navelem"><a class="el" href="spai_8cpp.html">spai.cpp</a></li>
    <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>