GNU Octave  4.2.1
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
svd.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 1996-2017 John W. Eaton
4 
5 This file is part of Octave.
6 
7 Octave is free software; you can redistribute it and/or modify it
8 under the terms of the GNU General Public License as published by the
9 Free Software Foundation; either version 3 of the License, or (at your
10 option) any later version.
11 
12 Octave is distributed in the hope that it will be useful, but WITHOUT
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
15 for more details.
16 
17 You should have received a copy of the GNU General Public License
18 along with Octave; see the file COPYING. If not, see
19 <http://www.gnu.org/licenses/>.
20 
21 */
22 
23 #if defined (HAVE_CONFIG_H)
24 # include "config.h"
25 #endif
26 
27 #include "svd.h"
28 
29 #include "defun.h"
30 #include "error.h"
31 #include "errwarn.h"
32 #include "ovl.h"
33 #include "pr-output.h"
34 #include "utils.h"
35 #include "variables.h"
36 
37 static std::string Vsvd_driver = "gesvd";
38 
39 template <typename T>
40 static typename octave::math::svd<T>::Type
42 {
43  return ((nargout == 0 || nargout == 1)
45  : ((nargin == 2)
48 }
49 
50 template <typename T>
51 static typename octave::math::svd<T>::Driver
52 svd_driver (void)
53 {
54  return (Vsvd_driver == "gesvd"
57 }
58 
60  doc: /* -*- texinfo -*-
61 @deftypefn {} {@var{s} =} svd (@var{A})
62 @deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A})
63 @deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A}, @var{econ})
64 @cindex singular value decomposition
65 Compute the singular value decomposition of @var{A}
66 @tex
67 $$
68  A = U S V^{\dagger}
69 $$
70 @end tex
71 @ifnottex
72 
73 @example
74 A = U*S*V'
75 @end example
76 
77 @end ifnottex
78 
79 The function @code{svd} normally returns only the vector of singular values.
80 When called with three return values, it computes
81 @tex
82 $U$, $S$, and $V$.
83 @end tex
84 @ifnottex
85 @var{U}, @var{S}, and @var{V}.
86 @end ifnottex
87 For example,
88 
89 @example
90 svd (hilb (3))
91 @end example
92 
93 @noindent
94 returns
95 
96 @example
97 @group
98 ans =
99 
100  1.4083189
101  0.1223271
102  0.0026873
103 @end group
104 @end example
105 
106 @noindent
107 and
108 
109 @example
110 [u, s, v] = svd (hilb (3))
111 @end example
112 
113 @noindent
114 returns
115 
116 @example
117 @group
118 u =
119 
120  -0.82704 0.54745 0.12766
121  -0.45986 -0.52829 -0.71375
122  -0.32330 -0.64901 0.68867
123 
124 s =
125 
126  1.40832 0.00000 0.00000
127  0.00000 0.12233 0.00000
128  0.00000 0.00000 0.00269
129 
130 v =
131 
132  -0.82704 0.54745 0.12766
133  -0.45986 -0.52829 -0.71375
134  -0.32330 -0.64901 0.68867
135 @end group
136 @end example
137 
138 If given a second argument, @code{svd} returns an economy-sized
139 decomposition, eliminating the unnecessary rows or columns of @var{U} or
140 @var{V}.
141 @seealso{svd_driver, svds, eig, lu, chol, hess, qr, qz}
142 @end deftypefn */)
143 {
144  int nargin = args.length ();
145 
146  if (nargin < 1 || nargin > 2 || nargout == 2 || nargout > 3)
147  print_usage ();
148 
149  octave_value arg = args(0);
150 
151  if (arg.ndims () != 2)
152  error ("svd: A must be a 2-D matrix");
153 
155 
156  bool isfloat = arg.is_single_type ();
157 
158  if (isfloat)
159  {
160  if (arg.is_real_type ())
161  {
163 
164  if (tmp.any_element_is_inf_or_nan ())
165  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
166 
168  (tmp, svd_type<FloatMatrix> (nargin, nargout),
169  svd_driver<FloatMatrix> ());
170 
171  FloatDiagMatrix sigma = result.singular_values ();
172 
173  if (nargout == 0 || nargout == 1)
174  retval(0) = sigma.extract_diag ();
175  else
176  retval = ovl (result.left_singular_matrix (),
177  sigma,
178  result.right_singular_matrix ());
179  }
180  else if (arg.is_complex_type ())
181  {
183 
184  if (ctmp.any_element_is_inf_or_nan ())
185  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
186 
188  (ctmp, svd_type<FloatComplexMatrix> (nargin, nargout),
189  svd_driver<FloatComplexMatrix> ());
190 
191  FloatDiagMatrix sigma = result.singular_values ();
192 
193  if (nargout == 0 || nargout == 1)
194  retval(0) = sigma.extract_diag ();
195  else
196  retval = ovl (result.left_singular_matrix (),
197  sigma,
198  result.right_singular_matrix ());
199  }
200  }
201  else
202  {
203  if (arg.is_real_type ())
204  {
205  Matrix tmp = arg.matrix_value ();
206 
207  if (tmp.any_element_is_inf_or_nan ())
208  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
209 
211  (tmp, svd_type<Matrix> (nargin, nargout),
212  svd_driver<Matrix> ());
213 
214  DiagMatrix sigma = result.singular_values ();
215 
216  if (nargout == 0 || nargout == 1)
217  retval(0) = sigma.extract_diag ();
218  else
219  retval = ovl (result.left_singular_matrix (),
220  sigma,
221  result.right_singular_matrix ());
222  }
223  else if (arg.is_complex_type ())
224  {
225  ComplexMatrix ctmp = arg.complex_matrix_value ();
226 
227  if (ctmp.any_element_is_inf_or_nan ())
228  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
229 
231  (ctmp, svd_type<ComplexMatrix> (nargin, nargout),
232  svd_driver<ComplexMatrix> ());
233 
234  DiagMatrix sigma = result.singular_values ();
235 
236  if (nargout == 0 || nargout == 1)
237  retval(0) = sigma.extract_diag ();
238  else
239  retval = ovl (result.left_singular_matrix (),
240  sigma,
241  result.right_singular_matrix ());
242  }
243  else
244  err_wrong_type_arg ("svd", arg);
245  }
246 
247  return retval;
248 }
249 
250 /*
251 %!assert (svd ([1, 2; 2, 1]), [3; 1], sqrt (eps))
252 
253 %!test
254 a = [1, 2; 3, 4] + [5, 6; 7, 8]*i;
255 [u,s,v] = svd (a);
256 assert (a, u * s * v', 128 * eps);
257 
258 
259 %!test
260 %! [u, s, v] = svd ([1, 2; 2, 1]);
261 %! x = 1 / sqrt (2);
262 %! assert (u, [-x, -x; -x, x], sqrt (eps));
263 %! assert (s, [3, 0; 0, 1], sqrt (eps));
264 %! assert (v, [-x, x; -x, -x], sqrt (eps));
265 
266 %!test
267 %! a = [1, 2, 3; 4, 5, 6];
268 %! [u, s, v] = svd (a);
269 %! assert (u * s * v', a, sqrt (eps));
270 
271 %!test
272 %! a = [1, 2; 3, 4; 5, 6];
273 %! [u, s, v] = svd (a);
274 %! assert (u * s * v', a, sqrt (eps));
275 
276 %!test
277 %! a = [1, 2, 3; 4, 5, 6];
278 %! [u, s, v] = svd (a, 1);
279 %! assert (u * s * v', a, sqrt (eps));
280 
281 %!test
282 %! a = [1, 2; 3, 4; 5, 6];
283 %! [u, s, v] = svd (a, 1);
284 %! assert (u * s * v', a, sqrt (eps));
285 
286 %!assert (svd (single ([1, 2; 2, 1])), single ([3; 1]), sqrt (eps ("single")))
287 
288 %!test
289 %! [u, s, v] = svd (single ([1, 2; 2, 1]));
290 %! x = single (1 / sqrt (2));
291 %! assert (u, [-x, -x; -x, x], sqrt (eps ("single")));
292 %! assert (s, single ([3, 0; 0, 1]), sqrt (eps ("single")));
293 %! assert (v, [-x, x; -x, -x], sqrt (eps ("single")));
294 
295 %!test
296 %! a = single ([1, 2, 3; 4, 5, 6]);
297 %! [u, s, v] = svd (a);
298 %! assert (u * s * v', a, sqrt (eps ("single")));
299 
300 %!test
301 %! a = single ([1, 2; 3, 4; 5, 6]);
302 %! [u, s, v] = svd (a);
303 %! assert (u * s * v', a, sqrt (eps ("single")));
304 
305 %!test
306 %! a = single ([1, 2, 3; 4, 5, 6]);
307 %! [u, s, v] = svd (a, 1);
308 %! assert (u * s * v', a, sqrt (eps ("single")));
309 
310 %!test
311 %! a = single ([1, 2; 3, 4; 5, 6]);
312 %! [u, s, v] = svd (a, 1);
313 %! assert (u * s * v', a, sqrt (eps ("single")));
314 
315 %!test
316 %! a = zeros (0, 5);
317 %! [u, s, v] = svd (a);
318 %! assert (size (u), [0, 0]);
319 %! assert (size (s), [0, 5]);
320 %! assert (size (v), [5, 5]);
321 
322 %!test
323 %! a = zeros (5, 0);
324 %! [u, s, v] = svd (a, 1);
325 %! assert (size (u), [5, 0]);
326 %! assert (size (s), [0, 0]);
327 %! assert (size (v), [0, 0]);
328 
329 %!error svd ()
330 %!error svd ([1, 2; 4, 5], 2, 3)
331 %!error [u, v] = svd ([1, 2; 3, 4])
332 */
333 
335  doc: /* -*- texinfo -*-
336 @deftypefn {} {@var{val} =} svd_driver ()
337 @deftypefnx {} {@var{old_val} =} svd_driver (@var{new_val})
338 @deftypefnx {} {} svd_driver (@var{new_val}, "local")
339 Query or set the underlying @sc{lapack} driver used by @code{svd}.
340 
341 Currently recognized values are @qcode{"gesvd"} and @qcode{"gesdd"}.
342 The default is @qcode{"gesvd"}.
343 
344 When called from inside a function with the @qcode{"local"} option, the
345 variable is changed locally for the function and any subroutines it calls.
346 The original variable value is restored when exiting the function.
347 @seealso{svd}
348 @end deftypefn */)
349 {
350  static const char *driver_names[] = { "gesvd", "gesdd", 0 };
351 
352  return SET_INTERNAL_VARIABLE_CHOICES (svd_driver, driver_names);
353 }
354 
355 /*
356 %!test
357 %! A = [1+1i, 1-1i, 0; 0, 2, 0; 1i, 1i, 1+2i];
358 %! old_driver = svd_driver ("gesvd");
359 %! [U1, S1, V1] = svd (A);
360 %! svd_driver ("gesdd");
361 %! [U2, S2, V2] = svd (A);
362 %! assert (U1, U2, 5*eps);
363 %! assert (S1, S2, 5*eps);
364 %! assert (V1, V2, 5*eps);
365 %! svd_driver (old_driver);
366 */
bool is_real_type(void) const
Definition: ov.h:667
int ndims(void) const
Definition: ov.h:495
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:809
OCTAVE_EXPORT octave_value_list isa nd deftypefn *return ovl(args(0).is_integer_type())
OCTINTERP_API void print_usage(void)
Definition: defun.cc:52
Definition: mx-defs.h:73
#define DEFUN(name, args_name, nargout_name, doc)
Definition: defun.h:46
void error(const char *fmt,...)
Definition: error.cc:570
T left_singular_matrix(void) const
Definition: svd.cc:49
bool any_element_is_inf_or_nan(void) const
Definition: fCNDArray.cc:502
ColumnVector extract_diag(octave_idx_type k=0) const
Definition: dDiagMatrix.h:104
octave_value arg
Definition: pr-output.cc:3440
T right_singular_matrix(void) const
Definition: svd.cc:60
JNIEnv void * args
Definition: ov-java.cc:67
static std::string Vsvd_driver
Definition: svd.cc:37
OCTAVE_EXPORT octave_value_list return the number of command line arguments passed to Octave If called with the optional argument the function xample nargout(@histc)
Definition: ov-usr-fcn.cc:935
OCTAVE_EXPORT octave_value_list isfloat
Definition: data.cc:3327
int nargin
Definition: graphics.cc:10115
bool is_complex_type(void) const
Definition: ov.h:670
DM_T singular_values(void) const
Definition: svd.h:86
double tmp
Definition: data.cc:6300
octave_value retval
Definition: data.cc:6294
#define SET_INTERNAL_VARIABLE_CHOICES(NM, CHOICES)
Definition: variables.h:136
FloatColumnVector extract_diag(octave_idx_type k=0) const
Definition: fDiagMatrix.h:105
Definition: dMatrix.h:37
bool any_element_is_inf_or_nan(void) const
Definition: CNDArray.cc:508
Matrix matrix_value(bool frc_str_conv=false) const
Definition: ov.h:787
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:156
With real return the complex result
Definition: data.cc:3375
static octave::math::svd< T >::Type svd_type(int nargin, int nargout)
Definition: svd.cc:41
bool any_element_is_inf_or_nan(void) const
Definition: fNDArray.cc:513
static octave::math::svd< T >::Driver svd_driver(void)
Definition: svd.cc:52
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:805
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:790
bool any_element_is_inf_or_nan(void) const
Definition: dNDArray.cc:561
bool is_single_type(void) const
Definition: ov.h:627
If this string is the system will ring the terminal sometimes it is useful to be able to print the original representation of the string
Definition: utils.cc:854