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