_svds.py 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563
  1. import os
  2. import numpy as np
  3. from .arpack import _arpack # type: ignore[attr-defined]
  4. from . import eigsh
  5. from scipy._lib._util import check_random_state
  6. from scipy.sparse.linalg._interface import LinearOperator, aslinearoperator
  7. from scipy.sparse.linalg._eigen.lobpcg import lobpcg # type: ignore[no-redef]
  8. if os.environ.get("SCIPY_USE_PROPACK"):
  9. from scipy.sparse.linalg._svdp import _svdp
  10. HAS_PROPACK = True
  11. else:
  12. HAS_PROPACK = False
  13. from scipy.linalg import svd
  14. arpack_int = _arpack.timing.nbx.dtype
  15. __all__ = ['svds']
  16. def _herm(x):
  17. return x.T.conj()
  18. def _iv(A, k, ncv, tol, which, v0, maxiter,
  19. return_singular, solver, random_state):
  20. # input validation/standardization for `solver`
  21. # out of order because it's needed for other parameters
  22. solver = str(solver).lower()
  23. solvers = {"arpack", "lobpcg", "propack"}
  24. if solver not in solvers:
  25. raise ValueError(f"solver must be one of {solvers}.")
  26. # input validation/standardization for `A`
  27. A = aslinearoperator(A) # this takes care of some input validation
  28. if not (np.issubdtype(A.dtype, np.complexfloating)
  29. or np.issubdtype(A.dtype, np.floating)):
  30. message = "`A` must be of floating or complex floating data type."
  31. raise ValueError(message)
  32. if np.prod(A.shape) == 0:
  33. message = "`A` must not be empty."
  34. raise ValueError(message)
  35. # input validation/standardization for `k`
  36. kmax = min(A.shape) if solver == 'propack' else min(A.shape) - 1
  37. if int(k) != k or not (0 < k <= kmax):
  38. message = "`k` must be an integer satisfying `0 < k < min(A.shape)`."
  39. raise ValueError(message)
  40. k = int(k)
  41. # input validation/standardization for `ncv`
  42. if solver == "arpack" and ncv is not None:
  43. if int(ncv) != ncv or not (k < ncv < min(A.shape)):
  44. message = ("`ncv` must be an integer satisfying "
  45. "`k < ncv < min(A.shape)`.")
  46. raise ValueError(message)
  47. ncv = int(ncv)
  48. # input validation/standardization for `tol`
  49. if tol < 0 or not np.isfinite(tol):
  50. message = "`tol` must be a non-negative floating point value."
  51. raise ValueError(message)
  52. tol = float(tol)
  53. # input validation/standardization for `which`
  54. which = str(which).upper()
  55. whichs = {'LM', 'SM'}
  56. if which not in whichs:
  57. raise ValueError(f"`which` must be in {whichs}.")
  58. # input validation/standardization for `v0`
  59. if v0 is not None:
  60. v0 = np.atleast_1d(v0)
  61. if not (np.issubdtype(v0.dtype, np.complexfloating)
  62. or np.issubdtype(v0.dtype, np.floating)):
  63. message = ("`v0` must be of floating or complex floating "
  64. "data type.")
  65. raise ValueError(message)
  66. shape = (A.shape[0],) if solver == 'propack' else (min(A.shape),)
  67. if v0.shape != shape:
  68. message = f"`v0` must have shape {shape}."
  69. raise ValueError(message)
  70. # input validation/standardization for `maxiter`
  71. if maxiter is not None and (int(maxiter) != maxiter or maxiter <= 0):
  72. message = "`maxiter` must be a positive integer."
  73. raise ValueError(message)
  74. maxiter = int(maxiter) if maxiter is not None else maxiter
  75. # input validation/standardization for `return_singular_vectors`
  76. # not going to be flexible with this; too complicated for little gain
  77. rs_options = {True, False, "vh", "u"}
  78. if return_singular not in rs_options:
  79. raise ValueError(f"`return_singular_vectors` must be in {rs_options}.")
  80. random_state = check_random_state(random_state)
  81. return (A, k, ncv, tol, which, v0, maxiter,
  82. return_singular, solver, random_state)
  83. def svds(A, k=6, ncv=None, tol=0, which='LM', v0=None,
  84. maxiter=None, return_singular_vectors=True,
  85. solver='arpack', random_state=None, options=None):
  86. """
  87. Partial singular value decomposition of a sparse matrix.
  88. Compute the largest or smallest `k` singular values and corresponding
  89. singular vectors of a sparse matrix `A`. The order in which the singular
  90. values are returned is not guaranteed.
  91. In the descriptions below, let ``M, N = A.shape``.
  92. Parameters
  93. ----------
  94. A : ndarray, sparse matrix, or LinearOperator
  95. Matrix to decompose of a floating point numeric dtype.
  96. k : int, default: 6
  97. Number of singular values and singular vectors to compute.
  98. Must satisfy ``1 <= k <= kmax``, where ``kmax=min(M, N)`` for
  99. ``solver='propack'`` and ``kmax=min(M, N) - 1`` otherwise.
  100. ncv : int, optional
  101. When ``solver='arpack'``, this is the number of Lanczos vectors
  102. generated. See :ref:`'arpack' <sparse.linalg.svds-arpack>` for details.
  103. When ``solver='lobpcg'`` or ``solver='propack'``, this parameter is
  104. ignored.
  105. tol : float, optional
  106. Tolerance for singular values. Zero (default) means machine precision.
  107. which : {'LM', 'SM'}
  108. Which `k` singular values to find: either the largest magnitude ('LM')
  109. or smallest magnitude ('SM') singular values.
  110. v0 : ndarray, optional
  111. The starting vector for iteration; see method-specific
  112. documentation (:ref:`'arpack' <sparse.linalg.svds-arpack>`,
  113. :ref:`'lobpcg' <sparse.linalg.svds-lobpcg>`), or
  114. :ref:`'propack' <sparse.linalg.svds-propack>` for details.
  115. maxiter : int, optional
  116. Maximum number of iterations; see method-specific
  117. documentation (:ref:`'arpack' <sparse.linalg.svds-arpack>`,
  118. :ref:`'lobpcg' <sparse.linalg.svds-lobpcg>`), or
  119. :ref:`'propack' <sparse.linalg.svds-propack>` for details.
  120. return_singular_vectors : {True, False, "u", "vh"}
  121. Singular values are always computed and returned; this parameter
  122. controls the computation and return of singular vectors.
  123. - ``True``: return singular vectors.
  124. - ``False``: do not return singular vectors.
  125. - ``"u"``: if ``M <= N``, compute only the left singular vectors and
  126. return ``None`` for the right singular vectors. Otherwise, compute
  127. all singular vectors.
  128. - ``"vh"``: if ``M > N``, compute only the right singular vectors and
  129. return ``None`` for the left singular vectors. Otherwise, compute
  130. all singular vectors.
  131. If ``solver='propack'``, the option is respected regardless of the
  132. matrix shape.
  133. solver : {'arpack', 'propack', 'lobpcg'}, optional
  134. The solver used.
  135. :ref:`'arpack' <sparse.linalg.svds-arpack>`,
  136. :ref:`'lobpcg' <sparse.linalg.svds-lobpcg>`, and
  137. :ref:`'propack' <sparse.linalg.svds-propack>` are supported.
  138. Default: `'arpack'`.
  139. random_state : {None, int, `numpy.random.Generator`,
  140. `numpy.random.RandomState`}, optional
  141. Pseudorandom number generator state used to generate resamples.
  142. If `random_state` is ``None`` (or `np.random`), the
  143. `numpy.random.RandomState` singleton is used.
  144. If `random_state` is an int, a new ``RandomState`` instance is used,
  145. seeded with `random_state`.
  146. If `random_state` is already a ``Generator`` or ``RandomState``
  147. instance then that instance is used.
  148. options : dict, optional
  149. A dictionary of solver-specific options. No solver-specific options
  150. are currently supported; this parameter is reserved for future use.
  151. Returns
  152. -------
  153. u : ndarray, shape=(M, k)
  154. Unitary matrix having left singular vectors as columns.
  155. s : ndarray, shape=(k,)
  156. The singular values.
  157. vh : ndarray, shape=(k, N)
  158. Unitary matrix having right singular vectors as rows.
  159. Notes
  160. -----
  161. This is a naive implementation using ARPACK or LOBPCG as an eigensolver
  162. on the matrix ``A.conj().T @ A`` or ``A @ A.conj().T``, depending on
  163. which one is smaller size, followed by the Rayleigh-Ritz method
  164. as postprocessing; see
  165. Using the normal matrix, in Rayleigh-Ritz method, (2022, Nov. 19),
  166. Wikipedia, https://w.wiki/4zms.
  167. Alternatively, the PROPACK solver can be called. ``form="array"``
  168. Choices of the input matrix ``A`` numeric dtype may be limited.
  169. Only ``solver="lobpcg"`` supports all floating point dtypes
  170. real: 'np.single', 'np.double', 'np.longdouble' and
  171. complex: 'np.csingle', 'np.cdouble', 'np.clongdouble'.
  172. The ``solver="arpack"`` supports only
  173. 'np.single', 'np.double', and 'np.cdouble'.
  174. Examples
  175. --------
  176. Construct a matrix ``A`` from singular values and vectors.
  177. >>> import numpy as np
  178. >>> from scipy.stats import ortho_group
  179. >>> from scipy.sparse.linalg import svds
  180. >>> from scipy.sparse import csr_matrix
  181. >>> rng = np.random.default_rng()
  182. Construct a dense matrix ``A`` from singular values and vectors.
  183. >>> orthogonal = ortho_group.rvs(10, random_state=rng)
  184. >>> s = [1e-3, 1, 2, 3, 4] # non-zero singular values
  185. >>> u = orthogonal[:, :5] # left singular vectors
  186. >>> vT = orthogonal[:, 5:].T # right singular vectors
  187. >>> A = u @ np.diag(s) @ vT
  188. With only four singular values/vectors, the SVD approximates the original
  189. matrix.
  190. >>> u4, s4, vT4 = svds(A, k=4)
  191. >>> A4 = u4 @ np.diag(s4) @ vT4
  192. >>> np.allclose(A4, A, atol=1e-3)
  193. True
  194. With all five non-zero singular values/vectors, we can reproduce
  195. the original matrix more accurately.
  196. >>> u5, s5, vT5 = svds(A, k=5)
  197. >>> A5 = u5 @ np.diag(s5) @ vT5
  198. >>> np.allclose(A5, A)
  199. True
  200. The singular values match the expected singular values.
  201. >>> np.allclose(s5, s)
  202. True
  203. Since the singular values are not close to each other in this example,
  204. every singular vector matches as expected up to a difference in sign.
  205. >>> (np.allclose(np.abs(u5), np.abs(u)) and
  206. ... np.allclose(np.abs(vT5), np.abs(vT)))
  207. True
  208. The singular vectors are also orthogonal.
  209. >>> (np.allclose(u5.T @ u5, np.eye(5)) and
  210. ... np.allclose(vT5 @ vT5.T, np.eye(5)))
  211. True
  212. If there are (nearly) multiple singular values, the corresponding
  213. individual singular vectors may be unstable, but the whole invariant
  214. subspace containing all such singular vectors is computed accurately
  215. as can be measured by angles between subspaces via 'subspace_angles'.
  216. >>> from scipy.linalg import subspace_angles as s_a
  217. >>> rng = np.random.default_rng()
  218. >>> s = [1, 1 + 1e-6] # non-zero singular values
  219. >>> u, _ = np.linalg.qr(rng.standard_normal((99, 2)))
  220. >>> v, _ = np.linalg.qr(rng.standard_normal((99, 2)))
  221. >>> vT = v.T
  222. >>> A = u @ np.diag(s) @ vT
  223. >>> A = A.astype(np.float32)
  224. >>> u2, s2, vT2 = svds(A, k=2)
  225. >>> np.allclose(s2, s)
  226. True
  227. The angles between the individual exact and computed singular vectors
  228. are not so small.
  229. >>> s_a(u2[:, :1], u[:, :1]) + s_a(u2[:, 1:], u[:, 1:]) > 1e-3
  230. True
  231. >>> (s_a(vT2[:1, :].T, vT[:1, :].T) +
  232. ... s_a(vT2[1:, :].T, vT[1:, :].T)) > 1e-3
  233. True
  234. As opposed to the angles between the 2-dimensional invariant subspaces
  235. that these vectors span, which are small for rights singular vectors
  236. >>> s_a(u2, u).sum() < 1e-6
  237. True
  238. as well as for left singular vectors.
  239. >>> s_a(vT2.T, vT.T).sum() < 1e-6
  240. True
  241. The next example follows that of 'sklearn.decomposition.TruncatedSVD'.
  242. >>> rng = np.random.RandomState(0)
  243. >>> X_dense = rng.random(size=(100, 100))
  244. >>> X_dense[:, 2 * np.arange(50)] = 0
  245. >>> X = csr_matrix(X_dense)
  246. >>> _, singular_values, _ = svds(X, k=5)
  247. >>> print(singular_values)
  248. [ 4.3293... 4.4491... 4.5420... 4.5987... 35.2410...]
  249. The function can be called without the transpose of the input matrix
  250. ever explicitly constructed.
  251. >>> from scipy.linalg import svd
  252. >>> from scipy.sparse import rand
  253. >>> from scipy.sparse.linalg import aslinearoperator
  254. >>> rng = np.random.RandomState(0)
  255. >>> G = rand(8, 9, density=0.5, random_state=rng)
  256. >>> Glo = aslinearoperator(G)
  257. >>> _, singular_values_svds, _ = svds(Glo, k=5)
  258. >>> _, singular_values_svd, _ = svd(G.toarray())
  259. >>> np.allclose(singular_values_svds, singular_values_svd[-4::-1])
  260. True
  261. The most memory efficient scenario is where neither
  262. the original matrix, nor its transpose, is explicitly constructed.
  263. Our example computes the smallest singular values and vectors
  264. of 'LinearOperator' constructed from the numpy function 'np.diff' used
  265. column-wise to be consistent with 'LinearOperator' operating on columns.
  266. >>> from scipy.sparse.linalg import LinearOperator, aslinearoperator
  267. >>> diff0 = lambda a: np.diff(a, axis=0)
  268. Let us create the matrix from 'diff0' to be used for validation only.
  269. >>> n = 5 # The dimension of the space.
  270. >>> M_from_diff0 = diff0(np.eye(n))
  271. >>> print(M_from_diff0.astype(int))
  272. [[-1 1 0 0 0]
  273. [ 0 -1 1 0 0]
  274. [ 0 0 -1 1 0]
  275. [ 0 0 0 -1 1]]
  276. The matrix 'M_from_diff0' is bi-diagonal and could be alternatively
  277. created directly by
  278. >>> M = - np.eye(n - 1, n, dtype=int)
  279. >>> np.fill_diagonal(M[:,1:], 1)
  280. >>> np.allclose(M, M_from_diff0)
  281. True
  282. Its transpose
  283. >>> print(M.T)
  284. [[-1 0 0 0]
  285. [ 1 -1 0 0]
  286. [ 0 1 -1 0]
  287. [ 0 0 1 -1]
  288. [ 0 0 0 1]]
  289. can be viewed as the incidence matrix; see
  290. Incidence matrix, (2022, Nov. 19), Wikipedia, https://w.wiki/5YXU,
  291. of a linear graph with 5 vertices and 4 edges. The 5x5 normal matrix
  292. 'M.T @ M' thus is
  293. >>> print(M.T @ M)
  294. [[ 1 -1 0 0 0]
  295. [-1 2 -1 0 0]
  296. [ 0 -1 2 -1 0]
  297. [ 0 0 -1 2 -1]
  298. [ 0 0 0 -1 1]]
  299. the graph Laplacian, while the actually used in 'svds' smaller size
  300. 4x4 normal matrix 'M @ M.T'
  301. >>> print(M @ M.T)
  302. [[ 2 -1 0 0]
  303. [-1 2 -1 0]
  304. [ 0 -1 2 -1]
  305. [ 0 0 -1 2]]
  306. is the so-called edge-based Laplacian; see
  307. Symmetric Laplacian via the incidence matrix, in Laplacian matrix,
  308. (2022, Nov. 19), Wikipedia, https://w.wiki/5YXW.
  309. The 'LinearOperator' setup needs the options 'rmatvec' and 'rmatmat'
  310. of multiplication by the matrix transpose 'M.T', but we want to be
  311. matrix-free to save memory, so knowing how 'M.T' looks like, we
  312. manually construct the following function to be used in 'rmatmat=diff0t'.
  313. >>> def diff0t(a):
  314. ... if a.ndim == 1:
  315. ... a = a[:,np.newaxis] # Turn 1D into 2D array
  316. ... d = np.zeros((a.shape[0] + 1, a.shape[1]), dtype=a.dtype)
  317. ... d[0, :] = - a[0, :]
  318. ... d[1:-1, :] = a[0:-1, :] - a[1:, :]
  319. ... d[-1, :] = a[-1, :]
  320. ... return d
  321. We check that our function 'diff0t' for the matrix transpose is valid.
  322. >>> np.allclose(M.T, diff0t(np.eye(n-1)))
  323. True
  324. Now we setup our matrix-free 'LinearOperator' called 'diff0_func_aslo'
  325. and for validation the matrix-based 'diff0_matrix_aslo'.
  326. >>> def diff0_func_aslo_def(n):
  327. ... return LinearOperator(matvec=diff0,
  328. ... matmat=diff0,
  329. ... rmatvec=diff0t,
  330. ... rmatmat=diff0t,
  331. ... shape=(n - 1, n))
  332. >>> diff0_func_aslo = diff0_func_aslo_def(n)
  333. >>> diff0_matrix_aslo = aslinearoperator(M_from_diff0)
  334. And validate both the matrix and its transpose in 'LinearOperator'.
  335. >>> np.allclose(diff0_func_aslo(np.eye(n)),
  336. ... diff0_matrix_aslo(np.eye(n)))
  337. True
  338. >>> np.allclose(diff0_func_aslo.T(np.eye(n-1)),
  339. ... diff0_matrix_aslo.T(np.eye(n-1)))
  340. True
  341. Having the 'LinearOperator' setup validated, we run the solver.
  342. >>> n = 100
  343. >>> diff0_func_aslo = diff0_func_aslo_def(n)
  344. >>> u, s, vT = svds(diff0_func_aslo, k=3, which='SM')
  345. The singular values squared and the singular vectors are known
  346. explicitly; see
  347. Pure Dirichlet boundary conditions, in
  348. Eigenvalues and eigenvectors of the second derivative,
  349. (2022, Nov. 19), Wikipedia, https://w.wiki/5YX6,
  350. since 'diff' corresponds to first
  351. derivative, and its smaller size n-1 x n-1 normal matrix
  352. 'M @ M.T' represent the discrete second derivative with the Dirichlet
  353. boundary conditions. We use these analytic expressions for validation.
  354. >>> se = 2. * np.sin(np.pi * np.arange(1, 4) / (2. * n))
  355. >>> ue = np.sqrt(2 / n) * np.sin(np.pi * np.outer(np.arange(1, n),
  356. ... np.arange(1, 4)) / n)
  357. >>> np.allclose(s, se, atol=1e-3)
  358. True
  359. >>> print(np.allclose(np.abs(u), np.abs(ue), atol=1e-6))
  360. True
  361. """
  362. args = _iv(A, k, ncv, tol, which, v0, maxiter, return_singular_vectors,
  363. solver, random_state)
  364. (A, k, ncv, tol, which, v0, maxiter,
  365. return_singular_vectors, solver, random_state) = args
  366. largest = (which == 'LM')
  367. n, m = A.shape
  368. if n >= m:
  369. X_dot = A.matvec
  370. X_matmat = A.matmat
  371. XH_dot = A.rmatvec
  372. XH_mat = A.rmatmat
  373. transpose = False
  374. else:
  375. X_dot = A.rmatvec
  376. X_matmat = A.rmatmat
  377. XH_dot = A.matvec
  378. XH_mat = A.matmat
  379. transpose = True
  380. dtype = getattr(A, 'dtype', None)
  381. if dtype is None:
  382. dtype = A.dot(np.zeros([m, 1])).dtype
  383. def matvec_XH_X(x):
  384. return XH_dot(X_dot(x))
  385. def matmat_XH_X(x):
  386. return XH_mat(X_matmat(x))
  387. XH_X = LinearOperator(matvec=matvec_XH_X, dtype=A.dtype,
  388. matmat=matmat_XH_X,
  389. shape=(min(A.shape), min(A.shape)))
  390. # Get a low rank approximation of the implicitly defined gramian matrix.
  391. # This is not a stable way to approach the problem.
  392. if solver == 'lobpcg':
  393. if k == 1 and v0 is not None:
  394. X = np.reshape(v0, (-1, 1))
  395. else:
  396. X = random_state.standard_normal(size=(min(A.shape), k))
  397. _, eigvec = lobpcg(XH_X, X, tol=tol ** 2, maxiter=maxiter,
  398. largest=largest)
  399. # lobpcg does not guarantee exactly orthonormal eigenvectors
  400. # until after gh-16320 is merged
  401. eigvec, _ = np.linalg.qr(eigvec)
  402. elif solver == 'propack':
  403. if not HAS_PROPACK:
  404. raise ValueError("`solver='propack'` is opt-in due "
  405. "to potential issues on Windows, "
  406. "it can be enabled by setting the "
  407. "`SCIPY_USE_PROPACK` environment "
  408. "variable before importing scipy")
  409. jobu = return_singular_vectors in {True, 'u'}
  410. jobv = return_singular_vectors in {True, 'vh'}
  411. irl_mode = (which == 'SM')
  412. res = _svdp(A, k=k, tol=tol**2, which=which, maxiter=None,
  413. compute_u=jobu, compute_v=jobv, irl_mode=irl_mode,
  414. kmax=maxiter, v0=v0, random_state=random_state)
  415. u, s, vh, _ = res # but we'll ignore bnd, the last output
  416. # PROPACK order appears to be largest first. `svds` output order is not
  417. # guaranteed, according to documentation, but for ARPACK and LOBPCG
  418. # they actually are ordered smallest to largest, so reverse for
  419. # consistency.
  420. s = s[::-1]
  421. u = u[:, ::-1]
  422. vh = vh[::-1]
  423. u = u if jobu else None
  424. vh = vh if jobv else None
  425. if return_singular_vectors:
  426. return u, s, vh
  427. else:
  428. return s
  429. elif solver == 'arpack' or solver is None:
  430. if v0 is None:
  431. v0 = random_state.standard_normal(size=(min(A.shape),))
  432. _, eigvec = eigsh(XH_X, k=k, tol=tol ** 2, maxiter=maxiter,
  433. ncv=ncv, which=which, v0=v0)
  434. # arpack do not guarantee exactly orthonormal eigenvectors
  435. # for clustered eigenvalues, especially in complex arithmetic
  436. eigvec, _ = np.linalg.qr(eigvec)
  437. # the eigenvectors eigvec must be orthonomal here; see gh-16712
  438. Av = X_matmat(eigvec)
  439. if not return_singular_vectors:
  440. s = svd(Av, compute_uv=False, overwrite_a=True)
  441. return s[::-1]
  442. # compute the left singular vectors of X and update the right ones
  443. # accordingly
  444. u, s, vh = svd(Av, full_matrices=False, overwrite_a=True)
  445. u = u[:, ::-1]
  446. s = s[::-1]
  447. vh = vh[::-1]
  448. jobu = return_singular_vectors in {True, 'u'}
  449. jobv = return_singular_vectors in {True, 'vh'}
  450. if transpose:
  451. u_tmp = eigvec @ _herm(vh) if jobu else None
  452. vh = _herm(u) if jobv else None
  453. u = u_tmp
  454. else:
  455. if not jobu:
  456. u = None
  457. vh = vh @ _herm(eigvec) if jobv else None
  458. return u, s, vh