Dense matrices over \(\ZZ/n\ZZ\) for \(n < 2^{11}\) using LinBox’s Modular<float>
¶
AUTHORS: - Burcin Erocal - Martin Albrecht
-
class
sage.matrix.matrix_modn_dense_float.
Matrix_modn_dense_float
¶ Bases:
sage.matrix.matrix_modn_dense_float.Matrix_modn_dense_template
Dense matrices over \(\ZZ/n\ZZ\) for \(n < 2^{11}\) using LinBox’s
Modular<float>
These are matrices with integer entries mod
n
represented as floating-point numbers in a 32-bit word for use with LinBox routines. This allows forn
up to \(2^{11}\). TheMatrix_modn_dense_double
class is used for larger moduli.Routines here are for the most basic access, see the \(matrix_modn_dense_template.pxi\) file for higher-level routines.
-
class
sage.matrix.matrix_modn_dense_float.
Matrix_modn_dense_template
¶ Bases:
sage.matrix.matrix_dense.Matrix_dense
Create a new matrix.
INPUT:
parent
– a matrix spaceentries
– seematrix()
copy
– ignored (for backwards compatibility)coerce
- perform modular reduction first?
EXAMPLES:
sage: A = random_matrix(GF(3),1000,1000) sage: type(A) <type 'sage.matrix.matrix_modn_dense_float.Matrix_modn_dense_float'> sage: A = random_matrix(Integers(10),1000,1000) sage: type(A) <type 'sage.matrix.matrix_modn_dense_float.Matrix_modn_dense_float'> sage: A = random_matrix(Integers(2^16),1000,1000) sage: type(A) <type 'sage.matrix.matrix_modn_dense_double.Matrix_modn_dense_double'>
-
charpoly
(var='x', algorithm='linbox')¶ Return the characteristic polynomial of
self
.INPUT:
var
- a variable namealgorithm
- ‘generic’, ‘linbox’ or ‘all’ (default: linbox)
EXAMPLES:
sage: A = random_matrix(GF(19), 10, 10); A [ 3 1 8 10 5 16 18 9 6 1] [ 5 14 4 4 14 15 5 11 3 0] [ 4 1 0 7 11 6 17 8 5 6] [ 4 6 9 4 8 1 18 17 8 18] [11 2 0 6 13 7 4 11 16 10] [12 6 12 3 15 10 5 11 3 8] [15 1 16 2 18 15 14 7 2 11] [16 16 17 7 14 12 7 7 0 5] [13 15 9 2 12 16 1 15 18 7] [10 8 16 18 9 18 2 13 5 10] sage: B = copy(A) sage: char_p = A.characteristic_polynomial(); char_p x^10 + 2*x^9 + 18*x^8 + 4*x^7 + 13*x^6 + 11*x^5 + 2*x^4 + 5*x^3 + 7*x^2 + 16*x + 6 sage: char_p(A) == 0 True sage: B == A # A is not modified True sage: min_p = A.minimal_polynomial(proof=True); min_p x^10 + 2*x^9 + 18*x^8 + 4*x^7 + 13*x^6 + 11*x^5 + 2*x^4 + 5*x^3 + 7*x^2 + 16*x + 6 sage: min_p.divides(char_p) True
sage: A = random_matrix(GF(2916337), 7, 7); A [ 514193 1196222 1242955 1040744 99523 2447069 40527] [ 930282 2685786 2892660 1347146 1126775 2131459 869381] [1853546 2266414 2897342 1342067 1054026 373002 84731] [1270068 2421818 569466 537440 572533 297105 1415002] [2079710 355705 2546914 2299052 2883413 1558788 1494309] [1027319 1572148 250822 522367 2516720 585897 2296292] [1797050 2128203 1161160 562535 2875615 1165768 286972] sage: B = copy(A) sage: char_p = A.characteristic_polynomial(); char_p x^7 + 1274305*x^6 + 1497602*x^5 + 12362*x^4 + 875330*x^3 + 31311*x^2 + 1858466*x + 700510 sage: char_p(A) == 0 True sage: B == A # A is not modified True sage: min_p = A.minimal_polynomial(proof=True); min_p x^7 + 1274305*x^6 + 1497602*x^5 + 12362*x^4 + 875330*x^3 + 31311*x^2 + 1858466*x + 700510 sage: min_p.divides(char_p) True sage: A = Mat(Integers(6),3,3)(range(9)) sage: A.charpoly() x^3
ALGORITHM: Uses LinBox if
self.base_ring()
is a field, otherwise use Hessenberg form algorithm.
-
determinant
()¶ Return the determinant of this matrix.
EXAMPLES:
sage: A = random_matrix(GF(7), 10, 10); A [3 1 6 6 4 4 2 2 3 5] [4 5 6 2 2 1 2 5 0 5] [3 2 0 5 0 1 5 4 2 3] [6 4 5 0 2 4 2 0 6 3] [2 2 4 2 4 5 3 4 4 4] [2 5 2 5 4 5 1 1 1 1] [0 6 3 4 2 2 3 5 1 1] [4 2 6 5 6 3 4 5 5 3] [5 2 4 3 6 2 3 6 2 1] [3 3 5 3 4 2 2 1 6 2] sage: A.determinant() 6
sage: A = random_matrix(GF(7), 100, 100) sage: A.determinant() 2 sage: A.transpose().determinant() 2 sage: B = random_matrix(GF(7), 100, 100) sage: B.determinant() 4 sage: (A*B).determinant() == A.determinant() * B.determinant() True
sage: A = random_matrix(GF(16007), 10, 10); A [ 5037 2388 4150 1400 345 5945 4240 14022 10514 700] [15552 8539 1927 3870 9867 3263 11637 609 15424 2443] [ 3761 15836 12246 15577 10178 13602 13183 15918 13942 2958] [ 4526 10817 6887 6678 1764 9964 6107 1705 5619 5811] [13537 15004 8307 11846 14779 550 14113 5477 7271 7091] [13338 4927 11406 13065 5437 12431 6318 5119 14198 496] [ 1044 179 12881 353 12975 12567 1092 10433 12304 954] [10072 8821 14118 13895 6543 13484 10685 14363 2612 11070] [15113 237 2612 14127 11589 5808 117 9656 15957 14118] [15233 11080 5716 9029 11402 9380 13045 13986 14544 5771] sage: A.determinant() 10207
sage: A = random_matrix(GF(16007), 100, 100) sage: A.determinant() 3576 sage: A.transpose().determinant() 3576 sage: B = random_matrix(GF(16007), 100, 100) sage: B.determinant() 4075 sage: (A*B).determinant() == A.determinant() * B.determinant() True
Parallel computation:
sage: A = random_matrix(GF(65521),200) sage: B = copy(A) sage: Parallelism().set('linbox', nproc=2) sage: d = A.determinant() sage: Parallelism().set('linbox', nproc=1) # switch off parallelization sage: e = B.determinant() sage: d==e True
-
echelonize
(algorithm='linbox_noefd', **kwds)¶ Put
self
in reduced row echelon form.INPUT:
self
- a mutable matrixalgorithm
linbox
- uses the LinBox library (wrapping fflas-ffpack)linbox_noefd
- uses the FFPACK directly, less memory and faster (default)gauss
- uses a custom slower \(O(n^3)\) Gauss elimination implemented in Sage.all
- compute using both algorithms and verify that the results are the same.
**kwds
- these are all ignored
OUTPUT:
self
is put in reduced row echelon form.the rank of self is computed and cached
the pivot columns of self are computed and cached.
the fact that self is now in echelon form is recorded and cached so future calls to echelonize return immediately.
EXAMPLES:
sage: A = random_matrix(GF(7), 10, 20); A [3 1 6 6 4 4 2 2 3 5 4 5 6 2 2 1 2 5 0 5] [3 2 0 5 0 1 5 4 2 3 6 4 5 0 2 4 2 0 6 3] [2 2 4 2 4 5 3 4 4 4 2 5 2 5 4 5 1 1 1 1] [0 6 3 4 2 2 3 5 1 1 4 2 6 5 6 3 4 5 5 3] [5 2 4 3 6 2 3 6 2 1 3 3 5 3 4 2 2 1 6 2] [0 5 6 3 2 5 6 6 3 2 1 4 5 0 2 6 5 2 5 1] [4 0 4 2 6 3 3 5 3 0 0 1 2 5 5 1 6 0 0 3] [2 0 1 0 0 3 0 2 4 2 2 4 4 4 5 4 1 2 3 4] [2 4 1 4 3 0 6 2 2 5 2 5 3 6 4 2 2 6 4 4] [0 0 2 2 1 6 2 0 5 0 4 3 1 6 0 6 0 4 6 5] sage: A.echelon_form() [1 0 0 0 0 0 0 0 0 0 6 2 6 0 1 1 2 5 6 2] [0 1 0 0 0 0 0 0 0 0 0 4 5 4 3 4 2 5 1 2] [0 0 1 0 0 0 0 0 0 0 6 3 4 6 1 0 3 6 5 6] [0 0 0 1 0 0 0 0 0 0 0 3 5 2 3 4 0 6 5 3] [0 0 0 0 1 0 0 0 0 0 0 6 3 4 5 3 0 4 3 2] [0 0 0 0 0 1 0 0 0 0 1 1 0 2 4 2 5 5 5 0] [0 0 0 0 0 0 1 0 0 0 1 0 1 3 2 0 0 0 5 3] [0 0 0 0 0 0 0 1 0 0 4 4 2 6 5 4 3 4 1 0] [0 0 0 0 0 0 0 0 1 0 1 0 4 2 3 5 4 6 4 0] [0 0 0 0 0 0 0 0 0 1 2 0 5 0 5 5 3 1 1 4]
sage: A = random_matrix(GF(13), 10, 10); A [ 8 3 11 11 9 4 8 7 9 9] [ 2 9 6 5 7 12 3 4 11 5] [12 6 11 12 4 3 3 8 9 5] [ 4 2 10 5 10 1 1 1 6 9] [12 8 5 5 11 4 1 2 8 11] [ 2 6 9 11 4 7 1 0 12 2] [ 8 9 0 7 7 7 10 4 1 4] [ 0 8 2 6 7 5 7 12 2 3] [ 2 11 12 3 4 7 2 9 6 1] [ 0 11 5 9 4 5 5 8 7 10] sage: MS = parent(A) sage: B = A.augment(MS(1)) sage: B.echelonize() sage: A.rank() 10 sage: C = B.submatrix(0,10,10,10); C [ 4 9 4 4 0 4 7 11 9 11] [11 7 6 8 2 8 6 11 9 5] [ 3 9 9 2 4 8 9 2 9 4] [ 7 0 11 4 0 9 6 11 8 1] [12 12 4 12 3 12 6 1 7 12] [12 2 11 6 6 6 7 0 10 6] [ 0 7 3 4 7 11 10 12 4 6] [ 5 11 0 5 3 11 4 12 5 12] [ 6 7 3 5 1 4 11 7 4 1] [ 4 9 6 7 11 1 2 12 6 7] sage: ~A == C True
sage: A = random_matrix(Integers(10), 10, 20) sage: A.echelon_form() Traceback (most recent call last): ... NotImplementedError: Echelon form not implemented over 'Ring of integers modulo 10'.
sage: A = random_matrix(GF(16007), 10, 20); A [15455 1177 10072 4693 3887 4102 10746 15265 6684 14559 4535 13921 9757 9525 9301 8566 2460 9609 3887 6205] [ 8602 10035 1242 9776 162 7893 12619 6660 13250 1988 14263 11377 2216 1247 7261 8446 15081 14412 7371 7948] [12634 7602 905 9617 13557 2694 13039 4936 12208 15480 3787 11229 593 12462 5123 14167 6460 3649 5821 6736] [10554 2511 11685 12325 12287 6534 11636 5004 6468 3180 3607 11627 13436 5106 3138 13376 8641 9093 2297 5893] [ 1025 11376 10288 609 12330 3021 908 13012 2112 11505 56 5971 338 2317 2396 8561 5593 3782 7986 13173] [ 7607 588 6099 12749 10378 111 2852 10375 8996 7969 774 13498 12720 4378 6817 6707 5299 9406 13318 2863] [15545 538 4840 1885 8471 1303 11086 14168 1853 14263 3995 12104 1294 7184 1188 11901 15971 2899 4632 711] [ 584 11745 7540 15826 15027 5953 7097 14329 10889 12532 13309 15041 6211 1749 10481 9999 2751 11068 21 2795] [ 761 11453 3435 10596 2173 7752 15941 14610 1072 8012 9458 5440 612 10581 10400 101 11472 13068 7758 7898] [10658 4035 6662 655 7546 4107 6987 1877 4072 4221 7679 14579 2474 8693 8127 12999 11141 605 9404 10003] sage: A.echelon_form() [ 1 0 0 0 0 0 0 0 0 0 8416 8364 10318 1782 13872 4566 14855 7678 11899 2652] [ 0 1 0 0 0 0 0 0 0 0 4782 15571 3133 10964 5581 10435 9989 14303 5951 8048] [ 0 0 1 0 0 0 0 0 0 0 15688 6716 13819 4144 257 5743 14865 15680 4179 10478] [ 0 0 0 1 0 0 0 0 0 0 4307 9488 2992 9925 13984 15754 8185 11598 14701 10784] [ 0 0 0 0 1 0 0 0 0 0 927 3404 15076 1040 2827 9317 14041 10566 5117 7452] [ 0 0 0 0 0 1 0 0 0 0 1144 10861 5241 6288 9282 5748 3715 13482 7258 9401] [ 0 0 0 0 0 0 1 0 0 0 769 1804 1879 4624 6170 7500 11883 9047 874 597] [ 0 0 0 0 0 0 0 1 0 0 15591 13686 5729 11259 10219 13222 15177 15727 5082 11211] [ 0 0 0 0 0 0 0 0 1 0 8375 14939 13471 12221 8103 4212 11744 10182 2492 11068] [ 0 0 0 0 0 0 0 0 0 1 6534 396 6780 14734 1206 3848 7712 9770 10755 410]
sage: A = random_matrix(Integers(10000), 10, 20) sage: A.echelon_form() Traceback (most recent call last): ... NotImplementedError: Echelon form not implemented over 'Ring of integers modulo 10000'.
Parallel computation:
sage: A = random_matrix(GF(65521),100,200) sage: Parallelism().set('linbox', nproc=2) sage: E = A.echelon_form() sage: Parallelism().set('linbox', nproc=1) # switch off parallelization sage: F = A.echelon_form() sage: E==F True
-
hessenbergize
()¶ Transforms self in place to its Hessenberg form.
EXAMPLES:
sage: A = random_matrix(GF(17), 10, 10, density=0.1); A [ 0 0 0 0 12 0 0 0 0 0] [ 0 0 0 4 0 0 0 0 0 0] [ 0 0 0 0 2 0 0 0 0 0] [ 0 14 0 0 0 0 0 0 0 0] [ 0 0 0 0 0 10 0 0 0 0] [ 0 0 0 0 0 16 0 0 0 0] [ 0 0 0 0 0 0 6 0 0 0] [15 0 0 0 0 0 0 0 0 0] [ 0 0 0 16 0 0 0 0 0 0] [ 0 5 0 0 0 0 0 0 0 0] sage: A.hessenbergize(); A [ 0 0 0 0 0 0 0 12 0 0] [15 0 0 0 0 0 0 0 0 0] [ 0 0 0 0 0 0 0 2 0 0] [ 0 0 0 0 14 0 0 0 0 0] [ 0 0 0 4 0 0 0 0 0 0] [ 0 0 0 0 5 0 0 0 0 0] [ 0 0 0 0 0 0 6 0 0 0] [ 0 0 0 0 0 0 0 0 0 10] [ 0 0 0 0 0 0 0 0 0 0] [ 0 0 0 0 0 0 0 0 0 16]
-
lift
()¶ Return the lift of this matrix to the integers.
EXAMPLES:
sage: A = matrix(GF(7),2,3,[1..6]) sage: A.lift() [1 2 3] [4 5 6] sage: A.lift().parent() Full MatrixSpace of 2 by 3 dense matrices over Integer Ring sage: A = matrix(GF(16007),2,3,[1..6]) sage: A.lift() [1 2 3] [4 5 6] sage: A.lift().parent() Full MatrixSpace of 2 by 3 dense matrices over Integer Ring
Subdivisions are preserved when lifting:
sage: A.subdivide([], [1,1]); A [1||2 3] [4||5 6] sage: A.lift() [1||2 3] [4||5 6]
-
minpoly
(var='x', algorithm='linbox', proof=None)¶ Returns the minimal polynomial of`` self``.
INPUT:
var
- a variable namealgorithm
-generic
orlinbox
(default:linbox
)proof
– (default:True
); whether to provably return the true minimal polynomial; ifFalse
, we only guarantee to return a divisor of the minimal polynomial. There are also certainly cases where the computed results is frequently not exactly equal to the minimal polynomial (but is instead merely a divisor of it).
Warning
If
proof=True
, minpoly is insanely slow compared toproof=False
. This matters since proof=True is the default, unless you first typeproof.linear_algebra(False)
.EXAMPLES:
sage: A = random_matrix(GF(17), 10, 10); A [ 2 14 0 15 11 10 16 2 9 4] [10 14 1 14 3 14 12 14 3 13] [10 1 14 6 2 14 13 7 6 14] [10 3 9 15 8 1 5 8 10 11] [ 5 12 4 9 15 2 6 11 2 12] [ 6 10 12 0 6 9 7 7 3 8] [ 2 9 1 5 12 13 7 16 7 11] [11 1 0 2 0 4 7 9 8 15] [ 5 3 16 2 11 10 12 14 0 7] [16 4 6 5 2 3 14 15 16 4] sage: B = copy(A) sage: min_p = A.minimal_polynomial(proof=True); min_p x^10 + 13*x^9 + 10*x^8 + 9*x^7 + 10*x^6 + 4*x^5 + 10*x^4 + 10*x^3 + 12*x^2 + 14*x + 7 sage: min_p(A) == 0 True sage: B == A True sage: char_p = A.characteristic_polynomial(); char_p x^10 + 13*x^9 + 10*x^8 + 9*x^7 + 10*x^6 + 4*x^5 + 10*x^4 + 10*x^3 + 12*x^2 + 14*x + 7 sage: min_p.divides(char_p) True
sage: A = random_matrix(GF(1214471), 10, 10); A [ 160562 831940 65852 173001 515930 714380 778254 844537 584888 392730] [ 502193 959391 614352 775603 240043 1156372 104118 1175992 612032 1049083] [ 660489 1066446 809624 15010 1002045 470722 314480 1155149 1173111 14213] [1190467 1079166 786442 429883 563611 625490 1015074 888047 1090092 892387] [ 4724 244901 696350 384684 254561 898612 44844 83752 1091581 349242] [ 130212 580087 253296 472569 913613 919150 38603 710029 438461 736442] [ 943501 792110 110470 850040 713428 668799 1122064 325250 1084368 520553] [1179743 791517 34060 1183757 1118938 642169 47513 73428 1076788 216479] [ 626571 105273 400489 1041378 1186801 158611 888598 1138220 1089631 56266] [1092400 890773 1060810 211135 719636 1011640 631366 427711 547497 1084281] sage: B = copy(A) sage: min_p = A.minimal_polynomial(proof=True); min_p x^10 + 384251*x^9 + 702437*x^8 + 960299*x^7 + 202699*x^6 + 409368*x^5 + 1109249*x^4 + 1163061*x^3 + 333802*x^2 + 273775*x + 55190 sage: min_p(A) == 0 True sage: B == A True sage: char_p = A.characteristic_polynomial(); char_p x^10 + 384251*x^9 + 702437*x^8 + 960299*x^7 + 202699*x^6 + 409368*x^5 + 1109249*x^4 + 1163061*x^3 + 333802*x^2 + 273775*x + 55190 sage: min_p.divides(char_p) True
EXAMPLES:
sage: R.<x>=GF(3)[] sage: A = matrix(GF(3),2,[0,0,1,2]) sage: A.minpoly() x^2 + x sage: A.minpoly(proof=False) in [x, x+1, x^2+x] True
-
randomize
(density=1, nonzero=False)¶ Randomize
density
proportion of the entries of this matrix, leaving the rest unchanged.INPUT:
density
- Integer; proportion (roughly) to be consideredfor changes
nonzero
- Bool (default:False
); whether the newentries are forced to be non-zero
OUTPUT:
None, the matrix is modified in-space
EXAMPLES:
sage: A = matrix(GF(5), 5, 5, 0) sage: A.randomize(0.5); A [0 0 0 2 0] [0 3 0 0 2] [4 0 0 0 0] [4 0 0 0 0] [0 1 0 0 0] sage: A.randomize(); A [3 3 2 1 2] [4 3 3 2 2] [0 3 3 3 3] [3 3 2 2 4] [2 2 2 1 4]
The matrix is updated instead of overwritten:
sage: A = random_matrix(GF(5), 100, 100, density=0.1) sage: A.density() 961/10000 sage: A.randomize(density=0.1) sage: A.density() 801/5000
-
rank
()¶ Return the rank of this matrix.
EXAMPLES:
sage: A = random_matrix(GF(3), 100, 100) sage: B = copy(A) sage: A.rank() 99 sage: B == A True sage: A = random_matrix(GF(3), 100, 100, density=0.01) sage: A.rank() 63 sage: A = matrix(GF(3), 100, 100) sage: A.rank() 0
Rank is not implemented over the integers modulo a composite yet.:
sage: M = matrix(Integers(4), 2, [2,2,2,2]) sage: M.rank() Traceback (most recent call last): ... NotImplementedError: Echelon form not implemented over 'Ring of integers modulo 4'.
sage: A = random_matrix(GF(16007), 100, 100) sage: B = copy(A) sage: A.rank() 100 sage: B == A True sage: MS = A.parent() sage: MS(1) == ~A*A True
-
right_kernel_matrix
(algorithm='linbox', basis='echelon')¶ Returns a matrix whose rows form a basis for the right kernel of
self
, whereself
is a matrix over a (small) finite field.INPUT:
algorithm
– (default:'linbox'
) a parameter that is passed on toself.echelon_form
, if computation of an echelon form is required; see that routine for allowable valuesbasis
– (default:'echelon'
) a keyword that describes the format of the basis returned, allowable values are:'echelon'
: the basis matrix is in echelon form'pivot'
: the basis matrix is such that the submatrix obtainedby taking the columns that in
self
contain no pivots, is the identity matrix
'computed'
: no work is done to transform the basis; inthe current implementation the result is the negative of that returned by
'pivot'
OUTPUT:
A matrix
X
whose rows are a basis for the right kernel ofself
. This means thatself * X.transpose()
is a zero matrix.The result is not cached, but the routine benefits when
self
is known to be in echelon form already.EXAMPLES:
sage: M = matrix(GF(5),6,6,range(36)) sage: M.right_kernel_matrix(basis='computed') [4 2 4 0 0 0] [3 3 0 4 0 0] [2 4 0 0 4 0] [1 0 0 0 0 4] sage: M.right_kernel_matrix(basis='pivot') [1 3 1 0 0 0] [2 2 0 1 0 0] [3 1 0 0 1 0] [4 0 0 0 0 1] sage: M.right_kernel_matrix() [1 0 0 0 0 4] [0 1 0 0 1 3] [0 0 1 0 2 2] [0 0 0 1 3 1] sage: M * M.right_kernel_matrix().transpose() [0 0 0 0] [0 0 0 0] [0 0 0 0] [0 0 0 0] [0 0 0 0] [0 0 0 0]
-
submatrix
(row=0, col=0, nrows=- 1, ncols=- 1)¶ Return the matrix constructed from self using the specified range of rows and columns.
INPUT:
row
,col
– index of the starting row and column. Indices start at zeronrows
,ncols
– (optional) number of rows and columns to take. If not provided, take all rows below and all columns to the right of the starting entry
See also
The functions
matrix_from_rows()
,matrix_from_columns()
, andmatrix_from_rows_and_columns()
allow one to select arbitrary subsets of rows and/or columns.EXAMPLES:
Take the \(3 \times 3\) submatrix starting from entry \((1,1)\) in a \(4 \times 4\) matrix:
sage: m = matrix(GF(17),4, [1..16]) sage: m.submatrix(1, 1) [ 6 7 8] [10 11 12] [14 15 16]
Same thing, except take only two rows:
sage: m.submatrix(1, 1, 2) [ 6 7 8] [10 11 12]
And now take only one column:
sage: m.submatrix(1, 1, 2, 1) [ 6] [10]
You can take zero rows or columns if you want:
sage: m.submatrix(0, 0, 0) [] sage: parent(m.submatrix(0, 0, 0)) Full MatrixSpace of 0 by 4 dense matrices over Finite Field of size 17
-
transpose
()¶ Return the transpose of
self
, without changingself
.EXAMPLES:
We create a matrix, compute its transpose, and note that the original matrix is not changed.
sage: M = MatrixSpace(GF(41), 2) sage: A = M([1,2,3,4]) sage: B = A.transpose() sage: B [1 3] [2 4] sage: A [1 2] [3 4]
.T
is a convenient shortcut for the transpose:sage: A.T [1 3] [2 4]
sage: A.subdivide(None, 1); A [1|2] [3|4] sage: A.transpose() [1 3] [---] [2 4]