# This is a cut-down and hacked version of the mpmath matrix and linear algebra code # See http://code.google.com/p/mpmath for the full source # The license is as follows: # Copyright (c) 2005-2010 Fredrik Johansson and mpmath contributors # # All rights reserved. # # Redistribution and use in source and binary forms, with or without # modification, are permitted provided that the following conditions are met: # # a. Redistributions of source code must retain the above copyright notice, # this list of conditions and the following disclaimer. # b. Redistributions in binary form must reproduce the above copyright # notice, this list of conditions and the following disclaimer in the # documentation and/or other materials provided with the distribution. # c. Neither the name of mpmath nor the names of its contributors # may be used to endorse or promote products derived from this software # without specific prior written permission. # # # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" # AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE # IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE # ARE DISCLAIMED. IN NO EVENT SHALL THE REGENTS OR CONTRIBUTORS BE LIABLE FOR # ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL # DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR # SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER # CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT # LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY # OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH # DAMAGE. import math from copy import copy rowsep = '\n' colsep = ' ' zero = 0.0 one = 1.0 eps = 2.2204460492503131e-16 inf = 1e300*1e300 ninf = -inf nan = inf-inf absmin = absmax = abs class matrix(object): """ Numerical matrix. Specify the dimensions or the data as a nested list. Elements default to zero. Use a flat list to create a column vector easily. By default, only mpf is used to store the data. You can specify another type using force_type=type. It's possible to specify None. Make sure force_type(force_type()) is fast. Creating matrices ----------------- Matrices in mpmath are implemented using dictionaries. Only non-zero values are stored, so it is cheap to represent sparse matrices. The most basic way to create one is to use the ``matrix`` class directly. You can create an empty matrix specifying the dimensions: >>> from mpmath import * >>> mp.dps = 15 >>> matrix(2) matrix( [['0.0', '0.0'], ['0.0', '0.0']]) >>> matrix(2, 3) matrix( [['0.0', '0.0', '0.0'], ['0.0', '0.0', '0.0']]) Calling ``matrix`` with one dimension will create a square matrix. To access the dimensions of a matrix, use the ``rows`` or ``cols`` keyword: >>> A = matrix(3, 2) >>> A matrix( [['0.0', '0.0'], ['0.0', '0.0'], ['0.0', '0.0']]) >>> A.rows 3 >>> A.cols 2 You can also change the dimension of an existing matrix. This will set the new elements to 0. If the new dimension is smaller than before, the concerning elements are discarded: >>> A.rows = 2 >>> A matrix( [['0.0', '0.0'], ['0.0', '0.0']]) Internally ``mpmathify`` is used every time an element is set. This is done using the syntax A[row,column], counting from 0: >>> A = matrix(2) >>> A[1,1] = 1 + 1j >>> A matrix( [['0.0', '0.0'], ['0.0', '(1.0 + 1.0j)']]) You can use the keyword ``force_type`` to change the function which is called on every new element: >>> matrix(2, 5, force_type=int) matrix( [[0, 0, 0, 0, 0], [0, 0, 0, 0, 0]]) A more comfortable way to create a matrix lets you use nested lists: >>> matrix([[1, 2], [3, 4]]) matrix( [['1.0', '2.0'], ['3.0', '4.0']]) If you want to preserve the type of the elements you can use ``force_type=None``: >>> matrix([[1, 2.5], [1j, mpf(2)]], force_type=None) matrix( [[1, 2.5], [1j, '2.0']]) Convenient advanced functions are available for creating various standard matrices, see ``zeros``, ``ones``, ``diag``, ``eye``, ``randmatrix`` and ``hilbert``. Vectors ....... Vectors may also be represented by the ``matrix`` class (with rows = 1 or cols = 1). For vectors there are some things which make life easier. A column vector can be created using a flat list, a row vectors using an almost flat nested list:: >>> matrix([1, 2, 3]) matrix( [['1.0'], ['2.0'], ['3.0']]) >>> matrix([[1, 2, 3]]) matrix( [['1.0', '2.0', '3.0']]) Optionally vectors can be accessed like lists, using only a single index:: >>> x = matrix([1, 2, 3]) >>> x[1] mpf('2.0') >>> x[1,0] mpf('2.0') Other ..... Like you probably expected, matrices can be printed:: >>> print randmatrix(3) # doctest:+SKIP [ 0.782963853573023 0.802057689719883 0.427895717335467] [0.0541876859348597 0.708243266653103 0.615134039977379] [ 0.856151514955773 0.544759264818486 0.686210904770947] Use ``nstr`` or ``nprint`` to specify the number of digits to print:: >>> nprint(randmatrix(5), 3) # doctest:+SKIP [2.07e-1 1.66e-1 5.06e-1 1.89e-1 8.29e-1] [6.62e-1 6.55e-1 4.47e-1 4.82e-1 2.06e-2] [4.33e-1 7.75e-1 6.93e-2 2.86e-1 5.71e-1] [1.01e-1 2.53e-1 6.13e-1 3.32e-1 2.59e-1] [1.56e-1 7.27e-2 6.05e-1 6.67e-2 2.79e-1] As matrices are mutable, you will need to copy them sometimes:: >>> A = matrix(2) >>> A matrix( [['0.0', '0.0'], ['0.0', '0.0']]) >>> B = A.copy() >>> B[0,0] = 1 >>> B matrix( [['1.0', '0.0'], ['0.0', '0.0']]) >>> A matrix( [['0.0', '0.0'], ['0.0', '0.0']]) Finally, it is possible to convert a matrix to a nested list. This is very useful, as most Python libraries involving matrices or arrays (namely NumPy or SymPy) support this format:: >>> B.tolist() [[mpf('1.0'), mpf('0.0')], [mpf('0.0'), mpf('0.0')]] Matrix operations ----------------- You can add and subtract matrices of compatible dimensions:: >>> A = matrix([[1, 2], [3, 4]]) >>> B = matrix([[-2, 4], [5, 9]]) >>> A + B matrix( [['-1.0', '6.0'], ['8.0', '13.0']]) >>> A - B matrix( [['3.0', '-2.0'], ['-2.0', '-5.0']]) >>> A + ones(3) # doctest:+ELLIPSIS Traceback (most recent call last): ... ValueError: incompatible dimensions for addition It is possible to multiply or add matrices and scalars. In the latter case the operation will be done element-wise:: >>> A * 2 matrix( [['2.0', '4.0'], ['6.0', '8.0']]) >>> A / 4 matrix( [['0.25', '0.5'], ['0.75', '1.0']]) >>> A - 1 matrix( [['0.0', '1.0'], ['2.0', '3.0']]) Of course you can perform matrix multiplication, if the dimensions are compatible:: >>> A * B matrix( [['8.0', '22.0'], ['14.0', '48.0']]) >>> matrix([[1, 2, 3]]) * matrix([[-6], [7], [-2]]) matrix( [['2.0']]) You can raise powers of square matrices:: >>> A**2 matrix( [['7.0', '10.0'], ['15.0', '22.0']]) Negative powers will calculate the inverse:: >>> A**-1 matrix( [['-2.0', '1.0'], ['1.5', '-0.5']]) >>> A * A**-1 matrix( [['1.0', '1.0842021724855e-19'], ['-2.16840434497101e-19', '1.0']]) Matrix transposition is straightforward:: >>> A = ones(2, 3) >>> A matrix( [['1.0', '1.0', '1.0'], ['1.0', '1.0', '1.0']]) >>> A.T matrix( [['1.0', '1.0'], ['1.0', '1.0'], ['1.0', '1.0']]) Norms ..... Sometimes you need to know how "large" a matrix or vector is. Due to their multidimensional nature it's not possible to compare them, but there are several functions to map a matrix or a vector to a positive real number, the so called norms. For vectors the p-norm is intended, usually the 1-, the 2- and the oo-norm are used. >>> x = matrix([-10, 2, 100]) >>> norm(x, 1) mpf('112.0') >>> norm(x, 2) mpf('100.5186549850325') >>> norm(x, inf) mpf('100.0') Please note that the 2-norm is the most used one, though it is more expensive to calculate than the 1- or oo-norm. It is possible to generalize some vector norms to matrix norm:: >>> A = matrix([[1, -1000], [100, 50]]) >>> mnorm(A, 1) mpf('1050.0') >>> mnorm(A, inf) mpf('1001.0') >>> mnorm(A, 'F') mpf('1006.2310867787777') The last norm (the "Frobenius-norm") is an approximation for the 2-norm, which is hard to calculate and not available. The Frobenius-norm lacks some mathematical properties you might expect from a norm. """ def __init__(self, *args, **kwargs): self.__data = {} # LU decompostion cache, this is useful when solving the same system # multiple times, when calculating the inverse and when calculating the # determinant self._LU = None if isinstance(args[0], (list, tuple)): if isinstance(args[0][0], (list, tuple)): # interpret nested list as matrix A = args[0] self.__rows = len(A) self.__cols = len(A[0]) for i, row in enumerate(A): for j, a in enumerate(row): self[i, j] = convert(a) else: # interpret list as row vector v = args[0] self.__rows = len(v) self.__cols = 1 for i, e in enumerate(v): self[i, 0] = e elif isinstance(args[0], int): # create empty matrix of given dimensions if len(args) == 1: self.__rows = self.__cols = args[0] else: assert isinstance(args[1], int), 'expected int' self.__rows = args[0] self.__cols = args[1] elif isinstance(args[0], matrix): A = args[0].copy() self.__data = A._matrix__data self.__rows = A._matrix__rows self.__cols = A._matrix__cols for i in xrange(A.__rows): for j in xrange(A.__cols): A[i,j] = convert(A[i,j]) elif hasattr(args[0], 'tolist'): A = matrix(args[0].tolist()) self.__data = A._matrix__data self.__rows = A._matrix__rows self.__cols = A._matrix__cols else: raise TypeError('could not interpret given arguments') def apply(self, f): """ Return a copy of self with the function `f` applied elementwise. """ new = matrix(self.__rows, self.__cols) for i in xrange(self.__rows): for j in xrange(self.__cols): new[i,j] = f(self[i,j]) return new def __nstr__(self, n=None, **kwargs): # Build table of string representations of the elements res = [] # Track per-column max lengths for pretty alignment maxlen = [0] * self.cols for i in range(self.rows): res.append([]) for j in range(self.cols): string = str(self[i,j]) res[-1].append(string) maxlen[j] = max(len(string), maxlen[j]) # Patch strings together for i, row in enumerate(res): for j, elem in enumerate(row): # Pad each element up to maxlen so the columns line up row[j] = elem.rjust(maxlen[j]) res[i] = "[" + colsep.join(row) + "]" return rowsep.join(res) def __str__(self): return self.__nstr__() def _toliststr(self, avoid_type=False): """ Create a list string from a matrix. If avoid_type: avoid multiple 'mpf's. """ s = '[' for i in xrange(self.__rows): s += '[' for j in xrange(self.__cols): a = "'" + str(self[i,j]) + "'" s += a + ', ' s = s[:-2] s += '],\n ' s = s[:-3] s += ']' return s def tolist(self): """ Convert the matrix to a nested list. """ return [[self[i,j] for j in range(self.__cols)] for i in range(self.__rows)] def __repr__(self): s = 'matrix(\n' s += self._toliststr(avoid_type=True) + ')' return s def __getitem__(self, key): if type(key) is int: # only sufficent for vectors if self.__rows == 1: key = (0, key) elif self.__cols == 1: key = (key, 0) else: raise IndexError('insufficient indices for matrix') if key in self.__data: return self.__data[key] else: if key[0] >= self.__rows or key[1] >= self.__cols: raise IndexError('matrix index out of range') return 0 def __setitem__(self, key, value): if type(key) is int: # only sufficent for vectors if self.__rows == 1: key = (0, key) elif self.__cols == 1: key = (key, 0) else: raise IndexError('insufficient indices for matrix') if key[0] >= self.__rows or key[1] >= self.__cols: raise IndexError('matrix index out of range') value = convert(value) if value: # only store non-zeros self.__data[key] = value elif key in self.__data: del self.__data[key] # TODO: maybe do this better, if the performance impact is significant if self._LU: self._LU = None def __iter__(self): for i in xrange(self.__rows): for j in xrange(self.__cols): yield self[i,j] def __mul__(self, other): if isinstance(other, matrix): # dot multiplication TODO: use Strassen's method? if self.__cols != other.__rows: raise ValueError('dimensions not compatible for multiplication') new = matrix(self.__rows, other.__cols) for i in xrange(self.__rows): for j in xrange(other.__cols): new[i, j] = fdot((self[i,k], other[k,j]) for k in xrange(other.__rows)) return new else: # try scalar multiplication new = matrix(self.__rows, self.__cols) for i in xrange(self.__rows): for j in xrange(self.__cols): new[i, j] = other * self[i, j] return new def __rmul__(self, other): # assume other is scalar and thus commutative assert not isinstance(other, matrix) return self.__mul__(other) def __pow__(self, other): # avoid cyclic import problems #from linalg import inverse if not isinstance(other, int): raise ValueError('only integer exponents are supported') if not self.__rows == self.__cols: raise ValueError('only powers of square matrices are defined') n = other if n == 0: return eye(self.__rows) if n < 0: n = -n neg = True else: neg = False i = n y = 1 z = self.copy() while i != 0: if i % 2 == 1: y = y * z z = z*z i = i // 2 if neg: y = inverse(y) return y def __div__(self, other): # assume other is scalar and do element-wise divison assert not isinstance(other, matrix) new = matrix(self.__rows, self.__cols) for i in xrange(self.__rows): for j in xrange(self.__cols): new[i,j] = self[i,j] / other return new __truediv__ = __div__ def __add__(self, other): if isinstance(other, matrix): if not (self.__rows == other.__rows and self.__cols == other.__cols): raise ValueError('incompatible dimensions for addition') new = matrix(self.__rows, self.__cols) for i in xrange(self.__rows): for j in xrange(self.__cols): new[i,j] = self[i,j] + other[i,j] return new else: # assume other is scalar and add element-wise new = matrix(self.__rows, self.__cols) for i in xrange(self.__rows): for j in xrange(self.__cols): new[i,j] += self[i,j] + other return new def __radd__(self, other): return self.__add__(other) def __sub__(self, other): if isinstance(other, matrix) and not (self.__rows == other.__rows and self.__cols == other.__cols): raise ValueError('incompatible dimensions for substraction') return self.__add__(other * (-1)) def __neg__(self): return (-1) * self def __rsub__(self, other): return -self + other def __eq__(self, other): return self.__rows == other.__rows and self.__cols == other.__cols \ and self.__data == other.__data def __len__(self): if self.rows == 1: return self.cols elif self.cols == 1: return self.rows else: return self.rows # do it like numpy def __getrows(self): return self.__rows def __setrows(self, value): for key in self.__data.copy().iterkeys(): if key[0] >= value: del self.__data[key] self.__rows = value rows = property(__getrows, __setrows, doc='number of rows') def __getcols(self): return self.__cols def __setcols(self, value): for key in self.__data.copy().iterkeys(): if key[1] >= value: del self.__data[key] self.__cols = value cols = property(__getcols, __setcols, doc='number of columns') def transpose(self): new = matrix(self.__cols, self.__rows) for i in xrange(self.__rows): for j in xrange(self.__cols): new[j,i] = self[i,j] return new T = property(transpose) def copy(self): new = matrix(self.__rows, self.__cols) new.__data = self.__data.copy() return new __copy__ = copy def column(self, n): m = matrix(self.rows, 1) for i in range(self.rows): m[i] = self[i,n] return m def convert(x): try: return float(x) except: return complex(x) def fdot(xs, ys=None): if ys is not None: xs = zip(xs, ys) return sum((x*y for (x,y) in xs), 0) def fsum(args, absolute=False, squared=False): if absolute: if squared: return sum((abs(x)**2 for x in args), zero) return sum((abs(x) for x in args), zero) if squared: return sum((x**2 for x in args), zero) return sum(args, zero) def fprod(args): prod = one for arg in args: prod *= arg return prod def eye(n, **kwargs): """ Create square identity matrix n x n. """ A = matrix(n, **kwargs) for i in xrange(n): A[i,i] = 1 return A def diag(diagonal, **kwargs): """ Create square diagonal matrix using given list. Example: >>> from mpmath import diag, mp >>> mp.pretty = False >>> diag([1, 2, 3]) matrix( [['1.0', '0.0', '0.0'], ['0.0', '2.0', '0.0'], ['0.0', '0.0', '3.0']]) """ A = matrix(len(diagonal), **kwargs) for i in xrange(len(diagonal)): A[i,i] = diagonal[i] return A def zeros(*args, **kwargs): """ Create matrix m x n filled with zeros. One given dimension will create square matrix n x n. Example: >>> from mpmath import zeros, mp >>> mp.pretty = False >>> zeros(2) matrix( [['0.0', '0.0'], ['0.0', '0.0']]) """ if len(args) == 1: m = n = args[0] elif len(args) == 2: m = args[0] n = args[1] else: raise TypeError('zeros expected at most 2 arguments, got %i' % len(args)) A = matrix(m, n, **kwargs) for i in xrange(m): for j in xrange(n): A[i,j] = 0 return A def ones(*args, **kwargs): """ Create matrix m x n filled with ones. One given dimension will create square matrix n x n. Example: >>> from mpmath import ones, mp >>> mp.pretty = False >>> ones(2) matrix( [['1.0', '1.0'], ['1.0', '1.0']]) """ if len(args) == 1: m = n = args[0] elif len(args) == 2: m = args[0] n = args[1] else: raise TypeError('ones expected at most 2 arguments, got %i' % len(args)) A = matrix(m, n, **kwargs) for i in xrange(m): for j in xrange(n): A[i,j] = 1 return A def hilbert(m, n=None): """ Create (pseudo) hilbert matrix m x n. One given dimension will create hilbert matrix n x n. The matrix is very ill-conditioned and symmetric, positive definite if square. """ if n is None: n = m A = matrix(m, n) for i in xrange(m): for j in xrange(n): A[i,j] = 1 / (i + j + 1) return A def swap_row(A, i, j): """ Swap row i with row j. """ if i == j: return if isinstance(A, matrix): for k in xrange(A.cols): A[i,k], A[j,k] = A[j,k], A[i,k] elif isinstance(A, list): A[i], A[j] = A[j], A[i] else: raise TypeError('could not interpret type') def extend(A, b): """ Extend matrix A with column b and return result. """ assert isinstance(A, matrix) assert A.rows == len(b) A = A.copy() A.cols += 1 for i in xrange(A.rows): A[i, A.cols-1] = b[i] return A def norm(x, p=2): r""" Gives the entrywise `p`-norm of an iterable *x*, i.e. the vector norm `\left(\sum_k |x_k|^p\right)^{1/p}`, for any given `1 \le p \le \infty`. Special cases: If *x* is not iterable, this just returns ``absmax(x)``. ``p=1`` gives the sum of absolute values. ``p=2`` is the standard Euclidean vector norm. ``p=inf`` gives the magnitude of the largest element. For *x* a matrix, ``p=2`` is the Frobenius norm. For operator matrix norms, use :func:`~mpmath.mnorm` instead. You can use the string 'inf' as well as float('inf') or mpf('inf') to specify the infinity norm. **Examples** >>> from mpmath import * >>> mp.dps = 15; mp.pretty = False >>> x = matrix([-10, 2, 100]) >>> norm(x, 1) mpf('112.0') >>> norm(x, 2) mpf('100.5186549850325') >>> norm(x, inf) mpf('100.0') """ try: iter(x) except TypeError: return absmax(x) if type(p) is not int: p = convert(p) if p == inf: return max(absmax(i) for i in x) elif p == 1: return fsum(x, absolute=1) elif p == 2: return math.sqrt(fsum(x, absolute=1, squared=1)) elif p > 1: return nthroot(fsum(abs(i)**p for i in x), p) else: raise ValueError('p has to be >= 1') def nthroot(x, n): r = 1./n x ** r def mnorm(A, p=1): r""" Gives the matrix (operator) `p`-norm of A. Currently ``p=1`` and ``p=inf`` are supported: ``p=1`` gives the 1-norm (maximal column sum) ``p=inf`` gives the `\infty`-norm (maximal row sum). You can use the string 'inf' as well as float('inf') or mpf('inf') ``p=2`` (not implemented) for a square matrix is the usual spectral matrix norm, i.e. the largest singular value. ``p='f'`` (or 'F', 'fro', 'Frobenius, 'frobenius') gives the Frobenius norm, which is the elementwise 2-norm. The Frobenius norm is an approximation of the spectral norm and satisfies .. math :: \frac{1}{\sqrt{\mathrm{rank}(A)}} \|A\|_F \le \|A\|_2 \le \|A\|_F The Frobenius norm lacks some mathematical properties that might be expected of a norm. For general elementwise `p`-norms, use :func:`~mpmath.norm` instead. **Examples** >>> from mpmath import * >>> mp.dps = 15; mp.pretty = False >>> A = matrix([[1, -1000], [100, 50]]) >>> mnorm(A, 1) mpf('1050.0') >>> mnorm(A, inf) mpf('1001.0') >>> mnorm(A, 'F') mpf('1006.2310867787777') """ A = matrix(A) if type(p) is not int: if type(p) is str and 'frobenius'.startswith(p.lower()): return norm(A, 2) p = convert(p) m, n = A.rows, A.cols if p == 1: return max(fsum((A[i,j] for i in xrange(m)), absolute=1) for j in xrange(n)) elif p == inf: return max(fsum((A[i,j] for j in xrange(n)), absolute=1) for i in xrange(m)) else: raise NotImplementedError("matrix p-norm for arbitrary p") def LU_decomp(A, overwrite=False, use_cache=True): """ LU-factorization of a n*n matrix using the Gauss algorithm. Returns L and U in one matrix and the pivot indices. Use overwrite to specify whether A will be overwritten with L and U. """ if not A.rows == A.cols: raise ValueError('need n*n matrix') # get from cache if possible if use_cache and isinstance(A, matrix) and A._LU: return A._LU if not overwrite: orig = A A = A.copy() tol = absmin(mnorm(A,1) * eps) # each pivot element has to be bigger n = A.rows p = [None]*(n - 1) for j in xrange(n - 1): # pivoting, choose max(abs(reciprocal row sum)*abs(pivot element)) biggest = 0 for k in xrange(j, n): s = fsum([absmin(A[k,l]) for l in xrange(j, n)]) if absmin(s) <= tol: raise ZeroDivisionError('matrix is numerically singular') current = 1/s * absmin(A[k,j]) if current > biggest: # TODO: what if equal? biggest = current p[j] = k # swap rows according to p swap_row(A, j, p[j]) if absmin(A[j,j]) <= tol: raise ZeroDivisionError('matrix is numerically singular') # calculate elimination factors and add rows for i in xrange(j + 1, n): A[i,j] /= A[j,j] for k in xrange(j + 1, n): A[i,k] -= A[i,j]*A[j,k] if absmin(A[n - 1,n - 1]) <= tol: raise ZeroDivisionError('matrix is numerically singular') # cache decomposition if not overwrite and isinstance(orig, matrix): orig._LU = (A, p) return A, p def L_solve(L, b, p=None): """ Solve the lower part of a LU factorized matrix for y. """ assert L.rows == L.cols, 'need n*n matrix' n = L.rows assert len(b) == n b = copy(b) if p: # swap b according to p for k in xrange(0, len(p)): swap_row(b, k, p[k]) # solve for i in xrange(1, n): for j in xrange(i): b[i] -= L[i,j] * b[j] return b def U_solve(U, y): """ Solve the upper part of a LU factorized matrix for x. """ assert U.rows == U.cols, 'need n*n matrix' n = U.rows assert len(y) == n x = copy(y) for i in xrange(n - 1, -1, -1): for j in xrange(i + 1, n): x[i] -= U[i,j] * x[j] x[i] /= U[i,i] return x def lu_solve(A, b, **kwargs): """ Ax = b => x Solve a determined or overdetermined linear equations system. Fast LU decomposition is used, which is less accurate than QR decomposition (especially for overdetermined systems), but it's twice as efficient. Use qr_solve if you want more precision or have to solve a very ill- conditioned system. If you specify real=True, it does not check for overdeterminded complex systems. """ global prec oldPrec = prec try: prec += 10 # do not overwrite A nor b A, b = matrix(A, **kwargs).copy(), matrix(b, **kwargs).copy() if A.rows < A.cols: raise ValueError('cannot solve underdetermined system') if A.rows > A.cols: # use least-squares method if overdetermined # (this increases errors) AH = A.H A = AH * A b = AH * b if (false): # TODO: necessary to check also b? x = cholesky_solve(A, b) else: x = lu_solve(A, b) else: # LU factorization A, p = LU_decomp(A) b = L_solve(A, b, p) x = U_solve(A, b) finally: prec = oldPrec return x def improve_solution(A, x, b, maxsteps=1): """ Improve a solution to a linear equation system iteratively. This re-uses the LU decomposition and is thus cheap. Usually 3 up to 4 iterations are giving the maximal improvement. """ assert A.rows == A.cols, 'need n*n matrix' # TODO: really? for _ in xrange(maxsteps): r = residual(A, x, b) if norm(r, 2) < 10*eps: break # this uses cached LU decomposition and is thus cheap dx = lu_solve(A, -r) x += dx return x def lu(A): """ A -> P, L, U LU factorisation of a square matrix A. L is the lower, U the upper part. P is the permutation matrix indicating the row swaps. P*A = L*U If you need efficiency, use the low-level method LU_decomp instead, it's much more memory efficient. """ # get factorization A, p = LU_decomp(A) n = A.rows L = matrix(n) U = matrix(n) for i in xrange(n): for j in xrange(n): if i > j: L[i,j] = A[i,j] elif i == j: L[i,j] = 1 U[i,j] = A[i,j] else: U[i,j] = A[i,j] # calculate permutation matrix P = eye(n) for k in xrange(len(p)): swap_row(P, k, p[k]) return P, L, U def unitvector(n, i): """ Return the i-th n-dimensional unit vector. """ assert 0 < i <= n, 'this unit vector does not exist' return [zero]*(i-1) + [one] + [zero]*(n-i) def inverse(A, **kwargs): """ Calculate the inverse of a matrix. If you want to solve an equation system Ax = b, it's recommended to use solve(A, b) instead, it's about 3 times more efficient. """ # do not overwrite A A = matrix(A, **kwargs).copy() n = A.rows # get LU factorisation A, p = LU_decomp(A) cols = [] # calculate unit vectors and solve corresponding system to get columns for i in xrange(1, n + 1): e = unitvector(n, i) y = L_solve(A, e, p) cols.append(U_solve(A, y)) # convert columns to matrix inv = [] for i in xrange(n): row = [] for j in xrange(n): row.append(cols[j][i]) inv.append(row) return matrix(inv, **kwargs) def householder(A): """ (A|b) -> H, p, x, res (A|b) is the coefficient matrix with left hand side of an optionally overdetermined linear equation system. H and p contain all information about the transformation matrices. x is the solution, res the residual. """ assert isinstance(A, matrix) m = A.rows n = A.cols assert m >= n - 1 # calculate Householder matrix p = [] for j in xrange(0, n - 1): s = fsum((A[i,j])**2 for i in xrange(j, m)) if not abs(s) > eps: raise ValueError('matrix is numerically singular') p.append(-sign(A[j,j]) * math.sqrt(s)) kappa = one / (s - p[j] * A[j,j]) A[j,j] -= p[j] for k in xrange(j+1, n): y = fsum(A[i,j] * A[i,k] for i in xrange(j, m)) * kappa for i in xrange(j, m): A[i,k] -= A[i,j] * y # solve Rx = c1 x = [A[i,n - 1] for i in xrange(n - 1)] for i in xrange(n - 2, -1, -1): x[i] -= fsum(A[i,j] * x[j] for j in xrange(i + 1, n - 1)) x[i] /= p[i] # calculate residual if not m == n - 1: r = [A[m-1-i, n-1] for i in xrange(m - n + 1)] else: # determined system, residual should be 0 r = [0]*m # maybe a bad idea, changing r[i] will change all elements return A, p, x, r def residual(A, x, b, **kwargs): """ Calculate the residual of a solution to a linear equation system. r = A*x - b for A*x = b """ A, x, b = matrix(A, **kwargs), matrix(x, **kwargs), matrix(b, **kwargs) return A*x - b def qr_solve(A, b, norm=None, **kwargs): """ Ax = b => x, ||Ax - b|| Solve a determined or overdetermined linear equations system and calculate the norm of the residual (error). QR decomposition using Householder factorization is applied, which gives very accurate results even for ill-conditioned matrices. qr_solve is twice as efficient. """ # do not overwrite A nor b A, b = matrix(A, **kwargs).copy(), matrix(b, **kwargs).copy() if A.rows < A.cols: raise ValueError('cannot solve underdetermined system') H, p, x, r = householder(extend(A, b)) res = norm(r) # calculate residual "manually" for determined systems if res == 0: res = norm(residual(A, x, b)) return matrix(x, **kwargs), res def cholesky(A, tol=None): r""" Cholesky decomposition of a symmetric positive-definite matrix `A`. Returns a lower triangular matrix `L` such that `A = L \times L^T`. More generally, for a complex Hermitian positive-definite matrix, a Cholesky decomposition satisfying `A = L \times L^H` is returned. The Cholesky decomposition can be used to solve linear equation systems twice as efficiently as LU decomposition, or to test whether `A` is positive-definite. The optional parameter ``tol`` determines the tolerance for verifying positive-definiteness. **Examples** Cholesky decomposition of a positive-definite symmetric matrix:: >>> from mpmath import * >>> mp.dps = 25; mp.pretty = True >>> A = eye(3) + hilbert(3) >>> nprint(A) [ 2.0 0.5 0.333333] [ 0.5 1.33333 0.25] [0.333333 0.25 1.2] >>> L = cholesky(A) >>> nprint(L) [ 1.41421 0.0 0.0] [0.353553 1.09924 0.0] [0.235702 0.15162 1.05899] >>> chop(A - L*L.T) [0.0 0.0 0.0] [0.0 0.0 0.0] [0.0 0.0 0.0] Cholesky decomposition of a Hermitian matrix:: >>> A = eye(3) + matrix([[0,0.25j,-0.5j],[-0.25j,0,0],[0.5j,0,0]]) >>> L = cholesky(A) >>> nprint(L) [ 1.0 0.0 0.0] [(0.0 - 0.25j) (0.968246 + 0.0j) 0.0] [ (0.0 + 0.5j) (0.129099 + 0.0j) (0.856349 + 0.0j)] >>> chop(A - L*L.H) [0.0 0.0 0.0] [0.0 0.0 0.0] [0.0 0.0 0.0] Attempted Cholesky decomposition of a matrix that is not positive definite:: >>> A = -eye(3) + hilbert(3) >>> L = cholesky(A) Traceback (most recent call last): ... ValueError: matrix is not positive-definite **References** 1. [Wikipedia]_ http://en.wikipedia.org/wiki/Cholesky_decomposition """ assert isinstance(A, matrix) if not A.rows == A.cols: raise ValueError('need n*n matrix') if tol is None: tol = +eps n = A.rows L = matrix(n) for j in xrange(n): c = re(A[j,j]) if abs(c-A[j,j]) > tol: raise ValueError('matrix is not Hermitian') s = c - fsum((L[j,k] for k in xrange(j)), absolute=True, squared=True) if s < tol: raise ValueError('matrix is not positive-definite') L[j,j] = math.sqrt(s) for i in xrange(j, n): it1 = (L[i,k] for k in xrange(j)) it2 = (L[j,k] for k in xrange(j)) t = fdot(it1, it2, conjugate=True) L[i,j] = (A[i,j] - t) / L[j,j] return L def cholesky_solve(A, b, **kwargs): """ Ax = b => x Solve a symmetric positive-definite linear equation system. This is twice as efficient as lu_solve. Typical use cases: * A.T*A * Hessian matrix * differential equations """ # do not overwrite A nor b A, b = matrix(A, **kwargs).copy(), matrix(b, **kwargs).copy() if A.rows != A.cols: raise ValueError('can only solve determined system') # Cholesky factorization L = cholesky(A) # solve n = L.rows assert len(b) == n for i in xrange(n): b[i] -= fsum(L[i,j] * b[j] for j in xrange(i)) b[i] /= L[i,i] x = U_solve(L.T, b) return x def det(A): """ Calculate the determinant of a matrix. """ # do not overwrite A A = matrix(A).copy() # use LU factorization to calculate determinant try: R, p = LU_decomp(A) except ZeroDivisionError: return 0 z = 1 for i, e in enumerate(p): if i != e: z *= -1 for i in xrange(A.rows): z *= R[i,i] return z def cond(A, norm=None): """ Calculate the condition number of a matrix using a specified matrix norm. The condition number estimates the sensitivity of a matrix to errors. Example: small input errors for ill-conditioned coefficient matrices alter the solution of the system dramatically. For ill-conditioned matrices it's recommended to use qr_solve() instead of lu_solve(). This does not help with input errors however, it just avoids to add additional errors. Definition: cond(A) = ||A|| * ||A**-1|| """ if norm is None: norm = lambda x: mnorm(x,1) return norm(A) * norm(inverse(A)) def lu_solve_mat(a, b): """Solve a * x = b where a and b are matrices.""" r = matrix(a.rows, b.cols) for i in range(b.cols): c = lu_solve(a, b.column(i)) for j in range(len(c)): r[j, i] = c[j] return r