GNU Octave  4.4.1
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
svd.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 1996-2018 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
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11 
12 Octave is distributed in the hope that it will be useful, but
13 WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License 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 <https://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 = "gesdd";
38 
39 template <typename T>
40 static typename octave::math::svd<T>::Type
41 svd_type (int nargin, int nargout, const octave_value_list & args, const T & A)
42 {
43  if (nargout == 0 || nargout == 1)
45  else if (nargin == 1)
47  else if (! args(1).is_real_scalar ())
49  else
50  {
51  if (A.rows () > A.columns ())
53  else
55  }
56 }
57 
58 template <typename T>
59 static typename octave::math::svd<T>::Driver
60 svd_driver (void)
61 {
62  return (Vsvd_driver == "gesvd"
65 }
66 
67 DEFUN (svd, args, nargout,
68  classes: double single
69  doc: /* -*- texinfo -*-
70 @deftypefn {} {@var{s} =} svd (@var{A})
71 @deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A})
72 @deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A}, "econ")
73 @deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A}, 0)
74 @cindex singular value decomposition
75 Compute the singular value decomposition of @var{A}.
76 
77 The singular value decomposition is defined by the relation
78 
79 @tex
80 $$
81  A = U S V^{\dagger}
82 $$
83 @end tex
84 @ifnottex
85 
86 @example
87 A = U*S*V'
88 @end example
89 
90 @end ifnottex
91 
92 The function @code{svd} normally returns only the vector of singular values.
93 When called with three return values, it computes
94 @tex
95 $U$, $S$, and $V$.
96 @end tex
97 @ifnottex
98 @var{U}, @var{S}, and @var{V}.
99 @end ifnottex
100 For example,
101 
102 @example
103 svd (hilb (3))
104 @end example
105 
106 @noindent
107 returns
108 
109 @example
110 @group
111 ans =
112 
113  1.4083189
114  0.1223271
115  0.0026873
116 @end group
117 @end example
118 
119 @noindent
120 and
121 
122 @example
123 [u, s, v] = svd (hilb (3))
124 @end example
125 
126 @noindent
127 returns
128 
129 @example
130 @group
131 u =
132 
133  -0.82704 0.54745 0.12766
134  -0.45986 -0.52829 -0.71375
135  -0.32330 -0.64901 0.68867
136 
137 s =
138 
139  1.40832 0.00000 0.00000
140  0.00000 0.12233 0.00000
141  0.00000 0.00000 0.00269
142 
143 v =
144 
145  -0.82704 0.54745 0.12766
146  -0.45986 -0.52829 -0.71375
147  -0.32330 -0.64901 0.68867
148 @end group
149 @end example
150 
151 When given a second argument that is not 0, @code{svd} returns an economy-sized
152 decomposition, eliminating the unnecessary rows or columns of @var{U} or
153 @var{V}.
154 
155 If the second argument is exactly 0, then the choice of decomposition is based
156 on the matrix @var{A}. If @var{A} has more rows than columns then an
157 economy-sized decomposition is returned, otherwise a regular decomposition
158 is calculated.
159 
160 Algorithm Notes: When calculating the full decomposition (left and right
161 singular matrices in addition to singular values) there is a choice of two
162 routines in @sc{lapack}. The default routine used by Octave is @code{gesdd}
163 which is 5X faster than the alternative @code{gesvd}, but may use more memory
164 and may be less accurate for some matrices. See the documentation for
165 @code{svd_driver} for more information.
166 @seealso{svd_driver, svds, eig, lu, chol, hess, qr, qz}
167 @end deftypefn */)
168 {
169  int nargin = args.length ();
170 
171  if (nargin < 1 || nargin > 2 || nargout == 2 || nargout > 3)
172  print_usage ();
173 
174  octave_value arg = args(0);
175 
176  if (arg.ndims () != 2)
177  error ("svd: A must be a 2-D matrix");
178 
180 
181  bool isfloat = arg.is_single_type ();
182 
183  if (isfloat)
184  {
185  if (arg.isreal ())
186  {
188 
189  if (tmp.any_element_is_inf_or_nan ())
190  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
191 
193  (tmp,
194  svd_type<FloatMatrix> (nargin, nargout, args, tmp),
195  svd_driver<FloatMatrix> ());
196 
197  FloatDiagMatrix sigma = result.singular_values ();
198 
199  if (nargout == 0 || nargout == 1)
200  retval(0) = sigma.extract_diag ();
201  else
202  retval = ovl (result.left_singular_matrix (),
203  sigma,
204  result.right_singular_matrix ());
205  }
206  else if (arg.iscomplex ())
207  {
209 
210  if (ctmp.any_element_is_inf_or_nan ())
211  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
212 
214  (ctmp,
215  svd_type<FloatComplexMatrix> (nargin, nargout, args, ctmp),
216  svd_driver<FloatComplexMatrix> ());
217 
218  FloatDiagMatrix sigma = result.singular_values ();
219 
220  if (nargout == 0 || nargout == 1)
221  retval(0) = sigma.extract_diag ();
222  else
223  retval = ovl (result.left_singular_matrix (),
224  sigma,
225  result.right_singular_matrix ());
226  }
227  }
228  else
229  {
230  if (arg.isreal ())
231  {
232  Matrix tmp = arg.matrix_value ();
233 
234  if (tmp.any_element_is_inf_or_nan ())
235  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
236 
238  (tmp,
239  svd_type<Matrix> (nargin, nargout, args, tmp),
240  svd_driver<Matrix> ());
241 
242  DiagMatrix sigma = result.singular_values ();
243 
244  if (nargout == 0 || nargout == 1)
245  retval(0) = sigma.extract_diag ();
246  else
247  retval = ovl (result.left_singular_matrix (),
248  sigma,
249  result.right_singular_matrix ());
250  }
251  else if (arg.iscomplex ())
252  {
254 
255  if (ctmp.any_element_is_inf_or_nan ())
256  error ("svd: cannot take SVD of matrix containing Inf or NaN values");
257 
259  (ctmp,
260  svd_type<ComplexMatrix> (nargin, nargout, args, ctmp),
261  svd_driver<ComplexMatrix> ());
262 
263  DiagMatrix sigma = result.singular_values ();
264 
265  if (nargout == 0 || nargout == 1)
266  retval(0) = sigma.extract_diag ();
267  else
268  retval = ovl (result.left_singular_matrix (),
269  sigma,
270  result.right_singular_matrix ());
271  }
272  else
273  err_wrong_type_arg ("svd", arg);
274  }
275 
276  return retval;
277 }
278 
279 /*
280 %!assert (svd ([1, 2; 2, 1]), [3; 1], sqrt (eps))
281 
282 %!test
283 a = [1, 2; 3, 4] + [5, 6; 7, 8]*i;
284 [u,s,v] = svd (a);
285 assert (a, u * s * v', 128 * eps);
286 
287 %!test
288 %! [u, s, v] = svd ([1, 2; 2, 1]);
289 %! x = 1 / sqrt (2);
290 %! assert (u, [-x, -x; -x, x], sqrt (eps));
291 %! assert (s, [3, 0; 0, 1], sqrt (eps));
292 %! assert (v, [-x, x; -x, -x], sqrt (eps));
293 
294 %!test
295 %! a = [1, 2, 3; 4, 5, 6];
296 %! [u, s, v] = svd (a);
297 %! assert (u * s * v', a, sqrt (eps));
298 
299 %!test
300 %! a = [1, 2; 3, 4; 5, 6];
301 %! [u, s, v] = svd (a);
302 %! assert (u * s * v', a, sqrt (eps));
303 
304 %!test
305 %! a = [1, 2, 3; 4, 5, 6];
306 %! [u, s, v] = svd (a, 1);
307 %! assert (u * s * v', a, sqrt (eps));
308 
309 %!test
310 %! a = [1, 2; 3, 4; 5, 6];
311 %! [u, s, v] = svd (a, 1);
312 %! assert (u * s * v', a, sqrt (eps));
313 
314 %!assert (svd (single ([1, 2; 2, 1])), single ([3; 1]), sqrt (eps ("single")))
315 
316 %!test
317 %! [u, s, v] = svd (single ([1, 2; 2, 1]));
318 %! x = single (1 / sqrt (2));
319 %! assert (u, [-x, -x; -x, x], sqrt (eps ("single")));
320 %! assert (s, single ([3, 0; 0, 1]), sqrt (eps ("single")));
321 %! assert (v, [-x, x; -x, -x], sqrt (eps ("single")));
322 
323 %!test
324 %! a = single ([1, 2, 3; 4, 5, 6]);
325 %! [u, s, v] = svd (a);
326 %! assert (u * s * v', a, sqrt (eps ("single")));
327 
328 %!test
329 %! a = single ([1, 2; 3, 4; 5, 6]);
330 %! [u, s, v] = svd (a);
331 %! assert (u * s * v', a, sqrt (eps ("single")));
332 
333 %!test
334 %! a = single ([1, 2, 3; 4, 5, 6]);
335 %! [u, s, v] = svd (a, 1);
336 %! assert (u * s * v', a, sqrt (eps ("single")));
337 
338 %!test
339 %! a = single ([1, 2; 3, 4; 5, 6]);
340 %! [u, s, v] = svd (a, 1);
341 %! assert (u * s * v', a, sqrt (eps ("single")));
342 
343 %!test
344 %! a = zeros (0, 5);
345 %! [u, s, v] = svd (a);
346 %! assert (size (u), [0, 0]);
347 %! assert (size (s), [0, 5]);
348 %! assert (size (v), [5, 5]);
349 
350 %!test
351 %! a = zeros (5, 0);
352 %! [u, s, v] = svd (a, 1);
353 %! assert (size (u), [5, 0]);
354 %! assert (size (s), [0, 0]);
355 %! assert (size (v), [0, 0]);
356 
357 %!test <*49309>
358 %! [~,~,v] = svd ([1, 1, 1], 0);
359 %! assert (size (v), [3 3]);
360 %! [~,~,v] = svd ([1, 1, 1], "econ");
361 %! assert (size (v), [3 1]);
362 
363 %!error svd ()
364 %!error svd ([1, 2; 4, 5], 2, 3)
365 %!error [u, v] = svd ([1, 2; 3, 4])
366 */
367 
368 DEFUN (svd_driver, args, nargout,
369  doc: /* -*- texinfo -*-
370 @deftypefn {} {@var{val} =} svd_driver ()
371 @deftypefnx {} {@var{old_val} =} svd_driver (@var{new_val})
372 @deftypefnx {} {} svd_driver (@var{new_val}, "local")
373 Query or set the underlying @sc{lapack} driver used by @code{svd}.
374 
375 Currently recognized values are @qcode{"gesdd"} and @qcode{"gesvd"}.
376 The default is @qcode{"gesdd"}.
377 
378 When called from inside a function with the @qcode{"local"} option, the
379 variable is changed locally for the function and any subroutines it calls.
380 The original variable value is restored when exiting the function.
381 
382 Algorithm Notes: The @sc{lapack} library provides two routines for calculating
383 the full singular value decomposition (left and right singular matrices as
384 well as singular values). When calculating just the singular values the
385 following discussion is not relevant.
386 
387 The default routine use by Octave is the newer @code{gesdd} which is based on a
388 Divide-and-Conquer algorithm that is 5X faster than the alternative
389 @code{gesvd}, which is based on QR factorization. However, the new algorithm
390 can use significantly more memory. For an @nospell{MxN} input matrix the
391 memory usage is of order O(min(M,N) ^ 2), whereas the alternative is of order
392 O(max(M,N)). In general, modern computers have abundant memory so Octave has
393 chosen to prioritize speed.
394 
395 In addition, there have been instances in the past where some input matrices
396 were not accurately decomposed by @code{gesdd}. This appears to have been
397 resolved with modern versions of @sc{lapack}. However, if certainty is
398 required the accuracy of the decomposition can always be tested after the fact
399 with
400 
401 @example
402 @group
403 [@var{u}, @var{s}, @var{v}] = svd (@var{x});
404 norm (@var{x} - @var{u}*@var{s}*@var{v'}, "fro")
405 @end group
406 @end example
407 
408 @seealso{svd}
409 @end deftypefn */)
410 {
411  static const char *driver_names[] = { "gesvd", "gesdd", nullptr };
412 
413  return SET_INTERNAL_VARIABLE_CHOICES (svd_driver, driver_names);
414 }
415 
416 /*
417 %!test
418 %! A = [1+1i, 1-1i, 0; 0, 2, 0; 1i, 1i, 1+2i];
419 %! old_driver = svd_driver ("gesvd");
420 %! [U1, S1, V1] = svd (A);
421 %! svd_driver ("gesdd");
422 %! [U2, S2, V2] = svd (A);
423 %! svd_driver (old_driver);
424 %! assert (U1, U2, 5*eps);
425 %! assert (S1, S2, 5*eps);
426 %! assert (V1, V2, 5*eps);
427 */
FloatColumnVector extract_diag(octave_idx_type k=0) const
Definition: fDiagMatrix.h:105
ColumnVector extract_diag(octave_idx_type k=0) const
Definition: dDiagMatrix.h:104
OCTINTERP_API void print_usage(void)
Definition: defun.cc:54
#define DEFUN(name, args_name, nargout_name, doc)
Macro to define a builtin function.
Definition: defun.h:53
void error(const char *fmt,...)
Definition: error.cc:578
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:837
int ndims(void) const
Definition: ov.h:478
octave_value arg
Definition: pr-output.cc:3244
bool any_element_is_inf_or_nan(void) const
Definition: fCNDArray.cc:500
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:856
static std::string Vsvd_driver
Definition: svd.cc:37
F77_RET_T const F77_INT F77_CMPLX * A
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:997
OCTAVE_EXPORT octave_value_list isfloat
Definition: data.cc:3212
bool is_single_type(void) const
Definition: ov.h:651
double tmp
Definition: data.cc:6252
octave_value retval
Definition: data.cc:6246
#define SET_INTERNAL_VARIABLE_CHOICES(NM, CHOICES)
Definition: variables.h:119
Definition: dMatrix.h:36
bool any_element_is_inf_or_nan(void) const
Definition: CNDArray.cc:506
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:162
With real return the complex result
Definition: data.cc:3260
bool isreal(void) const
Definition: ov.h:703
OCTAVE_EXPORT octave_value_list isa nd deftypefn *return ovl(args(0).isinteger())
static octave::math::svd< T >::Driver svd_driver(void)
Definition: svd.cc:60
static octave::math::svd< T >::Type svd_type(int nargin, int nargout, const octave_value_list &args, const T &A)
Definition: svd.cc:41
args.length() nargin
Definition: file-io.cc:589
bool iscomplex(void) const
Definition: ov.h:710
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition: ov.h:852
Definition: mx-defs.h:75
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:888
Matrix matrix_value(bool frc_str_conv=false) const
Definition: ov.h:834