annotate DEPENDENCIES/mingw32/Python27/Lib/site-packages/numpy/linalg/linalg.py @ 133:4acb5d8d80b6 tip

Don't fail environmental check if README.md exists (but .txt and no-suffix don't)
author Chris Cannam
date Tue, 30 Jul 2019 12:25:44 +0100
parents 2a2c65a20a8b
children
rev   line source
Chris@87 1 """Lite version of scipy.linalg.
Chris@87 2
Chris@87 3 Notes
Chris@87 4 -----
Chris@87 5 This module is a lite version of the linalg.py module in SciPy which
Chris@87 6 contains high-level Python interface to the LAPACK library. The lite
Chris@87 7 version only accesses the following LAPACK functions: dgesv, zgesv,
Chris@87 8 dgeev, zgeev, dgesdd, zgesdd, dgelsd, zgelsd, dsyevd, zheevd, dgetrf,
Chris@87 9 zgetrf, dpotrf, zpotrf, dgeqrf, zgeqrf, zungqr, dorgqr.
Chris@87 10 """
Chris@87 11 from __future__ import division, absolute_import, print_function
Chris@87 12
Chris@87 13
Chris@87 14 __all__ = ['matrix_power', 'solve', 'tensorsolve', 'tensorinv', 'inv',
Chris@87 15 'cholesky', 'eigvals', 'eigvalsh', 'pinv', 'slogdet', 'det',
Chris@87 16 'svd', 'eig', 'eigh', 'lstsq', 'norm', 'qr', 'cond', 'matrix_rank',
Chris@87 17 'LinAlgError']
Chris@87 18
Chris@87 19 import warnings
Chris@87 20
Chris@87 21 from numpy.core import (
Chris@87 22 array, asarray, zeros, empty, empty_like, transpose, intc, single, double,
Chris@87 23 csingle, cdouble, inexact, complexfloating, newaxis, ravel, all, Inf, dot,
Chris@87 24 add, multiply, sqrt, maximum, fastCopyAndTranspose, sum, isfinite, size,
Chris@87 25 finfo, errstate, geterrobj, longdouble, rollaxis, amin, amax, product, abs,
Chris@87 26 broadcast
Chris@87 27 )
Chris@87 28 from numpy.lib import triu, asfarray
Chris@87 29 from numpy.linalg import lapack_lite, _umath_linalg
Chris@87 30 from numpy.matrixlib.defmatrix import matrix_power
Chris@87 31 from numpy.compat import asbytes
Chris@87 32
Chris@87 33 # For Python2/3 compatibility
Chris@87 34 _N = asbytes('N')
Chris@87 35 _V = asbytes('V')
Chris@87 36 _A = asbytes('A')
Chris@87 37 _S = asbytes('S')
Chris@87 38 _L = asbytes('L')
Chris@87 39
Chris@87 40 fortran_int = intc
Chris@87 41
Chris@87 42 # Error object
Chris@87 43 class LinAlgError(Exception):
Chris@87 44 """
Chris@87 45 Generic Python-exception-derived object raised by linalg functions.
Chris@87 46
Chris@87 47 General purpose exception class, derived from Python's exception.Exception
Chris@87 48 class, programmatically raised in linalg functions when a Linear
Chris@87 49 Algebra-related condition would prevent further correct execution of the
Chris@87 50 function.
Chris@87 51
Chris@87 52 Parameters
Chris@87 53 ----------
Chris@87 54 None
Chris@87 55
Chris@87 56 Examples
Chris@87 57 --------
Chris@87 58 >>> from numpy import linalg as LA
Chris@87 59 >>> LA.inv(np.zeros((2,2)))
Chris@87 60 Traceback (most recent call last):
Chris@87 61 File "<stdin>", line 1, in <module>
Chris@87 62 File "...linalg.py", line 350,
Chris@87 63 in inv return wrap(solve(a, identity(a.shape[0], dtype=a.dtype)))
Chris@87 64 File "...linalg.py", line 249,
Chris@87 65 in solve
Chris@87 66 raise LinAlgError('Singular matrix')
Chris@87 67 numpy.linalg.LinAlgError: Singular matrix
Chris@87 68
Chris@87 69 """
Chris@87 70 pass
Chris@87 71
Chris@87 72 # Dealing with errors in _umath_linalg
Chris@87 73
Chris@87 74 _linalg_error_extobj = None
Chris@87 75
Chris@87 76 def _determine_error_states():
Chris@87 77 global _linalg_error_extobj
Chris@87 78 errobj = geterrobj()
Chris@87 79 bufsize = errobj[0]
Chris@87 80
Chris@87 81 with errstate(invalid='call', over='ignore',
Chris@87 82 divide='ignore', under='ignore'):
Chris@87 83 invalid_call_errmask = geterrobj()[1]
Chris@87 84
Chris@87 85 _linalg_error_extobj = [bufsize, invalid_call_errmask, None]
Chris@87 86
Chris@87 87 _determine_error_states()
Chris@87 88
Chris@87 89 def _raise_linalgerror_singular(err, flag):
Chris@87 90 raise LinAlgError("Singular matrix")
Chris@87 91
Chris@87 92 def _raise_linalgerror_nonposdef(err, flag):
Chris@87 93 raise LinAlgError("Matrix is not positive definite")
Chris@87 94
Chris@87 95 def _raise_linalgerror_eigenvalues_nonconvergence(err, flag):
Chris@87 96 raise LinAlgError("Eigenvalues did not converge")
Chris@87 97
Chris@87 98 def _raise_linalgerror_svd_nonconvergence(err, flag):
Chris@87 99 raise LinAlgError("SVD did not converge")
Chris@87 100
Chris@87 101 def get_linalg_error_extobj(callback):
Chris@87 102 extobj = list(_linalg_error_extobj)
Chris@87 103 extobj[2] = callback
Chris@87 104 return extobj
Chris@87 105
Chris@87 106 def _makearray(a):
Chris@87 107 new = asarray(a)
Chris@87 108 wrap = getattr(a, "__array_prepare__", new.__array_wrap__)
Chris@87 109 return new, wrap
Chris@87 110
Chris@87 111 def isComplexType(t):
Chris@87 112 return issubclass(t, complexfloating)
Chris@87 113
Chris@87 114 _real_types_map = {single : single,
Chris@87 115 double : double,
Chris@87 116 csingle : single,
Chris@87 117 cdouble : double}
Chris@87 118
Chris@87 119 _complex_types_map = {single : csingle,
Chris@87 120 double : cdouble,
Chris@87 121 csingle : csingle,
Chris@87 122 cdouble : cdouble}
Chris@87 123
Chris@87 124 def _realType(t, default=double):
Chris@87 125 return _real_types_map.get(t, default)
Chris@87 126
Chris@87 127 def _complexType(t, default=cdouble):
Chris@87 128 return _complex_types_map.get(t, default)
Chris@87 129
Chris@87 130 def _linalgRealType(t):
Chris@87 131 """Cast the type t to either double or cdouble."""
Chris@87 132 return double
Chris@87 133
Chris@87 134 _complex_types_map = {single : csingle,
Chris@87 135 double : cdouble,
Chris@87 136 csingle : csingle,
Chris@87 137 cdouble : cdouble}
Chris@87 138
Chris@87 139 def _commonType(*arrays):
Chris@87 140 # in lite version, use higher precision (always double or cdouble)
Chris@87 141 result_type = single
Chris@87 142 is_complex = False
Chris@87 143 for a in arrays:
Chris@87 144 if issubclass(a.dtype.type, inexact):
Chris@87 145 if isComplexType(a.dtype.type):
Chris@87 146 is_complex = True
Chris@87 147 rt = _realType(a.dtype.type, default=None)
Chris@87 148 if rt is None:
Chris@87 149 # unsupported inexact scalar
Chris@87 150 raise TypeError("array type %s is unsupported in linalg" %
Chris@87 151 (a.dtype.name,))
Chris@87 152 else:
Chris@87 153 rt = double
Chris@87 154 if rt is double:
Chris@87 155 result_type = double
Chris@87 156 if is_complex:
Chris@87 157 t = cdouble
Chris@87 158 result_type = _complex_types_map[result_type]
Chris@87 159 else:
Chris@87 160 t = double
Chris@87 161 return t, result_type
Chris@87 162
Chris@87 163
Chris@87 164 # _fastCopyAndTranpose assumes the input is 2D (as all the calls in here are).
Chris@87 165
Chris@87 166 _fastCT = fastCopyAndTranspose
Chris@87 167
Chris@87 168 def _to_native_byte_order(*arrays):
Chris@87 169 ret = []
Chris@87 170 for arr in arrays:
Chris@87 171 if arr.dtype.byteorder not in ('=', '|'):
Chris@87 172 ret.append(asarray(arr, dtype=arr.dtype.newbyteorder('=')))
Chris@87 173 else:
Chris@87 174 ret.append(arr)
Chris@87 175 if len(ret) == 1:
Chris@87 176 return ret[0]
Chris@87 177 else:
Chris@87 178 return ret
Chris@87 179
Chris@87 180 def _fastCopyAndTranspose(type, *arrays):
Chris@87 181 cast_arrays = ()
Chris@87 182 for a in arrays:
Chris@87 183 if a.dtype.type is type:
Chris@87 184 cast_arrays = cast_arrays + (_fastCT(a),)
Chris@87 185 else:
Chris@87 186 cast_arrays = cast_arrays + (_fastCT(a.astype(type)),)
Chris@87 187 if len(cast_arrays) == 1:
Chris@87 188 return cast_arrays[0]
Chris@87 189 else:
Chris@87 190 return cast_arrays
Chris@87 191
Chris@87 192 def _assertRank2(*arrays):
Chris@87 193 for a in arrays:
Chris@87 194 if len(a.shape) != 2:
Chris@87 195 raise LinAlgError('%d-dimensional array given. Array must be '
Chris@87 196 'two-dimensional' % len(a.shape))
Chris@87 197
Chris@87 198 def _assertRankAtLeast2(*arrays):
Chris@87 199 for a in arrays:
Chris@87 200 if len(a.shape) < 2:
Chris@87 201 raise LinAlgError('%d-dimensional array given. Array must be '
Chris@87 202 'at least two-dimensional' % len(a.shape))
Chris@87 203
Chris@87 204 def _assertSquareness(*arrays):
Chris@87 205 for a in arrays:
Chris@87 206 if max(a.shape) != min(a.shape):
Chris@87 207 raise LinAlgError('Array must be square')
Chris@87 208
Chris@87 209 def _assertNdSquareness(*arrays):
Chris@87 210 for a in arrays:
Chris@87 211 if max(a.shape[-2:]) != min(a.shape[-2:]):
Chris@87 212 raise LinAlgError('Last 2 dimensions of the array must be square')
Chris@87 213
Chris@87 214 def _assertFinite(*arrays):
Chris@87 215 for a in arrays:
Chris@87 216 if not (isfinite(a).all()):
Chris@87 217 raise LinAlgError("Array must not contain infs or NaNs")
Chris@87 218
Chris@87 219 def _assertNoEmpty2d(*arrays):
Chris@87 220 for a in arrays:
Chris@87 221 if a.size == 0 and product(a.shape[-2:]) == 0:
Chris@87 222 raise LinAlgError("Arrays cannot be empty")
Chris@87 223
Chris@87 224
Chris@87 225 # Linear equations
Chris@87 226
Chris@87 227 def tensorsolve(a, b, axes=None):
Chris@87 228 """
Chris@87 229 Solve the tensor equation ``a x = b`` for x.
Chris@87 230
Chris@87 231 It is assumed that all indices of `x` are summed over in the product,
Chris@87 232 together with the rightmost indices of `a`, as is done in, for example,
Chris@87 233 ``tensordot(a, x, axes=len(b.shape))``.
Chris@87 234
Chris@87 235 Parameters
Chris@87 236 ----------
Chris@87 237 a : array_like
Chris@87 238 Coefficient tensor, of shape ``b.shape + Q``. `Q`, a tuple, equals
Chris@87 239 the shape of that sub-tensor of `a` consisting of the appropriate
Chris@87 240 number of its rightmost indices, and must be such that
Chris@87 241 ``prod(Q) == prod(b.shape)`` (in which sense `a` is said to be
Chris@87 242 'square').
Chris@87 243 b : array_like
Chris@87 244 Right-hand tensor, which can be of any shape.
Chris@87 245 axes : tuple of ints, optional
Chris@87 246 Axes in `a` to reorder to the right, before inversion.
Chris@87 247 If None (default), no reordering is done.
Chris@87 248
Chris@87 249 Returns
Chris@87 250 -------
Chris@87 251 x : ndarray, shape Q
Chris@87 252
Chris@87 253 Raises
Chris@87 254 ------
Chris@87 255 LinAlgError
Chris@87 256 If `a` is singular or not 'square' (in the above sense).
Chris@87 257
Chris@87 258 See Also
Chris@87 259 --------
Chris@87 260 tensordot, tensorinv, einsum
Chris@87 261
Chris@87 262 Examples
Chris@87 263 --------
Chris@87 264 >>> a = np.eye(2*3*4)
Chris@87 265 >>> a.shape = (2*3, 4, 2, 3, 4)
Chris@87 266 >>> b = np.random.randn(2*3, 4)
Chris@87 267 >>> x = np.linalg.tensorsolve(a, b)
Chris@87 268 >>> x.shape
Chris@87 269 (2, 3, 4)
Chris@87 270 >>> np.allclose(np.tensordot(a, x, axes=3), b)
Chris@87 271 True
Chris@87 272
Chris@87 273 """
Chris@87 274 a, wrap = _makearray(a)
Chris@87 275 b = asarray(b)
Chris@87 276 an = a.ndim
Chris@87 277
Chris@87 278 if axes is not None:
Chris@87 279 allaxes = list(range(0, an))
Chris@87 280 for k in axes:
Chris@87 281 allaxes.remove(k)
Chris@87 282 allaxes.insert(an, k)
Chris@87 283 a = a.transpose(allaxes)
Chris@87 284
Chris@87 285 oldshape = a.shape[-(an-b.ndim):]
Chris@87 286 prod = 1
Chris@87 287 for k in oldshape:
Chris@87 288 prod *= k
Chris@87 289
Chris@87 290 a = a.reshape(-1, prod)
Chris@87 291 b = b.ravel()
Chris@87 292 res = wrap(solve(a, b))
Chris@87 293 res.shape = oldshape
Chris@87 294 return res
Chris@87 295
Chris@87 296 def solve(a, b):
Chris@87 297 """
Chris@87 298 Solve a linear matrix equation, or system of linear scalar equations.
Chris@87 299
Chris@87 300 Computes the "exact" solution, `x`, of the well-determined, i.e., full
Chris@87 301 rank, linear matrix equation `ax = b`.
Chris@87 302
Chris@87 303 Parameters
Chris@87 304 ----------
Chris@87 305 a : (..., M, M) array_like
Chris@87 306 Coefficient matrix.
Chris@87 307 b : {(..., M,), (..., M, K)}, array_like
Chris@87 308 Ordinate or "dependent variable" values.
Chris@87 309
Chris@87 310 Returns
Chris@87 311 -------
Chris@87 312 x : {(..., M,), (..., M, K)} ndarray
Chris@87 313 Solution to the system a x = b. Returned shape is identical to `b`.
Chris@87 314
Chris@87 315 Raises
Chris@87 316 ------
Chris@87 317 LinAlgError
Chris@87 318 If `a` is singular or not square.
Chris@87 319
Chris@87 320 Notes
Chris@87 321 -----
Chris@87 322 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 323 details.
Chris@87 324
Chris@87 325 The solutions are computed using LAPACK routine _gesv
Chris@87 326
Chris@87 327 `a` must be square and of full-rank, i.e., all rows (or, equivalently,
Chris@87 328 columns) must be linearly independent; if either is not true, use
Chris@87 329 `lstsq` for the least-squares best "solution" of the
Chris@87 330 system/equation.
Chris@87 331
Chris@87 332 References
Chris@87 333 ----------
Chris@87 334 .. [1] G. Strang, *Linear Algebra and Its Applications*, 2nd Ed., Orlando,
Chris@87 335 FL, Academic Press, Inc., 1980, pg. 22.
Chris@87 336
Chris@87 337 Examples
Chris@87 338 --------
Chris@87 339 Solve the system of equations ``3 * x0 + x1 = 9`` and ``x0 + 2 * x1 = 8``:
Chris@87 340
Chris@87 341 >>> a = np.array([[3,1], [1,2]])
Chris@87 342 >>> b = np.array([9,8])
Chris@87 343 >>> x = np.linalg.solve(a, b)
Chris@87 344 >>> x
Chris@87 345 array([ 2., 3.])
Chris@87 346
Chris@87 347 Check that the solution is correct:
Chris@87 348
Chris@87 349 >>> np.allclose(np.dot(a, x), b)
Chris@87 350 True
Chris@87 351
Chris@87 352 """
Chris@87 353 a, _ = _makearray(a)
Chris@87 354 _assertRankAtLeast2(a)
Chris@87 355 _assertNdSquareness(a)
Chris@87 356 b, wrap = _makearray(b)
Chris@87 357 t, result_t = _commonType(a, b)
Chris@87 358
Chris@87 359 # We use the b = (..., M,) logic, only if the number of extra dimensions
Chris@87 360 # match exactly
Chris@87 361 if b.ndim == a.ndim - 1:
Chris@87 362 if a.shape[-1] == 0 and b.shape[-1] == 0:
Chris@87 363 # Legal, but the ufunc cannot handle the 0-sized inner dims
Chris@87 364 # let the ufunc handle all wrong cases.
Chris@87 365 a = a.reshape(a.shape[:-1])
Chris@87 366 bc = broadcast(a, b)
Chris@87 367 return wrap(empty(bc.shape, dtype=result_t))
Chris@87 368
Chris@87 369 gufunc = _umath_linalg.solve1
Chris@87 370 else:
Chris@87 371 if b.size == 0:
Chris@87 372 if (a.shape[-1] == 0 and b.shape[-2] == 0) or b.shape[-1] == 0:
Chris@87 373 a = a[:,:1].reshape(a.shape[:-1] + (1,))
Chris@87 374 bc = broadcast(a, b)
Chris@87 375 return wrap(empty(bc.shape, dtype=result_t))
Chris@87 376
Chris@87 377 gufunc = _umath_linalg.solve
Chris@87 378
Chris@87 379 signature = 'DD->D' if isComplexType(t) else 'dd->d'
Chris@87 380 extobj = get_linalg_error_extobj(_raise_linalgerror_singular)
Chris@87 381 r = gufunc(a, b, signature=signature, extobj=extobj)
Chris@87 382
Chris@87 383 return wrap(r.astype(result_t))
Chris@87 384
Chris@87 385
Chris@87 386 def tensorinv(a, ind=2):
Chris@87 387 """
Chris@87 388 Compute the 'inverse' of an N-dimensional array.
Chris@87 389
Chris@87 390 The result is an inverse for `a` relative to the tensordot operation
Chris@87 391 ``tensordot(a, b, ind)``, i. e., up to floating-point accuracy,
Chris@87 392 ``tensordot(tensorinv(a), a, ind)`` is the "identity" tensor for the
Chris@87 393 tensordot operation.
Chris@87 394
Chris@87 395 Parameters
Chris@87 396 ----------
Chris@87 397 a : array_like
Chris@87 398 Tensor to 'invert'. Its shape must be 'square', i. e.,
Chris@87 399 ``prod(a.shape[:ind]) == prod(a.shape[ind:])``.
Chris@87 400 ind : int, optional
Chris@87 401 Number of first indices that are involved in the inverse sum.
Chris@87 402 Must be a positive integer, default is 2.
Chris@87 403
Chris@87 404 Returns
Chris@87 405 -------
Chris@87 406 b : ndarray
Chris@87 407 `a`'s tensordot inverse, shape ``a.shape[ind:] + a.shape[:ind]``.
Chris@87 408
Chris@87 409 Raises
Chris@87 410 ------
Chris@87 411 LinAlgError
Chris@87 412 If `a` is singular or not 'square' (in the above sense).
Chris@87 413
Chris@87 414 See Also
Chris@87 415 --------
Chris@87 416 tensordot, tensorsolve
Chris@87 417
Chris@87 418 Examples
Chris@87 419 --------
Chris@87 420 >>> a = np.eye(4*6)
Chris@87 421 >>> a.shape = (4, 6, 8, 3)
Chris@87 422 >>> ainv = np.linalg.tensorinv(a, ind=2)
Chris@87 423 >>> ainv.shape
Chris@87 424 (8, 3, 4, 6)
Chris@87 425 >>> b = np.random.randn(4, 6)
Chris@87 426 >>> np.allclose(np.tensordot(ainv, b), np.linalg.tensorsolve(a, b))
Chris@87 427 True
Chris@87 428
Chris@87 429 >>> a = np.eye(4*6)
Chris@87 430 >>> a.shape = (24, 8, 3)
Chris@87 431 >>> ainv = np.linalg.tensorinv(a, ind=1)
Chris@87 432 >>> ainv.shape
Chris@87 433 (8, 3, 24)
Chris@87 434 >>> b = np.random.randn(24)
Chris@87 435 >>> np.allclose(np.tensordot(ainv, b, 1), np.linalg.tensorsolve(a, b))
Chris@87 436 True
Chris@87 437
Chris@87 438 """
Chris@87 439 a = asarray(a)
Chris@87 440 oldshape = a.shape
Chris@87 441 prod = 1
Chris@87 442 if ind > 0:
Chris@87 443 invshape = oldshape[ind:] + oldshape[:ind]
Chris@87 444 for k in oldshape[ind:]:
Chris@87 445 prod *= k
Chris@87 446 else:
Chris@87 447 raise ValueError("Invalid ind argument.")
Chris@87 448 a = a.reshape(prod, -1)
Chris@87 449 ia = inv(a)
Chris@87 450 return ia.reshape(*invshape)
Chris@87 451
Chris@87 452
Chris@87 453 # Matrix inversion
Chris@87 454
Chris@87 455 def inv(a):
Chris@87 456 """
Chris@87 457 Compute the (multiplicative) inverse of a matrix.
Chris@87 458
Chris@87 459 Given a square matrix `a`, return the matrix `ainv` satisfying
Chris@87 460 ``dot(a, ainv) = dot(ainv, a) = eye(a.shape[0])``.
Chris@87 461
Chris@87 462 Parameters
Chris@87 463 ----------
Chris@87 464 a : (..., M, M) array_like
Chris@87 465 Matrix to be inverted.
Chris@87 466
Chris@87 467 Returns
Chris@87 468 -------
Chris@87 469 ainv : (..., M, M) ndarray or matrix
Chris@87 470 (Multiplicative) inverse of the matrix `a`.
Chris@87 471
Chris@87 472 Raises
Chris@87 473 ------
Chris@87 474 LinAlgError
Chris@87 475 If `a` is not square or inversion fails.
Chris@87 476
Chris@87 477 Notes
Chris@87 478 -----
Chris@87 479 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 480 details.
Chris@87 481
Chris@87 482 Examples
Chris@87 483 --------
Chris@87 484 >>> from numpy.linalg import inv
Chris@87 485 >>> a = np.array([[1., 2.], [3., 4.]])
Chris@87 486 >>> ainv = inv(a)
Chris@87 487 >>> np.allclose(np.dot(a, ainv), np.eye(2))
Chris@87 488 True
Chris@87 489 >>> np.allclose(np.dot(ainv, a), np.eye(2))
Chris@87 490 True
Chris@87 491
Chris@87 492 If a is a matrix object, then the return value is a matrix as well:
Chris@87 493
Chris@87 494 >>> ainv = inv(np.matrix(a))
Chris@87 495 >>> ainv
Chris@87 496 matrix([[-2. , 1. ],
Chris@87 497 [ 1.5, -0.5]])
Chris@87 498
Chris@87 499 Inverses of several matrices can be computed at once:
Chris@87 500
Chris@87 501 >>> a = np.array([[[1., 2.], [3., 4.]], [[1, 3], [3, 5]]])
Chris@87 502 >>> inv(a)
Chris@87 503 array([[[-2. , 1. ],
Chris@87 504 [ 1.5, -0.5]],
Chris@87 505 [[-5. , 2. ],
Chris@87 506 [ 3. , -1. ]]])
Chris@87 507
Chris@87 508 """
Chris@87 509 a, wrap = _makearray(a)
Chris@87 510 _assertRankAtLeast2(a)
Chris@87 511 _assertNdSquareness(a)
Chris@87 512 t, result_t = _commonType(a)
Chris@87 513
Chris@87 514 if a.shape[-1] == 0:
Chris@87 515 # The inner array is 0x0, the ufunc cannot handle this case
Chris@87 516 return wrap(empty_like(a, dtype=result_t))
Chris@87 517
Chris@87 518 signature = 'D->D' if isComplexType(t) else 'd->d'
Chris@87 519 extobj = get_linalg_error_extobj(_raise_linalgerror_singular)
Chris@87 520 ainv = _umath_linalg.inv(a, signature=signature, extobj=extobj)
Chris@87 521 return wrap(ainv.astype(result_t))
Chris@87 522
Chris@87 523
Chris@87 524 # Cholesky decomposition
Chris@87 525
Chris@87 526 def cholesky(a):
Chris@87 527 """
Chris@87 528 Cholesky decomposition.
Chris@87 529
Chris@87 530 Return the Cholesky decomposition, `L * L.H`, of the square matrix `a`,
Chris@87 531 where `L` is lower-triangular and .H is the conjugate transpose operator
Chris@87 532 (which is the ordinary transpose if `a` is real-valued). `a` must be
Chris@87 533 Hermitian (symmetric if real-valued) and positive-definite. Only `L` is
Chris@87 534 actually returned.
Chris@87 535
Chris@87 536 Parameters
Chris@87 537 ----------
Chris@87 538 a : (..., M, M) array_like
Chris@87 539 Hermitian (symmetric if all elements are real), positive-definite
Chris@87 540 input matrix.
Chris@87 541
Chris@87 542 Returns
Chris@87 543 -------
Chris@87 544 L : (..., M, M) array_like
Chris@87 545 Upper or lower-triangular Cholesky factor of `a`. Returns a
Chris@87 546 matrix object if `a` is a matrix object.
Chris@87 547
Chris@87 548 Raises
Chris@87 549 ------
Chris@87 550 LinAlgError
Chris@87 551 If the decomposition fails, for example, if `a` is not
Chris@87 552 positive-definite.
Chris@87 553
Chris@87 554 Notes
Chris@87 555 -----
Chris@87 556 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 557 details.
Chris@87 558
Chris@87 559 The Cholesky decomposition is often used as a fast way of solving
Chris@87 560
Chris@87 561 .. math:: A \\mathbf{x} = \\mathbf{b}
Chris@87 562
Chris@87 563 (when `A` is both Hermitian/symmetric and positive-definite).
Chris@87 564
Chris@87 565 First, we solve for :math:`\\mathbf{y}` in
Chris@87 566
Chris@87 567 .. math:: L \\mathbf{y} = \\mathbf{b},
Chris@87 568
Chris@87 569 and then for :math:`\\mathbf{x}` in
Chris@87 570
Chris@87 571 .. math:: L.H \\mathbf{x} = \\mathbf{y}.
Chris@87 572
Chris@87 573 Examples
Chris@87 574 --------
Chris@87 575 >>> A = np.array([[1,-2j],[2j,5]])
Chris@87 576 >>> A
Chris@87 577 array([[ 1.+0.j, 0.-2.j],
Chris@87 578 [ 0.+2.j, 5.+0.j]])
Chris@87 579 >>> L = np.linalg.cholesky(A)
Chris@87 580 >>> L
Chris@87 581 array([[ 1.+0.j, 0.+0.j],
Chris@87 582 [ 0.+2.j, 1.+0.j]])
Chris@87 583 >>> np.dot(L, L.T.conj()) # verify that L * L.H = A
Chris@87 584 array([[ 1.+0.j, 0.-2.j],
Chris@87 585 [ 0.+2.j, 5.+0.j]])
Chris@87 586 >>> A = [[1,-2j],[2j,5]] # what happens if A is only array_like?
Chris@87 587 >>> np.linalg.cholesky(A) # an ndarray object is returned
Chris@87 588 array([[ 1.+0.j, 0.+0.j],
Chris@87 589 [ 0.+2.j, 1.+0.j]])
Chris@87 590 >>> # But a matrix object is returned if A is a matrix object
Chris@87 591 >>> LA.cholesky(np.matrix(A))
Chris@87 592 matrix([[ 1.+0.j, 0.+0.j],
Chris@87 593 [ 0.+2.j, 1.+0.j]])
Chris@87 594
Chris@87 595 """
Chris@87 596 extobj = get_linalg_error_extobj(_raise_linalgerror_nonposdef)
Chris@87 597 gufunc = _umath_linalg.cholesky_lo
Chris@87 598 a, wrap = _makearray(a)
Chris@87 599 _assertRankAtLeast2(a)
Chris@87 600 _assertNdSquareness(a)
Chris@87 601 t, result_t = _commonType(a)
Chris@87 602 signature = 'D->D' if isComplexType(t) else 'd->d'
Chris@87 603 return wrap(gufunc(a, signature=signature, extobj=extobj).astype(result_t))
Chris@87 604
Chris@87 605 # QR decompostion
Chris@87 606
Chris@87 607 def qr(a, mode='reduced'):
Chris@87 608 """
Chris@87 609 Compute the qr factorization of a matrix.
Chris@87 610
Chris@87 611 Factor the matrix `a` as *qr*, where `q` is orthonormal and `r` is
Chris@87 612 upper-triangular.
Chris@87 613
Chris@87 614 Parameters
Chris@87 615 ----------
Chris@87 616 a : array_like, shape (M, N)
Chris@87 617 Matrix to be factored.
Chris@87 618 mode : {'reduced', 'complete', 'r', 'raw', 'full', 'economic'}, optional
Chris@87 619 If K = min(M, N), then
Chris@87 620
Chris@87 621 'reduced' : returns q, r with dimensions (M, K), (K, N) (default)
Chris@87 622 'complete' : returns q, r with dimensions (M, M), (M, N)
Chris@87 623 'r' : returns r only with dimensions (K, N)
Chris@87 624 'raw' : returns h, tau with dimensions (N, M), (K,)
Chris@87 625 'full' : alias of 'reduced', deprecated
Chris@87 626 'economic' : returns h from 'raw', deprecated.
Chris@87 627
Chris@87 628 The options 'reduced', 'complete, and 'raw' are new in numpy 1.8,
Chris@87 629 see the notes for more information. The default is 'reduced' and to
Chris@87 630 maintain backward compatibility with earlier versions of numpy both
Chris@87 631 it and the old default 'full' can be omitted. Note that array h
Chris@87 632 returned in 'raw' mode is transposed for calling Fortran. The
Chris@87 633 'economic' mode is deprecated. The modes 'full' and 'economic' may
Chris@87 634 be passed using only the first letter for backwards compatibility,
Chris@87 635 but all others must be spelled out. See the Notes for more
Chris@87 636 explanation.
Chris@87 637
Chris@87 638
Chris@87 639 Returns
Chris@87 640 -------
Chris@87 641 q : ndarray of float or complex, optional
Chris@87 642 A matrix with orthonormal columns. When mode = 'complete' the
Chris@87 643 result is an orthogonal/unitary matrix depending on whether or not
Chris@87 644 a is real/complex. The determinant may be either +/- 1 in that
Chris@87 645 case.
Chris@87 646 r : ndarray of float or complex, optional
Chris@87 647 The upper-triangular matrix.
Chris@87 648 (h, tau) : ndarrays of np.double or np.cdouble, optional
Chris@87 649 The array h contains the Householder reflectors that generate q
Chris@87 650 along with r. The tau array contains scaling factors for the
Chris@87 651 reflectors. In the deprecated 'economic' mode only h is returned.
Chris@87 652
Chris@87 653 Raises
Chris@87 654 ------
Chris@87 655 LinAlgError
Chris@87 656 If factoring fails.
Chris@87 657
Chris@87 658 Notes
Chris@87 659 -----
Chris@87 660 This is an interface to the LAPACK routines dgeqrf, zgeqrf,
Chris@87 661 dorgqr, and zungqr.
Chris@87 662
Chris@87 663 For more information on the qr factorization, see for example:
Chris@87 664 http://en.wikipedia.org/wiki/QR_factorization
Chris@87 665
Chris@87 666 Subclasses of `ndarray` are preserved except for the 'raw' mode. So if
Chris@87 667 `a` is of type `matrix`, all the return values will be matrices too.
Chris@87 668
Chris@87 669 New 'reduced', 'complete', and 'raw' options for mode were added in
Chris@87 670 Numpy 1.8 and the old option 'full' was made an alias of 'reduced'. In
Chris@87 671 addition the options 'full' and 'economic' were deprecated. Because
Chris@87 672 'full' was the previous default and 'reduced' is the new default,
Chris@87 673 backward compatibility can be maintained by letting `mode` default.
Chris@87 674 The 'raw' option was added so that LAPACK routines that can multiply
Chris@87 675 arrays by q using the Householder reflectors can be used. Note that in
Chris@87 676 this case the returned arrays are of type np.double or np.cdouble and
Chris@87 677 the h array is transposed to be FORTRAN compatible. No routines using
Chris@87 678 the 'raw' return are currently exposed by numpy, but some are available
Chris@87 679 in lapack_lite and just await the necessary work.
Chris@87 680
Chris@87 681 Examples
Chris@87 682 --------
Chris@87 683 >>> a = np.random.randn(9, 6)
Chris@87 684 >>> q, r = np.linalg.qr(a)
Chris@87 685 >>> np.allclose(a, np.dot(q, r)) # a does equal qr
Chris@87 686 True
Chris@87 687 >>> r2 = np.linalg.qr(a, mode='r')
Chris@87 688 >>> r3 = np.linalg.qr(a, mode='economic')
Chris@87 689 >>> np.allclose(r, r2) # mode='r' returns the same r as mode='full'
Chris@87 690 True
Chris@87 691 >>> # But only triu parts are guaranteed equal when mode='economic'
Chris@87 692 >>> np.allclose(r, np.triu(r3[:6,:6], k=0))
Chris@87 693 True
Chris@87 694
Chris@87 695 Example illustrating a common use of `qr`: solving of least squares
Chris@87 696 problems
Chris@87 697
Chris@87 698 What are the least-squares-best `m` and `y0` in ``y = y0 + mx`` for
Chris@87 699 the following data: {(0,1), (1,0), (1,2), (2,1)}. (Graph the points
Chris@87 700 and you'll see that it should be y0 = 0, m = 1.) The answer is provided
Chris@87 701 by solving the over-determined matrix equation ``Ax = b``, where::
Chris@87 702
Chris@87 703 A = array([[0, 1], [1, 1], [1, 1], [2, 1]])
Chris@87 704 x = array([[y0], [m]])
Chris@87 705 b = array([[1], [0], [2], [1]])
Chris@87 706
Chris@87 707 If A = qr such that q is orthonormal (which is always possible via
Chris@87 708 Gram-Schmidt), then ``x = inv(r) * (q.T) * b``. (In numpy practice,
Chris@87 709 however, we simply use `lstsq`.)
Chris@87 710
Chris@87 711 >>> A = np.array([[0, 1], [1, 1], [1, 1], [2, 1]])
Chris@87 712 >>> A
Chris@87 713 array([[0, 1],
Chris@87 714 [1, 1],
Chris@87 715 [1, 1],
Chris@87 716 [2, 1]])
Chris@87 717 >>> b = np.array([1, 0, 2, 1])
Chris@87 718 >>> q, r = LA.qr(A)
Chris@87 719 >>> p = np.dot(q.T, b)
Chris@87 720 >>> np.dot(LA.inv(r), p)
Chris@87 721 array([ 1.1e-16, 1.0e+00])
Chris@87 722
Chris@87 723 """
Chris@87 724 if mode not in ('reduced', 'complete', 'r', 'raw'):
Chris@87 725 if mode in ('f', 'full'):
Chris@87 726 msg = "".join((
Chris@87 727 "The 'full' option is deprecated in favor of 'reduced'.\n",
Chris@87 728 "For backward compatibility let mode default."))
Chris@87 729 warnings.warn(msg, DeprecationWarning)
Chris@87 730 mode = 'reduced'
Chris@87 731 elif mode in ('e', 'economic'):
Chris@87 732 msg = "The 'economic' option is deprecated.",
Chris@87 733 warnings.warn(msg, DeprecationWarning)
Chris@87 734 mode = 'economic'
Chris@87 735 else:
Chris@87 736 raise ValueError("Unrecognized mode '%s'" % mode)
Chris@87 737
Chris@87 738 a, wrap = _makearray(a)
Chris@87 739 _assertRank2(a)
Chris@87 740 _assertNoEmpty2d(a)
Chris@87 741 m, n = a.shape
Chris@87 742 t, result_t = _commonType(a)
Chris@87 743 a = _fastCopyAndTranspose(t, a)
Chris@87 744 a = _to_native_byte_order(a)
Chris@87 745 mn = min(m, n)
Chris@87 746 tau = zeros((mn,), t)
Chris@87 747 if isComplexType(t):
Chris@87 748 lapack_routine = lapack_lite.zgeqrf
Chris@87 749 routine_name = 'zgeqrf'
Chris@87 750 else:
Chris@87 751 lapack_routine = lapack_lite.dgeqrf
Chris@87 752 routine_name = 'dgeqrf'
Chris@87 753
Chris@87 754 # calculate optimal size of work data 'work'
Chris@87 755 lwork = 1
Chris@87 756 work = zeros((lwork,), t)
Chris@87 757 results = lapack_routine(m, n, a, m, tau, work, -1, 0)
Chris@87 758 if results['info'] != 0:
Chris@87 759 raise LinAlgError('%s returns %d' % (routine_name, results['info']))
Chris@87 760
Chris@87 761 # do qr decomposition
Chris@87 762 lwork = int(abs(work[0]))
Chris@87 763 work = zeros((lwork,), t)
Chris@87 764 results = lapack_routine(m, n, a, m, tau, work, lwork, 0)
Chris@87 765 if results['info'] != 0:
Chris@87 766 raise LinAlgError('%s returns %d' % (routine_name, results['info']))
Chris@87 767
Chris@87 768 # handle modes that don't return q
Chris@87 769 if mode == 'r':
Chris@87 770 r = _fastCopyAndTranspose(result_t, a[:, :mn])
Chris@87 771 return wrap(triu(r))
Chris@87 772
Chris@87 773 if mode == 'raw':
Chris@87 774 return a, tau
Chris@87 775
Chris@87 776 if mode == 'economic':
Chris@87 777 if t != result_t :
Chris@87 778 a = a.astype(result_t)
Chris@87 779 return wrap(a.T)
Chris@87 780
Chris@87 781 # generate q from a
Chris@87 782 if mode == 'complete' and m > n:
Chris@87 783 mc = m
Chris@87 784 q = empty((m, m), t)
Chris@87 785 else:
Chris@87 786 mc = mn
Chris@87 787 q = empty((n, m), t)
Chris@87 788 q[:n] = a
Chris@87 789
Chris@87 790 if isComplexType(t):
Chris@87 791 lapack_routine = lapack_lite.zungqr
Chris@87 792 routine_name = 'zungqr'
Chris@87 793 else:
Chris@87 794 lapack_routine = lapack_lite.dorgqr
Chris@87 795 routine_name = 'dorgqr'
Chris@87 796
Chris@87 797 # determine optimal lwork
Chris@87 798 lwork = 1
Chris@87 799 work = zeros((lwork,), t)
Chris@87 800 results = lapack_routine(m, mc, mn, q, m, tau, work, -1, 0)
Chris@87 801 if results['info'] != 0:
Chris@87 802 raise LinAlgError('%s returns %d' % (routine_name, results['info']))
Chris@87 803
Chris@87 804 # compute q
Chris@87 805 lwork = int(abs(work[0]))
Chris@87 806 work = zeros((lwork,), t)
Chris@87 807 results = lapack_routine(m, mc, mn, q, m, tau, work, lwork, 0)
Chris@87 808 if results['info'] != 0:
Chris@87 809 raise LinAlgError('%s returns %d' % (routine_name, results['info']))
Chris@87 810
Chris@87 811 q = _fastCopyAndTranspose(result_t, q[:mc])
Chris@87 812 r = _fastCopyAndTranspose(result_t, a[:, :mc])
Chris@87 813
Chris@87 814 return wrap(q), wrap(triu(r))
Chris@87 815
Chris@87 816
Chris@87 817 # Eigenvalues
Chris@87 818
Chris@87 819
Chris@87 820 def eigvals(a):
Chris@87 821 """
Chris@87 822 Compute the eigenvalues of a general matrix.
Chris@87 823
Chris@87 824 Main difference between `eigvals` and `eig`: the eigenvectors aren't
Chris@87 825 returned.
Chris@87 826
Chris@87 827 Parameters
Chris@87 828 ----------
Chris@87 829 a : (..., M, M) array_like
Chris@87 830 A complex- or real-valued matrix whose eigenvalues will be computed.
Chris@87 831
Chris@87 832 Returns
Chris@87 833 -------
Chris@87 834 w : (..., M,) ndarray
Chris@87 835 The eigenvalues, each repeated according to its multiplicity.
Chris@87 836 They are not necessarily ordered, nor are they necessarily
Chris@87 837 real for real matrices.
Chris@87 838
Chris@87 839 Raises
Chris@87 840 ------
Chris@87 841 LinAlgError
Chris@87 842 If the eigenvalue computation does not converge.
Chris@87 843
Chris@87 844 See Also
Chris@87 845 --------
Chris@87 846 eig : eigenvalues and right eigenvectors of general arrays
Chris@87 847 eigvalsh : eigenvalues of symmetric or Hermitian arrays.
Chris@87 848 eigh : eigenvalues and eigenvectors of symmetric/Hermitian arrays.
Chris@87 849
Chris@87 850 Notes
Chris@87 851 -----
Chris@87 852 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 853 details.
Chris@87 854
Chris@87 855 This is implemented using the _geev LAPACK routines which compute
Chris@87 856 the eigenvalues and eigenvectors of general square arrays.
Chris@87 857
Chris@87 858 Examples
Chris@87 859 --------
Chris@87 860 Illustration, using the fact that the eigenvalues of a diagonal matrix
Chris@87 861 are its diagonal elements, that multiplying a matrix on the left
Chris@87 862 by an orthogonal matrix, `Q`, and on the right by `Q.T` (the transpose
Chris@87 863 of `Q`), preserves the eigenvalues of the "middle" matrix. In other words,
Chris@87 864 if `Q` is orthogonal, then ``Q * A * Q.T`` has the same eigenvalues as
Chris@87 865 ``A``:
Chris@87 866
Chris@87 867 >>> from numpy import linalg as LA
Chris@87 868 >>> x = np.random.random()
Chris@87 869 >>> Q = np.array([[np.cos(x), -np.sin(x)], [np.sin(x), np.cos(x)]])
Chris@87 870 >>> LA.norm(Q[0, :]), LA.norm(Q[1, :]), np.dot(Q[0, :],Q[1, :])
Chris@87 871 (1.0, 1.0, 0.0)
Chris@87 872
Chris@87 873 Now multiply a diagonal matrix by Q on one side and by Q.T on the other:
Chris@87 874
Chris@87 875 >>> D = np.diag((-1,1))
Chris@87 876 >>> LA.eigvals(D)
Chris@87 877 array([-1., 1.])
Chris@87 878 >>> A = np.dot(Q, D)
Chris@87 879 >>> A = np.dot(A, Q.T)
Chris@87 880 >>> LA.eigvals(A)
Chris@87 881 array([ 1., -1.])
Chris@87 882
Chris@87 883 """
Chris@87 884 a, wrap = _makearray(a)
Chris@87 885 _assertNoEmpty2d(a)
Chris@87 886 _assertRankAtLeast2(a)
Chris@87 887 _assertNdSquareness(a)
Chris@87 888 _assertFinite(a)
Chris@87 889 t, result_t = _commonType(a)
Chris@87 890
Chris@87 891 extobj = get_linalg_error_extobj(
Chris@87 892 _raise_linalgerror_eigenvalues_nonconvergence)
Chris@87 893 signature = 'D->D' if isComplexType(t) else 'd->D'
Chris@87 894 w = _umath_linalg.eigvals(a, signature=signature, extobj=extobj)
Chris@87 895
Chris@87 896 if not isComplexType(t):
Chris@87 897 if all(w.imag == 0):
Chris@87 898 w = w.real
Chris@87 899 result_t = _realType(result_t)
Chris@87 900 else:
Chris@87 901 result_t = _complexType(result_t)
Chris@87 902
Chris@87 903 return w.astype(result_t)
Chris@87 904
Chris@87 905 def eigvalsh(a, UPLO='L'):
Chris@87 906 """
Chris@87 907 Compute the eigenvalues of a Hermitian or real symmetric matrix.
Chris@87 908
Chris@87 909 Main difference from eigh: the eigenvectors are not computed.
Chris@87 910
Chris@87 911 Parameters
Chris@87 912 ----------
Chris@87 913 a : (..., M, M) array_like
Chris@87 914 A complex- or real-valued matrix whose eigenvalues are to be
Chris@87 915 computed.
Chris@87 916 UPLO : {'L', 'U'}, optional
Chris@87 917 Same as `lower`, with 'L' for lower and 'U' for upper triangular.
Chris@87 918 Deprecated.
Chris@87 919
Chris@87 920 Returns
Chris@87 921 -------
Chris@87 922 w : (..., M,) ndarray
Chris@87 923 The eigenvalues, not necessarily ordered, each repeated according to
Chris@87 924 its multiplicity.
Chris@87 925
Chris@87 926 Raises
Chris@87 927 ------
Chris@87 928 LinAlgError
Chris@87 929 If the eigenvalue computation does not converge.
Chris@87 930
Chris@87 931 See Also
Chris@87 932 --------
Chris@87 933 eigh : eigenvalues and eigenvectors of symmetric/Hermitian arrays.
Chris@87 934 eigvals : eigenvalues of general real or complex arrays.
Chris@87 935 eig : eigenvalues and right eigenvectors of general real or complex
Chris@87 936 arrays.
Chris@87 937
Chris@87 938 Notes
Chris@87 939 -----
Chris@87 940 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 941 details.
Chris@87 942
Chris@87 943 The eigenvalues are computed using LAPACK routines _ssyevd, _heevd
Chris@87 944
Chris@87 945 Examples
Chris@87 946 --------
Chris@87 947 >>> from numpy import linalg as LA
Chris@87 948 >>> a = np.array([[1, -2j], [2j, 5]])
Chris@87 949 >>> LA.eigvalsh(a)
Chris@87 950 array([ 0.17157288+0.j, 5.82842712+0.j])
Chris@87 951
Chris@87 952 """
Chris@87 953 UPLO = UPLO.upper()
Chris@87 954 if UPLO not in ('L', 'U'):
Chris@87 955 raise ValueError("UPLO argument must be 'L' or 'U'")
Chris@87 956
Chris@87 957 extobj = get_linalg_error_extobj(
Chris@87 958 _raise_linalgerror_eigenvalues_nonconvergence)
Chris@87 959 if UPLO == 'L':
Chris@87 960 gufunc = _umath_linalg.eigvalsh_lo
Chris@87 961 else:
Chris@87 962 gufunc = _umath_linalg.eigvalsh_up
Chris@87 963
Chris@87 964 a, wrap = _makearray(a)
Chris@87 965 _assertNoEmpty2d(a)
Chris@87 966 _assertRankAtLeast2(a)
Chris@87 967 _assertNdSquareness(a)
Chris@87 968 t, result_t = _commonType(a)
Chris@87 969 signature = 'D->d' if isComplexType(t) else 'd->d'
Chris@87 970 w = gufunc(a, signature=signature, extobj=extobj)
Chris@87 971 return w.astype(_realType(result_t))
Chris@87 972
Chris@87 973 def _convertarray(a):
Chris@87 974 t, result_t = _commonType(a)
Chris@87 975 a = _fastCT(a.astype(t))
Chris@87 976 return a, t, result_t
Chris@87 977
Chris@87 978
Chris@87 979 # Eigenvectors
Chris@87 980
Chris@87 981
Chris@87 982 def eig(a):
Chris@87 983 """
Chris@87 984 Compute the eigenvalues and right eigenvectors of a square array.
Chris@87 985
Chris@87 986 Parameters
Chris@87 987 ----------
Chris@87 988 a : (..., M, M) array
Chris@87 989 Matrices for which the eigenvalues and right eigenvectors will
Chris@87 990 be computed
Chris@87 991
Chris@87 992 Returns
Chris@87 993 -------
Chris@87 994 w : (..., M) array
Chris@87 995 The eigenvalues, each repeated according to its multiplicity.
Chris@87 996 The eigenvalues are not necessarily ordered. The resulting
Chris@87 997 array will be always be of complex type. When `a` is real
Chris@87 998 the resulting eigenvalues will be real (0 imaginary part) or
Chris@87 999 occur in conjugate pairs
Chris@87 1000
Chris@87 1001 v : (..., M, M) array
Chris@87 1002 The normalized (unit "length") eigenvectors, such that the
Chris@87 1003 column ``v[:,i]`` is the eigenvector corresponding to the
Chris@87 1004 eigenvalue ``w[i]``.
Chris@87 1005
Chris@87 1006 Raises
Chris@87 1007 ------
Chris@87 1008 LinAlgError
Chris@87 1009 If the eigenvalue computation does not converge.
Chris@87 1010
Chris@87 1011 See Also
Chris@87 1012 --------
Chris@87 1013 eigvalsh : eigenvalues of a symmetric or Hermitian (conjugate symmetric)
Chris@87 1014 array.
Chris@87 1015
Chris@87 1016 eigvals : eigenvalues of a non-symmetric array.
Chris@87 1017
Chris@87 1018 Notes
Chris@87 1019 -----
Chris@87 1020 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 1021 details.
Chris@87 1022
Chris@87 1023 This is implemented using the _geev LAPACK routines which compute
Chris@87 1024 the eigenvalues and eigenvectors of general square arrays.
Chris@87 1025
Chris@87 1026 The number `w` is an eigenvalue of `a` if there exists a vector
Chris@87 1027 `v` such that ``dot(a,v) = w * v``. Thus, the arrays `a`, `w`, and
Chris@87 1028 `v` satisfy the equations ``dot(a[:,:], v[:,i]) = w[i] * v[:,i]``
Chris@87 1029 for :math:`i \\in \\{0,...,M-1\\}`.
Chris@87 1030
Chris@87 1031 The array `v` of eigenvectors may not be of maximum rank, that is, some
Chris@87 1032 of the columns may be linearly dependent, although round-off error may
Chris@87 1033 obscure that fact. If the eigenvalues are all different, then theoretically
Chris@87 1034 the eigenvectors are linearly independent. Likewise, the (complex-valued)
Chris@87 1035 matrix of eigenvectors `v` is unitary if the matrix `a` is normal, i.e.,
Chris@87 1036 if ``dot(a, a.H) = dot(a.H, a)``, where `a.H` denotes the conjugate
Chris@87 1037 transpose of `a`.
Chris@87 1038
Chris@87 1039 Finally, it is emphasized that `v` consists of the *right* (as in
Chris@87 1040 right-hand side) eigenvectors of `a`. A vector `y` satisfying
Chris@87 1041 ``dot(y.T, a) = z * y.T`` for some number `z` is called a *left*
Chris@87 1042 eigenvector of `a`, and, in general, the left and right eigenvectors
Chris@87 1043 of a matrix are not necessarily the (perhaps conjugate) transposes
Chris@87 1044 of each other.
Chris@87 1045
Chris@87 1046 References
Chris@87 1047 ----------
Chris@87 1048 G. Strang, *Linear Algebra and Its Applications*, 2nd Ed., Orlando, FL,
Chris@87 1049 Academic Press, Inc., 1980, Various pp.
Chris@87 1050
Chris@87 1051 Examples
Chris@87 1052 --------
Chris@87 1053 >>> from numpy import linalg as LA
Chris@87 1054
Chris@87 1055 (Almost) trivial example with real e-values and e-vectors.
Chris@87 1056
Chris@87 1057 >>> w, v = LA.eig(np.diag((1, 2, 3)))
Chris@87 1058 >>> w; v
Chris@87 1059 array([ 1., 2., 3.])
Chris@87 1060 array([[ 1., 0., 0.],
Chris@87 1061 [ 0., 1., 0.],
Chris@87 1062 [ 0., 0., 1.]])
Chris@87 1063
Chris@87 1064 Real matrix possessing complex e-values and e-vectors; note that the
Chris@87 1065 e-values are complex conjugates of each other.
Chris@87 1066
Chris@87 1067 >>> w, v = LA.eig(np.array([[1, -1], [1, 1]]))
Chris@87 1068 >>> w; v
Chris@87 1069 array([ 1. + 1.j, 1. - 1.j])
Chris@87 1070 array([[ 0.70710678+0.j , 0.70710678+0.j ],
Chris@87 1071 [ 0.00000000-0.70710678j, 0.00000000+0.70710678j]])
Chris@87 1072
Chris@87 1073 Complex-valued matrix with real e-values (but complex-valued e-vectors);
Chris@87 1074 note that a.conj().T = a, i.e., a is Hermitian.
Chris@87 1075
Chris@87 1076 >>> a = np.array([[1, 1j], [-1j, 1]])
Chris@87 1077 >>> w, v = LA.eig(a)
Chris@87 1078 >>> w; v
Chris@87 1079 array([ 2.00000000e+00+0.j, 5.98651912e-36+0.j]) # i.e., {2, 0}
Chris@87 1080 array([[ 0.00000000+0.70710678j, 0.70710678+0.j ],
Chris@87 1081 [ 0.70710678+0.j , 0.00000000+0.70710678j]])
Chris@87 1082
Chris@87 1083 Be careful about round-off error!
Chris@87 1084
Chris@87 1085 >>> a = np.array([[1 + 1e-9, 0], [0, 1 - 1e-9]])
Chris@87 1086 >>> # Theor. e-values are 1 +/- 1e-9
Chris@87 1087 >>> w, v = LA.eig(a)
Chris@87 1088 >>> w; v
Chris@87 1089 array([ 1., 1.])
Chris@87 1090 array([[ 1., 0.],
Chris@87 1091 [ 0., 1.]])
Chris@87 1092
Chris@87 1093 """
Chris@87 1094 a, wrap = _makearray(a)
Chris@87 1095 _assertRankAtLeast2(a)
Chris@87 1096 _assertNdSquareness(a)
Chris@87 1097 _assertFinite(a)
Chris@87 1098 t, result_t = _commonType(a)
Chris@87 1099
Chris@87 1100 extobj = get_linalg_error_extobj(
Chris@87 1101 _raise_linalgerror_eigenvalues_nonconvergence)
Chris@87 1102 signature = 'D->DD' if isComplexType(t) else 'd->DD'
Chris@87 1103 w, vt = _umath_linalg.eig(a, signature=signature, extobj=extobj)
Chris@87 1104
Chris@87 1105 if not isComplexType(t) and all(w.imag == 0.0):
Chris@87 1106 w = w.real
Chris@87 1107 vt = vt.real
Chris@87 1108 result_t = _realType(result_t)
Chris@87 1109 else:
Chris@87 1110 result_t = _complexType(result_t)
Chris@87 1111
Chris@87 1112 vt = vt.astype(result_t)
Chris@87 1113 return w.astype(result_t), wrap(vt)
Chris@87 1114
Chris@87 1115
Chris@87 1116 def eigh(a, UPLO='L'):
Chris@87 1117 """
Chris@87 1118 Return the eigenvalues and eigenvectors of a Hermitian or symmetric matrix.
Chris@87 1119
Chris@87 1120 Returns two objects, a 1-D array containing the eigenvalues of `a`, and
Chris@87 1121 a 2-D square array or matrix (depending on the input type) of the
Chris@87 1122 corresponding eigenvectors (in columns).
Chris@87 1123
Chris@87 1124 Parameters
Chris@87 1125 ----------
Chris@87 1126 A : (..., M, M) array
Chris@87 1127 Hermitian/Symmetric matrices whose eigenvalues and
Chris@87 1128 eigenvectors are to be computed.
Chris@87 1129 UPLO : {'L', 'U'}, optional
Chris@87 1130 Specifies whether the calculation is done with the lower triangular
Chris@87 1131 part of `a` ('L', default) or the upper triangular part ('U').
Chris@87 1132
Chris@87 1133 Returns
Chris@87 1134 -------
Chris@87 1135 w : (..., M) ndarray
Chris@87 1136 The eigenvalues, not necessarily ordered.
Chris@87 1137 v : {(..., M, M) ndarray, (..., M, M) matrix}
Chris@87 1138 The column ``v[:, i]`` is the normalized eigenvector corresponding
Chris@87 1139 to the eigenvalue ``w[i]``. Will return a matrix object if `a` is
Chris@87 1140 a matrix object.
Chris@87 1141
Chris@87 1142 Raises
Chris@87 1143 ------
Chris@87 1144 LinAlgError
Chris@87 1145 If the eigenvalue computation does not converge.
Chris@87 1146
Chris@87 1147 See Also
Chris@87 1148 --------
Chris@87 1149 eigvalsh : eigenvalues of symmetric or Hermitian arrays.
Chris@87 1150 eig : eigenvalues and right eigenvectors for non-symmetric arrays.
Chris@87 1151 eigvals : eigenvalues of non-symmetric arrays.
Chris@87 1152
Chris@87 1153 Notes
Chris@87 1154 -----
Chris@87 1155 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 1156 details.
Chris@87 1157
Chris@87 1158 The eigenvalues/eigenvectors are computed using LAPACK routines _ssyevd,
Chris@87 1159 _heevd
Chris@87 1160
Chris@87 1161 The eigenvalues of real symmetric or complex Hermitian matrices are
Chris@87 1162 always real. [1]_ The array `v` of (column) eigenvectors is unitary
Chris@87 1163 and `a`, `w`, and `v` satisfy the equations
Chris@87 1164 ``dot(a, v[:, i]) = w[i] * v[:, i]``.
Chris@87 1165
Chris@87 1166 References
Chris@87 1167 ----------
Chris@87 1168 .. [1] G. Strang, *Linear Algebra and Its Applications*, 2nd Ed., Orlando,
Chris@87 1169 FL, Academic Press, Inc., 1980, pg. 222.
Chris@87 1170
Chris@87 1171 Examples
Chris@87 1172 --------
Chris@87 1173 >>> from numpy import linalg as LA
Chris@87 1174 >>> a = np.array([[1, -2j], [2j, 5]])
Chris@87 1175 >>> a
Chris@87 1176 array([[ 1.+0.j, 0.-2.j],
Chris@87 1177 [ 0.+2.j, 5.+0.j]])
Chris@87 1178 >>> w, v = LA.eigh(a)
Chris@87 1179 >>> w; v
Chris@87 1180 array([ 0.17157288, 5.82842712])
Chris@87 1181 array([[-0.92387953+0.j , -0.38268343+0.j ],
Chris@87 1182 [ 0.00000000+0.38268343j, 0.00000000-0.92387953j]])
Chris@87 1183
Chris@87 1184 >>> np.dot(a, v[:, 0]) - w[0] * v[:, 0] # verify 1st e-val/vec pair
Chris@87 1185 array([2.77555756e-17 + 0.j, 0. + 1.38777878e-16j])
Chris@87 1186 >>> np.dot(a, v[:, 1]) - w[1] * v[:, 1] # verify 2nd e-val/vec pair
Chris@87 1187 array([ 0.+0.j, 0.+0.j])
Chris@87 1188
Chris@87 1189 >>> A = np.matrix(a) # what happens if input is a matrix object
Chris@87 1190 >>> A
Chris@87 1191 matrix([[ 1.+0.j, 0.-2.j],
Chris@87 1192 [ 0.+2.j, 5.+0.j]])
Chris@87 1193 >>> w, v = LA.eigh(A)
Chris@87 1194 >>> w; v
Chris@87 1195 array([ 0.17157288, 5.82842712])
Chris@87 1196 matrix([[-0.92387953+0.j , -0.38268343+0.j ],
Chris@87 1197 [ 0.00000000+0.38268343j, 0.00000000-0.92387953j]])
Chris@87 1198
Chris@87 1199 """
Chris@87 1200 UPLO = UPLO.upper()
Chris@87 1201 if UPLO not in ('L', 'U'):
Chris@87 1202 raise ValueError("UPLO argument must be 'L' or 'U'")
Chris@87 1203
Chris@87 1204 a, wrap = _makearray(a)
Chris@87 1205 _assertRankAtLeast2(a)
Chris@87 1206 _assertNdSquareness(a)
Chris@87 1207 t, result_t = _commonType(a)
Chris@87 1208
Chris@87 1209 extobj = get_linalg_error_extobj(
Chris@87 1210 _raise_linalgerror_eigenvalues_nonconvergence)
Chris@87 1211 if UPLO == 'L':
Chris@87 1212 gufunc = _umath_linalg.eigh_lo
Chris@87 1213 else:
Chris@87 1214 gufunc = _umath_linalg.eigh_up
Chris@87 1215
Chris@87 1216 signature = 'D->dD' if isComplexType(t) else 'd->dd'
Chris@87 1217 w, vt = gufunc(a, signature=signature, extobj=extobj)
Chris@87 1218 w = w.astype(_realType(result_t))
Chris@87 1219 vt = vt.astype(result_t)
Chris@87 1220 return w, wrap(vt)
Chris@87 1221
Chris@87 1222
Chris@87 1223 # Singular value decomposition
Chris@87 1224
Chris@87 1225 def svd(a, full_matrices=1, compute_uv=1):
Chris@87 1226 """
Chris@87 1227 Singular Value Decomposition.
Chris@87 1228
Chris@87 1229 Factors the matrix `a` as ``u * np.diag(s) * v``, where `u` and `v`
Chris@87 1230 are unitary and `s` is a 1-d array of `a`'s singular values.
Chris@87 1231
Chris@87 1232 Parameters
Chris@87 1233 ----------
Chris@87 1234 a : (..., M, N) array_like
Chris@87 1235 A real or complex matrix of shape (`M`, `N`) .
Chris@87 1236 full_matrices : bool, optional
Chris@87 1237 If True (default), `u` and `v` have the shapes (`M`, `M`) and
Chris@87 1238 (`N`, `N`), respectively. Otherwise, the shapes are (`M`, `K`)
Chris@87 1239 and (`K`, `N`), respectively, where `K` = min(`M`, `N`).
Chris@87 1240 compute_uv : bool, optional
Chris@87 1241 Whether or not to compute `u` and `v` in addition to `s`. True
Chris@87 1242 by default.
Chris@87 1243
Chris@87 1244 Returns
Chris@87 1245 -------
Chris@87 1246 u : { (..., M, M), (..., M, K) } array
Chris@87 1247 Unitary matrices. The actual shape depends on the value of
Chris@87 1248 ``full_matrices``. Only returned when ``compute_uv`` is True.
Chris@87 1249 s : (..., K) array
Chris@87 1250 The singular values for every matrix, sorted in descending order.
Chris@87 1251 v : { (..., N, N), (..., K, N) } array
Chris@87 1252 Unitary matrices. The actual shape depends on the value of
Chris@87 1253 ``full_matrices``. Only returned when ``compute_uv`` is True.
Chris@87 1254
Chris@87 1255 Raises
Chris@87 1256 ------
Chris@87 1257 LinAlgError
Chris@87 1258 If SVD computation does not converge.
Chris@87 1259
Chris@87 1260 Notes
Chris@87 1261 -----
Chris@87 1262 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 1263 details.
Chris@87 1264
Chris@87 1265 The decomposition is performed using LAPACK routine _gesdd
Chris@87 1266
Chris@87 1267 The SVD is commonly written as ``a = U S V.H``. The `v` returned
Chris@87 1268 by this function is ``V.H`` and ``u = U``.
Chris@87 1269
Chris@87 1270 If ``U`` is a unitary matrix, it means that it
Chris@87 1271 satisfies ``U.H = inv(U)``.
Chris@87 1272
Chris@87 1273 The rows of `v` are the eigenvectors of ``a.H a``. The columns
Chris@87 1274 of `u` are the eigenvectors of ``a a.H``. For row ``i`` in
Chris@87 1275 `v` and column ``i`` in `u`, the corresponding eigenvalue is
Chris@87 1276 ``s[i]**2``.
Chris@87 1277
Chris@87 1278 If `a` is a `matrix` object (as opposed to an `ndarray`), then so
Chris@87 1279 are all the return values.
Chris@87 1280
Chris@87 1281 Examples
Chris@87 1282 --------
Chris@87 1283 >>> a = np.random.randn(9, 6) + 1j*np.random.randn(9, 6)
Chris@87 1284
Chris@87 1285 Reconstruction based on full SVD:
Chris@87 1286
Chris@87 1287 >>> U, s, V = np.linalg.svd(a, full_matrices=True)
Chris@87 1288 >>> U.shape, V.shape, s.shape
Chris@87 1289 ((9, 9), (6, 6), (6,))
Chris@87 1290 >>> S = np.zeros((9, 6), dtype=complex)
Chris@87 1291 >>> S[:6, :6] = np.diag(s)
Chris@87 1292 >>> np.allclose(a, np.dot(U, np.dot(S, V)))
Chris@87 1293 True
Chris@87 1294
Chris@87 1295 Reconstruction based on reduced SVD:
Chris@87 1296
Chris@87 1297 >>> U, s, V = np.linalg.svd(a, full_matrices=False)
Chris@87 1298 >>> U.shape, V.shape, s.shape
Chris@87 1299 ((9, 6), (6, 6), (6,))
Chris@87 1300 >>> S = np.diag(s)
Chris@87 1301 >>> np.allclose(a, np.dot(U, np.dot(S, V)))
Chris@87 1302 True
Chris@87 1303
Chris@87 1304 """
Chris@87 1305 a, wrap = _makearray(a)
Chris@87 1306 _assertNoEmpty2d(a)
Chris@87 1307 _assertRankAtLeast2(a)
Chris@87 1308 t, result_t = _commonType(a)
Chris@87 1309
Chris@87 1310 extobj = get_linalg_error_extobj(_raise_linalgerror_svd_nonconvergence)
Chris@87 1311
Chris@87 1312 m = a.shape[-2]
Chris@87 1313 n = a.shape[-1]
Chris@87 1314 if compute_uv:
Chris@87 1315 if full_matrices:
Chris@87 1316 if m < n:
Chris@87 1317 gufunc = _umath_linalg.svd_m_f
Chris@87 1318 else:
Chris@87 1319 gufunc = _umath_linalg.svd_n_f
Chris@87 1320 else:
Chris@87 1321 if m < n:
Chris@87 1322 gufunc = _umath_linalg.svd_m_s
Chris@87 1323 else:
Chris@87 1324 gufunc = _umath_linalg.svd_n_s
Chris@87 1325
Chris@87 1326 signature = 'D->DdD' if isComplexType(t) else 'd->ddd'
Chris@87 1327 u, s, vt = gufunc(a, signature=signature, extobj=extobj)
Chris@87 1328 u = u.astype(result_t)
Chris@87 1329 s = s.astype(_realType(result_t))
Chris@87 1330 vt = vt.astype(result_t)
Chris@87 1331 return wrap(u), s, wrap(vt)
Chris@87 1332 else:
Chris@87 1333 if m < n:
Chris@87 1334 gufunc = _umath_linalg.svd_m
Chris@87 1335 else:
Chris@87 1336 gufunc = _umath_linalg.svd_n
Chris@87 1337
Chris@87 1338 signature = 'D->d' if isComplexType(t) else 'd->d'
Chris@87 1339 s = gufunc(a, signature=signature, extobj=extobj)
Chris@87 1340 s = s.astype(_realType(result_t))
Chris@87 1341 return s
Chris@87 1342
Chris@87 1343 def cond(x, p=None):
Chris@87 1344 """
Chris@87 1345 Compute the condition number of a matrix.
Chris@87 1346
Chris@87 1347 This function is capable of returning the condition number using
Chris@87 1348 one of seven different norms, depending on the value of `p` (see
Chris@87 1349 Parameters below).
Chris@87 1350
Chris@87 1351 Parameters
Chris@87 1352 ----------
Chris@87 1353 x : (M, N) array_like
Chris@87 1354 The matrix whose condition number is sought.
Chris@87 1355 p : {None, 1, -1, 2, -2, inf, -inf, 'fro'}, optional
Chris@87 1356 Order of the norm:
Chris@87 1357
Chris@87 1358 ===== ============================
Chris@87 1359 p norm for matrices
Chris@87 1360 ===== ============================
Chris@87 1361 None 2-norm, computed directly using the ``SVD``
Chris@87 1362 'fro' Frobenius norm
Chris@87 1363 inf max(sum(abs(x), axis=1))
Chris@87 1364 -inf min(sum(abs(x), axis=1))
Chris@87 1365 1 max(sum(abs(x), axis=0))
Chris@87 1366 -1 min(sum(abs(x), axis=0))
Chris@87 1367 2 2-norm (largest sing. value)
Chris@87 1368 -2 smallest singular value
Chris@87 1369 ===== ============================
Chris@87 1370
Chris@87 1371 inf means the numpy.inf object, and the Frobenius norm is
Chris@87 1372 the root-of-sum-of-squares norm.
Chris@87 1373
Chris@87 1374 Returns
Chris@87 1375 -------
Chris@87 1376 c : {float, inf}
Chris@87 1377 The condition number of the matrix. May be infinite.
Chris@87 1378
Chris@87 1379 See Also
Chris@87 1380 --------
Chris@87 1381 numpy.linalg.norm
Chris@87 1382
Chris@87 1383 Notes
Chris@87 1384 -----
Chris@87 1385 The condition number of `x` is defined as the norm of `x` times the
Chris@87 1386 norm of the inverse of `x` [1]_; the norm can be the usual L2-norm
Chris@87 1387 (root-of-sum-of-squares) or one of a number of other matrix norms.
Chris@87 1388
Chris@87 1389 References
Chris@87 1390 ----------
Chris@87 1391 .. [1] G. Strang, *Linear Algebra and Its Applications*, Orlando, FL,
Chris@87 1392 Academic Press, Inc., 1980, pg. 285.
Chris@87 1393
Chris@87 1394 Examples
Chris@87 1395 --------
Chris@87 1396 >>> from numpy import linalg as LA
Chris@87 1397 >>> a = np.array([[1, 0, -1], [0, 1, 0], [1, 0, 1]])
Chris@87 1398 >>> a
Chris@87 1399 array([[ 1, 0, -1],
Chris@87 1400 [ 0, 1, 0],
Chris@87 1401 [ 1, 0, 1]])
Chris@87 1402 >>> LA.cond(a)
Chris@87 1403 1.4142135623730951
Chris@87 1404 >>> LA.cond(a, 'fro')
Chris@87 1405 3.1622776601683795
Chris@87 1406 >>> LA.cond(a, np.inf)
Chris@87 1407 2.0
Chris@87 1408 >>> LA.cond(a, -np.inf)
Chris@87 1409 1.0
Chris@87 1410 >>> LA.cond(a, 1)
Chris@87 1411 2.0
Chris@87 1412 >>> LA.cond(a, -1)
Chris@87 1413 1.0
Chris@87 1414 >>> LA.cond(a, 2)
Chris@87 1415 1.4142135623730951
Chris@87 1416 >>> LA.cond(a, -2)
Chris@87 1417 0.70710678118654746
Chris@87 1418 >>> min(LA.svd(a, compute_uv=0))*min(LA.svd(LA.inv(a), compute_uv=0))
Chris@87 1419 0.70710678118654746
Chris@87 1420
Chris@87 1421 """
Chris@87 1422 x = asarray(x) # in case we have a matrix
Chris@87 1423 if p is None:
Chris@87 1424 s = svd(x, compute_uv=False)
Chris@87 1425 return s[0]/s[-1]
Chris@87 1426 else:
Chris@87 1427 return norm(x, p)*norm(inv(x), p)
Chris@87 1428
Chris@87 1429
Chris@87 1430 def matrix_rank(M, tol=None):
Chris@87 1431 """
Chris@87 1432 Return matrix rank of array using SVD method
Chris@87 1433
Chris@87 1434 Rank of the array is the number of SVD singular values of the array that are
Chris@87 1435 greater than `tol`.
Chris@87 1436
Chris@87 1437 Parameters
Chris@87 1438 ----------
Chris@87 1439 M : {(M,), (M, N)} array_like
Chris@87 1440 array of <=2 dimensions
Chris@87 1441 tol : {None, float}, optional
Chris@87 1442 threshold below which SVD values are considered zero. If `tol` is
Chris@87 1443 None, and ``S`` is an array with singular values for `M`, and
Chris@87 1444 ``eps`` is the epsilon value for datatype of ``S``, then `tol` is
Chris@87 1445 set to ``S.max() * max(M.shape) * eps``.
Chris@87 1446
Chris@87 1447 Notes
Chris@87 1448 -----
Chris@87 1449 The default threshold to detect rank deficiency is a test on the magnitude
Chris@87 1450 of the singular values of `M`. By default, we identify singular values less
Chris@87 1451 than ``S.max() * max(M.shape) * eps`` as indicating rank deficiency (with
Chris@87 1452 the symbols defined above). This is the algorithm MATLAB uses [1]. It also
Chris@87 1453 appears in *Numerical recipes* in the discussion of SVD solutions for linear
Chris@87 1454 least squares [2].
Chris@87 1455
Chris@87 1456 This default threshold is designed to detect rank deficiency accounting for
Chris@87 1457 the numerical errors of the SVD computation. Imagine that there is a column
Chris@87 1458 in `M` that is an exact (in floating point) linear combination of other
Chris@87 1459 columns in `M`. Computing the SVD on `M` will not produce a singular value
Chris@87 1460 exactly equal to 0 in general: any difference of the smallest SVD value from
Chris@87 1461 0 will be caused by numerical imprecision in the calculation of the SVD.
Chris@87 1462 Our threshold for small SVD values takes this numerical imprecision into
Chris@87 1463 account, and the default threshold will detect such numerical rank
Chris@87 1464 deficiency. The threshold may declare a matrix `M` rank deficient even if
Chris@87 1465 the linear combination of some columns of `M` is not exactly equal to
Chris@87 1466 another column of `M` but only numerically very close to another column of
Chris@87 1467 `M`.
Chris@87 1468
Chris@87 1469 We chose our default threshold because it is in wide use. Other thresholds
Chris@87 1470 are possible. For example, elsewhere in the 2007 edition of *Numerical
Chris@87 1471 recipes* there is an alternative threshold of ``S.max() *
Chris@87 1472 np.finfo(M.dtype).eps / 2. * np.sqrt(m + n + 1.)``. The authors describe
Chris@87 1473 this threshold as being based on "expected roundoff error" (p 71).
Chris@87 1474
Chris@87 1475 The thresholds above deal with floating point roundoff error in the
Chris@87 1476 calculation of the SVD. However, you may have more information about the
Chris@87 1477 sources of error in `M` that would make you consider other tolerance values
Chris@87 1478 to detect *effective* rank deficiency. The most useful measure of the
Chris@87 1479 tolerance depends on the operations you intend to use on your matrix. For
Chris@87 1480 example, if your data come from uncertain measurements with uncertainties
Chris@87 1481 greater than floating point epsilon, choosing a tolerance near that
Chris@87 1482 uncertainty may be preferable. The tolerance may be absolute if the
Chris@87 1483 uncertainties are absolute rather than relative.
Chris@87 1484
Chris@87 1485 References
Chris@87 1486 ----------
Chris@87 1487 .. [1] MATLAB reference documention, "Rank"
Chris@87 1488 http://www.mathworks.com/help/techdoc/ref/rank.html
Chris@87 1489 .. [2] W. H. Press, S. A. Teukolsky, W. T. Vetterling and B. P. Flannery,
Chris@87 1490 "Numerical Recipes (3rd edition)", Cambridge University Press, 2007,
Chris@87 1491 page 795.
Chris@87 1492
Chris@87 1493 Examples
Chris@87 1494 --------
Chris@87 1495 >>> from numpy.linalg import matrix_rank
Chris@87 1496 >>> matrix_rank(np.eye(4)) # Full rank matrix
Chris@87 1497 4
Chris@87 1498 >>> I=np.eye(4); I[-1,-1] = 0. # rank deficient matrix
Chris@87 1499 >>> matrix_rank(I)
Chris@87 1500 3
Chris@87 1501 >>> matrix_rank(np.ones((4,))) # 1 dimension - rank 1 unless all 0
Chris@87 1502 1
Chris@87 1503 >>> matrix_rank(np.zeros((4,)))
Chris@87 1504 0
Chris@87 1505 """
Chris@87 1506 M = asarray(M)
Chris@87 1507 if M.ndim > 2:
Chris@87 1508 raise TypeError('array should have 2 or fewer dimensions')
Chris@87 1509 if M.ndim < 2:
Chris@87 1510 return int(not all(M==0))
Chris@87 1511 S = svd(M, compute_uv=False)
Chris@87 1512 if tol is None:
Chris@87 1513 tol = S.max() * max(M.shape) * finfo(S.dtype).eps
Chris@87 1514 return sum(S > tol)
Chris@87 1515
Chris@87 1516
Chris@87 1517 # Generalized inverse
Chris@87 1518
Chris@87 1519 def pinv(a, rcond=1e-15 ):
Chris@87 1520 """
Chris@87 1521 Compute the (Moore-Penrose) pseudo-inverse of a matrix.
Chris@87 1522
Chris@87 1523 Calculate the generalized inverse of a matrix using its
Chris@87 1524 singular-value decomposition (SVD) and including all
Chris@87 1525 *large* singular values.
Chris@87 1526
Chris@87 1527 Parameters
Chris@87 1528 ----------
Chris@87 1529 a : (M, N) array_like
Chris@87 1530 Matrix to be pseudo-inverted.
Chris@87 1531 rcond : float
Chris@87 1532 Cutoff for small singular values.
Chris@87 1533 Singular values smaller (in modulus) than
Chris@87 1534 `rcond` * largest_singular_value (again, in modulus)
Chris@87 1535 are set to zero.
Chris@87 1536
Chris@87 1537 Returns
Chris@87 1538 -------
Chris@87 1539 B : (N, M) ndarray
Chris@87 1540 The pseudo-inverse of `a`. If `a` is a `matrix` instance, then so
Chris@87 1541 is `B`.
Chris@87 1542
Chris@87 1543 Raises
Chris@87 1544 ------
Chris@87 1545 LinAlgError
Chris@87 1546 If the SVD computation does not converge.
Chris@87 1547
Chris@87 1548 Notes
Chris@87 1549 -----
Chris@87 1550 The pseudo-inverse of a matrix A, denoted :math:`A^+`, is
Chris@87 1551 defined as: "the matrix that 'solves' [the least-squares problem]
Chris@87 1552 :math:`Ax = b`," i.e., if :math:`\\bar{x}` is said solution, then
Chris@87 1553 :math:`A^+` is that matrix such that :math:`\\bar{x} = A^+b`.
Chris@87 1554
Chris@87 1555 It can be shown that if :math:`Q_1 \\Sigma Q_2^T = A` is the singular
Chris@87 1556 value decomposition of A, then
Chris@87 1557 :math:`A^+ = Q_2 \\Sigma^+ Q_1^T`, where :math:`Q_{1,2}` are
Chris@87 1558 orthogonal matrices, :math:`\\Sigma` is a diagonal matrix consisting
Chris@87 1559 of A's so-called singular values, (followed, typically, by
Chris@87 1560 zeros), and then :math:`\\Sigma^+` is simply the diagonal matrix
Chris@87 1561 consisting of the reciprocals of A's singular values
Chris@87 1562 (again, followed by zeros). [1]_
Chris@87 1563
Chris@87 1564 References
Chris@87 1565 ----------
Chris@87 1566 .. [1] G. Strang, *Linear Algebra and Its Applications*, 2nd Ed., Orlando,
Chris@87 1567 FL, Academic Press, Inc., 1980, pp. 139-142.
Chris@87 1568
Chris@87 1569 Examples
Chris@87 1570 --------
Chris@87 1571 The following example checks that ``a * a+ * a == a`` and
Chris@87 1572 ``a+ * a * a+ == a+``:
Chris@87 1573
Chris@87 1574 >>> a = np.random.randn(9, 6)
Chris@87 1575 >>> B = np.linalg.pinv(a)
Chris@87 1576 >>> np.allclose(a, np.dot(a, np.dot(B, a)))
Chris@87 1577 True
Chris@87 1578 >>> np.allclose(B, np.dot(B, np.dot(a, B)))
Chris@87 1579 True
Chris@87 1580
Chris@87 1581 """
Chris@87 1582 a, wrap = _makearray(a)
Chris@87 1583 _assertNoEmpty2d(a)
Chris@87 1584 a = a.conjugate()
Chris@87 1585 u, s, vt = svd(a, 0)
Chris@87 1586 m = u.shape[0]
Chris@87 1587 n = vt.shape[1]
Chris@87 1588 cutoff = rcond*maximum.reduce(s)
Chris@87 1589 for i in range(min(n, m)):
Chris@87 1590 if s[i] > cutoff:
Chris@87 1591 s[i] = 1./s[i]
Chris@87 1592 else:
Chris@87 1593 s[i] = 0.;
Chris@87 1594 res = dot(transpose(vt), multiply(s[:, newaxis], transpose(u)))
Chris@87 1595 return wrap(res)
Chris@87 1596
Chris@87 1597 # Determinant
Chris@87 1598
Chris@87 1599 def slogdet(a):
Chris@87 1600 """
Chris@87 1601 Compute the sign and (natural) logarithm of the determinant of an array.
Chris@87 1602
Chris@87 1603 If an array has a very small or very large determinant, than a call to
Chris@87 1604 `det` may overflow or underflow. This routine is more robust against such
Chris@87 1605 issues, because it computes the logarithm of the determinant rather than
Chris@87 1606 the determinant itself.
Chris@87 1607
Chris@87 1608 Parameters
Chris@87 1609 ----------
Chris@87 1610 a : (..., M, M) array_like
Chris@87 1611 Input array, has to be a square 2-D array.
Chris@87 1612
Chris@87 1613 Returns
Chris@87 1614 -------
Chris@87 1615 sign : (...) array_like
Chris@87 1616 A number representing the sign of the determinant. For a real matrix,
Chris@87 1617 this is 1, 0, or -1. For a complex matrix, this is a complex number
Chris@87 1618 with absolute value 1 (i.e., it is on the unit circle), or else 0.
Chris@87 1619 logdet : (...) array_like
Chris@87 1620 The natural log of the absolute value of the determinant.
Chris@87 1621
Chris@87 1622 If the determinant is zero, then `sign` will be 0 and `logdet` will be
Chris@87 1623 -Inf. In all cases, the determinant is equal to ``sign * np.exp(logdet)``.
Chris@87 1624
Chris@87 1625 See Also
Chris@87 1626 --------
Chris@87 1627 det
Chris@87 1628
Chris@87 1629 Notes
Chris@87 1630 -----
Chris@87 1631 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 1632 details.
Chris@87 1633
Chris@87 1634 The determinant is computed via LU factorization using the LAPACK
Chris@87 1635 routine z/dgetrf.
Chris@87 1636
Chris@87 1637 .. versionadded:: 1.6.0.
Chris@87 1638
Chris@87 1639 Examples
Chris@87 1640 --------
Chris@87 1641 The determinant of a 2-D array ``[[a, b], [c, d]]`` is ``ad - bc``:
Chris@87 1642
Chris@87 1643 >>> a = np.array([[1, 2], [3, 4]])
Chris@87 1644 >>> (sign, logdet) = np.linalg.slogdet(a)
Chris@87 1645 >>> (sign, logdet)
Chris@87 1646 (-1, 0.69314718055994529)
Chris@87 1647 >>> sign * np.exp(logdet)
Chris@87 1648 -2.0
Chris@87 1649
Chris@87 1650 Computing log-determinants for a stack of matrices:
Chris@87 1651
Chris@87 1652 >>> a = np.array([ [[1, 2], [3, 4]], [[1, 2], [2, 1]], [[1, 3], [3, 1]] ])
Chris@87 1653 >>> a.shape
Chris@87 1654 (3, 2, 2)
Chris@87 1655 >>> sign, logdet = np.linalg.slogdet(a)
Chris@87 1656 >>> (sign, logdet)
Chris@87 1657 (array([-1., -1., -1.]), array([ 0.69314718, 1.09861229, 2.07944154]))
Chris@87 1658 >>> sign * np.exp(logdet)
Chris@87 1659 array([-2., -3., -8.])
Chris@87 1660
Chris@87 1661 This routine succeeds where ordinary `det` does not:
Chris@87 1662
Chris@87 1663 >>> np.linalg.det(np.eye(500) * 0.1)
Chris@87 1664 0.0
Chris@87 1665 >>> np.linalg.slogdet(np.eye(500) * 0.1)
Chris@87 1666 (1, -1151.2925464970228)
Chris@87 1667
Chris@87 1668 """
Chris@87 1669 a = asarray(a)
Chris@87 1670 _assertNoEmpty2d(a)
Chris@87 1671 _assertRankAtLeast2(a)
Chris@87 1672 _assertNdSquareness(a)
Chris@87 1673 t, result_t = _commonType(a)
Chris@87 1674 real_t = _realType(result_t)
Chris@87 1675 signature = 'D->Dd' if isComplexType(t) else 'd->dd'
Chris@87 1676 sign, logdet = _umath_linalg.slogdet(a, signature=signature)
Chris@87 1677 return sign.astype(result_t), logdet.astype(real_t)
Chris@87 1678
Chris@87 1679 def det(a):
Chris@87 1680 """
Chris@87 1681 Compute the determinant of an array.
Chris@87 1682
Chris@87 1683 Parameters
Chris@87 1684 ----------
Chris@87 1685 a : (..., M, M) array_like
Chris@87 1686 Input array to compute determinants for.
Chris@87 1687
Chris@87 1688 Returns
Chris@87 1689 -------
Chris@87 1690 det : (...) array_like
Chris@87 1691 Determinant of `a`.
Chris@87 1692
Chris@87 1693 See Also
Chris@87 1694 --------
Chris@87 1695 slogdet : Another way to representing the determinant, more suitable
Chris@87 1696 for large matrices where underflow/overflow may occur.
Chris@87 1697
Chris@87 1698 Notes
Chris@87 1699 -----
Chris@87 1700 Broadcasting rules apply, see the `numpy.linalg` documentation for
Chris@87 1701 details.
Chris@87 1702
Chris@87 1703 The determinant is computed via LU factorization using the LAPACK
Chris@87 1704 routine z/dgetrf.
Chris@87 1705
Chris@87 1706 Examples
Chris@87 1707 --------
Chris@87 1708 The determinant of a 2-D array [[a, b], [c, d]] is ad - bc:
Chris@87 1709
Chris@87 1710 >>> a = np.array([[1, 2], [3, 4]])
Chris@87 1711 >>> np.linalg.det(a)
Chris@87 1712 -2.0
Chris@87 1713
Chris@87 1714 Computing determinants for a stack of matrices:
Chris@87 1715
Chris@87 1716 >>> a = np.array([ [[1, 2], [3, 4]], [[1, 2], [2, 1]], [[1, 3], [3, 1]] ])
Chris@87 1717 >>> a.shape
Chris@87 1718 (2, 2, 2
Chris@87 1719 >>> np.linalg.det(a)
Chris@87 1720 array([-2., -3., -8.])
Chris@87 1721
Chris@87 1722 """
Chris@87 1723 a = asarray(a)
Chris@87 1724 _assertNoEmpty2d(a)
Chris@87 1725 _assertRankAtLeast2(a)
Chris@87 1726 _assertNdSquareness(a)
Chris@87 1727 t, result_t = _commonType(a)
Chris@87 1728 signature = 'D->D' if isComplexType(t) else 'd->d'
Chris@87 1729 return _umath_linalg.det(a, signature=signature).astype(result_t)
Chris@87 1730
Chris@87 1731 # Linear Least Squares
Chris@87 1732
Chris@87 1733 def lstsq(a, b, rcond=-1):
Chris@87 1734 """
Chris@87 1735 Return the least-squares solution to a linear matrix equation.
Chris@87 1736
Chris@87 1737 Solves the equation `a x = b` by computing a vector `x` that
Chris@87 1738 minimizes the Euclidean 2-norm `|| b - a x ||^2`. The equation may
Chris@87 1739 be under-, well-, or over- determined (i.e., the number of
Chris@87 1740 linearly independent rows of `a` can be less than, equal to, or
Chris@87 1741 greater than its number of linearly independent columns). If `a`
Chris@87 1742 is square and of full rank, then `x` (but for round-off error) is
Chris@87 1743 the "exact" solution of the equation.
Chris@87 1744
Chris@87 1745 Parameters
Chris@87 1746 ----------
Chris@87 1747 a : (M, N) array_like
Chris@87 1748 "Coefficient" matrix.
Chris@87 1749 b : {(M,), (M, K)} array_like
Chris@87 1750 Ordinate or "dependent variable" values. If `b` is two-dimensional,
Chris@87 1751 the least-squares solution is calculated for each of the `K` columns
Chris@87 1752 of `b`.
Chris@87 1753 rcond : float, optional
Chris@87 1754 Cut-off ratio for small singular values of `a`.
Chris@87 1755 Singular values are set to zero if they are smaller than `rcond`
Chris@87 1756 times the largest singular value of `a`.
Chris@87 1757
Chris@87 1758 Returns
Chris@87 1759 -------
Chris@87 1760 x : {(N,), (N, K)} ndarray
Chris@87 1761 Least-squares solution. If `b` is two-dimensional,
Chris@87 1762 the solutions are in the `K` columns of `x`.
Chris@87 1763 residuals : {(), (1,), (K,)} ndarray
Chris@87 1764 Sums of residuals; squared Euclidean 2-norm for each column in
Chris@87 1765 ``b - a*x``.
Chris@87 1766 If the rank of `a` is < N or M <= N, this is an empty array.
Chris@87 1767 If `b` is 1-dimensional, this is a (1,) shape array.
Chris@87 1768 Otherwise the shape is (K,).
Chris@87 1769 rank : int
Chris@87 1770 Rank of matrix `a`.
Chris@87 1771 s : (min(M, N),) ndarray
Chris@87 1772 Singular values of `a`.
Chris@87 1773
Chris@87 1774 Raises
Chris@87 1775 ------
Chris@87 1776 LinAlgError
Chris@87 1777 If computation does not converge.
Chris@87 1778
Chris@87 1779 Notes
Chris@87 1780 -----
Chris@87 1781 If `b` is a matrix, then all array results are returned as matrices.
Chris@87 1782
Chris@87 1783 Examples
Chris@87 1784 --------
Chris@87 1785 Fit a line, ``y = mx + c``, through some noisy data-points:
Chris@87 1786
Chris@87 1787 >>> x = np.array([0, 1, 2, 3])
Chris@87 1788 >>> y = np.array([-1, 0.2, 0.9, 2.1])
Chris@87 1789
Chris@87 1790 By examining the coefficients, we see that the line should have a
Chris@87 1791 gradient of roughly 1 and cut the y-axis at, more or less, -1.
Chris@87 1792
Chris@87 1793 We can rewrite the line equation as ``y = Ap``, where ``A = [[x 1]]``
Chris@87 1794 and ``p = [[m], [c]]``. Now use `lstsq` to solve for `p`:
Chris@87 1795
Chris@87 1796 >>> A = np.vstack([x, np.ones(len(x))]).T
Chris@87 1797 >>> A
Chris@87 1798 array([[ 0., 1.],
Chris@87 1799 [ 1., 1.],
Chris@87 1800 [ 2., 1.],
Chris@87 1801 [ 3., 1.]])
Chris@87 1802
Chris@87 1803 >>> m, c = np.linalg.lstsq(A, y)[0]
Chris@87 1804 >>> print m, c
Chris@87 1805 1.0 -0.95
Chris@87 1806
Chris@87 1807 Plot the data along with the fitted line:
Chris@87 1808
Chris@87 1809 >>> import matplotlib.pyplot as plt
Chris@87 1810 >>> plt.plot(x, y, 'o', label='Original data', markersize=10)
Chris@87 1811 >>> plt.plot(x, m*x + c, 'r', label='Fitted line')
Chris@87 1812 >>> plt.legend()
Chris@87 1813 >>> plt.show()
Chris@87 1814
Chris@87 1815 """
Chris@87 1816 import math
Chris@87 1817 a, _ = _makearray(a)
Chris@87 1818 b, wrap = _makearray(b)
Chris@87 1819 is_1d = len(b.shape) == 1
Chris@87 1820 if is_1d:
Chris@87 1821 b = b[:, newaxis]
Chris@87 1822 _assertRank2(a, b)
Chris@87 1823 m = a.shape[0]
Chris@87 1824 n = a.shape[1]
Chris@87 1825 n_rhs = b.shape[1]
Chris@87 1826 ldb = max(n, m)
Chris@87 1827 if m != b.shape[0]:
Chris@87 1828 raise LinAlgError('Incompatible dimensions')
Chris@87 1829 t, result_t = _commonType(a, b)
Chris@87 1830 result_real_t = _realType(result_t)
Chris@87 1831 real_t = _linalgRealType(t)
Chris@87 1832 bstar = zeros((ldb, n_rhs), t)
Chris@87 1833 bstar[:b.shape[0], :n_rhs] = b.copy()
Chris@87 1834 a, bstar = _fastCopyAndTranspose(t, a, bstar)
Chris@87 1835 a, bstar = _to_native_byte_order(a, bstar)
Chris@87 1836 s = zeros((min(m, n),), real_t)
Chris@87 1837 nlvl = max( 0, int( math.log( float(min(m, n))/2. ) ) + 1 )
Chris@87 1838 iwork = zeros((3*min(m, n)*nlvl+11*min(m, n),), fortran_int)
Chris@87 1839 if isComplexType(t):
Chris@87 1840 lapack_routine = lapack_lite.zgelsd
Chris@87 1841 lwork = 1
Chris@87 1842 rwork = zeros((lwork,), real_t)
Chris@87 1843 work = zeros((lwork,), t)
Chris@87 1844 results = lapack_routine(m, n, n_rhs, a, m, bstar, ldb, s, rcond,
Chris@87 1845 0, work, -1, rwork, iwork, 0)
Chris@87 1846 lwork = int(abs(work[0]))
Chris@87 1847 rwork = zeros((lwork,), real_t)
Chris@87 1848 a_real = zeros((m, n), real_t)
Chris@87 1849 bstar_real = zeros((ldb, n_rhs,), real_t)
Chris@87 1850 results = lapack_lite.dgelsd(m, n, n_rhs, a_real, m,
Chris@87 1851 bstar_real, ldb, s, rcond,
Chris@87 1852 0, rwork, -1, iwork, 0)
Chris@87 1853 lrwork = int(rwork[0])
Chris@87 1854 work = zeros((lwork,), t)
Chris@87 1855 rwork = zeros((lrwork,), real_t)
Chris@87 1856 results = lapack_routine(m, n, n_rhs, a, m, bstar, ldb, s, rcond,
Chris@87 1857 0, work, lwork, rwork, iwork, 0)
Chris@87 1858 else:
Chris@87 1859 lapack_routine = lapack_lite.dgelsd
Chris@87 1860 lwork = 1
Chris@87 1861 work = zeros((lwork,), t)
Chris@87 1862 results = lapack_routine(m, n, n_rhs, a, m, bstar, ldb, s, rcond,
Chris@87 1863 0, work, -1, iwork, 0)
Chris@87 1864 lwork = int(work[0])
Chris@87 1865 work = zeros((lwork,), t)
Chris@87 1866 results = lapack_routine(m, n, n_rhs, a, m, bstar, ldb, s, rcond,
Chris@87 1867 0, work, lwork, iwork, 0)
Chris@87 1868 if results['info'] > 0:
Chris@87 1869 raise LinAlgError('SVD did not converge in Linear Least Squares')
Chris@87 1870 resids = array([], result_real_t)
Chris@87 1871 if is_1d:
Chris@87 1872 x = array(ravel(bstar)[:n], dtype=result_t, copy=True)
Chris@87 1873 if results['rank'] == n and m > n:
Chris@87 1874 if isComplexType(t):
Chris@87 1875 resids = array([sum(abs(ravel(bstar)[n:])**2)],
Chris@87 1876 dtype=result_real_t)
Chris@87 1877 else:
Chris@87 1878 resids = array([sum((ravel(bstar)[n:])**2)],
Chris@87 1879 dtype=result_real_t)
Chris@87 1880 else:
Chris@87 1881 x = array(transpose(bstar)[:n,:], dtype=result_t, copy=True)
Chris@87 1882 if results['rank'] == n and m > n:
Chris@87 1883 if isComplexType(t):
Chris@87 1884 resids = sum(abs(transpose(bstar)[n:,:])**2, axis=0).astype(
Chris@87 1885 result_real_t)
Chris@87 1886 else:
Chris@87 1887 resids = sum((transpose(bstar)[n:,:])**2, axis=0).astype(
Chris@87 1888 result_real_t)
Chris@87 1889
Chris@87 1890 st = s[:min(n, m)].copy().astype(result_real_t)
Chris@87 1891 return wrap(x), wrap(resids), results['rank'], st
Chris@87 1892
Chris@87 1893
Chris@87 1894 def _multi_svd_norm(x, row_axis, col_axis, op):
Chris@87 1895 """Compute the extreme singular values of the 2-D matrices in `x`.
Chris@87 1896
Chris@87 1897 This is a private utility function used by numpy.linalg.norm().
Chris@87 1898
Chris@87 1899 Parameters
Chris@87 1900 ----------
Chris@87 1901 x : ndarray
Chris@87 1902 row_axis, col_axis : int
Chris@87 1903 The axes of `x` that hold the 2-D matrices.
Chris@87 1904 op : callable
Chris@87 1905 This should be either numpy.amin or numpy.amax.
Chris@87 1906
Chris@87 1907 Returns
Chris@87 1908 -------
Chris@87 1909 result : float or ndarray
Chris@87 1910 If `x` is 2-D, the return values is a float.
Chris@87 1911 Otherwise, it is an array with ``x.ndim - 2`` dimensions.
Chris@87 1912 The return values are either the minimum or maximum of the
Chris@87 1913 singular values of the matrices, depending on whether `op`
Chris@87 1914 is `numpy.amin` or `numpy.amax`.
Chris@87 1915
Chris@87 1916 """
Chris@87 1917 if row_axis > col_axis:
Chris@87 1918 row_axis -= 1
Chris@87 1919 y = rollaxis(rollaxis(x, col_axis, x.ndim), row_axis, -1)
Chris@87 1920 result = op(svd(y, compute_uv=0), axis=-1)
Chris@87 1921 return result
Chris@87 1922
Chris@87 1923
Chris@87 1924 def norm(x, ord=None, axis=None):
Chris@87 1925 """
Chris@87 1926 Matrix or vector norm.
Chris@87 1927
Chris@87 1928 This function is able to return one of seven different matrix norms,
Chris@87 1929 or one of an infinite number of vector norms (described below), depending
Chris@87 1930 on the value of the ``ord`` parameter.
Chris@87 1931
Chris@87 1932 Parameters
Chris@87 1933 ----------
Chris@87 1934 x : array_like
Chris@87 1935 Input array. If `axis` is None, `x` must be 1-D or 2-D.
Chris@87 1936 ord : {non-zero int, inf, -inf, 'fro'}, optional
Chris@87 1937 Order of the norm (see table under ``Notes``). inf means numpy's
Chris@87 1938 `inf` object.
Chris@87 1939 axis : {int, 2-tuple of ints, None}, optional
Chris@87 1940 If `axis` is an integer, it specifies the axis of `x` along which to
Chris@87 1941 compute the vector norms. If `axis` is a 2-tuple, it specifies the
Chris@87 1942 axes that hold 2-D matrices, and the matrix norms of these matrices
Chris@87 1943 are computed. If `axis` is None then either a vector norm (when `x`
Chris@87 1944 is 1-D) or a matrix norm (when `x` is 2-D) is returned.
Chris@87 1945
Chris@87 1946 Returns
Chris@87 1947 -------
Chris@87 1948 n : float or ndarray
Chris@87 1949 Norm of the matrix or vector(s).
Chris@87 1950
Chris@87 1951 Notes
Chris@87 1952 -----
Chris@87 1953 For values of ``ord <= 0``, the result is, strictly speaking, not a
Chris@87 1954 mathematical 'norm', but it may still be useful for various numerical
Chris@87 1955 purposes.
Chris@87 1956
Chris@87 1957 The following norms can be calculated:
Chris@87 1958
Chris@87 1959 ===== ============================ ==========================
Chris@87 1960 ord norm for matrices norm for vectors
Chris@87 1961 ===== ============================ ==========================
Chris@87 1962 None Frobenius norm 2-norm
Chris@87 1963 'fro' Frobenius norm --
Chris@87 1964 inf max(sum(abs(x), axis=1)) max(abs(x))
Chris@87 1965 -inf min(sum(abs(x), axis=1)) min(abs(x))
Chris@87 1966 0 -- sum(x != 0)
Chris@87 1967 1 max(sum(abs(x), axis=0)) as below
Chris@87 1968 -1 min(sum(abs(x), axis=0)) as below
Chris@87 1969 2 2-norm (largest sing. value) as below
Chris@87 1970 -2 smallest singular value as below
Chris@87 1971 other -- sum(abs(x)**ord)**(1./ord)
Chris@87 1972 ===== ============================ ==========================
Chris@87 1973
Chris@87 1974 The Frobenius norm is given by [1]_:
Chris@87 1975
Chris@87 1976 :math:`||A||_F = [\\sum_{i,j} abs(a_{i,j})^2]^{1/2}`
Chris@87 1977
Chris@87 1978 References
Chris@87 1979 ----------
Chris@87 1980 .. [1] G. H. Golub and C. F. Van Loan, *Matrix Computations*,
Chris@87 1981 Baltimore, MD, Johns Hopkins University Press, 1985, pg. 15
Chris@87 1982
Chris@87 1983 Examples
Chris@87 1984 --------
Chris@87 1985 >>> from numpy import linalg as LA
Chris@87 1986 >>> a = np.arange(9) - 4
Chris@87 1987 >>> a
Chris@87 1988 array([-4, -3, -2, -1, 0, 1, 2, 3, 4])
Chris@87 1989 >>> b = a.reshape((3, 3))
Chris@87 1990 >>> b
Chris@87 1991 array([[-4, -3, -2],
Chris@87 1992 [-1, 0, 1],
Chris@87 1993 [ 2, 3, 4]])
Chris@87 1994
Chris@87 1995 >>> LA.norm(a)
Chris@87 1996 7.745966692414834
Chris@87 1997 >>> LA.norm(b)
Chris@87 1998 7.745966692414834
Chris@87 1999 >>> LA.norm(b, 'fro')
Chris@87 2000 7.745966692414834
Chris@87 2001 >>> LA.norm(a, np.inf)
Chris@87 2002 4
Chris@87 2003 >>> LA.norm(b, np.inf)
Chris@87 2004 9
Chris@87 2005 >>> LA.norm(a, -np.inf)
Chris@87 2006 0
Chris@87 2007 >>> LA.norm(b, -np.inf)
Chris@87 2008 2
Chris@87 2009
Chris@87 2010 >>> LA.norm(a, 1)
Chris@87 2011 20
Chris@87 2012 >>> LA.norm(b, 1)
Chris@87 2013 7
Chris@87 2014 >>> LA.norm(a, -1)
Chris@87 2015 -4.6566128774142013e-010
Chris@87 2016 >>> LA.norm(b, -1)
Chris@87 2017 6
Chris@87 2018 >>> LA.norm(a, 2)
Chris@87 2019 7.745966692414834
Chris@87 2020 >>> LA.norm(b, 2)
Chris@87 2021 7.3484692283495345
Chris@87 2022
Chris@87 2023 >>> LA.norm(a, -2)
Chris@87 2024 nan
Chris@87 2025 >>> LA.norm(b, -2)
Chris@87 2026 1.8570331885190563e-016
Chris@87 2027 >>> LA.norm(a, 3)
Chris@87 2028 5.8480354764257312
Chris@87 2029 >>> LA.norm(a, -3)
Chris@87 2030 nan
Chris@87 2031
Chris@87 2032 Using the `axis` argument to compute vector norms:
Chris@87 2033
Chris@87 2034 >>> c = np.array([[ 1, 2, 3],
Chris@87 2035 ... [-1, 1, 4]])
Chris@87 2036 >>> LA.norm(c, axis=0)
Chris@87 2037 array([ 1.41421356, 2.23606798, 5. ])
Chris@87 2038 >>> LA.norm(c, axis=1)
Chris@87 2039 array([ 3.74165739, 4.24264069])
Chris@87 2040 >>> LA.norm(c, ord=1, axis=1)
Chris@87 2041 array([6, 6])
Chris@87 2042
Chris@87 2043 Using the `axis` argument to compute matrix norms:
Chris@87 2044
Chris@87 2045 >>> m = np.arange(8).reshape(2,2,2)
Chris@87 2046 >>> LA.norm(m, axis=(1,2))
Chris@87 2047 array([ 3.74165739, 11.22497216])
Chris@87 2048 >>> LA.norm(m[0, :, :]), LA.norm(m[1, :, :])
Chris@87 2049 (3.7416573867739413, 11.224972160321824)
Chris@87 2050
Chris@87 2051 """
Chris@87 2052 x = asarray(x)
Chris@87 2053
Chris@87 2054 # Check the default case first and handle it immediately.
Chris@87 2055 if ord is None and axis is None:
Chris@87 2056 x = x.ravel(order='K')
Chris@87 2057 if isComplexType(x.dtype.type):
Chris@87 2058 sqnorm = dot(x.real, x.real) + dot(x.imag, x.imag)
Chris@87 2059 else:
Chris@87 2060 sqnorm = dot(x, x)
Chris@87 2061 return sqrt(sqnorm)
Chris@87 2062
Chris@87 2063 # Normalize the `axis` argument to a tuple.
Chris@87 2064 nd = x.ndim
Chris@87 2065 if axis is None:
Chris@87 2066 axis = tuple(range(nd))
Chris@87 2067 elif not isinstance(axis, tuple):
Chris@87 2068 axis = (axis,)
Chris@87 2069
Chris@87 2070 if len(axis) == 1:
Chris@87 2071 if ord == Inf:
Chris@87 2072 return abs(x).max(axis=axis)
Chris@87 2073 elif ord == -Inf:
Chris@87 2074 return abs(x).min(axis=axis)
Chris@87 2075 elif ord == 0:
Chris@87 2076 # Zero norm
Chris@87 2077 return (x != 0).sum(axis=axis)
Chris@87 2078 elif ord == 1:
Chris@87 2079 # special case for speedup
Chris@87 2080 return add.reduce(abs(x), axis=axis)
Chris@87 2081 elif ord is None or ord == 2:
Chris@87 2082 # special case for speedup
Chris@87 2083 s = (x.conj() * x).real
Chris@87 2084 return sqrt(add.reduce(s, axis=axis))
Chris@87 2085 else:
Chris@87 2086 try:
Chris@87 2087 ord + 1
Chris@87 2088 except TypeError:
Chris@87 2089 raise ValueError("Invalid norm order for vectors.")
Chris@87 2090 if x.dtype.type is longdouble:
Chris@87 2091 # Convert to a float type, so integer arrays give
Chris@87 2092 # float results. Don't apply asfarray to longdouble arrays,
Chris@87 2093 # because it will downcast to float64.
Chris@87 2094 absx = abs(x)
Chris@87 2095 else:
Chris@87 2096 absx = x if isComplexType(x.dtype.type) else asfarray(x)
Chris@87 2097 if absx.dtype is x.dtype:
Chris@87 2098 absx = abs(absx)
Chris@87 2099 else:
Chris@87 2100 # if the type changed, we can safely overwrite absx
Chris@87 2101 abs(absx, out=absx)
Chris@87 2102 absx **= ord
Chris@87 2103 return add.reduce(absx, axis=axis) ** (1.0 / ord)
Chris@87 2104 elif len(axis) == 2:
Chris@87 2105 row_axis, col_axis = axis
Chris@87 2106 if not (-nd <= row_axis < nd and -nd <= col_axis < nd):
Chris@87 2107 raise ValueError('Invalid axis %r for an array with shape %r' %
Chris@87 2108 (axis, x.shape))
Chris@87 2109 if row_axis % nd == col_axis % nd:
Chris@87 2110 raise ValueError('Duplicate axes given.')
Chris@87 2111 if ord == 2:
Chris@87 2112 return _multi_svd_norm(x, row_axis, col_axis, amax)
Chris@87 2113 elif ord == -2:
Chris@87 2114 return _multi_svd_norm(x, row_axis, col_axis, amin)
Chris@87 2115 elif ord == 1:
Chris@87 2116 if col_axis > row_axis:
Chris@87 2117 col_axis -= 1
Chris@87 2118 return add.reduce(abs(x), axis=row_axis).max(axis=col_axis)
Chris@87 2119 elif ord == Inf:
Chris@87 2120 if row_axis > col_axis:
Chris@87 2121 row_axis -= 1
Chris@87 2122 return add.reduce(abs(x), axis=col_axis).max(axis=row_axis)
Chris@87 2123 elif ord == -1:
Chris@87 2124 if col_axis > row_axis:
Chris@87 2125 col_axis -= 1
Chris@87 2126 return add.reduce(abs(x), axis=row_axis).min(axis=col_axis)
Chris@87 2127 elif ord == -Inf:
Chris@87 2128 if row_axis > col_axis:
Chris@87 2129 row_axis -= 1
Chris@87 2130 return add.reduce(abs(x), axis=col_axis).min(axis=row_axis)
Chris@87 2131 elif ord in [None, 'fro', 'f']:
Chris@87 2132 return sqrt(add.reduce((x.conj() * x).real, axis=axis))
Chris@87 2133 else:
Chris@87 2134 raise ValueError("Invalid norm order for matrices.")
Chris@87 2135 else:
Chris@87 2136 raise ValueError("Improper number of dimensions to norm.")