Miscellaneous matrix functions¶
- sage.matrix.matrix_misc.permanental_minor_polynomial(A, permanent_only=False, var='t', prec=None)¶
Return the polynomial of the sums of permanental minors of
A
.INPUT:
A – a matrix
permanentonly – if True, return only the permanent of A
var – name of the polynomial variable
prec – if prec is not None, truncate the polynomial at precision prec
The polynomial of the sums of permanental minors is
min(nrows,ncols)∑i=0pi(A)xiwhere pi(A) is the i-th permanental minor of A (that can also be obtained through the method
permanental_minor()
viaA.permanental_minor(i)
).The algorithm implemented by that function has been developed by P. Butera and M. Pernici, see [BP2015]. Its complexity is O(2nm2n) where m and n are the number of rows and columns of A. Moreover, if A is a banded matrix with width w, that is Aij=0 for |i−j|>w and w<n/2, then the complexity of the algorithm is O(4w(w+1)n2).
INPUT:
A
– matrixpermanent_only
– optional boolean. IfTrue
, only the permanent is computed (might be faster).var
– a variable name
EXAMPLES:
sage: from sage.matrix.matrix_misc import permanental_minor_polynomial sage: m = matrix([[1,1],[1,2]]) sage: permanental_minor_polynomial(m) 3*t^2 + 5*t + 1 sage: permanental_minor_polynomial(m, permanent_only=True) 3 sage: permanental_minor_polynomial(m, prec=2) 5*t + 1
sage: M = MatrixSpace(ZZ,4,4) sage: A = M([1,0,1,0,1,0,1,0,1,0,10,10,1,0,1,1]) sage: permanental_minor_polynomial(A) 84*t^3 + 114*t^2 + 28*t + 1 sage: [A.permanental_minor(i) for i in range(5)] [1, 28, 114, 84, 0]
An example over Q:
sage: M = MatrixSpace(QQ,2,2) sage: A = M([1/5,2/7,3/2,4/5]) sage: permanental_minor_polynomial(A, True) 103/175
An example with polynomial coefficients:
sage: R.<a> = PolynomialRing(ZZ) sage: A = MatrixSpace(R,2)([[a,1], [a,a+1]]) sage: permanental_minor_polynomial(A, True) a^2 + 2*a
A usage of the
var
argument:sage: m = matrix(ZZ,4,[0,1,2,3,1,2,3,0,2,3,0,1,3,0,1,2]) sage: permanental_minor_polynomial(m, var='x') 164*x^4 + 384*x^3 + 172*x^2 + 24*x + 1
ALGORITHM:
The permanent perm(A) of a n×n matrix A is the coefficient of the x1x2…xn monomial in
n∏i=1(n∑j=1Aijxj)Evaluating this product one can neglect x2i, that is xi can be considered to be nilpotent of order 2.
To formalize this procedure, consider the algebra R=K[η1,η2,…,ηn] where the ηi are commuting, nilpotent of order 2 (i.e. η2i=0). Formally it is the quotient ring of the polynomial ring in η1,η2,…,ηn quotiented by the ideal generated by the η2i.
We will mostly consider the ring R[t] of polynomials over R. We denote a generic element of R[t] by p(η1,…,ηn) or p(ηi1,…,ηik) if we want to emphasize that some monomials in the ηi are missing.
Introduce an “integration” operation ⟨p⟩ over R and R[t] consisting in the sum of the coefficients of the non-vanishing monomials in ηi (i.e. the result of setting all variables ηi to 1). Let us emphasize that this is not a morphism of algebras as ⟨η1⟩2=1 while ⟨η21⟩=0!
Let us consider an example of computation. Let p1=1+tη1+tη2 and p2=1+tη1+tη3. Then
p1p2=1+2tη1+t(η2+η3)+t2(η1η2+η1η3+η2η3)and
⟨p1p2⟩=1+4t+3t2In this formalism, the permanent is just
perm(A)=⟨n∏i=1n∑j=1Aijηj⟩A useful property of ⟨.⟩ which makes this algorithm efficient for band matrices is the following: let p1(η1,…,ηn) and p2(ηj,…,ηn) be polynomials in R[t] where j≥1. Then one has
⟨p1(η1,…,ηn)p2⟩=⟨p1(1,…,1,ηj,…,ηn)p2⟩where η1,..,ηj−1 are replaced by 1 in p1. Informally, we can “integrate” these variables before performing the product. More generally, if a monomial ηi is missing in one of the terms of a product of two terms, then it can be integrated in the other term.
Now let us consider an m×n matrix with m≤n. The sum of permanental `k`-minors of `A` is
perm(A,k)=∑r,cperm(Ar,c)where the sum is over the k-subsets r of rows and k-subsets c of columns and Ar,c is the submatrix obtained from A by keeping only the rows r and columns c. Of course perm(A,min and note that perm(A,1) is just the sum of all entries of the matrix.
The generating function of these sums of permanental minors is
g(t) = \left\langle \prod_{i=1}^m \left(1 + t \sum_{j=1}^n A_{ij} \eta_j\right) \right\rangleIn fact the t^k coefficient of g(t) corresponds to choosing k rows of A; \eta_i is associated to the i-th column; nilpotency avoids having twice the same column in a product of A’s.
For more details, see the article [BP2015].
From a technical point of view, the product in K[\eta_1, \ldots, \eta_n][t] is implemented as a subroutine in
prm_mul()
. The indices of the rows and columns actually start at 0, so the variables are \eta_0, \ldots, \eta_{n-1}. Polynomials are represented in dictionary form: to a variable \eta_i is associated the key 2^i (or in Python1 << i
). The keys associated to products are obtained by considering the development in base 2: to the monomial \eta_{i_1} \ldots \eta_{i_k} is associated the key 2^{i_1} + \ldots + 2^{i_k}. So the product \eta_1 \eta_2 corresponds to the key 6 = (110)_2 while \eta_0 \eta_3 has key 9 = (1001)_2. In particular all operations on monomials are implemented via bitwise operations on the keys.
- sage.matrix.matrix_misc.prm_mul(p1, p2, mask_free, prec)¶
Return the product of
p1
andp2
, putting free variables inmask_free
to 1.This function is mainly use as a subroutine of
permanental_minor_polynomial()
.INPUT:
p1,p2 – polynomials as dictionaries
mask_free – an integer mask that give the list of free variables (the i-th variable is free if the i-th bit of
mask_free
is 1)prec – if prec is not None, truncate the product at precision prec
EXAMPLES:
sage: from sage.matrix.matrix_misc import prm_mul sage: t = polygen(ZZ, 't') sage: p1 = {0: 1, 1: t, 4: t} sage: p2 = {0: 1, 1: t, 2: t} sage: prm_mul(p1, p2, 1, None) {0: 2*t + 1, 2: t^2 + t, 4: t^2 + t, 6: t^2}
- sage.matrix.matrix_misc.row_iterator(A)¶