GNU Octave 11.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
svd.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1996-2026 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
41
42static std::string Vsvd_driver = "gesvd";
43
44template <typename T>
45static typename math::svd<T>::Type
46svd_type (int nargin, int nargout, const octave_value_list& args, const T& A)
47{
48 if (nargout == 0 || nargout == 1)
49 return math::svd<T>::Type::sigma_only;
50 else if (nargin == 1)
51 return math::svd<T>::Type::std;
52 else if (! args(1).is_real_scalar ())
53 return math::svd<T>::Type::economy;
54 else
55 {
56 if (A.rows () > A.columns ())
57 return math::svd<T>::Type::economy;
58 else
59 return math::svd<T>::Type::std;
60 }
61}
62
63template <typename T>
64static typename math::svd<T>::Driver
65svd_driver ()
66{
67 if (Vsvd_driver == "gejsv")
68 return math::svd<T>::Driver::GEJSV;
69 else if (Vsvd_driver == "gesdd")
70 return math::svd<T>::Driver::GESDD;
71 else
72 return math::svd<T>::Driver::GESVD; // default
73}
74
75DEFUN (svd, args, nargout,
76 classes: double single
77 doc: /* -*- texinfo -*-
78@deftypefn {} {@var{s} =} svd (@var{A})
79@deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A})
80@deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A}, "econ")
81@deftypefnx {} {[@var{U}, @var{S}, @var{V}] =} svd (@var{A}, 0)
82@cindex singular value decomposition
83Compute the singular value decomposition of @var{A}.
84
85The singular value decomposition is defined by the relation
86@tex
87$$
88 A = U S V^{\dagger}
89$$
90@end tex
91@ifnottex
92
93@example
94A = U*S*V'
95@end example
96
97@end ifnottex
98
99The function @code{svd} normally returns only the vector of singular values.
100When called with three return values, it computes
101@tex
102$U$, $S$, and $V$.
103@end tex
104@ifnottex
105@var{U}, @var{S}, and @var{V}.
106@end ifnottex
107For example,
108
109@example
110svd (hilb (3))
111@end example
112
113@noindent
114returns
115
116@example
117@group
118ans =
119
120 1.4083189
121 0.1223271
122 0.0026873
123@end group
124@end example
125
126@noindent
127and
128
129@example
130[u, s, v] = svd (hilb (3))
131@end example
132
133@noindent
134returns
135
136@example
137@group
138u =
139
140 -0.82704 0.54745 0.12766
141 -0.45986 -0.52829 -0.71375
142 -0.32330 -0.64901 0.68867
143
144s =
145
146 1.40832 0.00000 0.00000
147 0.00000 0.12233 0.00000
148 0.00000 0.00000 0.00269
149
150v =
151
152 -0.82704 0.54745 0.12766
153 -0.45986 -0.52829 -0.71375
154 -0.32330 -0.64901 0.68867
155@end group
156@end example
157
158When given a second argument that is not 0, @code{svd} returns an economy-sized
159decomposition, eliminating the unnecessary rows or columns of @var{U} or
160@var{V}.
161
162If the second argument is exactly 0, then the choice of decomposition is based
163on the matrix @var{A}. If @var{A} has more rows than columns then an
164economy-sized decomposition is returned, otherwise a regular decomposition
165is calculated.
166
167Algorithm Notes: When calculating the full decomposition (left and right
168singular matrices in addition to singular values) there is a choice of two
169routines in @sc{lapack}. The default routine used by Octave is @code{gesvd}.
170The alternative is @code{gesdd} which is 5X faster, but may use more memory
171and may be inaccurate for some input matrices. There is a third routine
172@code{gejsv}, suitable for better accuracy at extreme scale. See the
173documentation for @code{svd_driver} for more information on choosing a driver.
174@seealso{svd_driver, svds, eig, lu, chol, hess, qr, qz}
175@end deftypefn */)
176{
177 int nargin = args.length ();
178
179 if (nargin < 1 || nargin > 2 || nargout > 3)
180 print_usage ();
181
182 octave_value arg = args(0);
183
184 if (arg.ndims () != 2)
185 error ("svd: A must be a 2-D matrix");
186
187 octave_value_list retval;
188
189 bool isfloat = arg.is_single_type ();
190
191 if (isfloat)
192 {
193 if (arg.isreal ())
194 {
195 FloatMatrix tmp = arg.float_matrix_value ();
196
197 if (tmp.any_element_is_inf_or_nan ())
198 error ("svd: cannot take SVD of matrix containing Inf or NaN values");
199
200 math::svd<FloatMatrix> result
201 (tmp,
202 svd_type<FloatMatrix> (nargin, nargout, args, tmp),
203 svd_driver<FloatMatrix> ());
204
205 FloatDiagMatrix sigma = result.singular_values ();
206
207 if (nargout == 0 || nargout == 1)
208 retval(0) = sigma.extract_diag ();
209 else if (nargout == 2)
210 retval = ovl (result.left_singular_matrix (),
211 sigma);
212 else
213 retval = ovl (result.left_singular_matrix (),
214 sigma,
215 result.right_singular_matrix ());
216 }
217 else if (arg.iscomplex ())
218 {
220
221 if (ctmp.any_element_is_inf_or_nan ())
222 error ("svd: cannot take SVD of matrix containing Inf or NaN values");
223
224 math::svd<FloatComplexMatrix> result
225 (ctmp,
226 svd_type<FloatComplexMatrix> (nargin, nargout, args, ctmp),
227 svd_driver<FloatComplexMatrix> ());
228
229 FloatDiagMatrix sigma = result.singular_values ();
230
231 if (nargout == 0 || nargout == 1)
232 retval(0) = sigma.extract_diag ();
233 else if (nargout == 2)
234 retval = ovl (result.left_singular_matrix (),
235 sigma);
236 else
237 retval = ovl (result.left_singular_matrix (),
238 sigma,
239 result.right_singular_matrix ());
240 }
241 }
242 else
243 {
244 if (arg.isreal ())
245 {
246 Matrix tmp = arg.matrix_value ();
247
248 if (tmp.any_element_is_inf_or_nan ())
249 error ("svd: cannot take SVD of matrix containing Inf or NaN values");
250
251 math::svd<Matrix> result
252 (tmp,
253 svd_type<Matrix> (nargin, nargout, args, tmp),
254 svd_driver<Matrix> ());
255
256 DiagMatrix sigma = result.singular_values ();
257
258 if (nargout == 0 || nargout == 1)
259 retval(0) = sigma.extract_diag ();
260 else if (nargout == 2)
261 retval = ovl (result.left_singular_matrix (),
262 sigma);
263 else
264 retval = ovl (result.left_singular_matrix (),
265 sigma,
266 result.right_singular_matrix ());
267 }
268 else if (arg.iscomplex ())
269 {
271
272 if (ctmp.any_element_is_inf_or_nan ())
273 error ("svd: cannot take SVD of matrix containing Inf or NaN values");
274
275 math::svd<ComplexMatrix> result
276 (ctmp,
277 svd_type<ComplexMatrix> (nargin, nargout, args, ctmp),
278 svd_driver<ComplexMatrix> ());
279
280 DiagMatrix sigma = result.singular_values ();
281
282 if (nargout == 0 || nargout == 1)
283 retval(0) = sigma.extract_diag ();
284 else if (nargout == 2)
285 retval = ovl (result.left_singular_matrix (),
286 sigma);
287 else
288 retval = ovl (result.left_singular_matrix (),
289 sigma,
290 result.right_singular_matrix ());
291 }
292 else
293 err_wrong_type_arg ("svd", arg);
294 }
295
296 return retval;
297}
298
299/*
300%!assert (svd ([1, 2; 2, 1]), [3; 1], sqrt (eps))
301
302%!test
303%! a = [1, 2; 3, 4] + [5, 6; 7, 8]*i;
304%! [u,s,v] = svd (a);
305%! assert (a, u * s * v', 128 * eps);
306
307%!test
308%! [u, s, v] = svd ([1, 2; 2, 1]);
309%! x = 1 / sqrt (2);
310%! ## Because values are accurate only up to a phase factor, tests must check
311%! ## separately for magnitude and for phase.
312%! assert (abs (u), abs ([-x, -x; -x, x]), sqrt (eps));
313%! phase = prod (sign (u));
314%! assert (phase == [1, -1] || phase == [-1, 1]);
315%! assert (s, [3, 0; 0, 1], sqrt (eps));
316%! assert (abs (v), abs ([-x, x; -x, -x]), sqrt (eps));
317%! phase = prod (sign (v));
318%! assert (phase == [1, -1] || phase == [-1, 1]);
319
320%!test
321%! a = [1, 2, 3; 4, 5, 6];
322%! [u, s, v] = svd (a);
323%! assert (u * s * v', a, sqrt (eps));
324
325%!test
326%! a = [1, 2; 3, 4; 5, 6];
327%! [u, s, v] = svd (a);
328%! assert (u * s * v', a, sqrt (eps));
329
330%!test
331%! a = [1, 2, 3; 4, 5, 6];
332%! [u, s, v] = svd (a, 1);
333%! assert (u * s * v', a, sqrt (eps));
334
335%!test
336%! a = [1, 2; 3, 4; 5, 6];
337%! [u, s, v] = svd (a, 1);
338%! assert (u * s * v', a, sqrt (eps));
339
340%!assert (svd (single ([1, 2; 2, 1])), single ([3; 1]), sqrt (eps ("single")))
341
342%!test
343%! [u, s, v] = svd (single ([1, 2; 2, 1]));
344%! x = single (1 / sqrt (2));
345
346%! assert (abs (u), abs ([-x, -x; -x, x]), sqrt (eps ("single")));
347%! phase = prod (sign (u));
348%! assert (phase == [1, -1] || phase == [-1, 1]);
349%! assert (s, single ([3, 0; 0, 1]), sqrt (eps ("single")));
350%! assert (abs (v), abs ([-x, x; -x, -x]), sqrt (eps ("single")));
351%! phase = prod (sign (v));
352%! assert (phase == [1, -1] || phase == [-1, 1]);
353
354%!test
355%! a = single ([1, 2, 3; 4, 5, 6]);
356%! [u, s, v] = svd (a);
357%! assert (u * s * v', a, sqrt (eps ("single")));
358
359%!test
360%! a = single ([1, 2; 3, 4; 5, 6]);
361%! [u, s, v] = svd (a);
362%! assert (u * s * v', a, sqrt (eps ("single")));
363
364%!test
365%! a = single ([1, 2, 3; 4, 5, 6]);
366%! [u, s, v] = svd (a, 1);
367%! assert (u * s * v', a, sqrt (eps ("single")));
368
369%!test
370%! a = single ([1, 2; 3, 4; 5, 6]);
371%! [u, s, v] = svd (a, 1);
372%! assert (u * s * v', a, sqrt (eps ("single")));
373
374%!test
375%! a = zeros (0, 5);
376%! [u, s, v] = svd (a);
377%! assert (size (u), [0, 0]);
378%! assert (size (s), [0, 5]);
379%! assert (size (v), [5, 5]);
380
381%!test
382%! a = zeros (5, 0);
383%! [u, s, v] = svd (a, 1);
384%! assert (size (u), [5, 0]);
385%! assert (size (s), [0, 0]);
386%! assert (size (v), [0, 0]);
387
388%!test <*49309>
389%! [~,~,v] = svd ([1, 1, 1], 0);
390%! assert (size (v), [3 3]);
391%! [~,~,v] = svd ([1, 1, 1], "econ");
392%! assert (size (v), [3 1]);
393
394%!assert <*55710> (1 / svd (-0), Inf)
395
396%!test
397%! old_driver = svd_driver ("gejsv");
398%! s0 = [1e-20; 1e-10; 1]; # only gejsv can pass
399%! q = sqrt (0.5);
400%! a = s0 .* [q, 0, -q; -0.5, q, -0.5; 0.5, q, 0.5];
401%! s1 = svd (a);
402%! svd_driver (old_driver);
403%! assert (sort (s1), s0, -10 * eps);
404
405%!error svd ()
406%!error svd ([1, 2; 4, 5], 2, 3)
407*/
408
409DEFUN (svd_driver, args, nargout,
410 doc: /* -*- texinfo -*-
411@deftypefn {} {@var{val} =} svd_driver ()
412@deftypefnx {} {@var{old_val} =} svd_driver (@var{new_val})
413@deftypefnx {} {@var{old_val} =} svd_driver (@var{new_val}, "local")
414Query or set the underlying @sc{lapack} driver used by @code{svd}.
415
416Currently recognized values are @qcode{"gesdd"}, @qcode{"gesvd"}, and
417@qcode{"gejsv"}. The default is @qcode{"gesvd"}.
418
419When called from inside a function with the @qcode{"local"} option, the
420variable is changed locally for the function and any subroutines it calls.
421The original variable value is restored when exiting the function.
422
423Algorithm Notes: The @sc{lapack} library routines @code{gesvd} and @code{gesdd}
424are different only when calculating the full singular value decomposition (left
425and right singular matrices as well as singular values). When calculating just
426the singular values the following discussion is not relevant.
427
428The newer @code{gesdd} routine is based on a Divide-and-Conquer algorithm that
429is 5X faster than the alternative @code{gesvd}, which is based on QR
430factorization. However, the new algorithm can use significantly more memory.
431For an @nospell{MxN} input matrix the memory usage is of order O(min(M,N) ^ 2),
432whereas the alternative is of order O(max(M,N)).
433
434The routine @code{gejsv} uses a preconditioned Jacobi SVD algorithm. Unlike
435@code{gesvd} and @code{gesdd}, in @code{gejsv}, there is no bidiagonalization
436step that could contaminate accuracy in some extreme cases. Also, @code{gejsv}
437is known to be optimally accurate in some sense. However, the speed is slower
438(single threaded at its core) and uses more memory (O(min(M,N) ^ 2 + M + N)).
439
440Beyond speed and memory issues, there have been instances where some input
441matrices were not accurately decomposed by @code{gesdd}. See currently active
442bug @url{https://savannah.gnu.org/bugs/?55564}. Until these accuracy issues
443are resolved in a new version of the @sc{lapack} library, the default driver
444in Octave has been set to @qcode{"gesvd"}.
445
446@seealso{svd}
447@end deftypefn */)
448{
449 static const char *driver_names[] = { "gesvd", "gesdd", "gejsv", nullptr };
450
451 return set_internal_variable (Vsvd_driver, args, nargout,
452 "svd_driver", driver_names);
453}
454
455/*
456%!test
457%! A = [1+1i, 1-1i, 0; 0, 2, 0; 1i, 1i, 1+2i];
458%! old_driver = svd_driver ("gesvd");
459%! [U1, S1, V1] = svd (A);
460%! svd_driver ("gesdd");
461%! [U2, S2, V2] = svd (A);
462%! svd_driver ("gejsv");
463%! [U3, S3, V3] = svd (A);
464%! assert (svd_driver (), "gejsv");
465%! svd_driver (old_driver);
466%! assert (U1, U2, 6*eps);
467%! assert (S1, S2, 6*eps);
468%! assert (V1, V2, 6*eps);
469%! z = U1(1,:) ./ U3(1,:);
470%! assert (U1, U3 .* z, 100*eps);
471%! assert (S1, S3, 10*eps);
472%! assert (V1, V3 .* z, 100*eps);
473*/
474
475OCTAVE_END_NAMESPACE(octave)
bool any_element_is_inf_or_nan() const
Definition CNDArray.cc:271
ColumnVector extract_diag(octave_idx_type k=0) const
bool any_element_is_inf_or_nan() const
Definition fCNDArray.cc:271
FloatColumnVector extract_diag(octave_idx_type k=0) const
bool any_element_is_inf_or_nan() const
Definition fNDArray.cc:281
bool any_element_is_inf_or_nan() const
Definition dNDArray.cc:324
octave_idx_type length() const
Definition ovl.h:111
int ndims() const
Definition ov.h:549
bool isreal() const
Definition ov.h:736
ComplexMatrix complex_matrix_value(bool frc_str_conv=false) const
Definition ov.h:875
bool is_single_type() const
Definition ov.h:696
FloatMatrix float_matrix_value(bool frc_str_conv=false) const
Definition ov.h:860
bool iscomplex() const
Definition ov.h:739
Matrix matrix_value(bool frc_str_conv=false) const
Definition ov.h:857
FloatComplexMatrix float_complex_matrix_value(bool frc_str_conv=false) const
Definition ov.h:879
Definition svd.h:38
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
void print_usage()
Definition defun-int.h:72
#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:1008
void err_wrong_type_arg(const char *name, const char *s)
Definition errwarn.cc:166
F77_RET_T const F77_INT F77_CMPLX * A
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition ovl.h:217
octave_value set_internal_variable(bool &var, const octave_value_list &args, int nargout, const char *nm)
Definition variables.cc:583