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
|
<!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: viennacl/linalg/qr-method-common.hpp 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
 <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('qr-method-common_8hpp_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"> </span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark"> </span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark"> </span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark"> </span>Files</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark"> </span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark"> </span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark"> </span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark"> </span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark"> </span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark"> </span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark"> </span>Macros</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(11)"><span class="SelectionMark"> </span>Pages</a></div>
<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<iframe src="javascript:void(0)" frameborder="0"
name="MSearchResults" id="MSearchResults">
</iframe>
</div>
<div class="header">
<div class="headertitle">
<div class="title">qr-method-common.hpp</div> </div>
</div><!--header-->
<div class="contents">
<a href="qr-method-common_8hpp.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="preprocessor">#ifndef VIENNACL_LINALG_QR_METHOD_COMMON_HPP</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="preprocessor">#define VIENNACL_LINALG_QR_METHOD_COMMON_HPP</span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span> </div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment">/* =========================================================================</span></div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment"> Copyright (c) 2010-2016, Institute for Microelectronics,</span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment"> Institute for Analysis and Scientific Computing,</span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment"> TU Wien.</span></div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment"> Portions of this software are copyright by UChicago Argonne, LLC.</span></div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="comment"></span></div>
<div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="comment"> -----------------</span></div>
<div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="comment"> ViennaCL - The Vienna Computing Library</span></div>
<div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="comment"> -----------------</span></div>
<div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="comment"></span></div>
<div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="comment"> Project Head: Karl Rupp rupp@iue.tuwien.ac.at</span></div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span> <span class="comment"></span></div>
<div class="line"><a name="l00016"></a><span class="lineno"> 16</span> <span class="comment"> (A list of authors and contributors can be found in the manual)</span></div>
<div class="line"><a name="l00017"></a><span class="lineno"> 17</span> <span class="comment"></span></div>
<div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="comment"> License: MIT (X11), see file LICENSE in the base directory</span></div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="comment">============================================================================= */</span></div>
<div class="line"><a name="l00020"></a><span class="lineno"> 20</span> </div>
<div class="line"><a name="l00021"></a><span class="lineno"> 21</span> <span class="preprocessor">#include <cmath></span></div>
<div class="line"><a name="l00022"></a><span class="lineno"> 22</span> </div>
<div class="line"><a name="l00023"></a><span class="lineno"> 23</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00024"></a><span class="lineno"> 24</span> <span class="preprocessor">#include "<a class="code" href="device_8hpp.html">viennacl/ocl/device.hpp</a>"</span></div>
<div class="line"><a name="l00025"></a><span class="lineno"> 25</span> <span class="preprocessor">#include "<a class="code" href="ocl_2handle_8hpp.html">viennacl/ocl/handle.hpp</a>"</span></div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span> <span class="preprocessor">#include "<a class="code" href="kernel_8hpp.html">viennacl/ocl/kernel.hpp</a>"</span></div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span> <span class="preprocessor">#include "<a class="code" href="opencl_2kernels_2svd_8hpp.html">viennacl/linalg/opencl/kernels/svd.hpp</a>"</span></div>
<div class="line"><a name="l00028"></a><span class="lineno"> 28</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00029"></a><span class="lineno"> 29</span> </div>
<div class="line"><a name="l00030"></a><span class="lineno"> 30</span> <span class="preprocessor">#ifdef VIENNACL_WITH_CUDA</span></div>
<div class="line"><a name="l00031"></a><span class="lineno"> 31</span> <span class="preprocessor">#include "<a class="code" href="cuda_2matrix__operations_8hpp.html">viennacl/linalg/cuda/matrix_operations.hpp</a>"</span></div>
<div class="line"><a name="l00032"></a><span class="lineno"> 32</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span> <span class="preprocessor">#include "<a class="code" href="result__of_8hpp.html">viennacl/meta/result_of.hpp</a>"</span></div>
<div class="line"><a name="l00034"></a><span class="lineno"> 34</span> <span class="preprocessor">#include "<a class="code" href="vector_8hpp.html">viennacl/vector.hpp</a>"</span></div>
<div class="line"><a name="l00035"></a><span class="lineno"> 35</span> <span class="preprocessor">#include "<a class="code" href="matrix_8hpp.html">viennacl/matrix.hpp</a>"</span></div>
<div class="line"><a name="l00036"></a><span class="lineno"> 36</span> </div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span> <span class="comment">//#include <boost/numeric/ublas/vector.hpp></span></div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span> <span class="comment">//#include <boost/numeric/ublas/io.hpp></span></div>
<div class="line"><a name="l00039"></a><span class="lineno"> 39</span> </div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span> <span class="keyword">namespace </span>viennacl</div>
<div class="line"><a name="l00045"></a><span class="lineno"> 45</span> {</div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span> <span class="keyword">namespace </span>linalg</div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span> {</div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span> </div>
<div class="line"><a name="l00049"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a5b5cb81b46208cae772a5e0813f02517"> 49</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a5b5cb81b46208cae772a5e0813f02517">SVD_HOUSEHOLDER_UPDATE_QR_KERNEL</a> = <span class="stringliteral">"house_update_QR"</span>;</div>
<div class="line"><a name="l00050"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a349e6e78f95ef631a5f2736158766739"> 50</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a349e6e78f95ef631a5f2736158766739">SVD_MATRIX_TRANSPOSE_KERNEL</a> = <span class="stringliteral">"transpose_inplace"</span>;</div>
<div class="line"><a name="l00051"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a05717ddd22894376e780b075d59616ed"> 51</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a05717ddd22894376e780b075d59616ed">SVD_INVERSE_SIGNS_KERNEL</a> = <span class="stringliteral">"inverse_signs"</span>;</div>
<div class="line"><a name="l00052"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a15bedacb51ec9b425baf52ee7f68bdba"> 52</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a15bedacb51ec9b425baf52ee7f68bdba">SVD_GIVENS_PREV_KERNEL</a> = <span class="stringliteral">"givens_prev"</span>;</div>
<div class="line"><a name="l00053"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a8477d884683dafc7e6ca04f2e02a1d1f"> 53</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a8477d884683dafc7e6ca04f2e02a1d1f">SVD_FINAL_ITER_UPDATE_KERNEL</a> = <span class="stringliteral">"final_iter_update"</span>;</div>
<div class="line"><a name="l00054"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a3a439370108ddfa59c2fcb4136b13c8d"> 54</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a3a439370108ddfa59c2fcb4136b13c8d">SVD_UPDATE_QR_COLUMN_KERNEL</a> = <span class="stringliteral">"update_qr_column"</span>;</div>
<div class="line"><a name="l00055"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a8e62fa452895c889f7bad9700300994d"> 55</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a8e62fa452895c889f7bad9700300994d">SVD_HOUSEHOLDER_UPDATE_A_LEFT_KERNEL</a> = <span class="stringliteral">"house_update_A_left"</span>;</div>
<div class="line"><a name="l00056"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#ace85ea4095e9e7404b6f77b767ac000f"> 56</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#ace85ea4095e9e7404b6f77b767ac000f">SVD_HOUSEHOLDER_UPDATE_A_RIGHT_KERNEL</a> = <span class="stringliteral">"house_update_A_right"</span>;</div>
<div class="line"><a name="l00057"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg.html#a04a5ce1f7230502c9d1ee87dcc24c8b4"> 57</a></span> <span class="keyword">const</span> std::string <a class="code" href="namespaceviennacl_1_1linalg.html#a04a5ce1f7230502c9d1ee87dcc24c8b4">SVD_HOUSEHOLDER_UPDATE_QL_KERNEL</a> = <span class="stringliteral">"house_update_QL"</span>;</div>
<div class="line"><a name="l00058"></a><span class="lineno"> 58</span> </div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span> <span class="keyword">namespace </span>detail</div>
<div class="line"><a name="l00060"></a><span class="lineno"> 60</span> {</div>
<div class="line"><a name="l00061"></a><span class="lineno"> 61</span> <span class="keyword">static</span> <span class="keyword">const</span> <span class="keywordtype">double</span> <a class="code" href="tests_2src_2bisect_8cpp.html#a6ebf6899d6c1c8b7b9d09be872c05aae">EPS</a> = 1e-10;</div>
<div class="line"><a name="l00062"></a><span class="lineno"> 62</span> <span class="keyword">static</span> <span class="keyword">const</span> <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> ITER_MAX = 50;</div>
<div class="line"><a name="l00063"></a><span class="lineno"> 63</span> </div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span> <span class="keyword">template</span> <<span class="keyword">typename</span> SCALARTYPE></div>
<div class="line"><a name="l00065"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#a0737c94fac817b0d0c8444463be7ceb7"> 65</a></span> SCALARTYPE <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a0737c94fac817b0d0c8444463be7ceb7">pythag</a>(SCALARTYPE a, SCALARTYPE b)</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span> {</div>
<div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  <span class="keywordflow">return</span> std::sqrt(a*a + b*b);</div>
<div class="line"><a name="l00068"></a><span class="lineno"> 68</span> }</div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span> </div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span> <span class="keyword">template</span> <<span class="keyword">typename</span> SCALARTYPE></div>
<div class="line"><a name="l00071"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#a7bfa1c9cb0afb4ca0f975e690f4a8f07"> 71</a></span> SCALARTYPE <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a7bfa1c9cb0afb4ca0f975e690f4a8f07">sign</a>(SCALARTYPE val)</div>
<div class="line"><a name="l00072"></a><span class="lineno"> 72</span> {</div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <span class="keywordflow">return</span> (val >= 0) ? SCALARTYPE(1) : SCALARTYPE(-1);</div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span> }</div>
<div class="line"><a name="l00075"></a><span class="lineno"> 75</span> </div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span> <span class="comment">// DEPRECATED: Replace with viennacl::linalg::norm_2</span></div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span> <span class="keyword">template</span> <<span class="keyword">typename</span> VectorType></div>
<div class="line"><a name="l00078"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#ae059c4a3216ff6bdc44110377ce34e88"> 78</a></span> <span class="keyword">typename</span> VectorType::value_type <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#ae059c4a3216ff6bdc44110377ce34e88">norm_lcl</a>(VectorType <span class="keyword">const</span> & x, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">size</a>)</div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span> {</div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  <span class="keyword">typename</span> VectorType::value_type x_norm = 0.0;</div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  <span class="keywordflow">for</span>(<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> i = 0; i < <a class="code" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">size</a>; i++)</div>
<div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  x_norm += std::pow(x[i], 2);</div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  <span class="keywordflow">return</span> std::sqrt(x_norm);</div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span> }</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span> </div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span> <span class="keyword">template</span> <<span class="keyword">typename</span> VectorType></div>
<div class="line"><a name="l00087"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#ae152dd8fd5f272525bac1b843375d329"> 87</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#ae152dd8fd5f272525bac1b843375d329">normalize</a>(VectorType & x, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">size</a>)</div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span> {</div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  <span class="keyword">typename</span> VectorType::value_type x_norm = <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#ae059c4a3216ff6bdc44110377ce34e88">norm_lcl</a>(x, size);</div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="keywordflow">for</span>(<a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> i = 0; i < <a class="code" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">size</a>; i++)</div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  x[i] /= x_norm;</div>
<div class="line"><a name="l00092"></a><span class="lineno"> 92</span> }</div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span> </div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span> </div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span> </div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span> <span class="keyword">template</span> <<span class="keyword">typename</span> VectorType></div>
<div class="line"><a name="l00097"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#a366df7768bf52e820fee364334133db1"> 97</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a366df7768bf52e820fee364334133db1">householder_vector</a>(VectorType & v, <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="namespaceviennacl_1_1traits.html#a68603a1e1ca1bfb6d107831ba5096786">start</a>)</div>
<div class="line"><a name="l00098"></a><span class="lineno"> 98</span> {</div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> VectorType::value_type <a class="code" href="fft__1d_8cpp.html#ad5c19ca4f47d3f8ec21232a5af2624e5">ScalarType</a>;</div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  ScalarType x_norm = <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#ae059c4a3216ff6bdc44110377ce34e88">norm_lcl</a>(v, v.size());</div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  ScalarType alpha = -<a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a7bfa1c9cb0afb4ca0f975e690f4a8f07">sign</a>(v[start]) * x_norm;</div>
<div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  v[<a class="code" href="namespaceviennacl_1_1traits.html#a68603a1e1ca1bfb6d107831ba5096786">start</a>] += alpha;</div>
<div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#ae152dd8fd5f272525bac1b843375d329">normalize</a>(v, v.size());</div>
<div class="line"><a name="l00104"></a><span class="lineno"> 104</span> }</div>
<div class="line"><a name="l00105"></a><span class="lineno"> 105</span> </div>
<div class="line"><a name="l00106"></a><span class="lineno"> 106</span> <span class="keyword">template</span> <<span class="keyword">typename</span> SCALARTYPE></div>
<div class="line"><a name="l00107"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#aad4cd2951549618c93f74ff26444675f"> 107</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#aad4cd2951549618c93f74ff26444675f">transpose</a>(<a class="code" href="classviennacl_1_1matrix__base.html">matrix_base<SCALARTYPE></a> & A)</div>
<div class="line"><a name="l00108"></a><span class="lineno"> 108</span> {</div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  (void)A;</div>
<div class="line"><a name="l00110"></a><span class="lineno"> 110</span> <span class="preprocessor">#ifdef VIENNACL_WITH_OPENCL</span></div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span> </div>
<div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  <a class="code" href="classviennacl_1_1ocl_1_1context.html">viennacl::ocl::context</a> & ctx = <span class="keyword">const_cast<</span><a class="code" href="classviennacl_1_1ocl_1_1context.html">viennacl::ocl::context</a> &<span class="keyword">></span>(viennacl::traits::opencl_handle(A).context());</div>
<div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  <span class="keywordflow">if</span>(A.<a class="code" href="classviennacl_1_1matrix__base.html#a15ec879f9e85a6fc209454b0bf137084">row_major</a>())</div>
<div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  {</div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  <a class="code" href="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd.html#ab9e831dbb989de8ddc771c949658b514">viennacl::linalg::opencl::kernels::svd<SCALARTYPE, row_major>::init</a>(ctx);</div>
<div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  <a class="code" href="classviennacl_1_1ocl_1_1kernel.html">viennacl::ocl::kernel</a> & kernel = <a class="code" href="namespaceviennacl_1_1ocl.html#a61a73653d92f1eb2ef9c649ec253e29f">viennacl::ocl::get_kernel</a>(<a class="code" href="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd.html">viennacl::linalg::opencl::kernels::svd<SCALARTYPE, row_major>::program_name</a>(), <a class="code" href="namespaceviennacl_1_1linalg.html#a349e6e78f95ef631a5f2736158766739">SVD_MATRIX_TRANSPOSE_KERNEL</a>);</div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span> </div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  <a class="code" href="namespaceviennacl_1_1ocl.html#a5f2022f653ea1cf364d20e3ff84dcada">viennacl::ocl::enqueue</a>(kernel(A,</div>
<div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  static_cast<cl_uint>(A.<a class="code" href="classviennacl_1_1matrix__base.html#a82018b4e169973bdfb9d3be68ceb5be0">internal_size1</a>()),</div>
<div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  static_cast<cl_uint>(A.<a class="code" href="classviennacl_1_1matrix__base.html#a00e30fdc61081a73995b54175a91d220">internal_size2</a>())</div>
<div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  )</div>
<div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  );</div>
<div class="line"><a name="l00123"></a><span class="lineno"> 123</span>  }</div>
<div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  {</div>
<div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  <a class="code" href="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd.html#ab9e831dbb989de8ddc771c949658b514">viennacl::linalg::opencl::kernels::svd<SCALARTYPE, row_major>::init</a>(ctx);</div>
<div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  <a class="code" href="classviennacl_1_1ocl_1_1kernel.html">viennacl::ocl::kernel</a> & kernel = <a class="code" href="namespaceviennacl_1_1ocl.html#a61a73653d92f1eb2ef9c649ec253e29f">viennacl::ocl::get_kernel</a>(<a class="code" href="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd.html">viennacl::linalg::opencl::kernels::svd<SCALARTYPE, column_major>::program_name</a>(), <a class="code" href="namespaceviennacl_1_1linalg.html#a349e6e78f95ef631a5f2736158766739">SVD_MATRIX_TRANSPOSE_KERNEL</a>);</div>
<div class="line"><a name="l00128"></a><span class="lineno"> 128</span> </div>
<div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  <a class="code" href="namespaceviennacl_1_1ocl.html#a5f2022f653ea1cf364d20e3ff84dcada">viennacl::ocl::enqueue</a>(kernel(A,</div>
<div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  static_cast<cl_uint>(A.<a class="code" href="classviennacl_1_1matrix__base.html#a82018b4e169973bdfb9d3be68ceb5be0">internal_size1</a>()),</div>
<div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  static_cast<cl_uint>(A.<a class="code" href="classviennacl_1_1matrix__base.html#a00e30fdc61081a73995b54175a91d220">internal_size2</a>())</div>
<div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  )</div>
<div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  );</div>
<div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  }</div>
<div class="line"><a name="l00135"></a><span class="lineno"> 135</span> </div>
<div class="line"><a name="l00136"></a><span class="lineno"> 136</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00137"></a><span class="lineno"> 137</span> }</div>
<div class="line"><a name="l00138"></a><span class="lineno"> 138</span> </div>
<div class="line"><a name="l00139"></a><span class="lineno"> 139</span> </div>
<div class="line"><a name="l00140"></a><span class="lineno"> 140</span> </div>
<div class="line"><a name="l00141"></a><span class="lineno"> 141</span> <span class="keyword">template</span> <<span class="keyword">typename</span> T></div>
<div class="line"><a name="l00142"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#a807e733d7bf94d491377fff3905e25d2"> 142</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a807e733d7bf94d491377fff3905e25d2">cdiv</a>(T xr, T xi, T yr, T yi, T& cdivr, T& cdivi)</div>
<div class="line"><a name="l00143"></a><span class="lineno"> 143</span> {</div>
<div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  <span class="comment">// Complex scalar division.</span></div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  T r;</div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  T d;</div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  <span class="keywordflow">if</span> (std::fabs(yr) > std::fabs(yi))</div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  {</div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  r = yi / yr;</div>
<div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  d = yr + r * yi;</div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  cdivr = (xr + r * xi) / d;</div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  cdivi = (xi - r * xr) / d;</div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  }</div>
<div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  {</div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  r = yr / yi;</div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  d = yi + r * yr;</div>
<div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  cdivr = (r * xr + xi) / d;</div>
<div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  cdivi = (r * xi - xr) / d;</div>
<div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  }</div>
<div class="line"><a name="l00161"></a><span class="lineno"> 161</span> }</div>
<div class="line"><a name="l00162"></a><span class="lineno"> 162</span> </div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span> </div>
<div class="line"><a name="l00164"></a><span class="lineno"> 164</span> <span class="keyword">template</span><<span class="keyword">typename</span> SCALARTYPE></div>
<div class="line"><a name="l00165"></a><span class="lineno"><a class="line" href="namespaceviennacl_1_1linalg_1_1detail.html#a5852aa7527c0599681f1b0fbc753d357"> 165</a></span> <span class="keywordtype">void</span> <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a5852aa7527c0599681f1b0fbc753d357">prepare_householder_vector</a>(</div>
<div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  <a class="code" href="classviennacl_1_1matrix__base.html">matrix_base<SCALARTYPE></a>& A,</div>
<div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  <a class="code" href="classviennacl_1_1vector__base.html">vector_base<SCALARTYPE></a>& D,</div>
<div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">size</a>,</div>
<div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> row_start,</div>
<div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> col_start,</div>
<div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  <a class="code" href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">vcl_size_t</a> <a class="code" href="namespaceviennacl_1_1traits.html#a68603a1e1ca1bfb6d107831ba5096786">start</a>,</div>
<div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  <span class="keywordtype">bool</span> is_column</div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  )</div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span> {</div>
<div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <span class="comment">//boost::numeric::ublas::vector<SCALARTYPE> tmp = boost::numeric::ublas::scalar_vector<SCALARTYPE>(size, 0);</span></div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  std::vector<SCALARTYPE> tmp(size);</div>
<div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  <a class="code" href="namespaceviennacl_1_1linalg.html#a0aec751d37dbd7b6d45ea9340fa8a317">copy_vec</a>(A, D, row_start, col_start, is_column);</div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  <a class="code" href="namespaceviennacl.html#aaa5c8726b45bc89a523ca2fa8c42107a">fast_copy</a>(D.<a class="code" href="classviennacl_1_1vector__base.html#a762d98e2f912fc534951f25555b6077f">begin</a>(), D.<a class="code" href="classviennacl_1_1vector__base.html#a762d98e2f912fc534951f25555b6077f">begin</a>() + <a class="code" href="namespaceviennacl.html#afad612e26a65487b7e747c18c678ffd9">vcl_ptrdiff_t</a>(size - start), tmp.begin() + <a class="code" href="namespaceviennacl.html#afad612e26a65487b7e747c18c678ffd9">vcl_ptrdiff_t</a>(start));</div>
<div class="line"><a name="l00179"></a><span class="lineno"> 179</span> </div>
<div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  <a class="code" href="namespaceviennacl_1_1linalg_1_1detail.html#a366df7768bf52e820fee364334133db1">detail::householder_vector</a>(tmp, start);</div>
<div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  <a class="code" href="namespaceviennacl.html#aaa5c8726b45bc89a523ca2fa8c42107a">fast_copy</a>(tmp, D);</div>
<div class="line"><a name="l00182"></a><span class="lineno"> 182</span> }</div>
<div class="line"><a name="l00183"></a><span class="lineno"> 183</span> </div>
<div class="line"><a name="l00184"></a><span class="lineno"> 184</span> } <span class="comment">//detail</span></div>
<div class="line"><a name="l00185"></a><span class="lineno"> 185</span> }</div>
<div class="line"><a name="l00186"></a><span class="lineno"> 186</span> }</div>
<div class="line"><a name="l00187"></a><span class="lineno"> 187</span> </div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span> <span class="preprocessor">#endif</span></div>
<div class="ttc" id="tests_2src_2bisect_8cpp_html_a6ebf6899d6c1c8b7b9d09be872c05aae"><div class="ttname"><a href="tests_2src_2bisect_8cpp.html#a6ebf6899d6c1c8b7b9d09be872c05aae">EPS</a></div><div class="ttdeci">#define EPS</div><div class="ttdef"><b>Definition:</b> <a href="tests_2src_2bisect_8cpp_source.html#l00038">bisect.cpp:38</a></div></div>
<div class="ttc" id="device_8hpp_html"><div class="ttname"><a href="device_8hpp.html">device.hpp</a></div><div class="ttdoc">Represents an OpenCL device within ViennaCL. </div></div>
<div class="ttc" id="namespaceviennacl_1_1ocl_html_a61a73653d92f1eb2ef9c649ec253e29f"><div class="ttname"><a href="namespaceviennacl_1_1ocl.html#a61a73653d92f1eb2ef9c649ec253e29f">viennacl::ocl::get_kernel</a></div><div class="ttdeci">viennacl::ocl::kernel & get_kernel(std::string const &prog_name, std::string const &kernel_name)</div><div class="ttdoc">Convenience function for getting the kernel for a particular program from the current active context...</div><div class="ttdef"><b>Definition:</b> <a href="backend_8hpp_source.html#l00339">backend.hpp:339</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a05717ddd22894376e780b075d59616ed"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a05717ddd22894376e780b075d59616ed">viennacl::linalg::SVD_INVERSE_SIGNS_KERNEL</a></div><div class="ttdeci">const std::string SVD_INVERSE_SIGNS_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00051">qr-method-common.hpp:51</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a0aec751d37dbd7b6d45ea9340fa8a317"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a0aec751d37dbd7b6d45ea9340fa8a317">viennacl::linalg::copy_vec</a></div><div class="ttdeci">void copy_vec(matrix_base< SCALARTYPE > &A, vector_base< SCALARTYPE > &V, vcl_size_t row_start, vcl_size_t col_start, bool copy_col)</div><div class="ttdoc">This function copies a row or a column from a matrix to a vector. </div><div class="ttdef"><b>Definition:</b> <a href="matrix__operations_8hpp_source.html#l00942">matrix_operations.hpp:942</a></div></div>
<div class="ttc" id="classviennacl_1_1ocl_1_1kernel_html"><div class="ttname"><a href="classviennacl_1_1ocl_1_1kernel.html">viennacl::ocl::kernel</a></div><div class="ttdoc">Represents an OpenCL kernel within ViennaCL. </div><div class="ttdef"><b>Definition:</b> <a href="kernel_8hpp_source.html#l00058">kernel.hpp:58</a></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="classviennacl_1_1ocl_1_1context_html"><div class="ttname"><a href="classviennacl_1_1ocl_1_1context.html">viennacl::ocl::context</a></div><div class="ttdoc">Manages an OpenCL context and provides the respective convenience functions for creating buffers...</div><div class="ttdef"><b>Definition:</b> <a href="ocl_2context_8hpp_source.html#l00055">context.hpp:55</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a3a439370108ddfa59c2fcb4136b13c8d"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a3a439370108ddfa59c2fcb4136b13c8d">viennacl::linalg::SVD_UPDATE_QR_COLUMN_KERNEL</a></div><div class="ttdeci">const std::string SVD_UPDATE_QR_COLUMN_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00054">qr-method-common.hpp:54</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__base_html"><div class="ttname"><a href="classviennacl_1_1matrix__base.html">viennacl::matrix_base</a></div><div class="ttdef"><b>Definition:</b> <a href="matrix__def_8hpp_source.html#l00103">matrix_def.hpp:103</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_ace85ea4095e9e7404b6f77b767ac000f"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#ace85ea4095e9e7404b6f77b767ac000f">viennacl::linalg::SVD_HOUSEHOLDER_UPDATE_A_RIGHT_KERNEL</a></div><div class="ttdeci">const std::string SVD_HOUSEHOLDER_UPDATE_A_RIGHT_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00056">qr-method-common.hpp:56</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_a807e733d7bf94d491377fff3905e25d2"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a807e733d7bf94d491377fff3905e25d2">viennacl::linalg::detail::cdiv</a></div><div class="ttdeci">void cdiv(T xr, T xi, T yr, T yi, T &cdivr, T &cdivi)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00142">qr-method-common.hpp:142</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a5b5cb81b46208cae772a5e0813f02517"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a5b5cb81b46208cae772a5e0813f02517">viennacl::linalg::SVD_HOUSEHOLDER_UPDATE_QR_KERNEL</a></div><div class="ttdeci">const std::string SVD_HOUSEHOLDER_UPDATE_QR_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00049">qr-method-common.hpp:49</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_aad4cd2951549618c93f74ff26444675f"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#aad4cd2951549618c93f74ff26444675f">viennacl::linalg::detail::transpose</a></div><div class="ttdeci">void transpose(matrix_base< SCALARTYPE > &A)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00107">qr-method-common.hpp:107</a></div></div>
<div class="ttc" id="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd_html_ab9e831dbb989de8ddc771c949658b514"><div class="ttname"><a href="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd.html#ab9e831dbb989de8ddc771c949658b514">viennacl::linalg::opencl::kernels::svd::init</a></div><div class="ttdeci">static void init(viennacl::ocl::context &ctx)</div><div class="ttdef"><b>Definition:</b> <a href="opencl_2kernels_2svd_8hpp_source.html#l00652">svd.hpp:652</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_a5852aa7527c0599681f1b0fbc753d357"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a5852aa7527c0599681f1b0fbc753d357">viennacl::linalg::detail::prepare_householder_vector</a></div><div class="ttdeci">void prepare_householder_vector(matrix_base< SCALARTYPE > &A, vector_base< SCALARTYPE > &D, vcl_size_t size, vcl_size_t row_start, vcl_size_t col_start, vcl_size_t start, bool is_column)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00165">qr-method-common.hpp:165</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_ae152dd8fd5f272525bac1b843375d329"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#ae152dd8fd5f272525bac1b843375d329">viennacl::linalg::detail::normalize</a></div><div class="ttdeci">void normalize(VectorType &x, vcl_size_t size)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00087">qr-method-common.hpp:87</a></div></div>
<div class="ttc" id="opencl_2kernels_2svd_8hpp_html"><div class="ttname"><a href="opencl_2kernels_2svd_8hpp.html">svd.hpp</a></div><div class="ttdoc">OpenCL kernel file for singular value decomposition. </div></div>
<div class="ttc" id="namespaceviennacl_1_1traits_html_aa2344ea20469f55fbc15a8e9526494d0"><div class="ttname"><a href="namespaceviennacl_1_1traits.html#aa2344ea20469f55fbc15a8e9526494d0">viennacl::traits::size</a></div><div class="ttdeci">vcl_size_t size(VectorType const &vec)</div><div class="ttdoc">Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.) </div><div class="ttdef"><b>Definition:</b> <a href="size_8hpp_source.html#l00239">size.hpp:239</a></div></div>
<div class="ttc" id="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd_html"><div class="ttname"><a href="structviennacl_1_1linalg_1_1opencl_1_1kernels_1_1svd.html">viennacl::linalg::opencl::kernels::svd</a></div><div class="ttdoc">Main kernel class for generating OpenCL kernels for singular value decomposition of dense matrices...</div><div class="ttdef"><b>Definition:</b> <a href="opencl_2kernels_2svd_8hpp_source.html#l00644">svd.hpp:644</a></div></div>
<div class="ttc" id="classviennacl_1_1vector__base_html_a762d98e2f912fc534951f25555b6077f"><div class="ttname"><a href="classviennacl_1_1vector__base.html#a762d98e2f912fc534951f25555b6077f">viennacl::vector_base::begin</a></div><div class="ttdeci">iterator begin()</div><div class="ttdoc">Returns an iterator pointing to the beginning of the vector (STL like) </div><div class="ttdef"><b>Definition:</b> <a href="vector_8hpp_source.html#l00841">vector.hpp:841</a></div></div>
<div class="ttc" id="ocl_2handle_8hpp_html"><div class="ttname"><a href="ocl_2handle_8hpp.html">handle.hpp</a></div><div class="ttdoc">Implementation of a smart-pointer-like class for handling OpenCL handles. </div></div>
<div class="ttc" id="namespaceviennacl_1_1traits_html_a68603a1e1ca1bfb6d107831ba5096786"><div class="ttname"><a href="namespaceviennacl_1_1traits.html#a68603a1e1ca1bfb6d107831ba5096786">viennacl::traits::start</a></div><div class="ttdeci">result_of::size_type< T >::type start(T const &obj)</div><div class="ttdef"><b>Definition:</b> <a href="start_8hpp_source.html#l00044">start.hpp:44</a></div></div>
<div class="ttc" id="classviennacl_1_1vector__base_html"><div class="ttname"><a href="classviennacl_1_1vector__base.html">viennacl::vector_base</a></div><div class="ttdoc">Common base class for dense vectors, vector ranges, and vector slices. </div><div class="ttdef"><b>Definition:</b> <a href="vector__def_8hpp_source.html#l00104">vector_def.hpp:104</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_a366df7768bf52e820fee364334133db1"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a366df7768bf52e820fee364334133db1">viennacl::linalg::detail::householder_vector</a></div><div class="ttdeci">void householder_vector(VectorType &v, vcl_size_t start)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00097">qr-method-common.hpp:97</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a8e62fa452895c889f7bad9700300994d"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a8e62fa452895c889f7bad9700300994d">viennacl::linalg::SVD_HOUSEHOLDER_UPDATE_A_LEFT_KERNEL</a></div><div class="ttdeci">const std::string SVD_HOUSEHOLDER_UPDATE_A_LEFT_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00055">qr-method-common.hpp:55</a></div></div>
<div class="ttc" id="namespaceviennacl_html_a98a0afcc513111ffa0bd138f891930df"><div class="ttname"><a href="namespaceviennacl.html#a98a0afcc513111ffa0bd138f891930df">viennacl::vcl_size_t</a></div><div class="ttdeci">std::size_t vcl_size_t</div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00075">forwards.h:75</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_a0737c94fac817b0d0c8444463be7ceb7"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a0737c94fac817b0d0c8444463be7ceb7">viennacl::linalg::detail::pythag</a></div><div class="ttdeci">SCALARTYPE pythag(SCALARTYPE a, SCALARTYPE b)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00065">qr-method-common.hpp:65</a></div></div>
<div class="ttc" id="cuda_2matrix__operations_8hpp_html"><div class="ttname"><a href="cuda_2matrix__operations_8hpp.html">matrix_operations.hpp</a></div><div class="ttdoc">Implementations of dense matrix related operations, including matrix-vector products, using CUDA. </div></div>
<div class="ttc" id="namespaceviennacl_1_1ocl_html_a5f2022f653ea1cf364d20e3ff84dcada"><div class="ttname"><a href="namespaceviennacl_1_1ocl.html#a5f2022f653ea1cf364d20e3ff84dcada">viennacl::ocl::enqueue</a></div><div class="ttdeci">void enqueue(KernelType &k, viennacl::ocl::command_queue const &queue)</div><div class="ttdoc">Enqueues a kernel in the provided queue. </div><div class="ttdef"><b>Definition:</b> <a href="enqueue_8hpp_source.html#l00050">enqueue.hpp:50</a></div></div>
<div class="ttc" id="kernel_8hpp_html"><div class="ttname"><a href="kernel_8hpp.html">kernel.hpp</a></div><div class="ttdoc">Representation of an OpenCL kernel in ViennaCL. </div></div>
<div class="ttc" id="vector_8hpp_html"><div class="ttname"><a href="vector_8hpp.html">vector.hpp</a></div><div class="ttdoc">The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...</div></div>
<div class="ttc" id="classviennacl_1_1matrix__base_html_a15ec879f9e85a6fc209454b0bf137084"><div class="ttname"><a href="classviennacl_1_1matrix__base.html#a15ec879f9e85a6fc209454b0bf137084">viennacl::matrix_base::row_major</a></div><div class="ttdeci">bool row_major() const </div><div class="ttdef"><b>Definition:</b> <a href="matrix__def_8hpp_source.html#l00248">matrix_def.hpp:248</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a04a5ce1f7230502c9d1ee87dcc24c8b4"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a04a5ce1f7230502c9d1ee87dcc24c8b4">viennacl::linalg::SVD_HOUSEHOLDER_UPDATE_QL_KERNEL</a></div><div class="ttdeci">const std::string SVD_HOUSEHOLDER_UPDATE_QL_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00057">qr-method-common.hpp:57</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="namespaceviennacl_1_1linalg_html_a8477d884683dafc7e6ca04f2e02a1d1f"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a8477d884683dafc7e6ca04f2e02a1d1f">viennacl::linalg::SVD_FINAL_ITER_UPDATE_KERNEL</a></div><div class="ttdeci">const std::string SVD_FINAL_ITER_UPDATE_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00053">qr-method-common.hpp:53</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a349e6e78f95ef631a5f2736158766739"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a349e6e78f95ef631a5f2736158766739">viennacl::linalg::SVD_MATRIX_TRANSPOSE_KERNEL</a></div><div class="ttdeci">const std::string SVD_MATRIX_TRANSPOSE_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00050">qr-method-common.hpp:50</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__base_html_a00e30fdc61081a73995b54175a91d220"><div class="ttname"><a href="classviennacl_1_1matrix__base.html#a00e30fdc61081a73995b54175a91d220">viennacl::matrix_base::internal_size2</a></div><div class="ttdeci">size_type internal_size2() const </div><div class="ttdoc">Returns the internal number of columns. Usually required for launching OpenCL kernels only...</div><div class="ttdef"><b>Definition:</b> <a href="matrix__def_8hpp_source.html#l00240">matrix_def.hpp:240</a></div></div>
<div class="ttc" id="classviennacl_1_1matrix__base_html_a82018b4e169973bdfb9d3be68ceb5be0"><div class="ttname"><a href="classviennacl_1_1matrix__base.html#a82018b4e169973bdfb9d3be68ceb5be0">viennacl::matrix_base::internal_size1</a></div><div class="ttdeci">size_type internal_size1() const </div><div class="ttdoc">Returns the internal number of rows. Usually required for launching OpenCL kernels only...</div><div class="ttdef"><b>Definition:</b> <a href="matrix__def_8hpp_source.html#l00238">matrix_def.hpp:238</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_ae059c4a3216ff6bdc44110377ce34e88"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#ae059c4a3216ff6bdc44110377ce34e88">viennacl::linalg::detail::norm_lcl</a></div><div class="ttdeci">VectorType::value_type norm_lcl(VectorType const &x, vcl_size_t size)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00078">qr-method-common.hpp:78</a></div></div>
<div class="ttc" id="result__of_8hpp_html"><div class="ttname"><a href="result__of_8hpp.html">result_of.hpp</a></div><div class="ttdoc">A collection of compile time type deductions. </div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_html_a15bedacb51ec9b425baf52ee7f68bdba"><div class="ttname"><a href="namespaceviennacl_1_1linalg.html#a15bedacb51ec9b425baf52ee7f68bdba">viennacl::linalg::SVD_GIVENS_PREV_KERNEL</a></div><div class="ttdeci">const std::string SVD_GIVENS_PREV_KERNEL</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00052">qr-method-common.hpp:52</a></div></div>
<div class="ttc" id="namespaceviennacl_html_afad612e26a65487b7e747c18c678ffd9"><div class="ttname"><a href="namespaceviennacl.html#afad612e26a65487b7e747c18c678ffd9">viennacl::vcl_ptrdiff_t</a></div><div class="ttdeci">std::ptrdiff_t vcl_ptrdiff_t</div><div class="ttdef"><b>Definition:</b> <a href="forwards_8h_source.html#l00076">forwards.h:76</a></div></div>
<div class="ttc" id="namespaceviennacl_1_1linalg_1_1detail_html_a7bfa1c9cb0afb4ca0f975e690f4a8f07"><div class="ttname"><a href="namespaceviennacl_1_1linalg_1_1detail.html#a7bfa1c9cb0afb4ca0f975e690f4a8f07">viennacl::linalg::detail::sign</a></div><div class="ttdeci">SCALARTYPE sign(SCALARTYPE val)</div><div class="ttdef"><b>Definition:</b> <a href="qr-method-common_8hpp_source.html#l00071">qr-method-common.hpp:71</a></div></div>
<div class="ttc" id="namespaceviennacl_html_aaa5c8726b45bc89a523ca2fa8c42107a"><div class="ttname"><a href="namespaceviennacl.html#aaa5c8726b45bc89a523ca2fa8c42107a">viennacl::fast_copy</a></div><div class="ttdeci">void fast_copy(const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_begin, const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_end, CPU_ITERATOR cpu_begin)</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_c82e3d11dd171600f4a6e0cab1ec1e0d.html">viennacl</a></li><li class="navelem"><a class="el" href="dir_63cde087767c4ed65c7901ffc6e293fe.html">linalg</a></li><li class="navelem"><a class="el" href="qr-method-common_8hpp.html">qr-method-common.hpp</a></li>
<li class="footer">Generated on Wed Jan 20 2016 22:32:43 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>
|