from __future__ import print_function, division
import collections
from ..core.add import Add
from ..core.basic import Basic, Atom
from ..core.expr import Expr
from ..core.power import Pow
from ..core.symbol import Symbol, Dummy, symbols
from ..core.numbers import Integer, ilcm, Float
from ..core.singleton import S
from ..core.sympify import sympify
from ..core.compatibility import is_sequence, default_sort_key, range, NotIterable
from ..polys import PurePoly, roots, cancel, gcd
from ..simplify import simplify as _simplify, signsimp, nsimplify
from ..utilities.iterables import flatten, numbered_symbols
from ..functions.elementary.miscellaneous import sqrt, Max, Min
from ..functions import Abs, exp, factorial
from ..printing import sstr
from ..core.compatibility import reduce, as_int, string_types
from ..assumptions.refine import refine
from ..core.decorators import call_highest_priority
from types import FunctionType
from .common import (
a2idx,
classof,
MatrixError,
ShapeError,
NonSquareMatrixError,
MatrixCommon,
)
def _iszero(x):
"""Returns True if x is zero."""
try:
return x.is_zero
except AttributeError:
return None
[docs]class DeferredVector(Symbol, NotIterable):
"""A vector whose components are deferred (e.g. for use with lambdify)
Examples
========
>>> from .. import DeferredVector, lambdify
>>> X = DeferredVector( 'X' )
>>> X
X
>>> expr = (X[0] + 2, X[2] + 3)
>>> func = lambdify( X, expr)
>>> func( [1, 2, 3] )
(3, 6)
"""
def __getitem__(self, i):
if i == -0:
i = 0
if i < 0:
raise IndexError("DeferredVector index out of range")
component_name = "%s[%d]" % (self.name, i)
return Symbol(component_name)
def __str__(self):
return sstr(self)
def __repr__(self):
return "DeferredVector('%s')" % self.name
[docs]class MatrixDeterminant(MatrixCommon):
"""Provides basic matrix determinant operations.
Should not be instantiated directly."""
def _eval_berkowitz_toeplitz_matrix(self):
"""Return (A,T) where T the Toeplitz matrix used in the Berkowitz algorithm
corresponding to `self` and A is the first principal submatrix."""
# the 0 x 0 case is trivial
if self.rows == 0 and self.cols == 0:
return self._new(1, 1, [S.One])
#
# Partition self = [ a_11 R ]
# [ C A ]
#
a, R = self[0, 0], self[0, 1:]
C, A = self[1:, 0], self[1:, 1:]
#
# The Toeplitz matrix looks like
#
# [ 1 ]
# [ -a 1 ]
# [ -RC -a 1 ]
# [ -RAC -RC -a 1 ]
# [ -RA**2C -RAC -RC -a 1 ]
# etc.
# Compute the diagonal entries.
# Because multiplying matrix times vector is so much
# more efficient than matrix times matrix, recursively
# compute -R * A**n * C.
diags = [C]
for i in range(self.rows - 2):
diags.append(A * diags[i])
diags = [(-R * d)[0, 0] for d in diags]
diags = [S.One, -a] + diags
def entry(i, j):
if j > i:
return S.Zero
return diags[i - j]
toeplitz = self._new(self.cols + 1, self.rows, entry)
return (A, toeplitz)
def _eval_berkowitz_vector(self):
"""Run the Berkowitz algorithm and return a vector whose entries
are the coefficients of the characteristic polynomial of `self`.
Given N x N matrix, efficiently compute
coefficients of characteristic polynomials of 'self'
without division in the ground domain.
This method is particularly useful for computing determinant,
principal minors and characteristic polynomial when 'self'
has complicated coefficients e.g. polynomials. Semi-direct
usage of this algorithm is also important in computing
efficiently sub-resultant PRS.
Assuming that M is a square matrix of dimension N x N and
I is N x N identity matrix, then the Berkowitz vector is
an N x 1 vector whose entries are coefficients of the
polynomial
charpoly(M) = det(t*I - M)
As a consequence, all polynomials generated by Berkowitz
algorithm are monic.
For more information on the implemented algorithm refer to:
[1] S.J. Berkowitz, On computing the determinant in small
parallel time using a small number of processors, ACM,
Information Processing Letters 18, 1984, pp. 147-150
[2] M. Keber, Division-Free computation of sub-resultants
using Bezout matrices, Tech. Report MPI-I-2006-1-006,
Saarbrucken, 2006
"""
# handle the trivial cases
if self.rows == 0 and self.cols == 0:
return self._new(1, 1, [S.One])
elif self.rows == 1 and self.cols == 1:
return self._new(2, 1, [S.One, -self[0, 0]])
submat, toeplitz = self._eval_berkowitz_toeplitz_matrix()
return toeplitz * submat._eval_berkowitz_vector()
def _eval_det_bareiss(self):
"""Compute matrix determinant using Bareiss' fraction-free
algorithm which is an extension of the well known Gaussian
elimination method. This approach is best suited for dense
symbolic matrices and will result in a determinant with
minimal number of fractions. It means that less term
rewriting is needed on resulting formulae.
TODO: Implement algorithm for sparse matrices (SFF),
http://www.eecis.udel.edu/~saunders/papers/sffge/it5.ps.
"""
# XXX included as a workaround for issue #12362. Should use `_find_reasonable_pivot` instead
def _find_pivot(l):
for pos, val in enumerate(l):
if val:
return (pos, val, None, None)
return (None, None, None, None)
# Recursively implimented Bareiss' algorithm as per Deanna Richelle Leggett's
# thesis http://www.math.usm.edu/perry/Research/Thesis_DRL.pdf
def bareiss(mat, cumm=1):
if mat.rows == 0:
return S.One
elif mat.rows == 1:
return mat[0, 0]
# find a pivot and extract the remaining matrix
# XXX should use `_find_reasonable_pivot`. Blocked by issue #12362
pivot_pos, pivot_val, _, _ = _find_pivot(mat[:, 0])
if pivot_pos == None:
return S.Zero
# if we have a valid pivot, we'll do a "row swap", so keep the
# sign of the det
sign = (-1) ** (pivot_pos % 2)
# we want every row but the pivot row and every column
rows = list(i for i in range(mat.rows) if i != pivot_pos)
cols = list(range(mat.cols))
tmp_mat = mat.extract(rows, cols)
def entry(i, j):
ret = (
pivot_val * tmp_mat[i, j + 1]
- mat[pivot_pos, j + 1] * tmp_mat[i, 0]
) / cumm
if not ret.is_Atom:
cancel(ret)
return ret
return sign * bareiss(
self._new(mat.rows - 1, mat.cols - 1, entry), pivot_val
)
return cancel(bareiss(self))
def _eval_det_berkowitz(self):
"""Use the Berkowitz algorithm to compute the determinant."""
berk_vector = self._eval_berkowitz_vector()
return (-1) ** (len(berk_vector) - 1) * berk_vector[-1]
def _eval_det_lu(self, iszerofunc=_iszero, simpfunc=None):
"""Computes the determinant of a matrix from its LU decomposition.
This function uses the LU decomposition computed by
LUDecomposition_Simple().
The keyword arguments iszerofunc and simpfunc are passed to
LUDecomposition_Simple().
iszerofunc is a callable that returns a boolean indicating if its
input is zero, or None if it cannot make the determination.
simpfunc is a callable that simplifies its input.
The default is simpfunc=None, which indicate that the pivot search
algorithm should not attempt to simplify any candidate pivots.
If simpfunc fails to simplify its input, then it must return its input
instead of a copy."""
if self.rows == 0:
return S.One
# sympy/matrices/tests/test_matrices.py contains a test that
# suggests that the determinant of a 0 x 0 matrix is one, by
# convention.
lu, row_swaps = self.LUdecomposition_Simple(
iszerofunc=iszerofunc, simpfunc=None
)
# P*A = L*U => det(A) = det(L)*det(U)/det(P) = det(P)*det(U).
# Lower triangular factor L encoded in lu has unit diagonal => det(L) = 1.
# P is a permutation matrix => det(P) in {-1, 1} => 1/det(P) = det(P).
# LUdecomposition_Simple() returns a list of row exchange index pairs, rather
# than a permutation matrix, but det(P) = (-1)**len(row_swaps).
# Avoid forming the potentially time consuming product of U's diagonal entries
# if the product is zero.
# Bottom right entry of U is 0 => det(A) = 0.
# It may be impossible to determine if this entry of U is zero when it is symbolic.
if iszerofunc(lu[lu.rows - 1, lu.rows - 1]):
return S.Zero
# Compute det(P)
det = -S.One if len(row_swaps) % 2 else S.One
# Compute det(U) by calculating the product of U's diagonal entries.
# The upper triangular portion of lu is the upper triangular portion of the
# U factor in the LU decomposition.
for k in range(lu.rows):
det *= lu[k, k]
# return det(P)*det(U)
return det
def _eval_determinant(self):
"""Assumed to exist by matrix expressions; If we subclass
MatrixDeterminant, we can fully evaluate determinants."""
return self.det()
[docs] def adjugate(self, method="berkowitz"):
"""Returns the adjugate, or classical adjoint, of
a matrix. That is, the transpose of the matrix of cofactors.
http://en.wikipedia.org/wiki/Adjugate
See Also
========
cofactor_matrix
transpose
"""
return self.cofactor_matrix(method).transpose()
[docs] def charpoly(self, x=Dummy("lambda"), simplify=_simplify):
"""Computes characteristic polynomial det(x*I - self) where I is
the identity matrix.
A PurePoly is returned, so using different variables for ``x`` does
not affect the comparison or the polynomials:
Examples
========
>>> from .. import Matrix
>>> from ..abc import x, y
>>> A = Matrix([[1, 3], [2, 0]])
>>> A.charpoly(x) == A.charpoly(y)
True
Specifying ``x`` is optional; a Dummy with name ``lambda`` is used by
default (which looks good when pretty-printed in unicode):
>>> A.charpoly().as_expr()
_lambda**2 - _lambda - 6
No test is done to see that ``x`` doesn't clash with an existing
symbol, so using the default (``lambda``) or your own Dummy symbol is
the safest option:
>>> A = Matrix([[1, 2], [x, 0]])
>>> A.charpoly().as_expr()
_lambda**2 - _lambda - 2*x
>>> A.charpoly(x).as_expr()
x**2 - 3*x
Notes
=====
The Samuelson-Berkowitz algorithm is used to compute
the characteristic polynomial efficiently and without any
division operations. Thus the characteristic polynomial over any
commutative ring without zero divisors can be computed.
See Also
========
det
"""
if self.rows != self.cols:
raise NonSquareMatrixError()
berk_vector = self._eval_berkowitz_vector()
return PurePoly([simplify(a) for a in berk_vector], x)
[docs] def cofactor(self, i, j, method="berkowitz"):
"""Calculate the cofactor of an element.
See Also
========
cofactor_matrix
minor
minor_submatrix
"""
if self.rows != self.cols or self.rows < 1:
raise NonSquareMatrixError()
return (-1) ** ((i + j) % 2) * self.minor(i, j, method)
[docs] def cofactor_matrix(self, method="berkowitz"):
"""Return a matrix containing the cofactor of each element.
See Also
========
cofactor
minor
minor_submatrix
adjugate
"""
if self.rows != self.cols or self.rows < 1:
raise NonSquareMatrixError()
return self._new(self.rows, self.cols, lambda i, j: self.cofactor(i, j, method))
[docs] def det(self, method="bareiss"):
"""Computes the determinant of a matrix. If the matrix
is at most 3x3, a hard-coded formula is used.
Otherwise, the determinant using the method `method`.
Possible values for "method":
bareis
berkowitz
lu
"""
# sanitize `method`
method = method.lower()
if method == "bareis":
method = "bareiss"
if method == "det_lu":
method = "lu"
if method not in ("bareiss", "berkowitz", "lu"):
raise ValueError("Determinant method '%s' unrecognized" % method)
# if methods were made internal and all determinant calculations
# passed through here, then these lines could be factored out of
# the method routines
if self.rows != self.cols:
raise NonSquareMatrixError()
n = self.rows
if n == 0:
return S.One
elif n == 1:
return self[0, 0]
elif n == 2:
return self[0, 0] * self[1, 1] - self[0, 1] * self[1, 0]
elif n == 3:
return (
self[0, 0] * self[1, 1] * self[2, 2]
+ self[0, 1] * self[1, 2] * self[2, 0]
+ self[0, 2] * self[1, 0] * self[2, 1]
- self[0, 2] * self[1, 1] * self[2, 0]
- self[0, 0] * self[1, 2] * self[2, 1]
- self[0, 1] * self[1, 0] * self[2, 2]
)
if method == "bareiss":
return self._eval_det_bareiss()
elif method == "berkowitz":
return self._eval_det_berkowitz()
elif method == "lu":
return self._eval_det_lu()
[docs] def minor(self, i, j, method="berkowitz"):
"""Return the (i,j) minor of `self`. That is,
return the determinant of the matrix obtained by deleting
the `i`th row and `j`th column from `self`.
See Also
========
minor_submatrix
cofactor
det
"""
if self.rows != self.cols or self.rows < 1:
raise NonSquareMatrixError()
return self.minor_submatrix(i, j).det(method=method)
[docs] def minor_submatrix(self, i, j):
"""Return the submatrix obtained by removing the `i`th row
and `j`th column from `self`.
See Also
========
minor
cofactor
"""
if i < 0:
i += self.rows
if j < 0:
j += self.cols
if not 0 <= i < self.rows or not 0 <= j < self.cols:
raise ValueError(
"`i` and `j` must satisfy 0 <= i < `self.rows` "
"(%d)" % self.rows + "and 0 <= j < `self.cols` (%d)." % self.cols
)
rows = [a for a in range(self.rows) if a != i]
cols = [a for a in range(self.cols) if a != j]
return self.extract(rows, cols)
[docs]class MatrixReductions(MatrixDeterminant):
"""Provides basic matrix row/column operations.
Should not be instantiated directly."""
def _eval_col_op_swap(self, col1, col2):
def entry(i, j):
if j == col1:
return self[i, col2]
elif j == col2:
return self[i, col1]
return self[i, j]
return self._new(self.rows, self.cols, entry)
def _eval_col_op_multiply_col_by_const(self, col, k):
def entry(i, j):
if j == col:
return k * self[i, j]
return self[i, j]
return self._new(self.rows, self.cols, entry)
def _eval_col_op_add_multiple_to_other_col(self, col, k, col2):
def entry(i, j):
if j == col:
return self[i, j] + k * self[i, col2]
return self[i, j]
return self._new(self.rows, self.cols, entry)
def _eval_row_op_swap(self, row1, row2):
def entry(i, j):
if i == row1:
return self[row2, j]
elif i == row2:
return self[row1, j]
return self[i, j]
return self._new(self.rows, self.cols, entry)
def _eval_row_op_multiply_row_by_const(self, row, k):
def entry(i, j):
if i == row:
return k * self[i, j]
return self[i, j]
return self._new(self.rows, self.cols, entry)
def _eval_row_op_add_multiple_to_other_row(self, row, k, row2):
def entry(i, j):
if i == row:
return self[i, j] + k * self[row2, j]
return self[i, j]
return self._new(self.rows, self.cols, entry)
def _eval_echelon_form(self, iszerofunc, simpfunc):
"""Returns (mat, swaps) where `mat` is a row-equivalent matrix
in echelon form and `swaps` is a list of row-swaps performed."""
reduced, pivot_cols, swaps = self._row_reduce(
iszerofunc, simpfunc, normalize_last=True, normalize=False, zero_above=False
)
return reduced, pivot_cols, swaps
def _eval_is_echelon(self, iszerofunc):
if self.rows <= 0 or self.cols <= 0:
return True
zeros_below = all(iszerofunc(t) for t in self[1:, 0])
if iszerofunc(self[0, 0]):
return zeros_below and self[:, 1:]._eval_is_echelon(iszerofunc)
return zeros_below and self[1:, 1:]._eval_is_echelon(iszerofunc)
def _eval_rref(self, iszerofunc, simpfunc, normalize_last=True):
reduced, pivot_cols, swaps = self._row_reduce(
iszerofunc, simpfunc, normalize_last, normalize=True, zero_above=True
)
return reduced, pivot_cols
def _normalize_op_args(self, op, col, k, col1, col2, error_str="col"):
"""Validate the arguments for a row/column operation. `error_str`
can be one of "row" or "col" depending on the arguments being parsed."""
if op not in ["n->kn", "n<->m", "n->n+km"]:
raise ValueError(
"Unknown {} operation '{}'. Valid col operations "
"are 'n->kn', 'n<->m', 'n->n+km'".format(error_str, op)
)
# normalize and validate the arguments
if op == "n->kn":
col = col if col is not None else col1
if col is None or k is None:
raise ValueError(
"For a {0} operation 'n->kn' you must provide the "
"kwargs `{0}` and `k`".format(error_str)
)
if not 0 <= col <= self.cols:
raise ValueError(
"This matrix doesn't have a {} '{}'".format(error_str, col)
)
if op == "n<->m":
# we need two cols to swap. It doesn't matter
# how they were specified, so gather them together and
# remove `None`
cols = set((col, k, col1, col2)).difference([None])
if len(cols) > 2:
# maybe the user left `k` by mistake?
cols = set((col, col1, col2)).difference([None])
if len(cols) != 2:
raise ValueError(
"For a {0} operation 'n<->m' you must provide the "
"kwargs `{0}1` and `{0}2`".format(error_str)
)
col1, col2 = cols
if not 0 <= col1 <= self.cols:
raise ValueError(
"This matrix doesn't have a {} '{}'".format(error_str, col1)
)
if not 0 <= col2 <= self.cols:
raise ValueError(
"This matrix doesn't have a {} '{}'".format(error_str, col2)
)
if op == "n->n+km":
col = col1 if col is None else col
col2 = col1 if col2 is None else col2
if col is None or col2 is None or k is None:
raise ValueError(
"For a {0} operation 'n->n+km' you must provide the "
"kwargs `{0}`, `k`, and `{0}2`".format(error_str)
)
if col == col2:
raise ValueError(
"For a {0} operation 'n->n+km' `{0}` and `{0}2` must "
"be different.".format(error_str)
)
if not 0 <= col <= self.cols:
raise ValueError(
"This matrix doesn't have a {} '{}'".format(error_str, col)
)
if not 0 <= col2 <= self.cols:
raise ValueError(
"This matrix doesn't have a {} '{}'".format(error_str, col2)
)
return op, col, k, col1, col2
def _permute_complexity_right(self, iszerofunc):
"""Permute columns with complicated elements as
far right as they can go. Since the `sympy` row reduction
algorithms start on the left, having complexity right-shifted
speeds things up.
Returns a tuple (mat, perm) where perm is a permutation
of the columns to perform to shift the complex columns right, and mat
is the permuted matrix."""
def complexity(i):
# the complexity of a column will be judged by how many
# element's zero-ness cannot be determined
return sum(1 if iszerofunc(e) is None else 0 for e in self[:, i])
complex = [(complexity(i), i) for i in range(self.cols)]
perm = [j for (i, j) in sorted(complex)]
return (self.permute(perm, orientation="cols"), perm)
def _row_reduce(
self, iszerofunc, simpfunc, normalize_last=True, normalize=True, zero_above=True
):
"""Row reduce `self` and return a tuple (rref_matrix,
pivot_cols, swaps) where pivot_cols are the pivot columns
and swaps are any row swaps that were used in the process
of row reduction.
Parameters
==========
iszerofunc : determines if an entry can be used as a pivot
simpfunc : used to simplify elements and test if they are
zero if `iszerofunc` returns `None`
normalize_last : indicates where all row reduction should
happen in a fraction-free manner and then the rows are
normalized (so that the pivots are 1), or whether
rows should be normalized along the way (like the naive
row reduction algorithm)
normalize : whether pivot rows should be normalized so that
the pivot value is 1
zero_above : whether entries above the pivot should be zeroed.
If `zero_above=False`, an echelon matrix will be returned.
"""
rows, cols = self.rows, self.cols
mat = list(self)
def get_col(i):
return mat[i::cols]
def row_swap(i, j):
mat[i * cols : (i + 1) * cols], mat[j * cols : (j + 1) * cols] = (
mat[j * cols : (j + 1) * cols],
mat[i * cols : (i + 1) * cols],
)
def cross_cancel(a, i, b, j):
"""Does the row op row[i] = a*row[i] - b*row[j]"""
q = (j - i) * cols
for p in range(i * cols, (i + 1) * cols):
mat[p] = a * mat[p] - b * mat[p + q]
piv_row, piv_col = 0, 0
pivot_cols = []
swaps = []
# use a fraction free method to zero above and below each pivot
while piv_col < cols and piv_row < rows:
(
pivot_offset,
pivot_val,
assumed_nonzero,
newly_determined,
) = _find_reasonable_pivot(get_col(piv_col)[piv_row:], iszerofunc, simpfunc)
# _find_reasonable_pivot may have simplified some things
# in the process. Let's not let them go to waste
for (offset, val) in newly_determined:
offset += piv_row
mat[offset * cols + piv_col] = val
if pivot_offset is None:
piv_col += 1
continue
pivot_cols.append(piv_col)
if pivot_offset != 0:
row_swap(piv_row, pivot_offset + piv_row)
swaps.append((piv_row, pivot_offset + piv_row))
# if we aren't normalizing last, we normalize
# before we zero the other rows
if normalize_last is False:
i, j = piv_row, piv_col
mat[i * cols + j] = S.One
for p in range(i * cols + j + 1, (i + 1) * cols):
mat[p] = mat[p] / pivot_val
# after normalizing, the pivot value is 1
pivot_val = S.One
# zero above and below the pivot
for row in range(rows):
# don't zero our current row
if row == piv_row:
continue
# don't zero above the pivot unless we're told.
if zero_above is False and row < piv_row:
continue
# if we're already a zero, don't do anything
val = mat[row * cols + piv_col]
if iszerofunc(val):
continue
cross_cancel(pivot_val, row, val, piv_row)
piv_row += 1
# normalize each row
if normalize_last is True and normalize is True:
for piv_i, piv_j in enumerate(pivot_cols):
pivot_val = mat[piv_i * cols + piv_j]
mat[piv_i * cols + piv_j] = S.One
for p in range(piv_i * cols + piv_j + 1, (piv_i + 1) * cols):
mat[p] = mat[p] / pivot_val
return self._new(self.rows, self.cols, mat), tuple(pivot_cols), tuple(swaps)
[docs] def elementary_col_op(self, op="n->kn", col=None, k=None, col1=None, col2=None):
"""Perfoms the elementary column operation `op`.
`op` may be one of
* "n->kn" (column n goes to k*n)
* "n<->m" (swap column n and column m)
* "n->n+km" (column n goes to column n + k*column m)
Parameters
=========
op : string; the elementary row operation
col : the column to apply the column operation
k : the multiple to apply in the column operation
col1 : one column of a column swap
col2 : second column of a column swap or column "m" in the column operation
"n->n+km"
"""
op, col, k, col1, col2 = self._normalize_op_args(op, col, k, col1, col2, "col")
# now that we've validated, we're all good to dispatch
if op == "n->kn":
return self._eval_col_op_multiply_col_by_const(col, k)
if op == "n<->m":
return self._eval_col_op_swap(col1, col2)
if op == "n->n+km":
return self._eval_col_op_add_multiple_to_other_col(col, k, col2)
[docs] def elementary_row_op(self, op="n->kn", row=None, k=None, row1=None, row2=None):
"""Perfoms the elementary row operation `op`.
`op` may be one of
* "n->kn" (row n goes to k*n)
* "n<->m" (swap row n and row m)
* "n->n+km" (row n goes to row n + k*row m)
Parameters
==========
op : string; the elementary row operation
row : the row to apply the row operation
k : the multiple to apply in the row operation
row1 : one row of a row swap
row2 : second row of a row swap or row "m" in the row operation
"n->n+km"
"""
op, row, k, row1, row2 = self._normalize_op_args(op, row, k, row1, row2, "row")
# now that we've validated, we're all good to dispatch
if op == "n->kn":
return self._eval_row_op_multiply_row_by_const(row, k)
if op == "n<->m":
return self._eval_row_op_swap(row1, row2)
if op == "n->n+km":
return self._eval_row_op_add_multiple_to_other_row(row, k, row2)
@property
def is_echelon(self, iszerofunc=_iszero):
"""Returns `True` if he matrix is in echelon form.
That is, all rows of zeros are at the bottom, and below
each leading non-zero in a row are exclusively zeros."""
return self._eval_is_echelon(iszerofunc)
[docs] def rank(self, iszerofunc=_iszero, simplify=False):
"""
Returns the rank of a matrix
>>> from .. import Matrix
>>> from ..abc import x
>>> m = Matrix([[1, 2], [x, 1 - 1/x]])
>>> m.rank()
2
>>> n = Matrix(3, 3, range(1, 10))
>>> n.rank()
2
"""
simpfunc = simplify if isinstance(simplify, FunctionType) else _simplify
# for small matrices, we compute the rank explicitly
# if is_zero on elements doesn't answer the question
# for small matrices, we fall back to the full routine.
if self.rows <= 0 or self.cols <= 0:
return 0
if self.rows <= 1 or self.cols <= 1:
zeros = [iszerofunc(x) for x in self]
if False in zeros:
return 1
if self.rows == 2 and self.cols == 2:
zeros = [iszerofunc(x) for x in self]
if not False in zeros and not None in zeros:
return 0
det = self.det()
if iszerofunc(det) and False in zeros:
return 1
if iszerofunc(det) is False:
return 2
mat, _ = self._permute_complexity_right(iszerofunc=iszerofunc)
echelon_form, pivots, swaps = mat._eval_echelon_form(
iszerofunc=iszerofunc, simpfunc=simpfunc
)
return len(pivots)
[docs] def rref(
self, iszerofunc=_iszero, simplify=False, pivots=True, normalize_last=True
):
"""Return reduced row-echelon form of matrix and indices of pivot vars.
Parameters
==========
iszerofunc : Function
A function used for detecting whether an element can
act as a pivot. `lambda x: x.is_zero` is used by default.
simplify : Function
A function used to simplify elements when looking for a pivot.
By default SymPy's `simplify`is used.
pivots : True or False
If `True`, a tuple containing the row-reduced matrix and a tuple
of pivot columns is returned. If `False` just the row-reduced
matrix is returned.
normalize_last : True or False
If `True`, no pivots are normalized to `1` until after all entries
above and below each pivot are zeroed. This means the row
reduction algorithm is fraction free until the very last step.
If `False`, the naive row reduction procedure is used where
each pivot is normalized to be `1` before row operations are
used to zero above and below the pivot.
Notes
=====
The default value of `normalize_last=True` can provide significant
speedup to row reduction, especially on matrices with symbols. However,
if you depend on the form row reduction algorithm leaves entries
of the matrix, set `noramlize_last=False`
Examples
========
>>> from .. import Matrix
>>> from ..abc import x
>>> m = Matrix([[1, 2], [x, 1 - 1/x]])
>>> m.rref()
(Matrix([
[1, 0],
[0, 1]]), (0, 1))
>>> rref_matrix, rref_pivots = m.rref()
>>> rref_matrix
Matrix([
[1, 0],
[0, 1]])
>>> rref_pivots
(0, 1)
"""
simpfunc = simplify if isinstance(simplify, FunctionType) else _simplify
ret, pivot_cols = self._eval_rref(
iszerofunc=iszerofunc, simpfunc=simpfunc, normalize_last=normalize_last
)
if pivots:
ret = (ret, pivot_cols)
return ret
[docs]class MatrixSubspaces(MatrixReductions):
"""Provides methods relating to the fundamental subspaces
of a matrix. Should not be instantiated directly."""
[docs] def columnspace(self, simplify=False):
"""Returns a list of vectors (Matrix objects) that span columnspace of self
Examples
========
>>> from ..matrices import Matrix
>>> m = Matrix(3, 3, [1, 3, 0, -2, -6, 0, 3, 9, 6])
>>> m
Matrix([
[ 1, 3, 0],
[-2, -6, 0],
[ 3, 9, 6]])
>>> m.columnspace()
[Matrix([
[ 1],
[-2],
[ 3]]), Matrix([
[0],
[0],
[6]])]
See Also
========
nullspace
rowspace
"""
reduced, pivots = self.echelon_form(simplify=simplify, with_pivots=True)
return [self.col(i) for i in pivots]
[docs] def nullspace(self, simplify=False):
"""Returns list of vectors (Matrix objects) that span nullspace of self
Examples
========
>>> from ..matrices import Matrix
>>> m = Matrix(3, 3, [1, 3, 0, -2, -6, 0, 3, 9, 6])
>>> m
Matrix([
[ 1, 3, 0],
[-2, -6, 0],
[ 3, 9, 6]])
>>> m.nullspace()
[Matrix([
[-3],
[ 1],
[ 0]])]
See Also
========
columnspace
rowspace
"""
reduced, pivots = self.rref(simplify=simplify)
free_vars = [i for i in range(self.cols) if i not in pivots]
basis = []
for free_var in free_vars:
# for each free variable, we will set it to 1 and all others
# to 0. Then, we will use back substitution to solve the system
vec = [S.Zero] * self.cols
vec[free_var] = S.One
for piv_row, piv_col in enumerate(pivots):
for pos in pivots[piv_row + 1 :] + (free_var,):
vec[piv_col] -= reduced[piv_row, pos]
basis.append(vec)
return [self._new(self.cols, 1, b) for b in basis]
[docs] def rowspace(self, simplify=False):
"""Returns a list of vectors that span the row space of self."""
reduced, pivots = self.echelon_form(simplify=simplify, with_pivots=True)
return [reduced.row(i) for i in range(len(pivots))]
[docs] @classmethod
def orthogonalize(cls, *vecs, **kwargs):
"""Apply the Gram-Schmidt orthogonalization procedure
to vectors supplied in `vecs`.
Arguments
=========
vecs : vectors to be made orthogonal
normalize : bool. Whether the returned vectors
should be renormalized to be unit vectors.
"""
normalize = kwargs.get("normalize", False)
def project(a, b):
return b * (a.dot(b) / b.dot(b))
def perp_to_subspace(vec, basis):
"""projects vec onto the subspace given
by the orthogonal basis `basis`"""
components = [project(vec, b) for b in basis]
if len(basis) == 0:
return vec
return vec - reduce(lambda a, b: a + b, components)
ret = []
# make sure we start with a non-zero vector
while len(vecs) > 0 and vecs[0].is_zero:
del vecs[0]
for vec in vecs:
perp = perp_to_subspace(vec, ret)
if not perp.is_zero:
ret.append(perp)
if normalize:
ret = [vec / vec.norm() for vec in ret]
return ret
[docs]class MatrixEigen(MatrixSubspaces):
"""Provides basic matrix eigenvalue/vector operations.
Should not be instantiated directly."""
_cache_is_diagonalizable = None
_cache_eigenvects = None
[docs] def diagonalize(self, reals_only=False, sort=False, normalize=False):
"""
Return (P, D), where D is diagonal and
D = P^-1 * M * P
where M is current matrix.
Parameters
==========
reals_only : bool. Whether to throw an error if complex numbers are need
to diagonalize. (Default: False)
sort : bool. Sort the eigenvalues along the diagonal. (Default: False)
normalize : bool. If True, normalize the columns of P. (Default: False)
Examples
========
>>> from .. import Matrix
>>> m = Matrix(3, 3, [1, 2, 0, 0, 3, 0, 2, -4, 2])
>>> m
Matrix([
[1, 2, 0],
[0, 3, 0],
[2, -4, 2]])
>>> (P, D) = m.diagonalize()
>>> D
Matrix([
[1, 0, 0],
[0, 2, 0],
[0, 0, 3]])
>>> P
Matrix([
[-1, 0, -1],
[ 0, 0, -1],
[ 2, 1, 2]])
>>> P.inv() * m * P
Matrix([
[1, 0, 0],
[0, 2, 0],
[0, 0, 3]])
See Also
========
is_diagonal
is_diagonalizable
"""
if not self.is_square:
raise NonSquareMatrixError()
if not self.is_diagonalizable(reals_only=reals_only, clear_cache=False):
raise MatrixError("Matrix is not diagonalizable")
eigenvecs = self._cache_eigenvects
if eigenvecs is None:
eigenvecs = self.eigenvects(simplify=True)
if sort:
eigenvecs = sorted(eigenvecs, key=default_sort_key)
p_cols, diag = [], []
for val, mult, basis in eigenvecs:
diag += [val] * mult
p_cols += basis
if normalize:
p_cols = [v / v.norm() for v in p_cols]
return self.hstack(*p_cols), self.diag(*diag)
[docs] def eigenvals(self, error_when_incomplete=True, **flags):
"""Return eigenvalues using the Berkowitz agorithm to compute
the characteristic polynomial.
Parameters
==========
error_when_incomplete : bool
Raise an error when not all eigenvalues are computed. This is
caused by ``roots`` not returning a full list of eigenvalues.
Since the roots routine doesn't always work well with Floats,
they will be replaced with Rationals before calling that
routine. If this is not desired, set flag ``rational`` to False.
"""
mat = self
if not mat:
return {}
if flags.pop("rational", True):
if any(v.has(Float) for v in mat):
mat = mat.applyfunc(lambda x: nsimplify(x, rational=True))
flags.pop("simplify", None) # pop unsupported flag
eigs = roots(mat.charpoly(x=Dummy("x")), **flags)
# make sure the algebraic multiplicty sums to the
# size of the matrix
if error_when_incomplete and sum(m for m in eigs.values()) != self.cols:
raise MatrixError("Could not compute eigenvalues for {}".format(self))
return eigs
[docs] def eigenvects(self, error_when_incomplete=True, **flags):
"""Return list of triples (eigenval, multiplicity, basis).
The flag ``simplify`` has two effects:
1) if bool(simplify) is True, as_content_primitive()
will be used to tidy up normalization artifacts;
2) if nullspace needs simplification to compute the
basis, the simplify flag will be passed on to the
nullspace routine which will interpret it there.
Parameters
==========
error_when_incomplete : bool
Raise an error when not all eigenvalues are computed. This is
caused by ``roots`` not returning a full list of eigenvalues.
If the matrix contains any Floats, they will be changed to Rationals
for computation purposes, but the answers will be returned after being
evaluated with evalf. If it is desired to removed small imaginary
portions during the evalf step, pass a value for the ``chop`` flag.
"""
from ..matrices import eye
simplify = flags.get("simplify", True)
if not isinstance(simplify, FunctionType):
simpfunc = _simplify if simplify else lambda x: x
primitive = flags.get("simplify", False)
chop = flags.pop("chop", False)
flags.pop("multiple", None) # remove this if it's there
mat = self
# roots doesn't like Floats, so replace them with Rationals
has_floats = any(v.has(Float) for v in self)
if has_floats:
mat = mat.applyfunc(lambda x: nsimplify(x, rational=True))
def eigenspace(eigenval):
"""Get a basis for the eigenspace for a particular eigenvalue"""
m = mat - self.eye(mat.rows) * eigenval
ret = m.nullspace()
# the nullspace for a real eigenvalue should be
# non-trivial. If we didn't find an eigenvector, try once
# more a little harder
if len(ret) == 0 and simplify:
ret = m.nullspace(simplify=True)
if len(ret) == 0:
raise NotImplementedError(
"Can't evaluate eigenvector for eigenvalue %s" % eigenval
)
return ret
eigenvals = mat.eigenvals(
rational=False, error_when_incomplete=error_when_incomplete, **flags
)
ret = [
(val, mult, eigenspace(val))
for val, mult in sorted(eigenvals.items(), key=default_sort_key)
]
if primitive:
# if the primitive flag is set, get rid of any common
# integer denominators
def denom_clean(l):
from .. import gcd
return [(v / gcd(list(v))).applyfunc(simpfunc) for v in l]
ret = [(val, mult, denom_clean(es)) for val, mult, es in ret]
if has_floats:
# if we had floats to start with, turn the eigenvectors to floats
ret = [
(val.evalf(chop=chop), mult, [v.evalf(chop=chop) for v in es])
for val, mult, es in ret
]
return ret
[docs] def is_diagonalizable(self, reals_only=False, **kwargs):
"""Returns true if a matrix is diagonalizable.
Parameters
==========
reals_only : bool. If reals_only=True, determine whether the matrix can be
diagonalized without complex numbers. (Default: False)
kwargs
======
clear_cache : bool. If True, clear the result of any computations when finished.
(Default: True)
Examples
========
>>> from .. import Matrix
>>> m = Matrix(3, 3, [1, 2, 0, 0, 3, 0, 2, -4, 2])
>>> m
Matrix([
[1, 2, 0],
[0, 3, 0],
[2, -4, 2]])
>>> m.is_diagonalizable()
True
>>> m = Matrix(2, 2, [0, 1, 0, 0])
>>> m
Matrix([
[0, 1],
[0, 0]])
>>> m.is_diagonalizable()
False
>>> m = Matrix(2, 2, [0, 1, -1, 0])
>>> m
Matrix([
[ 0, 1],
[-1, 0]])
>>> m.is_diagonalizable()
True
>>> m.is_diagonalizable(reals_only=True)
False
See Also
========
is_diagonal
diagonalize
"""
clear_cache = kwargs.get("clear_cache", True)
if "clear_subproducts" in kwargs:
clear_cache = kwargs.get("clear_subproducts")
def cleanup():
"""Clears any cached values if requested"""
if clear_cache:
self._cache_eigenvects = None
self._cache_is_diagonalizable = None
if not self.is_square:
cleanup()
return False
# use the cached value if we have it
if self._cache_is_diagonalizable is not None:
ret = self._cache_is_diagonalizable
cleanup()
return ret
if all(e.is_real for e in self) and self.is_symmetric():
# every real symmetric matrix is real diagonalizable
self._cache_is_diagonalizable = True
cleanup()
return True
self._cache_eigenvects = self.eigenvects(simplify=True)
ret = True
for val, mult, basis in self._cache_eigenvects:
# if we have a complex eigenvalue
if reals_only and not val.is_real:
ret = False
# if the geometric multiplicity doesn't equal the algebraic
if mult != len(basis):
ret = False
cleanup()
return ret
[docs] def left_eigenvects(self, **flags):
"""Returns left eigenvectors and eigenvalues.
This function returns the list of triples (eigenval, multiplicity,
basis) for the left eigenvectors. Options are the same as for
eigenvects(), i.e. the ``**flags`` arguments gets passed directly to
eigenvects().
Examples
========
>>> from .. import Matrix
>>> M = Matrix([[0, 1, 1], [1, 0, 0], [1, 1, 1]])
>>> M.eigenvects()
[(-1, 1, [Matrix([
[-1],
[ 1],
[ 0]])]), (0, 1, [Matrix([
[ 0],
[-1],
[ 1]])]), (2, 1, [Matrix([
[2/3],
[1/3],
[ 1]])])]
>>> M.left_eigenvects()
[(-1, 1, [Matrix([[-2, 1, 1]])]), (0, 1, [Matrix([[-1, -1, 1]])]), (2,
1, [Matrix([[1, 1, 1]])])]
"""
eigs = self.transpose().eigenvects(**flags)
return [
(val, mult, [l.transpose() for l in basis]) for val, mult, basis in eigs
]
[docs] def singular_values(self):
"""Compute the singular values of a Matrix
Examples
========
>>> from .. import Matrix, Symbol
>>> x = Symbol('x', real=True)
>>> A = Matrix([[0, 1, 0], [0, x, 0], [-1, 0, 0]])
>>> A.singular_values()
[sqrt(x**2 + 1), 1, 0]
See Also
========
condition_number
"""
mat = self
# Compute eigenvalues of A.H A
valmultpairs = (mat.H * mat).eigenvals()
# Expands result from eigenvals into a simple list
vals = []
for k, v in valmultpairs.items():
vals += [sqrt(k)] * v # dangerous! same k in several spots!
# sort them in descending order
vals.sort(reverse=True, key=default_sort_key)
return vals
[docs]class MatrixCalculus(MatrixCommon):
"""Provides calculus-related matrix operations."""
[docs] def diff(self, *args):
"""Calculate the derivative of each element in the matrix.
``args`` will be passed to the ``integrate`` function.
Examples
========
>>> from ..matrices import Matrix
>>> from ..abc import x, y
>>> M = Matrix([[x, y], [1, 0]])
>>> M.diff(x)
Matrix([
[1, 0],
[0, 0]])
See Also
========
integrate
limit
"""
return self.applyfunc(lambda x: x.diff(*args))
[docs] def integrate(self, *args):
"""Integrate each element of the matrix. ``args`` will
be passed to the ``integrate`` function.
Examples
========
>>> from ..matrices import Matrix
>>> from ..abc import x, y
>>> M = Matrix([[x, y], [1, 0]])
>>> M.integrate((x, ))
Matrix([
[x**2/2, x*y],
[ x, 0]])
>>> M.integrate((x, 0, 2))
Matrix([
[2, 2*y],
[2, 0]])
See Also
========
limit
diff
"""
return self.applyfunc(lambda x: x.integrate(*args))
[docs] def jacobian(self, X):
"""Calculates the Jacobian matrix (derivative of a vector-valued function).
Parameters
==========
self : vector of expressions representing functions f_i(x_1, ..., x_n).
X : set of x_i's in order, it can be a list or a Matrix
Both self and X can be a row or a column matrix in any order
(i.e., jacobian() should always work).
Examples
========
>>> from .. import sin, cos, Matrix
>>> from ..abc import rho, phi
>>> X = Matrix([rho*cos(phi), rho*sin(phi), rho**2])
>>> Y = Matrix([rho, phi])
>>> X.jacobian(Y)
Matrix([
[cos(phi), -rho*sin(phi)],
[sin(phi), rho*cos(phi)],
[ 2*rho, 0]])
>>> X = Matrix([rho*cos(phi), rho*sin(phi)])
>>> X.jacobian(Y)
Matrix([
[cos(phi), -rho*sin(phi)],
[sin(phi), rho*cos(phi)]])
See Also
========
hessian
wronskian
"""
if not isinstance(X, MatrixBase):
X = self._new(X)
# Both X and self can be a row or a column matrix, so we need to make
# sure all valid combinations work, but everything else fails:
if self.shape[0] == 1:
m = self.shape[1]
elif self.shape[1] == 1:
m = self.shape[0]
else:
raise TypeError("self must be a row or a column matrix")
if X.shape[0] == 1:
n = X.shape[1]
elif X.shape[1] == 1:
n = X.shape[0]
else:
raise TypeError("X must be a row or a column matrix")
# m is the number of functions and n is the number of variables
# computing the Jacobian is now easy:
return self._new(m, n, lambda j, i: self[j].diff(X[i]))
[docs] def limit(self, *args):
"""Calculate the limit of each element in the matrix.
``args`` will be passed to the ``limit`` function.
Examples
========
>>> from ..matrices import Matrix
>>> from ..abc import x, y
>>> M = Matrix([[x, y], [1, 0]])
>>> M.limit(x, 2)
Matrix([
[2, y],
[1, 0]])
See Also
========
integrate
diff
"""
return self.applyfunc(lambda x: x.limit(*args))
# https://github.com/sympy/sympy/pull/12854
[docs]class MatrixDeprecated(MatrixCommon):
"""A class to house deprecated matrix methods."""
[docs] def berkowitz_charpoly(self, x=Dummy("lambda"), simplify=_simplify):
return self.charpoly(x=x)
[docs] def berkowitz_det(self):
"""Computes determinant using Berkowitz method.
See Also
========
det
berkowitz
"""
return self.det(method="berkowitz")
[docs] def berkowitz_eigenvals(self, **flags):
"""Computes eigenvalues of a Matrix using Berkowitz method.
See Also
========
berkowitz
"""
return self.eigenvals(**flags)
[docs] def berkowitz_minors(self):
"""Computes principal minors using Berkowitz method.
See Also
========
berkowitz
"""
sign, minors = S.One, []
for poly in self.berkowitz():
minors.append(sign * poly[-1])
sign = -sign
return tuple(minors)
[docs] def berkowitz(self):
from ..matrices import zeros
berk = ((1,),)
if not self:
return berk
if not self.is_square:
raise NonSquareMatrixError()
A, N = self, self.rows
transforms = [0] * (N - 1)
for n in range(N, 1, -1):
T, k = zeros(n + 1, n), n - 1
R, C = -A[k, :k], A[:k, k]
A, a = A[:k, :k], -A[k, k]
items = [C]
for i in range(0, n - 2):
items.append(A * items[i])
for i, B in enumerate(items):
items[i] = (R * B)[0, 0]
items = [S.One, a] + items
for i in range(n):
T[i:, i] = items[: n - i + 1]
transforms[k - 1] = T
polys = [self._new([S.One, -A[0, 0]])]
for i, T in enumerate(transforms):
polys.append(T * polys[i])
return berk + tuple(map(tuple, polys))
[docs] def cofactorMatrix(self, method="berkowitz"):
return self.cofactor_matrix(method=method)
[docs] def det_bareis(self):
return self.det(method="bareiss")
[docs] def det_bareiss(self):
"""Compute matrix determinant using Bareiss' fraction-free
algorithm which is an extension of the well known Gaussian
elimination method. This approach is best suited for dense
symbolic matrices and will result in a determinant with
minimal number of fractions. It means that less term
rewriting is needed on resulting formulae.
TODO: Implement algorithm for sparse matrices (SFF),
http://www.eecis.udel.edu/~saunders/papers/sffge/it5.ps.
See Also
========
det
berkowitz_det
"""
return self.det(method="bareiss")
[docs] def det_LU_decomposition(self):
"""Compute matrix determinant using LU decomposition
Note that this method fails if the LU decomposition itself
fails. In particular, if the matrix has no inverse this method
will fail.
TODO: Implement algorithm for sparse matrices (SFF),
http://www.eecis.udel.edu/~saunders/papers/sffge/it5.ps.
See Also
========
det
det_bareiss
berkowitz_det
"""
return self.det(method="lu")
[docs] def jordan_cell(self, eigenval, n):
return self.jordan_block(size=n, eigenvalue=eigenval)
[docs] def jordan_cells(self, calc_transformation=True):
P, J = self.jordan_form()
return P, J.get_diag_blocks()
[docs] def minorEntry(self, i, j, method="berkowitz"):
return self.minor(i, j, method=method)
[docs] def minorMatrix(self, i, j):
return self.minor_submatrix(i, j)
[docs] def permuteBkwd(self, perm):
"""Permute the rows of the matrix with the given permutation in reverse."""
return self.permute_rows(perm, direction="backward")
[docs] def permuteFwd(self, perm):
"""Permute the rows of the matrix with the given permutation."""
return self.permute_rows(perm, direction="forward")
[docs]class MatrixBase(MatrixDeprecated, MatrixCalculus, MatrixEigen, MatrixCommon):
"""Base class for matrix objects."""
# Added just for numpy compatibility
__array_priority__ = 11
is_Matrix = True
_class_priority = 3
_sympify = staticmethod(sympify)
__hash__ = None # Mutable
def __array__(self):
from .dense import matrix2numpy
return matrix2numpy(self)
def __getattr__(self, attr):
if attr in ("diff", "integrate", "limit"):
def doit(*args):
item_doit = lambda item: getattr(item, attr)(*args)
return self.applyfunc(item_doit)
return doit
else:
raise AttributeError(
"%s has no attribute %s." % (self.__class__.__name__, attr)
)
def __len__(self):
"""Return the number of elements of self.
Implemented mainly so bool(Matrix()) == False.
"""
return self.rows * self.cols
def __mathml__(self):
mml = ""
for i in range(self.rows):
mml += "<matrixrow>"
for j in range(self.cols):
mml += self[i, j].__mathml__()
mml += "</matrixrow>"
return "<matrix>" + mml + "</matrix>"
# needed for python 2 compatibility
def __ne__(self, other):
return not self == other
def _matrix_pow_by_jordan_blocks(self, num):
from ..matrices import diag, MutableMatrix
from .. import binomial
def jordan_cell_power(jc, n):
N = jc.shape[0]
l = jc[0, 0]
if l == 0 and (n < N - 1) != False:
raise ValueError("Matrix det == 0; not invertible")
elif l == 0 and N > 1 and n % 1 != 0:
raise ValueError("Non-integer power cannot be evaluated")
for i in range(N):
for j in range(N - i):
bn = binomial(n, i)
if isinstance(bn, binomial):
bn = bn._eval_expand_func()
jc[j, i + j] = l ** (n - i) * bn
P, J = self.jordan_form()
jordan_cells = J.get_diag_blocks()
# Make sure jordan_cells matrices are mutable:
jordan_cells = [MutableMatrix(j) for j in jordan_cells]
for j in jordan_cells:
jordan_cell_power(j, num)
return self._new(P * diag(*jordan_cells) * P.inv())
def __repr__(self):
return sstr(self)
def __str__(self):
if self.rows == 0 or self.cols == 0:
return "Matrix(%s, %s, [])" % (self.rows, self.cols)
return "Matrix(%s)" % str(self.tolist())
def _diagonalize_clear_subproducts(self):
del self._is_symbolic
del self._is_symmetric
del self._eigenvects
def _format_str(self, printer=None):
if not printer:
from ..printing.str import StrPrinter
printer = StrPrinter()
# Handle zero dimensions:
if self.rows == 0 or self.cols == 0:
return "Matrix(%s, %s, [])" % (self.rows, self.cols)
if self.rows == 1:
return "Matrix([%s])" % self.table(printer, rowsep=",\n")
return "Matrix([\n%s])" % self.table(printer, rowsep=",\n")
@classmethod
def _handle_creation_inputs(cls, *args, **kwargs):
"""Return the number of rows, cols and flat matrix elements.
Examples
========
>>> from .. import Matrix, I
Matrix can be constructed as follows:
* from a nested list of iterables
>>> Matrix( ((1, 2+I), (3, 4)) )
Matrix([
[1, 2 + I],
[3, 4]])
* from un-nested iterable (interpreted as a column)
>>> Matrix( [1, 2] )
Matrix([
[1],
[2]])
* from un-nested iterable with dimensions
>>> Matrix(1, 2, [1, 2] )
Matrix([[1, 2]])
* from no arguments (a 0 x 0 matrix)
>>> Matrix()
Matrix(0, 0, [])
* from a rule
>>> Matrix(2, 2, lambda i, j: i/(j + 1) )
Matrix([
[0, 0],
[1, 1/2]])
"""
from .sparse import SparseMatrix
flat_list = None
if len(args) == 1:
# Matrix(SparseMatrix(...))
if isinstance(args[0], SparseMatrix):
return args[0].rows, args[0].cols, flatten(args[0].tolist())
# Matrix(Matrix(...))
elif isinstance(args[0], MatrixBase):
return args[0].rows, args[0].cols, args[0]._mat
# Matrix(MatrixSymbol('X', 2, 2))
elif isinstance(args[0], Basic) and args[0].is_Matrix:
return args[0].rows, args[0].cols, args[0].as_explicit()._mat
# Matrix(numpy.ones((2, 2)))
elif hasattr(args[0], "__array__"):
# NumPy array or matrix or some other object that implements
# __array__. So let's first use this method to get a
# numpy.array() and then make a python list out of it.
arr = args[0].__array__()
if len(arr.shape) == 2:
rows, cols = arr.shape[0], arr.shape[1]
flat_list = [cls._sympify(i) for i in arr.ravel()]
return rows, cols, flat_list
elif len(arr.shape) == 1:
rows, cols = arr.shape[0], 1
flat_list = [S.Zero] * rows
for i in range(len(arr)):
flat_list[i] = cls._sympify(arr[i])
return rows, cols, flat_list
else:
raise NotImplementedError("SymPy supports just 1D and 2D matrices")
# Matrix([1, 2, 3]) or Matrix([[1, 2], [3, 4]])
elif is_sequence(args[0]) and not isinstance(args[0], DeferredVector):
in_mat = []
ncol = set()
for row in args[0]:
if isinstance(row, MatrixBase):
in_mat.extend(row.tolist())
if row.cols or row.rows: # only pay attention if it's not 0x0
ncol.add(row.cols)
else:
in_mat.append(row)
try:
ncol.add(len(row))
except TypeError:
ncol.add(1)
if len(ncol) > 1:
raise ValueError(
"Got rows of variable lengths: %s" % sorted(list(ncol))
)
cols = ncol.pop() if ncol else 0
rows = len(in_mat) if cols else 0
if rows:
if not is_sequence(in_mat[0]):
cols = 1
flat_list = [cls._sympify(i) for i in in_mat]
return rows, cols, flat_list
flat_list = []
for j in range(rows):
for i in range(cols):
flat_list.append(cls._sympify(in_mat[j][i]))
elif len(args) == 3:
rows = as_int(args[0])
cols = as_int(args[1])
if rows < 0 or cols < 0:
raise ValueError(
"Cannot create a {} x {} matrix. "
"Both dimensions must be positive".format(rows, cols)
)
# Matrix(2, 2, lambda i, j: i+j)
if len(args) == 3 and isinstance(args[2], collections.abc.Callable):
op = args[2]
flat_list = []
for i in range(rows):
flat_list.extend(
[
cls._sympify(op(cls._sympify(i), cls._sympify(j)))
for j in range(cols)
]
)
# Matrix(2, 2, [1, 2, 3, 4])
elif len(args) == 3 and is_sequence(args[2]):
flat_list = args[2]
if len(flat_list) != rows * cols:
raise ValueError("List length should be equal to rows*columns")
flat_list = [cls._sympify(i) for i in flat_list]
# Matrix()
elif len(args) == 0:
# Empty Matrix
rows = cols = 0
flat_list = []
if flat_list is None:
raise TypeError("Data type not understood")
return rows, cols, flat_list
def _setitem(self, key, value):
"""Helper to set value at location given by key.
Examples
========
>>> from .. import Matrix, I, zeros, ones
>>> m = Matrix(((1, 2+I), (3, 4)))
>>> m
Matrix([
[1, 2 + I],
[3, 4]])
>>> m[1, 0] = 9
>>> m
Matrix([
[1, 2 + I],
[9, 4]])
>>> m[1, 0] = [[0, 1]]
To replace row r you assign to position r*m where m
is the number of columns:
>>> M = zeros(4)
>>> m = M.cols
>>> M[3*m] = ones(1, m)*2; M
Matrix([
[0, 0, 0, 0],
[0, 0, 0, 0],
[0, 0, 0, 0],
[2, 2, 2, 2]])
And to replace column c you can assign to position c:
>>> M[2] = ones(m, 1)*4; M
Matrix([
[0, 0, 4, 0],
[0, 0, 4, 0],
[0, 0, 4, 0],
[2, 2, 4, 2]])
"""
from .dense import Matrix
is_slice = isinstance(key, slice)
i, j = key = self.key2ij(key)
is_mat = isinstance(value, MatrixBase)
if type(i) is slice or type(j) is slice:
if is_mat:
self.copyin_matrix(key, value)
return
if not isinstance(value, Expr) and is_sequence(value):
self.copyin_list(key, value)
return
raise ValueError("unexpected value: %s" % value)
else:
if not is_mat and not isinstance(value, Basic) and is_sequence(value):
value = Matrix(value)
is_mat = True
if is_mat:
if is_slice:
key = (slice(*divmod(i, self.cols)), slice(*divmod(j, self.cols)))
else:
key = (slice(i, i + value.rows), slice(j, j + value.cols))
self.copyin_matrix(key, value)
else:
return i, j, self._sympify(value)
return
[docs] def add(self, b):
"""Return self + b"""
return self + b
[docs] def cholesky_solve(self, rhs):
"""Solves Ax = B using Cholesky decomposition,
for a general square non-singular matrix.
For a non-square matrix with rows > cols,
the least squares solution is returned.
See Also
========
lower_triangular_solve
upper_triangular_solve
gauss_jordan_solve
diagonal_solve
LDLsolve
LUsolve
QRsolve
pinv_solve
"""
if self.is_symmetric():
L = self._cholesky()
elif self.rows >= self.cols:
L = (self.T * self)._cholesky()
rhs = self.T * rhs
else:
raise NotImplementedError(
"Under-determined System. " "Try M.gauss_jordan_solve(rhs)"
)
Y = L._lower_triangular_solve(rhs)
return (L.T)._upper_triangular_solve(Y)
[docs] def cholesky(self):
"""Returns the Cholesky decomposition L of a matrix A
such that L * L.T = A
A must be a square, symmetric, positive-definite
and non-singular matrix.
Examples
========
>>> from ..matrices import Matrix
>>> A = Matrix(((25, 15, -5), (15, 18, 0), (-5, 0, 11)))
>>> A.cholesky()
Matrix([
[ 5, 0, 0],
[ 3, 3, 0],
[-1, 1, 3]])
>>> A.cholesky() * A.cholesky().T
Matrix([
[25, 15, -5],
[15, 18, 0],
[-5, 0, 11]])
See Also
========
LDLdecomposition
LUdecomposition
QRdecomposition
"""
if not self.is_square:
raise NonSquareMatrixError("Matrix must be square.")
if not self.is_symmetric():
raise ValueError("Matrix must be symmetric.")
return self._cholesky()
[docs] def condition_number(self):
"""Returns the condition number of a matrix.
This is the maximum singular value divided by the minimum singular value
Examples
========
>>> from .. import Matrix, S
>>> A = Matrix([[1, 0, 0], [0, 10, 0], [0, 0, S.One/10]])
>>> A.condition_number()
100
See Also
========
singular_values
"""
if not self:
return S.Zero
singularvalues = self.singular_values()
return Max(*singularvalues) / Min(*singularvalues)
[docs] def copy(self):
"""
Returns the copy of a matrix.
Examples
========
>>> from .. import Matrix
>>> A = Matrix(2, 2, [1, 2, 3, 4])
>>> A.copy()
Matrix([
[1, 2],
[3, 4]])
"""
return self._new(self.rows, self.cols, self._mat)
[docs] def cross(self, b):
"""Return the cross product of `self` and `b` relaxing the condition
of compatible dimensions: if each has 3 elements, a matrix of the
same type and shape as `self` will be returned. If `b` has the same
shape as `self` then common identities for the cross product (like
`a x b = - b x a`) will hold.
See Also
========
dot
multiply
multiply_elementwise
"""
if not is_sequence(b):
raise TypeError(
"`b` must be an ordered iterable or Matrix, not %s." % type(b)
)
if not (self.rows * self.cols == b.rows * b.cols == 3):
raise ShapeError(
"Dimensions incorrect for cross product: %s x %s"
% ((self.rows, self.cols), (b.rows, b.cols))
)
else:
return self._new(
self.rows,
self.cols,
(
(self[1] * b[2] - self[2] * b[1]),
(self[2] * b[0] - self[0] * b[2]),
(self[0] * b[1] - self[1] * b[0]),
),
)
@property
def D(self):
"""Return Dirac conjugate (if self.rows == 4).
Examples
========
>>> from .. import Matrix, I, eye
>>> m = Matrix((0, 1 + I, 2, 3))
>>> m.D
Matrix([[0, 1 - I, -2, -3]])
>>> m = (eye(4) + I*eye(4))
>>> m[0, 3] = 2
>>> m.D
Matrix([
[1 - I, 0, 0, 0],
[ 0, 1 - I, 0, 0],
[ 0, 0, -1 + I, 0],
[ 2, 0, 0, -1 + I]])
If the matrix does not have 4 rows an AttributeError will be raised
because this property is only defined for matrices with 4 rows.
>>> Matrix(eye(2)).D
Traceback (most recent call last):
...
AttributeError: Matrix has no attribute D.
See Also
========
conjugate: By-element conjugation
H: Hermite conjugation
"""
from ..physics.matrices import mgamma
if self.rows != 4:
# In Python 3.2, properties can only return an AttributeError
# so we can't raise a ShapeError -- see commit which added the
# first line of this inline comment. Also, there is no need
# for a message since MatrixBase will raise the AttributeError
raise AttributeError
return self.H * mgamma(0)
[docs] def diagonal_solve(self, rhs):
"""Solves Ax = B efficiently, where A is a diagonal Matrix,
with non-zero diagonal entries.
Examples
========
>>> from ..matrices import Matrix, eye
>>> A = eye(2)*2
>>> B = Matrix([[1, 2], [3, 4]])
>>> A.diagonal_solve(B) == B/2
True
See Also
========
lower_triangular_solve
upper_triangular_solve
gauss_jordan_solve
cholesky_solve
LDLsolve
LUsolve
QRsolve
pinv_solve
"""
if not self.is_diagonal:
raise TypeError("Matrix should be diagonal")
if rhs.rows != self.rows:
raise TypeError("Size mis-match")
return self._diagonal_solve(rhs)
[docs] def dot(self, b):
"""Return the dot product of Matrix self and b relaxing the condition
of compatible dimensions: if either the number of rows or columns are
the same as the length of b then the dot product is returned. If self
is a row or column vector, a scalar is returned. Otherwise, a list
of results is returned (and in that case the number of columns in self
must match the length of b).
Examples
========
>>> from .. import Matrix
>>> M = Matrix([[1, 2, 3], [4, 5, 6], [7, 8, 9]])
>>> v = [1, 1, 1]
>>> M.row(0).dot(v)
6
>>> M.col(0).dot(v)
12
>>> M.dot(v)
[6, 15, 24]
See Also
========
cross
multiply
multiply_elementwise
"""
from .dense import Matrix
if not isinstance(b, MatrixBase):
if is_sequence(b):
if len(b) != self.cols and len(b) != self.rows:
raise ShapeError(
"Dimensions incorrect for dot product: %s, %s"
% (self.shape, len(b))
)
return self.dot(Matrix(b))
else:
raise TypeError(
"`b` must be an ordered iterable or Matrix, not %s." % type(b)
)
mat = self
if mat.cols == b.rows:
if b.cols != 1:
mat = mat.T
b = b.T
prod = flatten((mat * b).tolist())
if len(prod) == 1:
return prod[0]
return prod
if mat.cols == b.cols:
return mat.dot(b.T)
elif mat.rows == b.rows:
return mat.T.dot(b)
else:
raise ShapeError(
"Dimensions incorrect for dot product: %s, %s" % (self.shape, b.shape)
)
[docs] def dual(self):
"""Returns the dual of a matrix, which is:
`(1/2)*levicivita(i, j, k, l)*M(k, l)` summed over indices `k` and `l`
Since the levicivita method is anti_symmetric for any pairwise
exchange of indices, the dual of a symmetric matrix is the zero
matrix. Strictly speaking the dual defined here assumes that the
'matrix' `M` is a contravariant anti_symmetric second rank tensor,
so that the dual is a covariant second rank tensor.
"""
from .. import LeviCivita
from ..matrices import zeros
M, n = self[:, :], self.rows
work = zeros(n)
if self.is_symmetric():
return work
for i in range(1, n):
for j in range(1, n):
acum = 0
for k in range(1, n):
acum += LeviCivita(i, j, 0, k) * M[0, k]
work[i, j] = acum
work[j, i] = -acum
for l in range(1, n):
acum = 0
for a in range(1, n):
for b in range(1, n):
acum += LeviCivita(0, l, a, b) * M[a, b]
acum /= 2
work[0, l] = -acum
work[l, 0] = acum
return work
[docs] def exp(self):
"""Return the exponentiation of a square matrix."""
if not self.is_square:
raise NonSquareMatrixError(
"Exponentiation is valid only for square matrices"
)
try:
P, J = self.jordan_form()
cells = J.get_diag_blocks()
except MatrixError:
raise NotImplementedError(
"Exponentiation is implemented only for matrices for which the Jordan normal form can be computed"
)
def _jblock_exponential(b):
# This function computes the matrix exponential for one single Jordan block
nr = b.rows
l = b[0, 0]
if nr == 1:
res = exp(l)
else:
from .. import eye
# extract the diagonal part
d = b[0, 0] * eye(nr)
# and the nilpotent part
n = b - d
# compute its exponential
nex = eye(nr)
for i in range(1, nr):
nex = nex + n ** i / factorial(i)
# combine the two parts
res = exp(b[0, 0]) * nex
return res
blocks = list(map(_jblock_exponential, cells))
from ..matrices import diag
eJ = diag(*blocks)
# n = self.rows
ret = P * eJ * P.inv()
return type(self)(ret)
[docs] def gauss_jordan_solve(self, b, freevar=False):
"""
Solves Ax = b using Gauss Jordan elimination.
There may be zero, one, or infinite solutions. If one solution
exists, it will be returned. If infinite solutions exist, it will
be returned parametrically. If no solutions exist, It will throw
ValueError.
Parameters
==========
b : Matrix
The right hand side of the equation to be solved for. Must have
the same number of rows as matrix A.
freevar : List
If the system is underdetermined (e.g. A has more columns than
rows), infinite solutions are possible, in terms of an arbitrary
values of free variables. Then the index of the free variables
in the solutions (column Matrix) will be returned by freevar, if
the flag `freevar` is set to `True`.
Returns
=======
x : Matrix
The matrix that will satisfy Ax = B. Will have as many rows as
matrix A has columns, and as many columns as matrix B.
params : Matrix
If the system is underdetermined (e.g. A has more columns than
rows), infinite solutions are possible, in terms of an arbitrary
parameters. These arbitrary parameters are returned as params
Matrix.
Examples
========
>>> from .. import Matrix
>>> A = Matrix([[1, 2, 1, 1], [1, 2, 2, -1], [2, 4, 0, 6]])
>>> b = Matrix([7, 12, 4])
>>> sol, params = A.gauss_jordan_solve(b)
>>> sol
Matrix([
[-2*_tau0 - 3*_tau1 + 2],
[ _tau0],
[ 2*_tau1 + 5],
[ _tau1]])
>>> params
Matrix([
[_tau0],
[_tau1]])
>>> A = Matrix([[1, 2, 3], [4, 5, 6], [7, 8, 10]])
>>> b = Matrix([3, 6, 9])
>>> sol, params = A.gauss_jordan_solve(b)
>>> sol
Matrix([
[-1],
[ 2],
[ 0]])
>>> params
Matrix(0, 1, [])
See Also
========
lower_triangular_solve
upper_triangular_solve
cholesky_solve
diagonal_solve
LDLsolve
LUsolve
QRsolve
pinv
References
==========
.. [1] http://en.wikipedia.org/wiki/Gaussian_elimination
"""
from ..matrices import Matrix, zeros
aug = self.hstack(self.copy(), b.copy())
row, col = aug[:, :-1].shape
# solve by reduced row echelon form
A, pivots = aug.rref(simplify=True)
A, v = A[:, :-1], A[:, -1]
pivots = list(filter(lambda p: p < col, pivots))
rank = len(pivots)
# Bring to block form
permutation = Matrix(range(col)).T
A = A.vstack(A, permutation)
for i, c in enumerate(pivots):
A.col_swap(i, c)
A, permutation = A[:-1, :], A[-1, :]
# check for existence of solutions
# rank of aug Matrix should be equal to rank of coefficient matrix
if not v[rank:, 0].is_zero:
raise ValueError("Linear system has no solution")
# Get index of free symbols (free parameters)
free_var_index = permutation[
len(pivots) :
] # non-pivots columns are free variables
# Free parameters
dummygen = numbered_symbols("tau", Dummy)
tau = Matrix([next(dummygen) for k in range(col - rank)]).reshape(col - rank, 1)
# Full parametric solution
V = A[:rank, rank:]
vt = v[:rank, 0]
free_sol = tau.vstack(vt - V * tau, tau)
# Undo permutation
sol = zeros(col, 1)
for k, v in enumerate(free_sol):
sol[permutation[k], 0] = v
if freevar:
return sol, tau, free_var_index
else:
return sol, tau
[docs] def inv_mod(self, m):
r"""
Returns the inverse of the matrix `K` (mod `m`), if it exists.
Method to find the matrix inverse of `K` (mod `m`) implemented in this function:
* Compute `\mathrm{adj}(K) = \mathrm{cof}(K)^t`, the adjoint matrix of `K`.
* Compute `r = 1/\mathrm{det}(K) \pmod m`.
* `K^{-1} = r\cdot \mathrm{adj}(K) \pmod m`.
Examples
========
>>> from .. import Matrix
>>> A = Matrix(2, 2, [1, 2, 3, 4])
>>> A.inv_mod(5)
Matrix([
[3, 1],
[4, 2]])
>>> A.inv_mod(3)
Matrix([
[1, 1],
[0, 1]])
"""
from ..ntheory import totient
if not self.is_square:
raise NonSquareMatrixError()
N = self.cols
phi = totient(m)
det_K = self.det()
if gcd(det_K, m) != 1:
raise ValueError("Matrix is not invertible (mod %d)" % m)
det_inv = pow(int(det_K), int(phi - 1), int(m))
K_adj = self.adjugate()
K_inv = self.__class__(
N, N, [det_inv * K_adj[i, j] % m for i in range(N) for j in range(N)]
)
return K_inv
[docs] def inverse_ADJ(self, iszerofunc=_iszero):
"""Calculates the inverse using the adjugate matrix and a determinant.
See Also
========
inv
inverse_LU
inverse_GE
"""
if not self.is_square:
raise NonSquareMatrixError("A Matrix must be square to invert.")
d = self.det(method="berkowitz")
zero = d.equals(0)
if zero is None:
# if equals() can't decide, will rref be able to?
ok = self.rref(simplify=True)[0]
zero = any(iszerofunc(ok[j, j]) for j in range(ok.rows))
if zero:
raise ValueError("Matrix det == 0; not invertible.")
return self.adjugate() / d
[docs] def inverse_GE(self, iszerofunc=_iszero):
"""Calculates the inverse using Gaussian elimination.
See Also
========
inv
inverse_LU
inverse_ADJ
"""
from .dense import Matrix
if not self.is_square:
raise NonSquareMatrixError("A Matrix must be square to invert.")
big = Matrix.hstack(self.as_mutable(), Matrix.eye(self.rows))
red = big.rref(iszerofunc=iszerofunc, simplify=True)[0]
if any(iszerofunc(red[j, j]) for j in range(red.rows)):
raise ValueError("Matrix det == 0; not invertible.")
return self._new(red[:, big.rows :])
[docs] def inverse_LU(self, iszerofunc=_iszero):
"""Calculates the inverse using LU decomposition.
See Also
========
inv
inverse_GE
inverse_ADJ
"""
if not self.is_square:
raise NonSquareMatrixError()
ok = self.rref(simplify=True)[0]
if any(iszerofunc(ok[j, j]) for j in range(ok.rows)):
raise ValueError("Matrix det == 0; not invertible.")
return self.LUsolve(self.eye(self.rows), iszerofunc=_iszero)
[docs] def inv(self, method=None, **kwargs):
"""
Return the inverse of a matrix.
CASE 1: If the matrix is a dense matrix.
Return the matrix inverse using the method indicated (default
is Gauss elimination).
Parameters
==========
method : ('GE', 'LU', or 'ADJ')
Notes
=====
According to the ``method`` keyword, it calls the appropriate method:
GE .... inverse_GE(); default
LU .... inverse_LU()
ADJ ... inverse_ADJ()
See Also
========
inverse_LU
inverse_GE
inverse_ADJ
Raises
------
ValueError
If the determinant of the matrix is zero.
CASE 2: If the matrix is a sparse matrix.
Return the matrix inverse using Cholesky or LDL (default).
kwargs
======
method : ('CH', 'LDL')
Notes
=====
According to the ``method`` keyword, it calls the appropriate method:
LDL ... inverse_LDL(); default
CH .... inverse_CH()
Raises
------
ValueError
If the determinant of the matrix is zero.
"""
if not self.is_square:
raise NonSquareMatrixError()
if method is not None:
kwargs["method"] = method
return self._eval_inverse(**kwargs)
[docs] def is_nilpotent(self):
"""Checks if a matrix is nilpotent.
A matrix B is nilpotent if for some integer k, B**k is
a zero matrix.
Examples
========
>>> from .. import Matrix
>>> a = Matrix([[0, 0, 0], [1, 0, 0], [1, 1, 0]])
>>> a.is_nilpotent()
True
>>> a = Matrix([[1, 0, 1], [1, 0, 0], [1, 1, 0]])
>>> a.is_nilpotent()
False
"""
if not self:
return True
if not self.is_square:
raise NonSquareMatrixError("Nilpotency is valid only for square matrices")
x = Dummy("x")
if self.charpoly(x).args[0] == x ** self.rows:
return True
return False
[docs] def key2bounds(self, keys):
"""Converts a key with potentially mixed types of keys (integer and slice)
into a tuple of ranges and raises an error if any index is out of self's
range.
See Also
========
key2ij
"""
islice, jslice = [isinstance(k, slice) for k in keys]
if islice:
if not self.rows:
rlo = rhi = 0
else:
rlo, rhi = keys[0].indices(self.rows)[:2]
else:
rlo = a2idx(keys[0], self.rows)
rhi = rlo + 1
if jslice:
if not self.cols:
clo = chi = 0
else:
clo, chi = keys[1].indices(self.cols)[:2]
else:
clo = a2idx(keys[1], self.cols)
chi = clo + 1
return rlo, rhi, clo, chi
[docs] def key2ij(self, key):
"""Converts key into canonical form, converting integers or indexable
items into valid integers for self's range or returning slices
unchanged.
See Also
========
key2bounds
"""
if is_sequence(key):
if not len(key) == 2:
raise TypeError("key must be a sequence of length 2")
return [
a2idx(i, n) if not isinstance(i, slice) else i
for i, n in zip(key, self.shape)
]
elif isinstance(key, slice):
return key.indices(len(self))[:2]
else:
return divmod(a2idx(key, len(self)), self.cols)
[docs] def LDLdecomposition(self):
"""Returns the LDL Decomposition (L, D) of matrix A,
such that L * D * L.T == A
This method eliminates the use of square root.
Further this ensures that all the diagonal entries of L are 1.
A must be a square, symmetric, positive-definite
and non-singular matrix.
Examples
========
>>> from ..matrices import Matrix, eye
>>> A = Matrix(((25, 15, -5), (15, 18, 0), (-5, 0, 11)))
>>> L, D = A.LDLdecomposition()
>>> L
Matrix([
[ 1, 0, 0],
[ 3/5, 1, 0],
[-1/5, 1/3, 1]])
>>> D
Matrix([
[25, 0, 0],
[ 0, 9, 0],
[ 0, 0, 9]])
>>> L * D * L.T * A.inv() == eye(A.rows)
True
See Also
========
cholesky
LUdecomposition
QRdecomposition
"""
if not self.is_square:
raise NonSquareMatrixError("Matrix must be square.")
if not self.is_symmetric():
raise ValueError("Matrix must be symmetric.")
return self._LDLdecomposition()
[docs] def LDLsolve(self, rhs):
"""Solves Ax = B using LDL decomposition,
for a general square and non-singular matrix.
For a non-square matrix with rows > cols,
the least squares solution is returned.
Examples
========
>>> from ..matrices import Matrix, eye
>>> A = eye(2)*2
>>> B = Matrix([[1, 2], [3, 4]])
>>> A.LDLsolve(B) == B/2
True
See Also
========
LDLdecomposition
lower_triangular_solve
upper_triangular_solve
gauss_jordan_solve
cholesky_solve
diagonal_solve
LUsolve
QRsolve
pinv_solve
"""
if self.is_symmetric():
L, D = self.LDLdecomposition()
elif self.rows >= self.cols:
L, D = (self.T * self).LDLdecomposition()
rhs = self.T * rhs
else:
raise NotImplementedError(
"Under-determined System. " "Try M.gauss_jordan_solve(rhs)"
)
Y = L._lower_triangular_solve(rhs)
Z = D._diagonal_solve(Y)
return (L.T)._upper_triangular_solve(Z)
[docs] def lower_triangular_solve(self, rhs):
"""Solves Ax = B, where A is a lower triangular matrix.
See Also
========
upper_triangular_solve
gauss_jordan_solve
cholesky_solve
diagonal_solve
LDLsolve
LUsolve
QRsolve
pinv_solve
"""
if not self.is_square:
raise NonSquareMatrixError("Matrix must be square.")
if rhs.rows != self.rows:
raise ShapeError("Matrices size mismatch.")
if not self.is_lower:
raise ValueError("Matrix must be lower triangular.")
return self._lower_triangular_solve(rhs)
[docs] def LUdecomposition(self, iszerofunc=_iszero, simpfunc=None, rankcheck=False):
"""Returns (L, U, perm) where L is a lower triangular matrix with unit
diagonal, U is an upper triangular matrix, and perm is a list of row
swap index pairs. If A is the original matrix, then
A = (L*U).permuteBkwd(perm), and the row permutation matrix P such
that P*A = L*U can be computed by P=eye(A.row).permuteFwd(perm).
See documentation for LUCombined for details about the keyword argument
rankcheck, iszerofunc, and simpfunc.
Examples
========
>>> from .. import Matrix
>>> a = Matrix([[4, 3], [6, 3]])
>>> L, U, _ = a.LUdecomposition()
>>> L
Matrix([
[ 1, 0],
[3/2, 1]])
>>> U
Matrix([
[4, 3],
[0, -3/2]])
See Also
========
cholesky
LDLdecomposition
QRdecomposition
LUdecomposition_Simple
LUdecompositionFF
LUsolve
"""
combined, p = self.LUdecomposition_Simple(
iszerofunc=iszerofunc, simpfunc=simpfunc, rankcheck=rankcheck
)
# L is lower triangular self.rows x self.rows
# U is upper triangular self.rows x self.cols
# L has unit diagonal. For each column in combined, the subcolumn
# below the diagonal of combined is shared by L.
# If L has more columns than combined, then the remaining subcolumns
# below the diagonal of L are zero.
# The upper triangular portion of L and combined are equal.
def entry_L(i, j):
if i < j:
# Super diagonal entry
return S.Zero
elif i == j:
return S.One
elif j < combined.cols:
return combined[i, j]
# Subdiagonal entry of L with no corresponding
# entry in combined
return S.Zero
def entry_U(i, j):
return S.Zero if i > j else combined[i, j]
L = self._new(combined.rows, combined.rows, entry_L)
U = self._new(combined.rows, combined.cols, entry_U)
return L, U, p
[docs] def LUdecomposition_Simple(
self, iszerofunc=_iszero, simpfunc=None, rankcheck=False
):
"""Compute an lu decomposition of m x n matrix A, where P*A = L*U
* L is m x m lower triangular with unit diagonal
* U is m x n upper triangular
* P is an m x m permutation matrix
Returns an m x n matrix lu, and an m element list perm where each
element of perm is a pair of row exchange indices.
The factors L and U are stored in lu as follows:
The subdiagonal elements of L are stored in the subdiagonal elements
of lu, that is lu[i, j] = L[i, j] whenever i > j.
The elements on the diagonal of L are all 1, and are not explicitly
stored.
U is stored in the upper triangular portion of lu, that is
lu[i ,j] = U[i, j] whenever i <= j.
The output matrix can be visualized as:
Matrix([
[u, u, u, u],
[l, u, u, u],
[l, l, u, u],
[l, l, l, u]])
where l represents a subdiagonal entry of the L factor, and u
represents an entry from the upper triangular entry of the U
factor.
perm is a list row swap index pairs such that if A is the original
matrix, then A = (L*U).permuteBkwd(perm), and the row permutation
matrix P such that P*A = L*U can be computed by
soP=eye(A.row).permuteFwd(perm).
The keyword argument rankcheck determines if this function raises a
ValueError when passed a matrix whose rank is strictly less than
min(num rows, num cols). The default behavior is to decompose a rank
deficient matrix. Pass rankcheck=True to raise a
ValueError instead. (This mimics the previous behavior of this function).
The keyword arguments iszerofunc and simpfunc are used by the pivot
search algorithm.
iszerofunc is a callable that returns a boolean indicating if its
input is zero, or None if it cannot make the determination.
simpfunc is a callable that simplifies its input.
The default is simpfunc=None, which indicate that the pivot search
algorithm should not attempt to simplify any candidate pivots.
If simpfunc fails to simplify its input, then it must return its input
instead of a copy.
When a matrix contains symbolic entries, the pivot search algorithm
differs from the case where every entry can be categorized as zero or
nonzero.
The algorithm searches column by column through the submatrix whose
top left entry coincides with the pivot position.
If it exists, the pivot is the first entry in the current search
column that iszerofunc guarantees is nonzero.
If no such candidate exists, then each candidate pivot is simplified
if simpfunc is not None.
The search is repeated, with the difference that a candidate may be
the pivot if `iszerofunc()` cannot guarantee that it is nonzero.
In the second search the pivot is the first candidate that
iszerofunc can guarantee is nonzero.
If no such candidate exists, then the pivot is the first candidate
for which iszerofunc returns None.
If no such candidate exists, then the search is repeated in the next
column to the right.
The pivot search algorithm differs from the one in `rref()`, which
relies on `_find_reasonable_pivot()`.
Future versions of `LUdecomposition_simple()` may use
`_find_reasonable_pivot()`.
See Also
========
LUdecomposition
LUdecompositionFF
LUsolve
"""
if rankcheck:
# https://github.com/sympy/sympy/issues/9796
pass
if self.rows == 0 or self.cols == 0:
# Define LU decomposition of a matrix with no entries as a matrix
# of the same dimensions with all zero entries.
return self.zeros(self.rows, self.cols), []
lu = self.as_mutable()
row_swaps = []
pivot_col = 0
for pivot_row in range(0, lu.rows - 1):
# Search for pivot. Prefer entry that iszeropivot determines
# is nonzero, over entry that iszeropivot cannot guarantee
# is zero.
# XXX `_find_reasonable_pivot` uses slow zero testing. Blocked by bug #10279
# Future versions of LUdecomposition_simple can pass iszerofunc and simpfunc
# to _find_reasonable_pivot().
# In pass 3 of _find_reasonable_pivot(), the predicate in `if x.equals(S.Zero):`
# calls sympy.simplify(), and not the simplification function passed in via
# the keyword argument simpfunc.
iszeropivot = True
while pivot_col != self.cols and iszeropivot:
sub_col = (lu[r, pivot_col] for r in range(pivot_row, self.rows))
(
pivot_row_offset,
pivot_value,
is_assumed_non_zero,
ind_simplified_pairs,
) = _find_reasonable_pivot_naive(sub_col, iszerofunc, simpfunc)
iszeropivot = pivot_value is None
if iszeropivot:
# All candidate pivots in this column are zero.
# Proceed to next column.
pivot_col += 1
if rankcheck and pivot_col != pivot_row:
# All entries including and below the pivot position are
# zero, which indicates that the rank of the matrix is
# strictly less than min(num rows, num cols)
# Mimic behavior of previous implementation, by throwing a
# ValueError.
raise ValueError(
"Rank of matrix is strictly less than"
" number of rows or columns."
" Pass keyword argument"
" rankcheck=False to compute"
" the LU decomposition of this matrix."
)
candidate_pivot_row = (
None if pivot_row_offset is None else pivot_row + pivot_row_offset
)
if candidate_pivot_row is None and iszeropivot:
# If candidate_pivot_row is None and iszeropivot is True
# after pivot search has completed, then the submatrix
# below and to the right of (pivot_row, pivot_col) is
# all zeros, indicating that Gaussian elimination is
# complete.
return lu, row_swaps
# Update entries simplified during pivot search.
for offset, val in ind_simplified_pairs:
lu[pivot_row + offset, pivot_col] = val
if pivot_row != candidate_pivot_row:
# Row swap book keeping:
# Record which rows were swapped.
# Update stored portion of L factor by multiplying L on the
# left and right with the current permutation.
# Swap rows of U.
row_swaps.append([pivot_row, candidate_pivot_row])
# Update L.
lu[pivot_row, 0:pivot_row], lu[candidate_pivot_row, 0:pivot_row] = (
lu[candidate_pivot_row, 0:pivot_row],
lu[pivot_row, 0:pivot_row],
)
# Swap pivot row of U with candidate pivot row.
(
lu[pivot_row, pivot_col : lu.cols],
lu[candidate_pivot_row, pivot_col : lu.cols],
) = (
lu[candidate_pivot_row, pivot_col : lu.cols],
lu[pivot_row, pivot_col : lu.cols],
)
# Introduce zeros below the pivot by adding a multiple of the
# pivot row to a row under it, and store the result in the
# row under it.
# Only entries in the target row whose index is greater than
# start_col may be nonzero.
start_col = pivot_col + 1
for row in range(pivot_row + 1, lu.rows):
# Store factors of L in the subcolumn below
# (pivot_row, pivot_row).
lu[row, pivot_row] = lu[row, pivot_col] / lu[pivot_row, pivot_col]
# Form the linear combination of the pivot row and the current
# row below the pivot row that zeros the entries below the pivot.
# Employing slicing instead of a loop here raises
# NotImplementedError: Cannot add Zero to MutableSparseMatrix
# in sympy/matrices/tests/test_sparse.py.
# c = pivot_row + 1 if pivot_row == pivot_col else pivot_col
for c in range(start_col, lu.cols):
lu[row, c] = lu[row, c] - lu[row, pivot_row] * lu[pivot_row, c]
if pivot_row != pivot_col:
# matrix rank < min(num rows, num cols),
# so factors of L are not stored directly below the pivot.
# These entries are zero by construction, so don't bother
# computing them.
for row in range(pivot_row + 1, lu.rows):
lu[row, pivot_col] = S.Zero
pivot_col += 1
if pivot_col == lu.cols:
# All candidate pivots are zero implies that Gaussian
# elimination is complete.
return lu, row_swaps
return lu, row_swaps
[docs] def LUdecompositionFF(self):
"""Compute a fraction-free LU decomposition.
Returns 4 matrices P, L, D, U such that PA = L D**-1 U.
If the elements of the matrix belong to some integral domain I, then all
elements of L, D and U are guaranteed to belong to I.
**Reference**
- W. Zhou & D.J. Jeffrey, "Fraction-free matrix factors: new forms
for LU and QR factors". Frontiers in Computer Science in China,
Vol 2, no. 1, pp. 67-80, 2008.
See Also
========
LUdecomposition
LUdecomposition_Simple
LUsolve
"""
from ..matrices import SparseMatrix
zeros = SparseMatrix.zeros
eye = SparseMatrix.eye
n, m = self.rows, self.cols
U, L, P = self.as_mutable(), eye(n), eye(n)
DD = zeros(n, n)
oldpivot = 1
for k in range(n - 1):
if U[k, k] == 0:
for kpivot in range(k + 1, n):
if U[kpivot, k]:
break
else:
raise ValueError("Matrix is not full rank")
U[k, k:], U[kpivot, k:] = U[kpivot, k:], U[k, k:]
L[k, :k], L[kpivot, :k] = L[kpivot, :k], L[k, :k]
P[k, :], P[kpivot, :] = P[kpivot, :], P[k, :]
L[k, k] = Ukk = U[k, k]
DD[k, k] = oldpivot * Ukk
for i in range(k + 1, n):
L[i, k] = Uik = U[i, k]
for j in range(k + 1, m):
U[i, j] = (Ukk * U[i, j] - U[k, j] * Uik) / oldpivot
U[i, k] = 0
oldpivot = Ukk
DD[n - 1, n - 1] = oldpivot
return P, L, DD, U
[docs] def LUsolve(self, rhs, iszerofunc=_iszero):
"""Solve the linear system Ax = rhs for x where A = self.
This is for symbolic matrices, for real or complex ones use
mpmath.lu_solve or mpmath.qr_solve.
See Also
========
lower_triangular_solve
upper_triangular_solve
gauss_jordan_solve
cholesky_solve
diagonal_solve
LDLsolve
QRsolve
pinv_solve
LUdecomposition
"""
if rhs.rows != self.rows:
raise ShapeError("`self` and `rhs` must have the same number of rows.")
A, perm = self.LUdecomposition_Simple(iszerofunc=_iszero)
n = self.rows
b = rhs.permute_rows(perm).as_mutable()
# forward substitution, all diag entries are scaled to 1
for i in range(n):
for j in range(i):
scale = A[i, j]
b.zip_row_op(i, j, lambda x, y: x - y * scale)
# backward substitution
for i in range(n - 1, -1, -1):
for j in range(i + 1, n):
scale = A[i, j]
b.zip_row_op(i, j, lambda x, y: x - y * scale)
scale = A[i, i]
b.row_op(i, lambda x, _: x / scale)
return rhs.__class__(b)
[docs] def multiply(self, b):
"""Returns self*b
See Also
========
dot
cross
multiply_elementwise
"""
return self * b
[docs] def normalized(self):
"""Return the normalized version of ``self``.
See Also
========
norm
"""
if self.rows != 1 and self.cols != 1:
raise ShapeError("A Matrix must be a vector to normalize.")
norm = self.norm()
out = self.applyfunc(lambda i: i / norm)
return out
[docs] def norm(self, ord=None):
"""Return the Norm of a Matrix or Vector.
In the simplest case this is the geometric size of the vector
Other norms can be specified by the ord parameter
===== ============================ ==========================
ord norm for matrices norm for vectors
===== ============================ ==========================
None Frobenius norm 2-norm
'fro' Frobenius norm - does not exist
inf -- max(abs(x))
-inf -- min(abs(x))
1 -- as below
-1 -- as below
2 2-norm (largest sing. value) as below
-2 smallest singular value as below
other - does not exist sum(abs(x)**ord)**(1./ord)
===== ============================ ==========================
Examples
========
>>> from .. import Matrix, Symbol, trigsimp, cos, sin, oo
>>> x = Symbol('x', real=True)
>>> v = Matrix([cos(x), sin(x)])
>>> trigsimp( v.norm() )
1
>>> v.norm(10)
(sin(x)**10 + cos(x)**10)**(1/10)
>>> A = Matrix([[1, 1], [1, 1]])
>>> A.norm(2)# Spectral norm (max of |Ax|/|x| under 2-vector-norm)
2
>>> A.norm(-2) # Inverse spectral norm (smallest singular value)
0
>>> A.norm() # Frobenius Norm
2
>>> Matrix([1, -2]).norm(oo)
2
>>> Matrix([-1, 2]).norm(-oo)
1
See Also
========
normalized
"""
# Row or Column Vector Norms
vals = list(self.values()) or [0]
if self.rows == 1 or self.cols == 1:
if ord == 2 or ord is None: # Common case sqrt(<x, x>)
return sqrt(Add(*(abs(i) ** 2 for i in vals)))
elif ord == 1: # sum(abs(x))
return Add(*(abs(i) for i in vals))
elif ord == S.Infinity: # max(abs(x))
return Max(*[abs(i) for i in vals])
elif ord == S.NegativeInfinity: # min(abs(x))
return Min(*[abs(i) for i in vals])
# Otherwise generalize the 2-norm, Sum(x_i**ord)**(1/ord)
# Note that while useful this is not mathematically a norm
try:
return Pow(Add(*(abs(i) ** ord for i in vals)), S(1) / ord)
except (NotImplementedError, TypeError):
raise ValueError("Expected order to be Number, Symbol, oo")
# Matrix Norms
else:
if ord == 2: # Spectral Norm
# Maximum singular value
return Max(*self.singular_values())
elif ord == -2:
# Minimum singular value
return Min(*self.singular_values())
elif (
ord is None
or isinstance(ord, string_types)
and ord.lower() in ["f", "fro", "frobenius", "vector"]
):
# Reshape as vector and send back to norm function
return self.vec().norm(ord=2)
else:
raise NotImplementedError("Matrix Norms under development")
[docs] def pinv_solve(self, B, arbitrary_matrix=None):
"""Solve Ax = B using the Moore-Penrose pseudoinverse.
There may be zero, one, or infinite solutions. If one solution
exists, it will be returned. If infinite solutions exist, one will
be returned based on the value of arbitrary_matrix. If no solutions
exist, the least-squares solution is returned.
Parameters
==========
B : Matrix
The right hand side of the equation to be solved for. Must have
the same number of rows as matrix A.
arbitrary_matrix : Matrix
If the system is underdetermined (e.g. A has more columns than
rows), infinite solutions are possible, in terms of an arbitrary
matrix. This parameter may be set to a specific matrix to use
for that purpose; if so, it must be the same shape as x, with as
many rows as matrix A has columns, and as many columns as matrix
B. If left as None, an appropriate matrix containing dummy
symbols in the form of ``wn_m`` will be used, with n and m being
row and column position of each symbol.
Returns
=======
x : Matrix
The matrix that will satisfy Ax = B. Will have as many rows as
matrix A has columns, and as many columns as matrix B.
Examples
========
>>> from .. import Matrix
>>> A = Matrix([[1, 2, 3], [4, 5, 6]])
>>> B = Matrix([7, 8])
>>> A.pinv_solve(B)
Matrix([
[ _w0_0/6 - _w1_0/3 + _w2_0/6 - 55/18],
[-_w0_0/3 + 2*_w1_0/3 - _w2_0/3 + 1/9],
[ _w0_0/6 - _w1_0/3 + _w2_0/6 + 59/18]])
>>> A.pinv_solve(B, arbitrary_matrix=Matrix([0, 0, 0]))
Matrix([
[-55/18],
[ 1/9],
[ 59/18]])
See Also
========
lower_triangular_solve
upper_triangular_solve
gauss_jordan_solve
cholesky_solve
diagonal_solve
LDLsolve
LUsolve
QRsolve
pinv
Notes
=====
This may return either exact solutions or least squares solutions.
To determine which, check ``A * A.pinv() * B == B``. It will be
True if exact solutions exist, and False if only a least-squares
solution exists. Be aware that the left hand side of that equation
may need to be simplified to correctly compare to the right hand
side.
References
==========
.. [1] https://en.wikipedia.org/wiki/Moore-Penrose_pseudoinverse#Obtaining_all_solutions_of_a_linear_system
"""
from ..matrices import eye
A = self
A_pinv = self.pinv()
if arbitrary_matrix is None:
rows, cols = A.cols, B.cols
w = symbols("w:{0}_:{1}".format(rows, cols), cls=Dummy)
arbitrary_matrix = self.__class__(cols, rows, w).T
return A_pinv * B + (eye(A.cols) - A_pinv * A) * arbitrary_matrix
[docs] def pinv(self):
"""Calculate the Moore-Penrose pseudoinverse of the matrix.
The Moore-Penrose pseudoinverse exists and is unique for any matrix.
If the matrix is invertible, the pseudoinverse is the same as the
inverse.
Examples
========
>>> from .. import Matrix
>>> Matrix([[1, 2, 3], [4, 5, 6]]).pinv()
Matrix([
[-17/18, 4/9],
[ -1/9, 1/9],
[ 13/18, -2/9]])
See Also
========
inv
pinv_solve
References
==========
.. [1] https://en.wikipedia.org/wiki/Moore-Penrose_pseudoinverse
"""
A = self
AH = self.H
# Trivial case: pseudoinverse of all-zero matrix is its transpose.
if A.is_zero:
return AH
try:
if self.rows >= self.cols:
return (AH * A).inv() * AH
else:
return AH * (A * AH).inv()
except ValueError:
# Matrix is not full rank, so A*AH cannot be inverted.
raise NotImplementedError(
"Rank-deficient matrices are not yet " "supported."
)
[docs] def print_nonzero(self, symb="X"):
"""Shows location of non-zero entries for fast shape lookup.
Examples
========
>>> from ..matrices import Matrix, eye
>>> m = Matrix(2, 3, lambda i, j: i*3+j)
>>> m
Matrix([
[0, 1, 2],
[3, 4, 5]])
>>> m.print_nonzero()
[ XX]
[XXX]
>>> m = eye(4)
>>> m.print_nonzero("x")
[x ]
[ x ]
[ x ]
[ x]
"""
s = []
for i in range(self.rows):
line = []
for j in range(self.cols):
if self[i, j] == 0:
line.append(" ")
else:
line.append(str(symb))
s.append("[%s]" % "".join(line))
print("\n".join(s))
[docs] def project(self, v):
"""Return the projection of ``self`` onto the line containing ``v``.
Examples
========
>>> from .. import Matrix, S, sqrt
>>> V = Matrix([sqrt(3)/2, S.Half])
>>> x = Matrix([[1, 0]])
>>> V.project(x)
Matrix([[sqrt(3)/2, 0]])
>>> V.project(-x)
Matrix([[sqrt(3)/2, 0]])
"""
return v * (self.dot(v) / v.dot(v))
[docs] def QRdecomposition(self):
"""Return Q, R where A = Q*R, Q is orthogonal and R is upper triangular.
Examples
========
This is the example from wikipedia:
>>> from .. import Matrix
>>> A = Matrix([[12, -51, 4], [6, 167, -68], [-4, 24, -41]])
>>> Q, R = A.QRdecomposition()
>>> Q
Matrix([
[ 6/7, -69/175, -58/175],
[ 3/7, 158/175, 6/175],
[-2/7, 6/35, -33/35]])
>>> R
Matrix([
[14, 21, -14],
[ 0, 175, -70],
[ 0, 0, 35]])
>>> A == Q*R
True
QR factorization of an identity matrix:
>>> A = Matrix([[1, 0, 0], [0, 1, 0], [0, 0, 1]])
>>> Q, R = A.QRdecomposition()
>>> Q
Matrix([
[1, 0, 0],
[0, 1, 0],
[0, 0, 1]])
>>> R
Matrix([
[1, 0, 0],
[0, 1, 0],
[0, 0, 1]])
See Also
========
cholesky
LDLdecomposition
LUdecomposition
QRsolve
"""
cls = self.__class__
mat = self.as_mutable()
if not mat.rows >= mat.cols:
raise MatrixError("The number of rows must be greater than columns")
n = mat.rows
m = mat.cols
rank = n
row_reduced = mat.rref()[0]
for i in range(row_reduced.rows):
if row_reduced.row(i).norm() == 0:
rank -= 1
if not rank == mat.cols:
raise MatrixError("The rank of the matrix must match the columns")
Q, R = mat.zeros(n, m), mat.zeros(m)
for j in range(m): # for each column vector
tmp = mat[:, j] # take original v
for i in range(j):
# subtract the project of mat on new vector
tmp -= Q[:, i] * mat[:, j].dot(Q[:, i])
tmp.expand()
# normalize it
R[j, j] = tmp.norm()
Q[:, j] = tmp / R[j, j]
if Q[:, j].norm() != 1:
raise NotImplementedError("Could not normalize the vector %d." % j)
for i in range(j):
R[i, j] = Q[:, i].dot(mat[:, j])
return cls(Q), cls(R)
[docs] def QRsolve(self, b):
"""Solve the linear system 'Ax = b'.
'self' is the matrix 'A', the method argument is the vector
'b'. The method returns the solution vector 'x'. If 'b' is a
matrix, the system is solved for each column of 'b' and the
return value is a matrix of the same shape as 'b'.
This method is slower (approximately by a factor of 2) but
more stable for floating-point arithmetic than the LUsolve method.
However, LUsolve usually uses an exact arithmetic, so you don't need
to use QRsolve.
This is mainly for educational purposes and symbolic matrices, for real
(or complex) matrices use mpmath.qr_solve.
See Also
========
lower_triangular_solve
upper_triangular_solve
gauss_jordan_solve
cholesky_solve
diagonal_solve
LDLsolve
LUsolve
pinv_solve
QRdecomposition
"""
Q, R = self.as_mutable().QRdecomposition()
y = Q.T * b
# back substitution to solve R*x = y:
# We build up the result "backwards" in the vector 'x' and reverse it
# only in the end.
x = []
n = R.rows
for j in range(n - 1, -1, -1):
tmp = y[j, :]
for k in range(j + 1, n):
tmp -= R[j, k] * x[n - 1 - k]
x.append(tmp / R[j, j])
return self._new([row._mat for row in reversed(x)])
[docs] def solve_least_squares(self, rhs, method="CH"):
"""Return the least-square fit to the data.
By default the cholesky_solve routine is used (method='CH'); other
methods of matrix inversion can be used. To find out which are
available, see the docstring of the .inv() method.
Examples
========
>>> from ..matrices import Matrix, ones
>>> A = Matrix([1, 2, 3])
>>> B = Matrix([2, 3, 4])
>>> S = Matrix(A.row_join(B))
>>> S
Matrix([
[1, 2],
[2, 3],
[3, 4]])
If each line of S represent coefficients of Ax + By
and x and y are [2, 3] then S*xy is:
>>> r = S*Matrix([2, 3]); r
Matrix([
[ 8],
[13],
[18]])
But let's add 1 to the middle value and then solve for the
least-squares value of xy:
>>> xy = S.solve_least_squares(Matrix([8, 14, 18])); xy
Matrix([
[ 5/3],
[10/3]])
The error is given by S*xy - r:
>>> S*xy - r
Matrix([
[1/3],
[1/3],
[1/3]])
>>> _.norm().n(2)
0.58
If a different xy is used, the norm will be higher:
>>> xy += ones(2, 1)/10
>>> (S*xy - r).norm().n(2)
1.5
"""
if method == "CH":
return self.cholesky_solve(rhs)
t = self.T
return (t * self).inv(method=method) * t * rhs
[docs] def solve(self, rhs, method="GE"):
"""Return solution to self*soln = rhs using given inversion method.
For a list of possible inversion methods, see the .inv() docstring.
"""
if not self.is_square:
if self.rows < self.cols:
raise ValueError(
"Under-determined system. " "Try M.gauss_jordan_solve(rhs)"
)
elif self.rows > self.cols:
raise ValueError(
"For over-determined system, M, having "
"more rows than columns, try M.solve_least_squares(rhs)."
)
else:
return self.inv(method=method) * rhs
[docs] def table(
self, printer, rowstart="[", rowend="]", rowsep="\n", colsep=", ", align="right"
):
r"""
String form of Matrix as a table.
``printer`` is the printer to use for on the elements (generally
something like StrPrinter())
``rowstart`` is the string used to start each row (by default '[').
``rowend`` is the string used to end each row (by default ']').
``rowsep`` is the string used to separate rows (by default a newline).
``colsep`` is the string used to separate columns (by default ', ').
``align`` defines how the elements are aligned. Must be one of 'left',
'right', or 'center'. You can also use '<', '>', and '^' to mean the
same thing, respectively.
This is used by the string printer for Matrix.
Examples
========
>>> from .. import Matrix
>>> from ..printing.str import StrPrinter
>>> M = Matrix([[1, 2], [-33, 4]])
>>> printer = StrPrinter()
>>> M.table(printer)
'[ 1, 2]\n[-33, 4]'
>>> print(M.table(printer))
[ 1, 2]
[-33, 4]
>>> print(M.table(printer, rowsep=',\n'))
[ 1, 2],
[-33, 4]
>>> print('[%s]' % M.table(printer, rowsep=',\n'))
[[ 1, 2],
[-33, 4]]
>>> print(M.table(printer, colsep=' '))
[ 1 2]
[-33 4]
>>> print(M.table(printer, align='center'))
[ 1 , 2]
[-33, 4]
>>> print(M.table(printer, rowstart='{', rowend='}'))
{ 1, 2}
{-33, 4}
"""
# Handle zero dimensions:
if self.rows == 0 or self.cols == 0:
return "[]"
# 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):
s = printer._print(self[i, j])
res[-1].append(s)
maxlen[j] = max(len(s), maxlen[j])
# Patch strings together
align = {
"left": "ljust",
"right": "rjust",
"center": "center",
"<": "ljust",
">": "rjust",
"^": "center",
}[align]
for i, row in enumerate(res):
for j, elem in enumerate(row):
row[j] = getattr(elem, align)(maxlen[j])
res[i] = rowstart + colsep.join(row) + rowend
return rowsep.join(res)
[docs] def upper_triangular_solve(self, rhs):
"""Solves Ax = B, where A is an upper triangular matrix.
See Also
========
lower_triangular_solve
gauss_jordan_solve
cholesky_solve
diagonal_solve
LDLsolve
LUsolve
QRsolve
pinv_solve
"""
if not self.is_square:
raise NonSquareMatrixError("Matrix must be square.")
if rhs.rows != self.rows:
raise TypeError("Matrix size mismatch.")
if not self.is_upper:
raise TypeError("Matrix is not upper triangular.")
return self._upper_triangular_solve(rhs)
[docs] def vech(self, diagonal=True, check_symmetry=True):
"""Return the unique elements of a symmetric Matrix as a one column matrix
by stacking the elements in the lower triangle.
Arguments:
diagonal -- include the diagonal cells of self or not
check_symmetry -- checks symmetry of self but not completely reliably
Examples
========
>>> from .. import Matrix
>>> m=Matrix([[1, 2], [2, 3]])
>>> m
Matrix([
[1, 2],
[2, 3]])
>>> m.vech()
Matrix([
[1],
[2],
[3]])
>>> m.vech(diagonal=False)
Matrix([[2]])
See Also
========
vec
"""
from ..matrices import zeros
c = self.cols
if c != self.rows:
raise ShapeError("Matrix must be square")
if check_symmetry:
self.simplify()
if self != self.transpose():
raise ValueError(
"Matrix appears to be asymmetric; consider check_symmetry=False"
)
count = 0
if diagonal:
v = zeros(c * (c + 1) // 2, 1)
for j in range(c):
for i in range(j, c):
v[count] = self[i, j]
count += 1
else:
v = zeros(c * (c - 1) // 2, 1)
for j in range(c):
for i in range(j + 1, c):
v[count] = self[i, j]
count += 1
return v
[docs]def classof(A, B):
"""
Get the type of the result when combining matrices of different types.
Currently the strategy is that immutability is contagious.
Examples
========
>>> from .. import Matrix, ImmutableMatrix
>>> from .matrices import classof
>>> M = Matrix([[1, 2], [3, 4]]) # a Mutable Matrix
>>> IM = ImmutableMatrix([[1, 2], [3, 4]])
>>> classof(M, IM)
<class 'sympy.matrices.immutable.ImmutableDenseMatrix'>
"""
try:
if A._class_priority > B._class_priority:
return A.__class__
else:
return B.__class__
except Exception:
pass
try:
import numpy
if isinstance(A, numpy.ndarray):
return B.__class__
if isinstance(B, numpy.ndarray):
return A.__class__
except Exception:
pass
raise TypeError("Incompatible classes %s, %s" % (A.__class__, B.__class__))
[docs]def a2idx(j, n=None):
"""Return integer after making positive and validating against n."""
if type(j) is not int:
try:
j = j.__index__()
except AttributeError:
raise IndexError("Invalid index a[%r]" % (j,))
if n is not None:
if j < 0:
j += n
if not (j >= 0 and j < n):
raise IndexError("Index out of range: a[%s]" % j)
return int(j)
def _find_reasonable_pivot(col, iszerofunc=_iszero, simpfunc=_simplify):
"""Find the lowest index of an item in `col` that is
suitable for a pivot. If `col` consists only of
Floats, the pivot with the largest norm is returned.
Otherwise, the first element where `iszerofunc` returns
False is used. If `iszerofunc` doesn't return false,
items are simplified and retested until a suitable
pivot is found.
Returns a 4-tuple
(pivot_offset, pivot_val, assumed_nonzero, newly_determined)
where pivot_offset is the index of the pivot, pivot_val is
the (possibly simplified) value of the pivot, assumed_nonzero
is True if an assumption that the pivot was non-zero
was made without being proved, and newly_determined are
elements that were simplified during the process of pivot
finding."""
newly_determined = []
col = list(col)
# a column that contains a mix of floats and integers
# but at least one float is considered a numerical
# column, and so we do partial pivoting
if all(isinstance(x, (Float, Integer)) for x in col) and any(
isinstance(x, Float) for x in col
):
col_abs = [abs(x) for x in col]
max_value = max(col_abs)
if iszerofunc(max_value):
# just because iszerofunc returned True, doesn't
# mean the value is numerically zero. Make sure
# to replace all entries with numerical zeros
if max_value != 0:
newly_determined = [(i, 0) for i, x in enumerate(col) if x != 0]
return (None, None, False, newly_determined)
index = col_abs.index(max_value)
return (index, col[index], False, newly_determined)
# PASS 1 (iszerofunc directly)
possible_zeros = []
for i, x in enumerate(col):
is_zero = iszerofunc(x)
# is someone wrote a custom iszerofunc, it may return
# BooleanFalse or BooleanTrue instead of True or False,
# so use == for comparison instead of `is`
if is_zero == False:
# we found something that is definitely not zero
return (i, x, False, newly_determined)
possible_zeros.append(is_zero)
# by this point, we've found no certain non-zeros
if all(possible_zeros):
# if everything is definitely zero, we have
# no pivot
return (None, None, False, newly_determined)
# PASS 2 (iszerofunc after simplify)
# we haven't found any for-sure non-zeros, so
# go through the elements iszerofunc couldn't
# make a determination about and opportunistically
# simplify to see if we find something
for i, x in enumerate(col):
if possible_zeros[i] is not None:
continue
simped = simpfunc(x)
is_zero = iszerofunc(simped)
if is_zero == True or is_zero == False:
newly_determined.append((i, simped))
if is_zero == False:
return (i, simped, False, newly_determined)
possible_zeros[i] = is_zero
# after simplifying, some things that were recognized
# as zeros might be zeros
if all(possible_zeros):
# if everything is definitely zero, we have
# no pivot
return (None, None, False, newly_determined)
# PASS 3 (.equals(0))
# some expressions fail to simplify to zero, but
# `.equals(0)` evaluates to True. As a last-ditch
# attempt, apply `.equals` to these expressions
for i, x in enumerate(col):
if possible_zeros[i] is not None:
continue
if x.equals(S.Zero):
# `.iszero` may return False with
# an implicit assumption (e.g., `x.equals(0)`
# when `x` is a symbol), so only treat it
# as proved when `.equals(0)` returns True
possible_zeros[i] = True
newly_determined.append((i, S.Zero))
if all(possible_zeros):
return (None, None, False, newly_determined)
# at this point there is nothing that could definitely
# be a pivot. To maintain compatibility with existing
# behavior, we'll assume that an illdetermined thing is
# non-zero. We should probably raise a warning in this case
i = possible_zeros.index(None)
return (i, col[i], True, newly_determined)
def _find_reasonable_pivot_naive(col, iszerofunc=_iszero, simpfunc=None):
"""
Helper that computes the pivot value and location from a
sequence of contiguous matrix column elements. As a side effect
of the pivot search, this function may simplify some of the elements
of the input column. A list of these simplified entries and their
indices are also returned.
This function mimics the behavior of _find_reasonable_pivot(),
but does less work trying to determine if an indeterminate candidate
pivot simplifies to zero. This more naive approach can be much faster,
with the trade-off that it may erroneously return a pivot that is zero.
`col` is a sequence of contiguous column entries to be searched for
a suitable pivot.
`iszerofunc` is a callable that returns a Boolean that indicates
if its input is zero, or None if no such determination can be made.
`simpfunc` is a callable that simplifies its input. It must return
its input if it does not simplify its input. Passing in
`simpfunc=None` indicates that the pivot search should not attempt
to simplify any candidate pivots.
Returns a 4-tuple:
(pivot_offset, pivot_val, assumed_nonzero, newly_determined)
`pivot_offset` is the sequence index of the pivot.
`pivot_val` is the value of the pivot.
pivot_val and col[pivot_index] are equivalent, but will be different
when col[pivot_index] was simplified during the pivot search.
`assumed_nonzero` is a boolean indicating if the pivot cannot be
guaranteed to be zero. If assumed_nonzero is true, then the pivot
may or may not be non-zero. If assumed_nonzero is false, then
the pivot is non-zero.
`newly_determined` is a list of index-value pairs of pivot candidates
that were simplified during the pivot search.
"""
# indeterminates holds the index-value pairs of each pivot candidate
# that is neither zero or non-zero, as determined by iszerofunc().
# If iszerofunc() indicates that a candidate pivot is guaranteed
# non-zero, or that every candidate pivot is zero then the contents
# of indeterminates are unused.
# Otherwise, the only viable candidate pivots are symbolic.
# In this case, indeterminates will have at least one entry,
# and all but the first entry are ignored when simpfunc is None.
indeterminates = []
for i, col_val in enumerate(col):
col_val_is_zero = iszerofunc(col_val)
if col_val_is_zero == False:
# This pivot candidate is non-zero.
return i, col_val, False, []
elif col_val_is_zero is None:
# The candidate pivot's comparison with zero
# is indeterminate.
indeterminates.append((i, col_val))
if len(indeterminates) == 0:
# All candidate pivots are guaranteed to be zero, i.e. there is
# no pivot.
return None, None, False, []
if simpfunc is None:
# Caller did not pass in a simplification function that might
# determine if an indeterminate pivot candidate is guaranteed
# to be nonzero, so assume the first indeterminate candidate
# is non-zero.
return indeterminates[0][0], indeterminates[0][1], True, []
# newly_determined holds index-value pairs of candidate pivots
# that were simplified during the search for a non-zero pivot.
newly_determined = []
for i, col_val in indeterminates:
tmp_col_val = simpfunc(col_val)
if id(col_val) != id(tmp_col_val):
# simpfunc() simplified this candidate pivot.
newly_determined.append((i, tmp_col_val))
if iszerofunc(tmp_col_val) == False:
# Candidate pivot simplified to a guaranteed non-zero value.
return i, tmp_col_val, False, newly_determined
return indeterminates[0][0], indeterminates[0][1], True, newly_determined
class _MinimalMatrix(object):
"""Class providing the minimum functionality
for a matrix-like object and implementing every method
required for a `MatrixRequired`. This class does not have everything
needed to become a full-fledged sympy object, but it will satisfy the
requirements of anything inheriting from `MatrixRequired`. If you wish
to make a specialized matrix type, make sure to implement these
methods and properties with the exception of `__init__` and `__repr__`
which are included for convenience."""
is_MatrixLike = True
_sympify = staticmethod(sympify)
_class_priority = 3
is_Matrix = True
is_MatrixExpr = False
@classmethod
def _new(cls, *args, **kwargs):
return cls(*args, **kwargs)
def __init__(self, rows, cols=None, mat=None):
if isinstance(mat, FunctionType):
# if we passed in a function, use that to populate the indices
mat = list(mat(i, j) for i in range(rows) for j in range(cols))
try:
if cols is None and mat is None:
mat = rows
rows, cols = mat.shape
except AttributeError:
pass
try:
# if we passed in a list of lists, flatten it and set the size
if cols is None and mat is None:
mat = rows
cols = len(mat[0])
rows = len(mat)
mat = [x for l in mat for x in l]
except (IndexError, TypeError):
pass
self.mat = tuple(self._sympify(x) for x in mat)
self.rows, self.cols = rows, cols
if self.rows is None or self.cols is None:
raise NotImplementedError("Cannot initialize matrix with given parameters")
def __getitem__(self, key):
def _normalize_slices(row_slice, col_slice):
"""Ensure that row_slice and col_slice don't have
`None` in their arguments. Any integers are converted
to slices of length 1"""
if not isinstance(row_slice, slice):
row_slice = slice(row_slice, row_slice + 1, None)
row_slice = slice(*row_slice.indices(self.rows))
if not isinstance(col_slice, slice):
col_slice = slice(col_slice, col_slice + 1, None)
col_slice = slice(*col_slice.indices(self.cols))
return (row_slice, col_slice)
def _coord_to_index(i, j):
"""Return the index in _mat corresponding
to the (i,j) position in the matrix."""
return i * self.cols + j
if isinstance(key, tuple):
i, j = key
if isinstance(i, slice) or isinstance(j, slice):
# if the coordinates are not slices, make them so
# and expand the slices so they don't contain `None`
i, j = _normalize_slices(i, j)
rowsList, colsList = (
list(range(self.rows))[i],
list(range(self.cols))[j],
)
indices = (i * self.cols + j for i in rowsList for j in colsList)
return self._new(
len(rowsList), len(colsList), list(self.mat[i] for i in indices)
)
# if the key is a tuple of ints, change
# it to an array index
key = _coord_to_index(i, j)
return self.mat[key]
def __eq__(self, other):
return self.shape == other.shape and list(self) == list(other)
def __len__(self):
return self.rows * self.cols
def __repr__(self):
return "_MinimalMatrix({}, {}, {})".format(self.rows, self.cols, self.mat)
@property
def shape(self):
return (self.rows, self.cols)
class _MatrixWrapper(object):
"""Wrapper class providing the minimum functionality
for a matrix-like object: .rows, .cols, .shape, indexability,
and iterability. CommonMatrix math operations should work
on matrix-like objects. For example, wrapping a numpy
matrix in a MatrixWrapper allows it to be passed to CommonMatrix.
"""
is_MatrixLike = True
def __init__(self, mat, shape=None):
self.mat = mat
self.rows, self.cols = mat.shape if shape is None else shape
def __getattr__(self, attr):
"""Most attribute access is passed straight through
to the stored matrix"""
return getattr(self.mat, attr)
def __getitem__(self, key):
return self.mat.__getitem__(key)
def _matrixify(mat):
"""If `mat` is a Matrix or is matrix-like,
return a Matrix or MatrixWrapper object. Otherwise
`mat` is passed through without modification."""
if getattr(mat, "is_Matrix", False):
return mat
if hasattr(mat, "shape"):
if len(mat.shape) == 2:
return _MatrixWrapper(mat)
return mat