File: least-squares_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 (284 lines) | stat: -rw-r--r-- 25,155 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
<!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: least-squares.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('least-squares_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">least-squares.cpp</div>  </div>
</div><!--header-->
<div class="contents">
<p>This tutorial shows how least Squares problems for matrices from ViennaCL or Boost.uBLAS can be solved solved.</p>
<p>We start with including the respective header files: </p>
<div class="fragment"><div class="line"><span class="comment">// activate ublas support in ViennaCL</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS</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">// Boost headers</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/vector.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/vector_proxy.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/lu.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/io.hpp&gt;</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// ViennaCL headers</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="matrix_8hpp.html">viennacl/matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__proxy_8hpp.html">viennacl/matrix_proxy.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="qr_8hpp.html">viennacl/linalg/qr.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="lu_8hpp.html">viennacl/linalg/lu.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="direct__solve_8hpp.html">viennacl/linalg/direct_solve.hpp</a>&quot;</span></div>
</div><!-- fragment --><p> The minimization problem of finding x such that <img class="formulaInl" alt="$ \Vert Ax - b \Vert $" src="form_112.png"/> is solved as follows:</p>
<ul>
<li>Compute the QR-factorization of A = QR.</li>
<li>Compute <img class="formulaInl" alt="$ b' = Q^{\mathrm{T}} b $" src="form_113.png"/> for the equivalent minimization problem <img class="formulaInl" alt="$ \Vert Rx - Q^{\mathrm{T}} b $" src="form_114.png"/>.</li>
<li>Solve the triangular system <img class="formulaInl" alt="$ \tilde{R} x = b' $" src="form_115.png"/>, where <img class="formulaInl" alt="$ \tilde{R} $" src="form_116.png"/> is the upper square matrix of R.</li>
</ul>
<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> (<span class="keywordtype">int</span>, <span class="keyword">const</span> <span class="keywordtype">char</span> **)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">typedef</span> <span class="keywordtype">float</span>               <a name="a1"></a><a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>;     <span class="comment">//feel free to change this to &#39;double&#39; if supported by your hardware</span></div>
<div class="line"></div>
<div class="line">  <span class="keyword">typedef</span> boost::numeric::ublas::matrix&lt;ScalarType&gt;              MatrixType;</div>
<div class="line">  <span class="keyword">typedef</span> boost::numeric::ublas::vector&lt;ScalarType&gt;              VectorType;</div>
<div class="line">  <span class="keyword">typedef</span> <a name="_a2"></a><a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix&lt;ScalarType, viennacl::column_major&gt;</a>   VCLMatrixType;</div>
<div class="line">  <span class="keyword">typedef</span> <a name="_a3"></a><a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a>                           VCLVectorType;</div>
</div><!-- fragment --><p> Create vectors and matrices with data: </p>
<div class="fragment"><div class="line">VectorType ublas_b(4);</div>
<div class="line">ublas_b(0) = -4;</div>
<div class="line">ublas_b(1) =  2;</div>
<div class="line">ublas_b(2) =  5;</div>
<div class="line">ublas_b(3) = -1;</div>
<div class="line"></div>
<div class="line">MatrixType ublas_A(4, 3);</div>
<div class="line"></div>
<div class="line">ublas_A(0, 0) =  2; ublas_A(0, 1) = -1; ublas_A(0, 2) =  1;</div>
<div class="line">ublas_A(1, 0) =  1; ublas_A(1, 1) = -5; ublas_A(1, 2) =  2;</div>
<div class="line">ublas_A(2, 0) = -3; ublas_A(2, 1) =  1; ublas_A(2, 2) = -4;</div>
<div class="line">ublas_A(3, 0) =  1; ublas_A(3, 1) = -1; ublas_A(3, 2) =  1;</div>
</div><!-- fragment --><p> Setup the matrix and vector with ViennaCL objects and copy the data from the uBLAS objects: </p>
<div class="fragment"><div class="line">VCLVectorType <a name="a4"></a><a class="code" href="least-squares_8cpp.html#a0aff2af0ec74ea3ddb0129b72d9b921f">vcl_b</a>(ublas_b.size());</div>
<div class="line">VCLMatrixType vcl_A(ublas_A.size1(), ublas_A.size2());</div>
<div class="line"></div>
<div class="line"><a name="a5"></a><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_b, <a class="code" href="least-squares_8cpp.html#a0aff2af0ec74ea3ddb0129b72d9b921f">vcl_b</a>);</div>
<div class="line"><a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_A, vcl_A);</div>
</div><!-- fragment --> <h2>Option 1: Using Boost.uBLAS</h2>
<p>The implementation in ViennaCL accepts both uBLAS and ViennaCL types. We start with a single-threaded implementation using Boost.uBLAS. </p>
<div class="fragment"><div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;--- Boost.uBLAS ---&quot;</span> &lt;&lt; std::endl;</div>
</div><!-- fragment --><p> The first (and computationally most expensive) step is to compute the QR factorization of A. Since we do not need A later, we directly overwrite A with the householder reflectors and the upper triangular matrix R. The returned vector holds the scalar coefficients (betas) for the Householder reflections <img class="formulaInl" alt="$ I - \beta v v^{\mathrm{T}} $" src="form_117.png"/> </p>
<div class="fragment"><div class="line">std::vector&lt;ScalarType&gt; ublas_betas = <a name="a6"></a><a class="code" href="namespaceviennacl_1_1linalg.html#a83f573386903bae3a5379ab4cbfc5e2d">viennacl::linalg::inplace_qr</a>(ublas_A);</div>
</div><!-- fragment --><p> Compute the modified RHS of the minimization problem from the QR factorization, but do not form <img class="formulaInl" alt="$ Q^{\mathrm{T}} $" src="form_118.png"/> explicitly: b' := Q^T b </p>
<div class="fragment"><div class="line"><a name="a7"></a><a class="code" href="namespaceviennacl_1_1linalg.html#adf77bb0fe92ef09522d0458ee9df03db">viennacl::linalg::inplace_qr_apply_trans_Q</a>(ublas_A, ublas_betas, ublas_b);</div>
</div><!-- fragment --><p> Final step: triangular solve: Rx = b'', where b'' are the first three entries in b' We only need the upper left square part of A, which defines the upper triangular matrix R </p>
<div class="fragment"><div class="line"><a name="a8"></a><a class="code" href="namespaceviennacl.html#ae92c62d9fd59870c1f6b881e391d32aa">boost::numeric::ublas::range</a> ublas_range(0, 3);</div>
<div class="line">boost::numeric::ublas::matrix_range&lt;MatrixType&gt; ublas_R(ublas_A, ublas_range, ublas_range);</div>
<div class="line">boost::numeric::ublas::vector_range&lt;VectorType&gt; ublas_b2(ublas_b, ublas_range);</div>
<div class="line"><a name="a9"></a><a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#acb0b34e3ac22c72081eceaabeab80f22">boost::numeric::ublas::inplace_solve</a>(ublas_R, ublas_b2, boost::numeric::ublas::upper_tag());</div>
<div class="line"></div>
<div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;Result: &quot;</span> &lt;&lt; ublas_b2 &lt;&lt; std::endl;</div>
</div><!-- fragment --> <h2>Option 2: Use ViennaCL types</h2>
<p>ViennaCL is used for the computationally intensive BLAS 3 computations. Boost.uBLAS is used for the panel factorization on the host (CPU).</p>
<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">// activate ublas support in ViennaCL</span></div>
<div class="line"><span class="preprocessor">#define VIENNACL_WITH_UBLAS</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">// Boost headers</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/vector.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/vector_proxy.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/lu.hpp&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;boost/numeric/ublas/io.hpp&gt;</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="comment">// ViennaCL headers</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="matrix_8hpp.html">viennacl/matrix.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="matrix__proxy_8hpp.html">viennacl/matrix_proxy.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="qr_8hpp.html">viennacl/linalg/qr.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="lu_8hpp.html">viennacl/linalg/lu.hpp</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="direct__solve_8hpp.html">viennacl/linalg/direct_solve.hpp</a>&quot;</span></div>
<div class="line"></div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> <a class="code" href="tests_2src_2bisect_8cpp.html#ae66f6b31b5ad750f1fe042a706a4e3d4">main</a> (<span class="keywordtype">int</span>, <span class="keyword">const</span> <span class="keywordtype">char</span> **)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">typedef</span> <span class="keywordtype">float</span>               <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>;     <span class="comment">//feel free to change this to &#39;double&#39; if supported by your hardware</span></div>
<div class="line"></div>
<div class="line">  <span class="keyword">typedef</span> boost::numeric::ublas::matrix&lt;ScalarType&gt;              MatrixType;</div>
<div class="line">  <span class="keyword">typedef</span> boost::numeric::ublas::vector&lt;ScalarType&gt;              VectorType;</div>
<div class="line">  <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1matrix.html">viennacl::matrix&lt;ScalarType, viennacl::column_major&gt;</a>   VCLMatrixType;</div>
<div class="line">  <span class="keyword">typedef</span> <a class="code" href="classviennacl_1_1vector.html">viennacl::vector&lt;ScalarType&gt;</a>                           VCLVectorType;</div>
<div class="line"></div>
<div class="line">  VectorType ublas_b(4);</div>
<div class="line">  ublas_b(0) = -4;</div>
<div class="line">  ublas_b(1) =  2;</div>
<div class="line">  ublas_b(2) =  5;</div>
<div class="line">  ublas_b(3) = -1;</div>
<div class="line"></div>
<div class="line">  MatrixType ublas_A(4, 3);</div>
<div class="line"></div>
<div class="line">  ublas_A(0, 0) =  2; ublas_A(0, 1) = -1; ublas_A(0, 2) =  1;</div>
<div class="line">  ublas_A(1, 0) =  1; ublas_A(1, 1) = -5; ublas_A(1, 2) =  2;</div>
<div class="line">  ublas_A(2, 0) = -3; ublas_A(2, 1) =  1; ublas_A(2, 2) = -4;</div>
<div class="line">  ublas_A(3, 0) =  1; ublas_A(3, 1) = -1; ublas_A(3, 2) =  1;</div>
<div class="line"></div>
<div class="line">  VCLVectorType <a class="code" href="least-squares_8cpp.html#a0aff2af0ec74ea3ddb0129b72d9b921f">vcl_b</a>(ublas_b.size());</div>
<div class="line">  VCLMatrixType vcl_A(ublas_A.size1(), ublas_A.size2());</div>
<div class="line"></div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_b, <a class="code" href="least-squares_8cpp.html#a0aff2af0ec74ea3ddb0129b72d9b921f">vcl_b</a>);</div>
<div class="line">  <a class="code" href="namespaceviennacl.html#a10b7f8cf6b8864a7aa196d670481a453">viennacl::copy</a>(ublas_A, vcl_A);</div>
<div class="line"></div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;--- Boost.uBLAS ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  std::vector&lt;ScalarType&gt; ublas_betas = <a class="code" href="namespaceviennacl_1_1linalg.html#a83f573386903bae3a5379ab4cbfc5e2d">viennacl::linalg::inplace_qr</a>(ublas_A);</div>
<div class="line"></div>
<div class="line">  <a class="code" href="namespaceviennacl_1_1linalg.html#adf77bb0fe92ef09522d0458ee9df03db">viennacl::linalg::inplace_qr_apply_trans_Q</a>(ublas_A, ublas_betas, ublas_b);</div>
<div class="line"></div>
<div class="line">  <a class="code" href="namespaceviennacl.html#ae92c62d9fd59870c1f6b881e391d32aa">boost::numeric::ublas::range</a> ublas_range(0, 3);</div>
<div class="line">  boost::numeric::ublas::matrix_range&lt;MatrixType&gt; ublas_R(ublas_A, ublas_range, ublas_range);</div>
<div class="line">  boost::numeric::ublas::vector_range&lt;VectorType&gt; ublas_b2(ublas_b, ublas_range);</div>
<div class="line">  <a class="code" href="namespaceviennacl_1_1linalg_1_1cuda.html#acb0b34e3ac22c72081eceaabeab80f22">boost::numeric::ublas::inplace_solve</a>(ublas_R, ublas_b2, boost::numeric::ublas::upper_tag());</div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Result: &quot;</span> &lt;&lt; ublas_b2 &lt;&lt; std::endl;</div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;--- ViennaCL (hybrid implementation)  ---&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">  std::vector&lt;ScalarType&gt; <a name="a10"></a><a class="code" href="least-squares_8cpp.html#ac887df28e24b43a0bda14bd218659c78">hybrid_betas</a> = <a class="code" href="namespaceviennacl_1_1linalg.html#a83f573386903bae3a5379ab4cbfc5e2d">viennacl::linalg::inplace_qr</a>(vcl_A);</div>
<div class="line"></div>
<div class="line">  <a class="code" href="namespaceviennacl_1_1linalg.html#adf77bb0fe92ef09522d0458ee9df03db">viennacl::linalg::inplace_qr_apply_trans_Q</a>(vcl_A, hybrid_betas, <a class="code" href="least-squares_8cpp.html#a0aff2af0ec74ea3ddb0129b72d9b921f">vcl_b</a>);</div>
<div class="line"></div>
<div class="line">  <a name="_a11"></a><a class="code" href="classviennacl_1_1basic__range.html">viennacl::range</a> vcl_range(0, 3);</div>
<div class="line">  <a name="_a12"></a><a class="code" href="classviennacl_1_1matrix__range.html">viennacl::matrix_range&lt;VCLMatrixType&gt;</a> <a name="a13"></a><a class="code" href="least-squares_8cpp.html#aba586554adcf743c2cd6b09d207957f5">vcl_R</a>(vcl_A, vcl_range, vcl_range);</div>
<div class="line">  <a name="_a14"></a><a class="code" href="classviennacl_1_1vector__range.html">viennacl::vector_range&lt;VCLVectorType&gt;</a> <a name="a15"></a><a class="code" href="least-squares_8cpp.html#aa87c3ccaa3cee96a4db6f62dff634562">vcl_b2</a>(<a class="code" href="least-squares_8cpp.html#a0aff2af0ec74ea3ddb0129b72d9b921f">vcl_b</a>, vcl_range);</div>
<div class="line">  <a name="a16"></a><a class="code" href="namespaceviennacl_1_1linalg.html#adf06f61be9418e5817eed66004f9dd2b">viennacl::linalg::inplace_solve</a>(vcl_R, vcl_b2, <a name="_a17"></a><a class="code" href="structviennacl_1_1linalg_1_1upper__tag.html">viennacl::linalg::upper_tag</a>());</div>
<div class="line"></div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;Result: &quot;</span> &lt;&lt; vcl_b2 &lt;&lt; std::endl;</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>