source: sage/matrix/matrix_mod2_dense.pyx @ 3931:2e3dbe94ad80

Revision 3931:2e3dbe94ad80, 25.7 KB checked in by 'Martin Albrecht <malb@…, 6 years ago (diff)

fixed bug where rank wasn't computed correctly

Line 
1"""
2Dense matrice over GF(2) based on code by Gregory Bard.
3
4This implementation uses a packed representation of boolean matrices
5and provides a quite fast echelon form implemenation (M4RI).
6
7#For some solutions LinBox is used.
8
9AUTHOR: Martin Albrecht <malb@informatik.uni-bremen.de>
10
11EXAMPLES:
12    sage: a = matrix(GF(2),3,range(9),sparse=False); a
13    [0 1 0]
14    [1 0 1]
15    [0 1 0]
16    sage: a.rank()
17    2
18    sage: type(a)
19    <type 'sage.matrix.matrix_mod2_dense.Matrix_mod2_dense'>
20    sage: a[0,0] = 1
21    sage: a.rank()
22    3
23    sage: parent(a)
24    Full MatrixSpace of 3 by 3 dense matrices over Finite Field of size 2
25
26    sage: a^2
27    [0 1 1]
28    [1 0 0]
29    [1 0 1]
30    sage: a+a
31    [0 0 0]
32    [0 0 0]
33    [0 0 0]
34
35    sage: b = a.new_matrix(2,3,range(6)); b
36    [0 1 0]
37    [1 0 1]
38
39    sage: a*b
40    Traceback (most recent call last):
41    ...
42    TypeError: incompatible dimensions
43    sage: b*a
44    [1 0 1]
45    [1 0 0]
46
47    sage: a == loads(dumps(a))
48    True
49    sage: b == loads(dumps(b))
50    True
51
52    sage: a.echelonize(); a
53    [1 0 0]
54    [0 1 0]
55    [0 0 1]
56    sage: b.echelonize(); b
57    [1 0 1]
58    [0 1 0]
59
60TESTS:
61    sage: FF = FiniteField(2)
62    sage: V = VectorSpace(FF,2)
63    sage: v = V([0,1]); v
64    (0, 1)
65    sage: W = V.subspace([v])
66    sage: W
67    Vector space of degree 2 and dimension 1 over Finite Field of size 2
68    Basis matrix:
69    [0 1]
70    sage: v in W
71    True
72
73    sage: M = Matrix(GF(2), [[1,1,0],[0,1,0]])
74    sage: M.row_space()
75    Vector space of degree 3 and dimension 2 over Finite Field of size 2
76    Basis matrix:
77    [1 0 0]
78    [0 1 0]
79
80    sage: M = Matrix(GF(2), [[1,1,0],[0,0,1]])
81    sage: M.row_space()
82    Vector space of degree 3 and dimension 2 over Finite Field of size 2
83    Basis matrix:
84    [1 1 0]
85    [0 0 1]
86
87TODO:
88   - make linbox frontend and use it
89     - charpoly ?
90     - minpoly ?
91     - rank ?
92   - make Matrix_modn_frontend and use it (?)
93"""
94
95##############################################################################
96#       Copyright (C) 2004,2005,2006 William Stein <wstein@gmail.com>
97#  Distributed under the terms of the GNU General Public License (GPL)
98#  The full text of the GPL is available at:
99#                  http://www.gnu.org/licenses/
100##############################################################################
101
102include "../ext/interrupt.pxi"
103include "../ext/cdefs.pxi"
104include '../ext/stdsage.pxi'
105include '../ext/random.pxi'
106
107cimport matrix_dense
108from sage.structure.element cimport Matrix
109from sage.structure.element cimport ModuleElement, Element
110
111from sage.misc.functional import log
112
113from sage.misc.misc import verbose, get_verbose, cputime
114
115## from sage.libs.linbox.linbox cimport Linbox_mod2_dense
116## cdef Linbox_mod2_dense linbox
117## linbox = Linbox_mod2_dense()
118
119cdef object called
120
121cdef void init_m4ri():
122    global called
123    if called is None:
124        # TODO: remove packingmask
125        setupPackingMasks()
126        buildAllCodes()
127        called = True
128
129init_m4ri()
130
131cdef class Matrix_mod2_dense(matrix_dense.Matrix_dense):   # dense or sparse
132    """
133   
134    """
135    ########################################################################
136    # LEVEL 1 functionality
137    ########################################################################
138    def __new__(self, parent, entries, copy, coerce, alloc=True):
139        """
140        Creates a new dense matrix over GF(2).
141
142        INPUT:
143            parent -- MatrixSpace (always)
144            entries -- ignored
145            copy -- ignored
146            coerce -- ignored
147            alloc -- if True a zero matrix is allocated (default:True)
148       
149        """
150        matrix_dense.Matrix_dense.__init__(self, parent)
151
152        if alloc:
153            self._entries = createPackedMatrix(self._nrows, self._ncols)
154
155        # cache elements
156        self._zero = self._base_ring(0)
157        self._one = self._base_ring(1)
158
159    def __dealloc__(self):
160        if self._entries:
161            destroyPackedMatrix(self._entries)
162            self._entries = NULL
163
164    def __init__(self, parent, entries, copy, coerce):
165        """
166        Dense matrix over GF(2) constructor.
167
168        INPUT:
169            parent -- MatrixSpace.
170            entries -- may be list or 0 or 1
171            copy -- ignored, elements are always copied
172            coerce -- ignored, elements are always coerced to ints % 2
173
174        EXAMPLES:
175            sage: type(random_matrix(GF(2),2,2))
176            <type 'sage.matrix.matrix_mod2_dense.Matrix_mod2_dense'>
177
178            sage: Matrix(GF(2),3,3,1)
179            [1 0 0]
180            [0 1 0]
181            [0 0 1]
182
183            sage: Matrix(GF(2),2,2,[1,1,1,0])
184            [1 1]
185            [1 0]
186
187            sage: Matrix(GF(2),2,2,4)
188            [0 0]
189            [0 0]
190        """
191        cdef int i,j,e
192
193        # scalar ? 
194        if not isinstance(entries, list):
195            if int(entries) % 2 == 1:
196                makeIdentityPacked(self._entries)
197            return
198
199        # all entries are given as a long list
200        if len(entries) != self._nrows * self._ncols:
201            raise IndexError, "The vector of entries has the wrong length."
202       
203        k = 0
204        R = self.base_ring()
205
206        cdef PyObject** w
207        w = FAST_SEQ_UNSAFE(entries)
208       
209        for i from 0 <= i < self._nrows:
210            if PyErr_CheckSignals(): raise KeyboardInterrupt
211            for j from 0 <= j < self._ncols:
212                writePackedCell(self._entries,i,j, int(<object> w[k]) % 2)
213                k = k + 1
214           
215    def __richcmp__(Matrix self, right, int op):  # always need for mysterious reasons.
216        return self._richcmp(right, op)
217
218    def __hash__(self):
219        return self._hash()
220   
221    cdef set_unsafe(self, Py_ssize_t i, Py_ssize_t j, value):
222        writePackedCell(self._entries, i, j, int(value))
223
224    cdef get_unsafe(self, Py_ssize_t i, Py_ssize_t j):
225        if readPackedCell(self._entries, i, j):
226            return self._one
227        else:
228            return self._zero
229
230
231    def str(self):
232        cdef int i,j
233        s = []
234        for i from 0 <= i < self._nrows:
235            rl = []
236            for j from 0 <= j < self._ncols:
237                rl.append(str(readPackedCell(self._entries,i,j)))
238            s.append( " ".join(rl) )
239        return "[" + "]\n[".join(s) + "]"
240
241    ########################################################################
242    # LEVEL 2 functionality
243    #   * def _pickle
244    #   * def _unpickle
245    #   * cdef _mul_c_impl
246    #   * cdef _cmp_c_impl
247    #   * _list -- list of underlying elements (need not be a copy)
248    #   * _dict -- sparse dictionary of underlying elements (need not be a copy)
249    ########################################################################
250    # def _pickle(self):
251    # def _unpickle(self, data, int version):   # use version >= 0
252
253    cdef ModuleElement _add_c_impl(self, ModuleElement right):
254        """
255        Matrix addition.
256
257        INPUT:
258            right -- matrix of dimension self.nrows() x self.ncols()
259
260        EXAMPLES:
261            sage: A = random_matrix(GF(2),10,10)
262            sage: A + A == Matrix(GF(2),10,10,0)
263            True
264       
265        """
266        cdef Matrix_mod2_dense A
267        A = Matrix_mod2_dense.__new__(Matrix_mod2_dense, self._parent, 0, 0, 0, alloc=False)
268       
269        A._entries = addPacked(self._entries,(<Matrix_mod2_dense>right)._entries)
270
271        return A
272
273    cdef ModuleElement _sub_c_impl(self, ModuleElement right):
274        """
275        Matrix addition.
276
277        INPUT:
278            right -- matrix of dimension self.nrows() x self.ncols()
279
280        EXAMPLES:
281            sage: A = random_matrix(GF(2),10,10)
282            sage: A - A == Matrix(GF(2),10,10,0)
283            True
284        """
285        return self._add_c_impl(right)
286       
287    cdef Matrix _matrix_times_matrix_c_impl(self, Matrix right):
288        """
289        Matrix multiplication.
290
291        EXAMPLES:
292            sage: A = random_matrix(GF(2),200,200)
293            sage: A*A == A._multiply_linbox(A) # optional
294            True
295
296        ALGORITHM: Uses the 'Method of the Four Russians
297        Multiplication', see self._multiply_m4rm.
298
299        """
300        if get_verbose() >= 2:
301            verbose('matrix multiply of %s x %s matrix by %s x %s matrix'%(
302                self._nrows, self._ncols, right._nrows, right._ncols))
303
304##         if self._ncols < 1000 and right.nrows() < 1000:
305##              return self._multiply_classical(right)
306##         else:
307##             # uses way more RAM but is faster
308##              return self._multiply_linbox(right)
309
310       
311        cdef int n = self._ncols
312        cdef int k = round(min(0.75 * log(n,2), 16))
313
314        if k < 1:
315            k = 1
316
317##         if ( self.nrows() < right.ncols() ):
318##             return self._multiply_m4rm_c(right,k,1) # transpose
319##         else:
320        return self._multiply_m4rm_c(right,k,0)
321
322##     def _multiply_linbox(Matrix_mod2_dense self, Matrix_mod2_dense right):
323##         """
324##         Multiply matrices using LinBox.
325
326##         INPUT:
327##             right -- Matrix
328
329##         EXAMPLE:
330##               sage: A = Matrix(GF(2), 4, 3, [0, 0, 0, 0, 1, 0, 0, 1, 1, 0, 0, 1] )
331##               sage: B = Matrix(GF(2), 3, 4, [0, 0, 1, 0, 1, 0, 0, 1, 1, 1, 0, 0] )
332##               sage: A
333##               [0 0 0]
334##               [0 1 0]
335##               [0 1 1]
336##               [0 0 1]
337##               sage: B
338##               [0 0 1 0]
339##               [1 0 0 1]
340##               [1 1 0 0]
341##               sage: A._multiply_linbox(B)
342##               [0 0 0 0]
343##               [1 0 0 1]
344##               [0 1 0 1]
345##               [1 1 0 0]
346
347##         ALGORITHM: Uses LinBox
348
349##         """
350##         if get_verbose() >= 2:
351##             verbose('linbox multiply of %s x %s matrix by %s x %s matrix'%(
352##                 self._nrows, self._ncols, right._nrows, right._ncols))
353##         cdef int e
354##         cdef Matrix_mod2_dense ans
355       
356##         ans = self.new_matrix(nrows = self.nrows(), ncols = right.ncols())
357
358##         linbox.set(self._entries)
359
360##         _sig_on
361##         linbox.matrix_matrix_multiply(ans._entries, right._entries)
362##         _sig_off
363##         return ans
364
365    def _multiply_m4rm(Matrix_mod2_dense self, Matrix_mod2_dense right, k=0, transpose=False):
366        """
367        Multiply matrices using the 'Method of the Four Russians Multiplication' (M4RM).
368
369        The algorithm is based on an algorithm by Arlazarov, Dinic,
370        Kronrod, and Faradzev [ADKF70] and appeared in [AHU]. This
371        implementation is based on a description given in Gregory
372        Bard's 'Method of the Four Russians Inversion' paper [B06].
373
374        INPUT:
375            right     -- Matrix
376            k         -- parameter k for the Gray Code table size. If k=0 a suitable value is
377                         chosen by the function. (0<= k <= 16, default: 0)
378            transpose -- Calculate C = [B^T * A^T]^T which equals A * B. This might safe memory for large k
379                         if B.ncols() > A.nrows() and might also be faster (default:False)
380
381        EXAMPLE:
382              sage: A = Matrix(GF(2), 4, 3, [0, 0, 0, 0, 1, 0, 0, 1, 1, 0, 0, 1] )
383              sage: B = Matrix(GF(2), 3, 4, [0, 0, 1, 0, 1, 0, 0, 1, 1, 1, 0, 0] )
384              sage: A
385              [0 0 0]
386              [0 1 0]
387              [0 1 1]
388              [0 0 1]
389              sage: B
390              [0 0 1 0]
391              [1 0 0 1]
392              [1 1 0 0]
393              sage: A._multiply_m4rm(B)
394              [0 0 0 0]
395              [1 0 0 1]
396              [0 1 0 1]
397              [1 1 0 0]
398
399        ALGORITHM: Uses M4RM
400
401        REFERENCES:
402            [AHU]    A. Aho, J. Hopcroft, and J. Ullman. 'Chapter 6: Matrix Multiplication and Related Oper-
403                     ations.' The Design and Analysis of Computer Algorithms. Addison-Wesley, 1974.
404
405            [ADKF70] V. Arlazarov, E. Dinic, M. Kronrod, and I. Faradzev. 'On Economical Construction of
406                     the Transitive Closure of a Directed Graph.' Dokl. Akad. Nauk. SSSR No. 194 (in Russian),
407                     English Translation in Soviet Math Dokl. No. 11, 1970.
408
409            [Bard06] G. Bard. 'Accelerating Cryptanalysis with the Method of Four Russians'. Cryptography
410                     E-Print Archive (http://eprint.iacr.org/2006/251.pdf), 2006.
411
412        """
413
414
415        if k == 0:
416            n = self._ncols
417            k = round(min(0.75 * log(n,2), 16))
418       
419        if k<1 or k>16:
420            raise RuntimeError,"k must be between 1 and 16 or 0"
421        k = round(k)
422
423        if self._ncols != right._nrows:
424            raise ArithmeticError, "left ncols must match right nrows"
425
426        return self._multiply_m4rm_c(right, k, transpose)
427       
428    cdef Matrix_mod2_dense _multiply_m4rm_c(Matrix_mod2_dense self, Matrix_mod2_dense right, int k, int transpose):
429        if get_verbose() >= 2:
430            verbose('m4rm multiply of %s x %s matrix by %s x %s matrix'%(
431                self._nrows, self._ncols, right._nrows, right._ncols))
432
433        cdef Matrix_mod2_dense ans
434       
435        ans = self.new_matrix(nrows = self.nrows(), ncols = right.ncols())
436        destroyPackedMatrix(ans._entries)
437
438        if not transpose:
439            ans._entries = m4rmPacked(self._entries, right._entries, k)
440        else:
441            ans._entries = m4rmTransposePacked(self._entries, right._entries, k)
442        return ans
443
444
445    def _multiply_classical(Matrix_mod2_dense self, Matrix_mod2_dense right):
446        """
447        Classical n^3 multiplication.
448
449
450        EXAMPLE:
451              sage: A = Matrix(GF(2), 4, 3, [0, 0, 0, 0, 1, 0, 0, 1, 1, 0, 0, 1] )
452              sage: B = Matrix(GF(2), 3, 4, [0, 0, 1, 0, 1, 0, 0, 1, 1, 1, 0, 0] )
453              sage: A
454              [0 0 0]
455              [0 1 0]
456              [0 1 1]
457              [0 0 1]
458              sage: B
459              [0 0 1 0]
460              [1 0 0 1]
461              [1 1 0 0]
462              sage: A._multiply_classical(B)
463              [0 0 0 0]
464              [1 0 0 1]
465              [0 1 0 1]
466              [1 1 0 0]
467
468        """
469        cdef Matrix_mod2_dense A
470        A = self.new_matrix(nrows = self.nrows(), ncols = right.ncols())
471        destroyPackedMatrix(A._entries)
472
473        A._entries = matrixTimesMatrixPacked(self._entries,(<Matrix_mod2_dense>right)._entries)
474        return A
475
476    # cdef int _cmp_c_impl(self, Matrix right) except -2:
477
478    def __neg__(self):
479        """
480
481        EXAMPLES:
482            sage: A = random_matrix(GF(2),100,100)
483            sage: A - A == A - -A
484            True
485        """
486        return self.copy()
487
488    def __invert__(self):
489        """
490        Inverts self using M4RI.
491
492        If self is not invertible a ZeroDivisionError is raised.
493
494        EXAMPLE:
495              sage: A = Matrix(GF(2),3,3, [0, 0, 1, 0, 1, 1, 1, 0, 1])
496              sage: MS = A.parent()
497              sage: A
498              [0 0 1]
499              [0 1 1]
500              [1 0 1]
501              sage: ~A
502              [1 0 1]
503              [1 1 0]
504              [1 0 0]
505              sage: A * ~A == ~A * A == MS(1)
506              True
507
508        ALGORITHM: Uses M4RI.
509        """
510        cdef int k
511        cdef packedmatrix *I
512        cdef Matrix_mod2_dense A
513        k = 8
514
515        if self._nrows != self._ncols:
516            raise ArithmeticError, "self must be a square matrix"
517
518        I = createPackedMatrix(self._nrows,self._ncols)
519        makeIdentityPacked(I)
520
521        A = Matrix_mod2_dense.__new__(Matrix_mod2_dense, self._parent, 0, 0, 0, alloc = False)
522        A._entries = invertPackedFlexRussian(self._entries, I, k)
523        destroyPackedMatrix(I)
524
525        if A._entries==NULL:
526            raise ZeroDivisionError, "self is not invertible"
527        else:
528            return A
529
530    def __copy__(self):
531        """
532        Returns a copy of self.
533
534        EXAMPLES:
535             sage: MS = MatrixSpace(GF(2),3,3)
536             sage: A = MS(1)
537             sage: A.copy() == A
538             True
539             sage: A.copy() is A
540             False
541
542             sage: A = random_matrix(GF(2),100,100)
543             sage: A.copy() == A
544             True
545             sage: A.copy() is A
546             False
547
548             sage: A.echelonize()
549             sage: A.copy() == A
550             True
551       
552        """
553        cdef Matrix_mod2_dense A
554        cdef int width
555       
556        A = Matrix_mod2_dense.__new__(Matrix_mod2_dense, self._parent, 0, 0, 0)
557
558        memcpy(A._entries.values, self._entries.values, (RADIX>>3) * self._entries.width * self._nrows)
559        memcpy(A._entries.rowswap, self._entries.rowswap, self._nrows * sizeof(int))
560
561        return A
562       
563    def _list(self):
564        """
565        Returns list of the elements of self in row major ordering.
566
567        EXAMPLE:
568            sage: A = Matrix(GF(2),2,2,[1,0,1,1])
569            sage: A
570            [1 0]
571            [1 1]
572            sage: A.list()
573            [1, 0, 1, 1]
574       
575        """
576        cdef int i,j
577        l = []
578        for i from 0 <= i < self._nrows:
579            for j from 0 <= j < self._ncols:
580                if readPackedCell(self._entries,i,j):
581                    l.append(self._one)
582                else:
583                    l.append(self._zero)
584        return l
585
586    # def _dict(self):   
587
588    ########################################################################
589    # LEVEL 3 functionality (Optional)
590    #    * __deepcopy__
591    #    * Matrix windows -- only if you need strassen for that base
592    ########################################################################
593
594    def echelonize(self, algorithm='m4ri', cutoff=0, **kwds):
595        """
596        Puts self in reduced row echelon form.
597
598        INPUT:
599            self -- a mutable matrix
600            algorithm -- 'm4ri' -- uses M4RI (default)
601                      -- 'linbox' -- uses LinBox's Strassen algorithm
602                                     (needs more RAM, slower for many examples, but
603                                     asymtotically faster)
604            k --  the parameter 'k' of the M4RI algorithm. It MUST be between
605                  1 and 16 (inclusive). If it is not specified it will be calculated as
606                  3/4 * log_2( min(nrows, ncols) ) as suggested in the M4RI paper.
607
608        EXAMPLE:
609             sage: A = random_matrix(GF(2), 10, 10)
610             sage: B = A.copy(); B.echelonize() # fastest
611             sage: C = A.copy(); C.echelonize(k=2) # force k
612             sage: D = A.copy(); D.echelonize(algorithm='linbox') # optional for LinBox
613             sage: E = A.copy(); E.echelonize(algorithm='classical') # force Gaussian elimination
614             sage: B == C == E
615             True
616
617        TESTS:
618             sage: VF2 = VectorSpace(GF(2),2)
619             sage: WF2 = VF2.submodule([VF2([1,1])])
620             sage: WF2
621             Vector space of degree 2 and dimension 1 over Finite Field of size 2
622             Basis matrix:
623             [1 1]
624
625        ALGORITHM: Uses Gregory Bard's M4RI algorithm and implementation or
626                   LinBox.
627
628        REFERENCES:
629            [Bard06] G. Bard. 'Accelerating Cryptanalysis with the Method of Four Russians'. Cryptography
630                     E-Print Archive (http://eprint.iacr.org/2006/251.pdf), 2006.
631        """
632        cdef int k, n
633
634        x = self.fetch('in_echelon_form')
635        if not x is None: return  # already known to be in echelon form
636
637        if algorithm == 'm4ri':
638
639            self.check_mutability()
640            self.clear_cache()       
641
642            if 'k' in kwds:
643                k = int(kwds['k'])
644
645                if k<1 or k>16:
646                    raise RuntimeError,"k must be between 1 and 16"
647                k = round(k)
648            else:
649                n = min(self._nrows, self._ncols)
650                k = round(min(0.75 * log(n,2), 16))
651                if k<1:
652                    k = 1
653
654            _sig_on
655            r =  simpleFourRussiansPackedFlex(self._entries, 1, k)
656            _sig_off
657           
658            self.cache('in_echelon_form',True)
659            self.cache('rank', r)
660            self.cache('pivots', self._pivots())
661
662        elif algorithm == 'linbox':
663
664            #self._echelonize_linbox()
665            raise NotImplementedError
666
667        elif algorithm == 'classical':
668
669            # for debugging purposes only, it is slow
670            self._echelon_in_place_classical()
671        else:
672            raise ValueError, "no algorithm '%s'"%algorithm
673       
674##     def _echelonize_linbox(self):
675##         """
676##         Puts self in row echelon form using LinBox.
677##         """
678##         self.check_mutability()
679##         self.clear_cache()       
680       
681##         t = verbose('calling linbox echelonize')
682##         _sig_on
683##         linbox.set(self._entries)
684##         r = linbox.echelonize()
685##         _sig_off
686##         verbose('done with linbox echelonize',t)
687
688##         self.cache('in_echelon_form',True)
689##         self.cache('rank', r)
690##         #self.cache('pivots', self._pivots())
691
692    def _pivots(self):
693        """
694        Returns the pivot columns of self if self is in row echelon form.
695        """
696        if not self.fetch('in_echelon_form'):
697            raise RuntimeError, "self must be in reduced row echelon form first."
698        pivots = []
699        cdef Py_ssize_t i, j, nc
700        nc = self._ncols
701        i = 0
702        while i < self._nrows:
703            for j from i <= j < nc:
704                if readPackedCell(self._entries, i, j):
705                    pivots.append(j)
706                    i += 1
707                    break
708            if j == nc:
709                break
710        return pivots
711
712    def randomize(self, density=1):
713        """
714        Randomize density proportion of the entries of this matrix,
715        leaving the rest unchanged.
716
717        NOTE: The random() function doesn't seem to be as random as
718        expected. The lower-order bits seem to have a strong bias
719        towards zero. Even stranger, if you create two 1000x1000
720        matrices over GF(2) they always appear to have the same
721        reduced row echelon form, i.e. they span the same space. The
722        higher-order bits seem to be much more random and thus we
723        shift first and mod p then.
724        """
725       
726        density = float(density)
727        if density == 0:
728            return
729
730        self.check_mutability()
731        self.clear_cache()
732
733        cdef int nc
734        if density == 1:
735            fillRandomlyPacked(self._entries)
736        else:
737            nc = self._ncols
738            num_per_row = int(density * nc)
739            _sig_on
740            for i from 0 <= i < self._nrows:
741                for j from 0 <= j < num_per_row:
742                    k = ((random()>>16)+random())%nc # is this safe?
743                    # 16-bit seems safe
744                    writePackedCell(self._entries, i, j, (random()>>16) % 2)
745            _sig_off
746
747
748    cdef rescale_row_c(self, Py_ssize_t row, multiple, Py_ssize_t start_col):
749        if (int(multiple)%2) == 0:
750            rowClearPackedOffset(self._entries, row, start_col);
751           
752##     cdef rescale_col_c(self, Py_ssize_t col, multiple, Py_ssize_t start_row):
753##         pass
754
755    cdef add_multiple_of_row_c(self,  Py_ssize_t row_to, Py_ssize_t row_from, multiple,
756                               Py_ssize_t start_col):
757        if (int(multiple)%2) != 0:
758            rowAddPackedOffset(self._entries, row_from, row_to, start_col)
759   
760##     cdef add_multiple_of_column_c(self, Py_ssize_t col_to, Py_ssize_t col_from, s,
761##                                    Py_ssize_t start_row):
762##         pass
763
764    cdef swap_rows_c(self, Py_ssize_t row1, Py_ssize_t row2):
765        rowSwapPacked(self._entries, row1, row2)
766   
767##     cdef swap_columns_c(self, Py_ssize_t col1, Py_ssize_t col2):
768##         pass
769
770    def _magma_init_(self):
771        """
772        Returns a string of self in MAGMA form. Does not return MAGMA
773        object but string.
774        """
775        cdef int i,j
776        K = self._base_ring._magma_init_()
777        if self._nrows == self._ncols:
778            s = 'MatrixAlgebra(%s, %s)'%(K, self.nrows())
779        else:
780            s = 'RMatrixSpace(%s, %s, %s)'%(K, self.nrows(), self.ncols())
781        v = []
782        for i from 0 <= i < self._nrows:
783            for j from 0 <= j < self._ncols:
784                v.append(str(readPackedCell(self._entries,i,j)))
785        return s + '![%s]'%(','.join(v))
786
787    def transpose(self):
788        """
789        Returns transpose of self and leaves self untouched.
790
791        EXAMPLE:
792            sage: A = Matrix(GF(2),3,5,[1, 0, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 1, 0])
793            sage: A
794            [1 0 1 0 0]
795            [0 1 1 0 0]
796            [1 1 0 1 0]
797            sage: B = A.transpose(); B
798            [1 0 1]
799            [0 1 1]
800            [1 1 0]
801            [0 0 1]
802            [0 0 0]
803            sage: B.transpose() == A
804            True
805       
806        """
807        cdef Matrix_mod2_dense A = self.new_matrix(ncols = self._nrows,  nrows = self._ncols)
808        destroyPackedMatrix(A._entries)
809
810        A._entries = transposePacked(self._entries)
811        return A
812
813    cdef int _cmp_c_impl(self, Element right) except -2:
814        return comparePackedMatrix(self._entries, (<Matrix_mod2_dense>right)._entries)
815
816
817    def augment(self, Matrix_mod2_dense right):
818        """
819        Augements self with right.
820
821        EXAMPLE:
822            sage: MS = MatrixSpace(GF(2),3,3)
823            sage: A = MS([0, 1, 0, 1, 1, 0, 1, 1, 1]); A
824            [0 1 0]
825            [1 1 0]
826            [1 1 1]
827            sage: B = A.augment(MS(1)); B
828            [0 1 0 1 0 0]
829            [1 1 0 0 1 0]
830            [1 1 1 0 0 1]
831            sage: B.echelonize(); B
832            [1 0 0 1 1 0]
833            [0 1 0 1 0 0]
834            [0 0 1 0 1 1]
835            sage: C = B.matrix_from_columns([3,4,5]); C
836            [1 1 0]
837            [1 0 0]
838            [0 1 1]
839            sage: C == ~A
840            True
841            sage: C*A == MS(1)
842            True
843           
844        """
845        cdef Matrix_mod2_dense A
846
847        A = self.new_matrix(ncols = self._ncols + right._ncols)
848        destroyPackedMatrix(A._entries)
849       
850        A._entries = concatPacked(self._entries, right._entries)
851        return A
852
853       
Note: See TracBrowser for help on using the repository browser.