Dense matrices over univariate polynomials over fields¶
The implementation inherits from Matrix_generic_dense but some algorithms are optimized for polynomial matrices.
AUTHORS:
- Kwankyu Lee (2016-12-15): initial version with code moved from other files.
- Johan Rosenkilde (2017-02-07): added weak_popov_form()
-
class
sage.matrix.matrix_polynomial_dense.
Matrix_polynomial_dense
¶ Bases:
sage.matrix.matrix_generic_dense.Matrix_generic_dense
Dense matrix over a univariate polynomial ring over a field.
-
hermite_form
(include_zero_rows=True, transformation=False)¶ Return the Hermite form of this matrix.
The Hermite form is also normalized, i.e., the pivot polynomials are monic.
INPUT:
include_zero_rows
– boolean (default:True
); ifFalse
, the zero rows in the output matrix are deletedtransformation
– boolean (default:False
); ifTrue
, return the transformation matrix
OUTPUT:
- the Hermite normal form \(H\) of this matrix \(A\)
- (optional) transformation matrix \(U\) such that \(UA = H\)
EXAMPLES:
sage: M.<x> = GF(7)[] sage: A = matrix(M, 2, 3, [x, 1, 2*x, x, 1+x, 2]) sage: A.hermite_form() [ x 1 2*x] [ 0 x 5*x + 2] sage: A.hermite_form(transformation=True) ( [ x 1 2*x] [1 0] [ 0 x 5*x + 2], [6 1] ) sage: A = matrix(M, 2, 3, [x, 1, 2*x, 2*x, 2, 4*x]) sage: A.hermite_form(transformation=True, include_zero_rows=False) ([ x 1 2*x], [0 4]) sage: H, U = A.hermite_form(transformation=True, include_zero_rows=True); H, U ( [ x 1 2*x] [0 4] [ 0 0 0], [5 1] ) sage: U * A == H True sage: H, U = A.hermite_form(transformation=True, include_zero_rows=False) sage: U * A [ x 1 2*x] sage: U * A == H True
-
is_weak_popov
()¶ Return
True
if the matrix is in weak Popov form.OUTPUT:
A matrix over a polynomial ring is in weak Popov form if all leading positions are different [MS2003]. A leading position is the position \(i\) in a row with the highest degree; in case of tie, the maximal \(i\) is used (i.e. furthest to the right).
EXAMPLES:
A matrix with the same leading position in two rows is not in weak Popov form:
sage: PF = PolynomialRing(GF(2^12,'a'),'x') sage: A = matrix(PF,3,[x, x^2, x^3,\ x^2, x^2, x^2,\ x^3, x^2, x ]) sage: A.is_weak_popov() False
If a matrix has different leading positions, it is in weak Popov form:
sage: B = matrix(PF,3,[1, 1, x^3,\ x^2, 1, 1,\ 1,x^ 2, 1 ]) sage: B.is_weak_popov() True
Weak Popov form is not restricted to square matrices:
sage: PF = PolynomialRing(GF(7),'x') sage: D = matrix(PF,2,4,[x^2+1, 1, 2, x,\ 3*x+2, 0, 0, 0 ]) sage: D.is_weak_popov() False
Even a matrix with more rows than columns can still be in weak Popov form:
sage: E = matrix(PF,4,2,[4*x^3+x, x^2+5*x+2,\ 0, 0,\ 4, x,\ 0, 0 ]) sage: E.is_weak_popov() True
A matrix with fewer columns than non-zero rows is never in weak Popov form:
sage: F = matrix(PF,3,2,[x^2, x,\ x^3+2, x,\ 4, 5]) sage: F.is_weak_popov() False
See also
weak_popov_form
AUTHOR:
- David Moedinger (2014-07-30)
-
row_reduced_form
(transformation=None, shifts=None)¶ Return a row reduced form of this matrix.
A matrix \(M\) is row reduced if the (row-wise) leading term matrix has the same rank as \(M\). The (row-wise) leading term matrix of a polynomial matrix \(M\) is the matrix over \(k\) whose \((i,j)\)‘th entry is the \(x^{d_i}\) coefficient of \(M[i,j]\), where \(d_i\) is the greatest degree among polynomials in the \(i\)‘th row of \(M_0\).
A row reduced form is non-canonical so a given matrix has many row reduced forms.
INPUT:
transformation
– (default:False
). If thisTrue
, the transformation matrix \(U\) will be returned as well: this is an invertible matrix over \(k[x]\) such thatself
equals \(UW\), where \(W\) is the output matrix.shifts
– (default:None
) A tuple or list of integers \(s_1, \ldots, s_n\), where \(n\) is the number of columns of the matrix. If given, a “shifted row reduced form” is computed, i.e. such that the matrix \(A\,\mathrm{diag}(x^{s_1}, \ldots, x^{s_n})\) is row reduced, where \(\mathrm{diag}\) denotes a diagonal matrix.
OUTPUT:
- \(W\) – a row reduced form of this matrix.
EXAMPLES:
sage: R.<t> = GF(3)['t'] sage: M = matrix([[(t-1)^2],[(t-1)]]) sage: M.row_reduced_form() [ 0] [t + 2]
If the matrix is an \(n \times 1\) matrix with at least one non-zero entry, \(W\) has a single non-zero entry and that entry is a scalar multiple of the greatest-common-divisor of the entries of the matrix:
sage: M1 = matrix([[t*(t-1)*(t+1)],[t*(t-2)*(t+2)],[t]]) sage: output1 = M1.row_reduced_form() sage: output1 [0] [0] [t]
The following is the first half of example 5 in [Hes2002] except that we have transposed the matrix; [Hes2002] uses column operations and we use row:
sage: R.<t> = QQ['t'] sage: M = matrix([[t^3 - t,t^2 - 2],[0,t]]).transpose() sage: M.row_reduced_form() [ t -t^2] [t^2 - 2 t]
The next example demonstrates what happens when the matrix is a zero matrix:
sage: R.<t> = GF(5)['t'] sage: M = matrix(R, 2, [0,0,0,0]) sage: M.row_reduced_form() [0 0] [0 0]
In the following example, the original matrix is already row reduced, but the output is a different matrix. This is because currently this method simply computes a weak Popov form, which is always also a row reduced matrix (see
weak_popov_form()
). This behavior is likely to change when a faster algorithm designed specifically for row reduced form is implemented in Sage:sage: R.<t> = QQ['t'] sage: M = matrix([[t,t,t],[0,0,t]]); M [t t t] [0 0 t] sage: M.row_reduced_form() [ t t t] [-t -t 0]
The last example shows the usage of the transformation parameter:
sage: Fq.<a> = GF(2^3) sage: Fx.<x> = Fq[] sage: A = matrix(Fx,[[x^2+a,x^4+a],[x^3,a*x^4]]) sage: W,U = A.row_reduced_form(transformation=True); sage: W,U ( [ x^2 + a x^4 + a] [1 0] [x^3 + a*x^2 + a^2 a^2], [a 1] ) sage: U*W == A True sage: U.is_invertible() True
-
weak_popov_form
(transformation=False, shifts=None)¶ Return a weak Popov form of the matrix.
A matrix is in weak Popov form if the leading positions of the nonzero rows are all different. The leading position of a row is the right-most position whose entry has the maximal degree in the row.
The weak Popov form is non-canonical, so an input matrix have many weak Popov forms.
INPUT:
transformation
– boolean (default:False
) IfTrue
, the transformation matrix is returned together with the weak Popov form.shifts
– (default:None
) A tuple or list of integers \(s_1, \ldots, s_n\), where \(n\) is the number of columns of the matrix. If given, a “shifted weak Popov form” is computed, i.e. such that the matrix \(A\,\mathrm{diag}(x^{s_1}, \ldots, x^{s_n})\) is in weak Popov form, where \(\mathrm{diag}\) denotes a diagonal matrix.
ALGORITHM:
This method implements the Mulders-Storjohann algorithm of [MS2003].
EXAMPLES:
sage: F.<a> = GF(2^4,'a') sage: PF.<x> = F[] sage: A = matrix(PF,[[1, a*x^17 + 1 ],\ [0, a*x^11 + a^2*x^7 + 1 ]]) sage: M, U = A.weak_popov_form(transformation=True) sage: U * A == M True sage: M.is_weak_popov() True sage: U.is_invertible() True
A zero matrix will return itself:
sage: Z = matrix(PF,5,3) sage: Z.weak_popov_form() [0 0 0] [0 0 0] [0 0 0] [0 0 0] [0 0 0]
Shifted weak popov form is computed if
shifts
is given:sage: PF.<x> = QQ[] sage: A = matrix(PF,3,[x, x^2, x^3,\ x^2, x^1, 0,\ x^3, x^3, x^3]) sage: A.weak_popov_form() [ x x^2 x^3] [ x^2 x 0] [ x^3 - x x^3 - x^2 0] sage: H,U = A.weak_popov_form(transformation=True, shifts=[16,8,0]) sage: H [ x x^2 x^3] [ 0 -x^2 + x -x^4 + x^3] [ 0 0 -x^5 + x^4 + x^3] sage: U * A == H True
See also
is_weak_popov
-