source: sage/matrix/matrix_space.py @ 12287:d7533ae4895e

Revision 12287:d7533ae4895e, 42.7 KB checked in by Mike Hansen <mhansen@…>, 4 years ago (diff)

Updates for Pynac-0.1.7, main symbolics switch (#5930)

Line 
1r"""
2Matrix Spaces.
3
4You can create any space `\text{Mat}_{n\times m}(R)` of
5either dense or sparse matrices with given number of rows and
6columns over any commutative or noncommutative ring.
7
8EXAMPLES::
9
10    sage: MS = MatrixSpace(QQ,6,6,sparse=True); MS
11    Full MatrixSpace of 6 by 6 sparse matrices over Rational Field
12    sage: MS.base_ring()
13    Rational Field
14    sage: MS = MatrixSpace(ZZ,3,5,sparse=False); MS
15    Full MatrixSpace of 3 by 5 dense matrices over Integer Ring
16
17TESTS::
18
19    sage: matrix(RR,2,2,sparse=True)
20    [0.000000000000000 0.000000000000000]
21    [0.000000000000000 0.000000000000000]
22    sage: matrix(GF(11),2,2,sparse=True)
23    [0 0]
24    [0 0]
25"""
26
27# System imports
28import types
29import weakref
30import operator
31
32# SAGE matrix imports
33import matrix
34import matrix_generic_dense
35import matrix_generic_sparse
36
37import matrix_modn_dense
38import matrix_modn_sparse
39
40import matrix_mod2_dense
41#import matrix_mod2_sparse
42
43import matrix_integer_dense
44import matrix_integer_sparse
45
46import matrix_rational_dense
47import matrix_rational_sparse
48
49import matrix_mpolynomial_dense
50
51#import padics.matrix_padic_capped_relative_dense
52
53## import matrix_cyclo_dense
54## import matrix_cyclo_sparse
55
56
57# IMPORTANT - these two guys get imported below only later
58# since they currently force numpy to import, which takes
59# a *long* time.
60#import matrix_real_double_dense
61#import matrix_complex_double_dense
62
63import sage.groups.matrix_gps.matrix_group_element
64
65
66# SAGE imports
67import sage.structure.coerce
68import sage.structure.parent_gens as parent_gens
69import sage.rings.ring as ring
70import sage.rings.rational_field as rational_field
71import sage.rings.integer_ring as integer_ring
72import sage.rings.integer as integer
73import sage.rings.field as field
74import sage.rings.finite_field as finite_field
75import sage.rings.principal_ideal_domain as principal_ideal_domain
76import sage.rings.integral_domain as integral_domain
77import sage.rings.number_field.all
78import sage.rings.integer_mod_ring
79import sage.rings.polynomial.multi_polynomial_ring_generic
80import sage.misc.latex as latex
81#import sage.rings.real_double as real_double
82import sage.misc.mrange
83import sage.modules.free_module_element
84import sage.modules.free_module
85from sage.structure.sequence import Sequence
86
87def is_MatrixSpace(x):
88    """
89    Returns True if self is an instance of MatrixSpace returns false if
90    self is not an instance of MatrixSpace
91   
92    EXAMPLES::
93   
94        sage: from sage.matrix.matrix_space import is_MatrixSpace
95        sage: MS = MatrixSpace(QQ,2)
96        sage: A = MS.random_element()
97        sage: is_MatrixSpace(MS)
98        True
99        sage: is_MatrixSpace(A)
100        False
101        sage: is_MatrixSpace(5)
102        False
103    """
104
105    return isinstance(x, MatrixSpace_generic)
106
107_cache = {}
108def MatrixSpace(base_ring, nrows, ncols=None, sparse=False):
109    """
110    Create with the command
111   
112    MatrixSpace(base_ring , nrows [, ncols] [, sparse])
113   
114    The default value of the optional argument sparse is False. The
115    default value of the optional argument ncols is nrows.
116   
117    INPUT:
118   
119   
120    -  ``base_ring`` - a ring
121   
122    -  ``nrows`` - int, the number of rows
123   
124    -  ``ncols`` - (default nrows) int, the number of
125       columns
126   
127    -  ``sparse`` - (default false) whether or not matrices
128       are given a sparse representation
129   
130   
131    OUTPUT: The unique space of all nrows x ncols matrices over
132    base_ring.
133   
134    EXAMPLES::
135   
136        sage: MS = MatrixSpace(RationalField(),2)
137        sage: MS.base_ring()
138        Rational Field
139        sage: MS.dimension()
140        4
141        sage: MS.dims()
142        (2, 2)
143        sage: B = MS.basis()
144        sage: B
145        [
146        [1 0]
147        [0 0],
148        [0 1]
149        [0 0],
150        [0 0]
151        [1 0],
152        [0 0]
153        [0 1]
154        ]
155        sage: B[0]
156        [1 0]
157        [0 0]
158        sage: B[1]
159        [0 1]
160        [0 0]
161        sage: B[2]
162        [0 0]
163        [1 0]
164        sage: B[3]
165        [0 0]
166        [0 1]
167        sage: A = MS.matrix([1,2,3,4])
168        sage: A
169        [1 2]
170        [3 4]
171        sage: MS2 = MatrixSpace(RationalField(),2,3)
172        sage: B = MS2.matrix([1,2,3,4,5,6])
173        sage: A*B
174        [ 9 12 15]
175        [19 26 33]
176   
177    ::
178   
179        sage: M = MatrixSpace(ZZ, 10)
180        sage: M
181        Full MatrixSpace of 10 by 10 dense matrices over Integer Ring
182        sage: loads(M.dumps()) == M
183        True
184    """
185    if not sage.rings.ring.is_Ring(base_ring):
186        raise TypeError, "base_ring (=%s) must be a ring"%base_ring
187   
188    if ncols is None: ncols = nrows
189    nrows = int(nrows); ncols = int(ncols); sparse=bool(sparse)
190    key = (base_ring, nrows, ncols, sparse)
191    if _cache.has_key(key):
192        M = _cache[key]()
193        if not M is None: return M
194
195    M = MatrixSpace_generic(base_ring, nrows, ncols, sparse)
196
197    _cache[key] = weakref.ref(M)
198    return M
199
200
201   
202class MatrixSpace_generic(parent_gens.ParentWithGens):
203    """
204    The space of all nrows x ncols matrices over base_ring.
205   
206    EXAMPLES::
207   
208        sage: MatrixSpace(ZZ,10,5)
209        Full MatrixSpace of 10 by 5 dense matrices over Integer Ring
210        sage: MatrixSpace(ZZ,10,2^31)
211        Traceback (most recent call last):                                   # 32-bit
212        ...                                                                  # 32-bit
213        ValueError: number of rows and columns must be less than 2^31 (on a 32-bit computer -- use a 64-bit computer for matrices with up to 2^63-1 rows and columns)           # 32-bit
214        Full MatrixSpace of 10 by 2147483648 dense matrices over Integer Ring   # 64-bit
215    """
216    def __init__(self,  base_ring,             
217                        nrows, 
218                        ncols=None, 
219                        sparse=False):
220        parent_gens.ParentWithGens.__init__(self, base_ring)
221        if not isinstance(base_ring, ring.Ring):
222            raise TypeError, "base_ring must be a ring"
223        if ncols == None: ncols = nrows
224        nrows = int(nrows)
225        ncols = int(ncols)
226        if nrows < 0:
227            raise ArithmeticError, "nrows must be nonnegative"
228        if ncols < 0:
229            raise ArithmeticError, "ncols must be nonnegative"
230
231        if nrows >= 2**63 or ncols >= 2**63:
232            raise ValueError, "number of rows and columns must be less than 2^63"
233        elif nrows >= 2**31 or ncols >= 2**31 and not sage.misc.misc.is_64_bit:
234            raise ValueError, "number of rows and columns must be less than 2^31 (on a 32-bit computer -- use a 64-bit computer for matrices with up to 2^63-1 rows and columns)"
235           
236        self.__nrows = nrows
237        self.__is_sparse = sparse
238        if ncols == None:
239            self.__ncols = nrows
240        else:
241            self.__ncols = ncols
242        self.__matrix_class = self._get_matrix_class()
243
244    def __reduce__(self):
245        """
246        EXAMPLES::
247       
248            sage: A = Mat(ZZ,5,7,sparse=True)
249            sage: A
250            Full MatrixSpace of 5 by 7 sparse matrices over Integer Ring
251            sage: loads(dumps(A)) == A
252            True
253        """
254        return MatrixSpace, (self.base_ring(), self.__nrows, self.__ncols, self.__is_sparse)
255
256    def __call__(self, entries=0, coerce=True, copy=True, rows=None):
257        """
258        EXAMPLES::
259       
260            sage: k = GF(7); G = MatrixGroup([matrix(k,2,[1,1,0,1]), matrix(k,2,[1,0,0,2])])
261            sage: g = G.0
262            sage: MatrixSpace(k,2)(g)
263            [1 1]
264            [0 1]
265       
266        ::
267       
268            sage: MS = MatrixSpace(ZZ,2,4)
269            sage: M2 = MS(range(8)); M2
270            [0 1 2 3]
271            [4 5 6 7]
272            sage: M2.columns()
273            [(0, 4), (1, 5), (2, 6), (3, 7)]
274            sage: MS(M2.columns())
275            [0 1 2 3]
276            [4 5 6 7]
277            sage: M2 == MS(M2.columns())
278            True
279            sage: M2 == MS(M2.rows())
280            True
281       
282        ::
283       
284            sage: MS = MatrixSpace(ZZ,2,4, sparse=True)
285            sage: M2 = MS(range(8)); M2
286            [0 1 2 3]
287            [4 5 6 7]
288            sage: M2.columns()
289            [(0, 4), (1, 5), (2, 6), (3, 7)]
290            sage: MS(M2.columns())
291            [0 1 2 3]
292            [4 5 6 7]
293            sage: M2 == MS(M2.columns())
294            True
295            sage: M2 == MS(M2.rows())
296            True
297       
298        ::
299       
300            sage: MS = MatrixSpace(ZZ,2,2)
301            sage: MS([1,2,3,4])
302            [1 2]
303            [3 4]
304            sage: MS([1,2,3,4], rows=True)
305            [1 2]
306            [3 4]
307            sage: MS([1,2,3,4], rows=False)
308            [1 3]
309            [2 4]
310       
311        ::
312       
313            sage: MS = MatrixSpace(ZZ,2,2, sparse=True)
314            sage: MS([1,2,3,4])
315            [1 2]
316            [3 4]
317            sage: MS([1,2,3,4], rows=True)
318            [1 2]
319            [3 4]
320            sage: MS([1,2,3,4], rows=False)
321            [1 3]
322            [2 4]
323
324            sage: MS = MatrixSpace(ZZ, 2)
325            sage: g = Gamma0(5)([1,1,0,1])
326            sage: MS(g)
327            [1 1]
328            [0 1]
329        """
330        if entries is None:
331            entries = 0
332
333        if entries == 0 and hasattr(self, '__zero_matrix'):
334            return self.zero_matrix()
335
336        if isinstance(entries, (list, tuple)) and len(entries) > 0 and \
337           sage.modules.free_module_element.is_FreeModuleElement(entries[0]):
338            #Try to determine whether or not the entries should
339            #be rows or columns
340            if rows is None:
341                #If the matrix is square, default to rows
342                if self.__ncols == self.__nrows:
343                    rows = True
344                elif len(entries[0]) == self.__ncols:
345                    rows = True
346                elif len(entries[0]) == self.__nrows:
347                    rows = False
348                else:
349                    raise ValueError, "incorrect dimensions"
350           
351            if self.__is_sparse:
352                e = {}
353                zero = self.base_ring()(0)
354                for i in xrange(len(entries)):
355                    for j, x in entries[i].iteritems():
356                        if x != zero:
357                            if rows:
358                                e[(i,j)] = x
359                            else:
360                                e[(j,i)] = x
361                entries = e
362            else:
363                entries = sum([v.list() for v in entries],[])
364
365        if rows is None:
366            rows = True
367
368        if not self.__is_sparse and isinstance(entries, dict):
369            entries = dict_to_list(entries, self.__nrows, self.__ncols)
370            coerce = True
371            copy = False
372        elif self.__is_sparse and isinstance(entries, (list, tuple)):
373            entries = list_to_dict(entries, self.__nrows, self.__ncols, rows=rows)
374            coerce = True
375            copy = False
376        elif sage.groups.matrix_gps.matrix_group_element.is_MatrixGroupElement(entries) \
377             or isinstance(entries, sage.modular.arithgroup.arithgroup_element.ArithmeticSubgroupElement):
378            return self(entries.matrix(), copy=False)
379
380        return self.matrix(entries, copy=copy, coerce=coerce, rows=rows)
381
382    def change_ring(self, R):
383        """
384        Return matrix space over R with otherwise same parameters as self.
385       
386        INPUT:
387       
388       
389        -  ``R`` - ring
390       
391       
392        OUTPUT: a matrix space
393       
394        EXAMPLES::
395       
396            sage: Mat(QQ,3,5).change_ring(GF(7))
397            Full MatrixSpace of 3 by 5 dense matrices over Finite Field of size 7
398        """
399        try:
400            return self.__change_ring[R]
401        except AttributeError:
402            self.__change_ring = {}
403        except KeyError:
404            pass
405        M = MatrixSpace(R, self.__nrows, self.__ncols, self.__is_sparse)
406        self.__change_ring[R] = M
407        return M
408
409    def base_extend(self, R):
410        """
411        Return base extension of this matrix space to R.
412       
413        INPUT:
414       
415       
416        -  ``R`` - ring
417       
418       
419        OUTPUT: a matrix space
420       
421        EXAMPLES::
422       
423            sage: Mat(ZZ,3,5).base_extend(QQ)
424            Full MatrixSpace of 3 by 5 dense matrices over Rational Field
425            sage: Mat(QQ,3,5).base_extend(GF(7))
426            Traceback (most recent call last):
427            ...
428            TypeError: no base extension defined
429        """
430        if R.has_coerce_map_from(self.base_ring()):
431            return self.change_ring(R)
432        raise TypeError, "no base extension defined"
433       
434    def construction(self):
435        """
436        EXAMPLES::
437       
438            sage: A = matrix(ZZ, 2, [1..4], sparse=True)
439            sage: A.parent().construction()
440            (MatrixFunctor, Integer Ring)
441            sage: A.parent().construction()[0](QQ['x'])
442            Full MatrixSpace of 2 by 2 sparse matrices over Univariate Polynomial Ring in x over Rational Field
443            sage: parent(A/2)
444            Full MatrixSpace of 2 by 2 sparse matrices over Rational Field
445        """
446        from sage.categories.pushout import MatrixFunctor
447        return MatrixFunctor(self.__nrows, self.__ncols, is_sparse=self.is_sparse()), self.base_ring()
448
449    def get_action_impl(self, S, op, self_on_left):
450        try:
451            if op is operator.mul:
452                import action as matrix_action
453                if self_on_left:
454                    if is_MatrixSpace(S):
455                        return matrix_action.MatrixMatrixAction(self, S)
456                    elif sage.modules.free_module.is_FreeModule(S):
457                        return matrix_action.MatrixVectorAction(self, S)
458                    else:
459                        # action of basering
460                        return sage.structure.coerce.RightModuleAction(S, self)
461                else:
462                    if sage.modules.free_module.is_FreeModule(S):
463                        return matrix_action.VectorMatrixAction(self, S)
464                    else:
465                        # action of basering
466                        return sage.structure.coerce.LeftModuleAction(S, self)
467            else:
468                return None
469        except TypeError:
470            return None
471
472    def _coerce_impl(self, x):
473        """
474        EXAMPLES::
475       
476            sage: MS1 = MatrixSpace(QQ,3)
477            sage: MS2 = MatrixSpace(ZZ,4,5,true)
478            sage: A = MS1(range(9))
479            sage: D = MS2(range(20))
480            sage: MS1._coerce_(A)
481            [0 1 2]
482            [3 4 5]
483            [6 7 8]
484            sage: MS2._coerce_(D)
485            [ 0  1  2  3  4]
486            [ 5  6  7  8  9]
487            [10 11 12 13 14]
488            [15 16 17 18 19]
489        """
490        if isinstance(x, matrix.Matrix):
491            if self.is_sparse() and x.is_dense():
492                raise TypeError, "cannot coerce dense matrix into sparse space for arithmetic"
493            if x.nrows() == self.nrows() and x.ncols() == self.ncols():
494                if self.base_ring().has_coerce_map_from(x.base_ring()):
495                    return self(x)
496                raise TypeError, "no canonical coercion"
497        return self._coerce_try(x, self.base_ring())
498
499    def __cmp__(self, other):
500        """
501        Compare this matrix space with other. Sparse and dense matrix
502        spaces with otherwise the same parameters are considered equal.
503       
504        If other is not a matrix space, return something arbitrary but
505        deterministic. Otherwise, compare based on base ring, then on
506        number of rows and columns.
507       
508        EXAMPLES::
509       
510            sage: Mat(ZZ,1000) == Mat(QQ,1000)
511            False
512            sage: Mat(ZZ,10) == Mat(ZZ,10)
513            True
514            sage: Mat(ZZ,10, sparse=False) == Mat(ZZ,10, sparse=True)
515            True
516        """
517        if isinstance(other, MatrixSpace_generic):
518            return cmp((self.base_ring(), self.__nrows, self.__ncols),
519                       (other.base_ring(), other.__nrows, other.__ncols))
520        return cmp(type(self), type(other))
521       
522    def _repr_(self):
523        """
524        Returns the string representation of a MatrixSpace
525       
526        EXAMPLES::
527       
528            sage: MS = MatrixSpace(ZZ,2,4,true)
529            sage: repr(MS)
530            'Full MatrixSpace of 2 by 4 sparse matrices over Integer Ring'
531            sage: MS
532            Full MatrixSpace of 2 by 4 sparse matrices over Integer Ring
533        """
534        if self.is_sparse():
535            s = "sparse"
536        else:
537            s = "dense"
538        return "Full MatrixSpace of %s by %s %s matrices over %s"%(
539                    self.__nrows, self.__ncols, s, self.base_ring())
540
541    def _latex_(self):
542        r"""
543        Returns the latex representation of a MatrixSpace
544       
545        EXAMPLES::
546       
547            sage: MS3 = MatrixSpace(QQ,6,6,true)
548            sage: latex(MS3)
549            \mathrm{Mat}_{6\times 6}(\Bold{Q})
550        """
551        return "\\mathrm{Mat}_{%s\\times %s}(%s)"%(self.nrows(), self.ncols(),
552                                                      latex.latex(self.base_ring()))
553
554    def __iter__(self):
555        r"""
556        Returns a generator object which iterates through the elements of
557        self. The order in which the elements are generated is based on a
558        'weight' of a matrix which is the number of iterations on the base
559        ring that are required to reach that matrix.
560       
561        The ordering is similar to a degree negative lexicographic order in
562        monomials in a multivariate polynomial ring.
563       
564        EXAMPLES: Consider the case of 2 x 2 matrices over GF(5).
565       
566        ::
567       
568            sage: list( GF(5) )
569            [0, 1, 2, 3, 4]
570            sage: MS = MatrixSpace(GF(5), 2, 2)
571            sage: l = list(MS)
572       
573        Then, consider the following matrices::
574       
575            sage: A = MS([2,1,0,1]); A
576            [2 1]
577            [0 1]
578            sage: B = MS([1,2,1,0]); B
579            [1 2]
580            [1 0]
581            sage: C = MS([1,2,0,0]); C
582            [1 2]
583            [0 0]
584       
585        A appears before B since the weight of one of A's entries exceeds
586        the weight of the corresponding entry in B earliest in the list.
587       
588        ::
589       
590            sage: l.index(A)
591            41
592            sage: l.index(B)
593            46
594       
595        However, A would come after the matrix C since C has a lower weight
596        than A.
597       
598        ::
599       
600            sage: l.index(A)
601            41
602            sage: l.index(C)
603            19
604       
605        The weights of matrices over other base rings are not as obvious.
606        For example, the weight of
607       
608        ::
609       
610            sage: MS = MatrixSpace(ZZ, 2, 2)
611            sage: MS([-1,0,0,0])
612            [-1  0]
613            [ 0  0]
614       
615        is 2 since
616       
617        ::
618       
619            sage: i = iter(ZZ)
620            sage: i.next()
621            0
622            sage: i.next()
623            1
624            sage: i.next()
625            -1
626       
627        Some more examples::
628       
629            sage: MS = MatrixSpace(GF(2),2)
630            sage: a = list(MS)
631            sage: len(a)
632            16
633            sage: for m in a: print m, '\n-'
634            [0 0]
635            [0 0]
636            -
637            [1 0]
638            [0 0]
639            -
640            [0 1]
641            [0 0]
642            -
643            [0 0]
644            [1 0]
645            -
646            [0 0]
647            [0 1]
648            -
649            [1 1]
650            [0 0]
651            -           
652            [1 0]
653            [1 0]
654            -
655            [1 0]
656            [0 1]
657            -
658            [0 1]
659            [1 0]
660            -
661            [0 1]
662            [0 1]
663            -
664            [0 0]
665            [1 1]
666            -
667            [1 1]
668            [1 0]
669            -
670            [1 1]
671            [0 1]
672            -
673            [1 0]
674            [1 1]
675            -
676            [0 1]
677            [1 1]
678            -
679            [1 1]
680            [1 1]
681            -
682       
683        ::
684       
685            sage: MS = MatrixSpace(GF(2),2, 3)
686            sage: a = list(MS)
687            sage: len(a)
688            64
689            sage: a[0]
690            [0 0 0]
691            [0 0 0]
692       
693        ::
694       
695            sage: MS = MatrixSpace(ZZ, 2, 3)
696            sage: i = iter(MS)
697            sage: a = [ i.next() for _ in range(6) ]
698            sage: a[0]
699            [0 0 0]
700            [0 0 0]
701            sage: a[4]
702            [0 0 0]
703            [1 0 0]
704       
705        For degenerate cases, where either the number of rows or columns
706        (or both) are zero, then the single element of the space is
707        returned.
708       
709        ::
710       
711            sage: list( MatrixSpace(GF(2), 2, 0) )
712            [[]]
713            sage: list( MatrixSpace(GF(2), 0, 2) )
714            [[]]
715            sage: list( MatrixSpace(GF(2), 0, 0) )
716            [[]]
717       
718        If the base ring does not support iteration (for example, with the
719        reals), then the matrix space over that ring does not support
720        iteration either.
721       
722        ::
723       
724            sage: MS = MatrixSpace(RR, 2)
725            sage: a = list(MS)
726            Traceback (most recent call last):
727            ...
728            NotImplementedError: object does not support iteration
729        """
730        #Make sure that we can interate over the base ring
731        base_ring = self.base_ring()
732        base_iter = iter(base_ring)
733
734        number_of_entries = (self.__nrows*self.__ncols)
735
736        #If the number of entries is zero, then just
737        #yield the empty matrix in that case and return
738        if number_of_entries == 0:
739            yield self(0)
740            raise StopIteration
741       
742        import sage.combinat.integer_vector
743
744        if not base_ring.is_finite():
745            #When the base ring is not finite, then we should go
746            #through and yield the matrices by "weight", which is
747            #the total number of iterations that need to be done
748            #on the base ring to reach the matrix.
749            base_elements = [ base_iter.next() ]
750            weight = 0
751            while True:
752                for iv in sage.combinat.integer_vector.IntegerVectors(weight, number_of_entries):
753                    yield self(entries=[base_elements[i] for i in iv], rows=True)       
754                       
755                weight += 1
756                base_elements.append( base_iter.next() )
757        else:
758            #In the finite case, we do a similar thing except that
759            #the "weight" of each entry is bounded by the number
760            #of elements in the base ring
761            order = base_ring.order()
762            done = False
763            base_elements = list(base_ring)
764            for weight in range((order-1)*number_of_entries+1):
765                for iv in sage.combinat.integer_vector.IntegerVectors(weight, number_of_entries, max_part=(order-1)):
766                   yield self(entries=[base_elements[i] for i in iv], rows=True)
767
768         
769    def _get_matrix_class(self):
770        """
771        Returns the class of self
772       
773        EXAMPLES::
774       
775            sage: MS1 = MatrixSpace(QQ,4)
776            sage: MS2 = MatrixSpace(ZZ,4,5,true)
777            sage: MS1._get_matrix_class()
778            <type 'sage.matrix.matrix_rational_dense.Matrix_rational_dense'>
779            sage: MS2._get_matrix_class()
780            <type 'sage.matrix.matrix_integer_sparse.Matrix_integer_sparse'>
781            sage: type(matrix(SR, 2, 2, 0))
782            <type 'sage.matrix.matrix_symbolic_dense.Matrix_symbolic_dense'>
783        """
784        R = self.base_ring()
785        if self.is_dense():
786            if sage.rings.integer_ring.is_IntegerRing(R):
787                return matrix_integer_dense.Matrix_integer_dense
788            elif sage.rings.rational_field.is_RationalField(R):
789                return matrix_rational_dense.Matrix_rational_dense
790            elif sage.rings.number_field.number_field.is_CyclotomicField(R):
791                import matrix_cyclo_dense
792                return matrix_cyclo_dense.Matrix_cyclo_dense
793            elif R==sage.rings.real_double.RDF:
794                import matrix_real_double_dense
795                return matrix_real_double_dense.Matrix_real_double_dense
796            elif R==sage.rings.complex_double.CDF:
797                import matrix_complex_double_dense
798                return matrix_complex_double_dense.Matrix_complex_double_dense
799            elif sage.rings.integer_mod_ring.is_IntegerModRing(R) and R.order() < matrix_modn_dense.MAX_MODULUS:
800                if R.order() == 2:
801                    return matrix_mod2_dense.Matrix_mod2_dense
802                return matrix_modn_dense.Matrix_modn_dense
803            elif sage.rings.polynomial.multi_polynomial_ring_generic.is_MPolynomialRing(R) and R.base_ring().is_field():
804                return matrix_mpolynomial_dense.Matrix_mpolynomial_dense
805            #elif isinstance(R, sage.rings.padics.padic_ring_capped_relative.pAdicRingCappedRelative):
806            #    return padics.matrix_padic_capped_relative_dense
807            # the default     
808            else:
809                from sage.symbolic.ring import SR   # causes circular imports
810                if R is SR:
811                    import matrix_symbolic_dense
812                    return matrix_symbolic_dense.Matrix_symbolic_dense
813                return matrix_generic_dense.Matrix_generic_dense
814       
815        else:           
816            if sage.rings.integer_mod_ring.is_IntegerModRing(R) and R.order() < matrix_modn_sparse.MAX_MODULUS:
817                return matrix_modn_sparse.Matrix_modn_sparse
818            elif sage.rings.rational_field.is_RationalField(R):
819                return matrix_rational_sparse.Matrix_rational_sparse
820            elif sage.rings.integer_ring.is_IntegerRing(R):
821                return matrix_integer_sparse.Matrix_integer_sparse
822            # the default
823            return matrix_generic_sparse.Matrix_generic_sparse
824           
825       
826    def basis(self):
827        """
828        Returns a basis for this matrix space.
829       
830        .. warning::
831
832           This will of course compute every generator of this matrix
833           space. So for large matrices, this could take a long time,
834           waste a massive amount of memory (for dense matrices), and
835           is likely not very useful. Don't use this on large matrix
836           spaces.
837       
838        EXAMPLES::
839       
840            sage: Mat(ZZ,2,2).basis()
841            [
842            [1 0]
843            [0 0],
844            [0 1]
845            [0 0],
846            [0 0]
847            [1 0],
848            [0 0]
849            [0 1]
850            ]
851        """
852        v = [self.zero_matrix() for _ in range(self.dimension())]
853        one = self.base_ring()(1)
854        i = 0
855        for r in range(self.__nrows):
856            for c in range(self.__ncols):
857                v[i][r,c] = one
858                v[i].set_immutable()
859                i += 1
860        return Sequence(v, universe=self, check=False, immutable=True, cr=True)
861
862    def dimension(self):
863        """
864        Returns (m rows) \* (n cols) of self as Integer
865       
866        EXAMPLES::
867       
868            sage: MS = MatrixSpace(ZZ,4,6)
869            sage: u = MS.dimension()
870            sage: u - 24 == 0
871            True
872        """
873        return self.__nrows * self.__ncols
874   
875    def dims(self):
876        """
877        Returns (m row, n col) representation of self dimension
878       
879        EXAMPLES::
880       
881            sage: MS = MatrixSpace(ZZ,4,6)
882            sage: MS.dims()
883            (4, 6)
884        """
885        return (self.__nrows, self.__ncols)
886       
887    def identity_matrix(self):
888        """
889        Create an identity matrix in self. (Must be a space of square
890        matrices).
891       
892        EXAMPLES::
893       
894            sage: MS1 = MatrixSpace(ZZ,4)
895            sage: MS2 = MatrixSpace(QQ,3,4)
896            sage: I = MS1.identity_matrix()
897            sage: I
898            [1 0 0 0]
899            [0 1 0 0]
900            [0 0 1 0]
901            [0 0 0 1]
902            sage: Er = MS2.identity_matrix()
903            Traceback (most recent call last):
904            ...
905            TypeError: self must be a space of square matrices
906        """
907        if self.__nrows != self.__ncols:
908            raise TypeError, "self must be a space of square matrices"
909        A = self(0)
910        for i in xrange(self.__nrows):
911            A[i,i] = 1
912        return A
913       
914    def is_dense(self):
915        """
916        Returns True if matrices in self are dense and False otherwise.
917       
918        EXAMPLES::
919       
920            sage: Mat(RDF,2,3).is_sparse()
921            False
922            sage: Mat(RR,123456,22,sparse=True).is_sparse()
923            True
924        """
925        return not self.__is_sparse
926       
927    def is_sparse(self):
928        """
929        Returns True if matrices in self are sparse and False otherwise.
930       
931        EXAMPLES::
932       
933            sage: Mat(GF(2011),10000).is_sparse()
934            False
935            sage: Mat(GF(2011),10000,sparse=True).is_sparse()
936            True
937        """
938        return self.__is_sparse
939       
940    def is_finite(self):
941        """
942        EXAMPLES::
943       
944            sage: MatrixSpace(GF(101), 10000).is_finite()
945            True
946            sage: MatrixSpace(QQ, 2).is_finite()
947            False
948        """
949        return self.base_ring().is_finite()
950
951    def gen(self, n):
952        """
953        Return the n-th generator of this matrix space.
954       
955        This doesn't compute all basis matrices, so it is reasonably
956        intelligent.
957       
958        EXAMPLES::
959       
960            sage: M = Mat(GF(7),10000,5); M.ngens()
961            50000
962            sage: a = M.10
963            sage: a[:4]
964            [0 0 0 0 0]
965            [0 0 0 0 0]
966            [1 0 0 0 0]
967            [0 0 0 0 0]
968        """
969        if hasattr(self, '__basis'):
970            return self.__basis[n]
971        r = n // self.__ncols
972        c = n - (r * self.__ncols)
973        z = self.zero_matrix()
974        z[r,c] = 1
975        return z
976
977    def zero_matrix(self):
978        """
979        Return the zero matrix.
980        """
981        try:
982            z = self.__zero_matrix
983        except AttributeError:
984            z = self(0)
985            self.__zero_matrix = z
986        return z.__copy__()
987   
988    def ngens(self):
989        """
990        Return the number of generators of this matrix space, which is the
991        number of entries in the matrices in this space.
992       
993        EXAMPLES::
994       
995            sage: M = Mat(GF(7),100,200); M.ngens()
996            20000
997        """
998        return self.dimension()
999 
1000    def matrix(self, x=0, coerce=True, copy=True, rows=True):
1001        """
1002        Create a matrix in self. The entries can be specified either as a
1003        single list of length nrows\*ncols, or as a list of lists.
1004       
1005        EXAMPLES::
1006       
1007            sage: M = MatrixSpace(ZZ, 2)
1008            sage: M.matrix([[1,0],[0,-1]])
1009            [ 1  0]
1010            [ 0 -1]
1011            sage: M.matrix([1,0,0,-1])
1012            [ 1  0]
1013            [ 0 -1]
1014            sage: M.matrix([1,2,3,4])
1015            [1 2]
1016            [3 4]
1017            sage: M.matrix([1,2,3,4],rows=False)
1018            [1 3]
1019            [2 4]
1020        """
1021        if isinstance(x, (types.GeneratorType, xrange)):
1022            x = list(x)
1023        elif isinstance(x, (int, integer.Integer)) and x==1:
1024            return self.identity_matrix()
1025        if matrix.is_Matrix(x):
1026            if x.parent() is self:
1027                if x.is_immutable():
1028                    return x
1029                else:
1030                    return x.copy()
1031            x = x.list()
1032        if isinstance(x, list) and len(x) > 0:
1033            if isinstance(x[0], list):
1034                x = sum(x,[])
1035            elif hasattr(x[0], "is_vector"): # TODO: is this the best way to test that?
1036                e = []
1037                for v in x:
1038                    e = e + v.list()
1039                copy = False # deep copy?
1040                x = e
1041            elif isinstance(x[0], tuple):
1042                x = list(sum(x,()))
1043
1044            if not rows:
1045                new_x = []
1046                for k in range(len(x)):
1047                    i = k % self.__ncols
1048                    j = k // self.__ncols
1049                    new_x.append( x[ i*self.__nrows + j ] )
1050                x = new_x
1051           
1052        return self.__matrix_class(self, entries=x, copy=copy, coerce=coerce) 
1053     
1054    def matrix_space(self, nrows=None, ncols=None, sparse=False):
1055        """
1056        Return the matrix space with given number of rows, columns and
1057        sparcity over the same base ring as self, and defaults the same as
1058        self.
1059       
1060        EXAMPLES::
1061       
1062            sage: M = Mat(GF(7),100,200)
1063            sage: M.matrix_space(5000)
1064            Full MatrixSpace of 5000 by 200 dense matrices over Finite Field of size 7
1065            sage: M.matrix_space(ncols=5000)
1066            Full MatrixSpace of 100 by 5000 dense matrices over Finite Field of size 7
1067            sage: M.matrix_space(sparse=True)
1068            Full MatrixSpace of 100 by 200 sparse matrices over Finite Field of size 7
1069        """
1070        if nrows is None:
1071            nrows = self.__nrows
1072        if ncols is None:
1073            ncols = self.__ncols
1074        return MatrixSpace(self.base_ring(), nrows, ncols, 
1075                        sparse=sparse)
1076 
1077    def ncols(self):
1078        """
1079        Return the number of columns of matrices in this space.
1080       
1081        EXAMPLES::
1082       
1083            sage: M = Mat(ZZ['x'],200000,500000,sparse=True)
1084            sage: M.ncols()
1085            500000
1086        """
1087        return self.__ncols
1088       
1089    def nrows(self):
1090        """
1091        Return the number of rows of matrices in this space.
1092       
1093        EXAMPLES::
1094       
1095            sage: M = Mat(ZZ,200000,500000)
1096            sage: M.nrows()
1097            200000
1098        """
1099        return self.__nrows
1100
1101    def row_space(self):
1102        """
1103        Return the module spanned by all rows of matrices in this matrix
1104        space. This is a free module of rank the number of rows. It will be
1105        sparse or dense as this matrix space is sparse or dense.
1106       
1107        EXAMPLES::
1108       
1109            sage: M = Mat(ZZ,20,5,sparse=False); M.row_space()
1110            Ambient free module of rank 5 over the principal ideal domain Integer Ring
1111        """
1112        try:
1113            return self.__row_space
1114        except AttributeError:
1115            self.__row_space = sage.modules.free_module.FreeModule(self.base_ring(),
1116                                                self.ncols(), sparse=self.is_sparse())
1117            return self.__row_space
1118   
1119    def column_space(self):
1120        """
1121        Return the module spanned by all columns of matrices in this matrix
1122        space. This is a free module of rank the number of columns. It will
1123        be sparse or dense as this matrix space is sparse or dense.
1124       
1125        EXAMPLES::
1126       
1127            sage: M = Mat(GF(9,'a'),20,5,sparse=True); M.column_space()
1128            Sparse vector space of dimension 20 over Finite Field in a of size 3^2
1129        """
1130        try:
1131            return self.__column_space
1132        except AttributeError:
1133            self.__column_space = sage.modules.free_module.FreeModule(self.base_ring(), self.nrows(),
1134                                                                   sparse=self.is_sparse())
1135            return self.__column_space
1136
1137    def random_element(self, density=1, *args, **kwds):
1138        """
1139        INPUT:
1140       
1141       
1142        -  ``density`` - integer (default: 1) rough measure of
1143           the proportion of nonzero entries in the random matrix
1144       
1145        -  ``*args, **kwds`` - rest of parameters may be
1146           passed to the random_element function of the base ring. ("may be",
1147           since this function calls the randomize function on the zero
1148           matrix, which need not call the random_element function of the
1149           base ring at all in general.)
1150       
1151       
1152        EXAMPLES::
1153       
1154            sage: Mat(ZZ,2,5).random_element()
1155            [ -8   2   0   0   1]
1156            [ -1   2   1 -95  -1]
1157            sage: Mat(QQ,2,5).random_element(density=0.5)
1158            [  2   0   0   0   1]
1159            [  0   0   0 1/2   0]
1160            sage: Mat(QQ,3,sparse=True).random_element()
1161            [  -1  1/3    1]
1162            [   0   -1    0]
1163            [  -1    1 -1/4]
1164            sage: Mat(GF(9,'a'),3,sparse=True).random_element()
1165            [      1       2       1]
1166            [2*a + 1       a       2]
1167            [      2 2*a + 2       1]
1168        """
1169        Z = self.zero_matrix()
1170        Z.randomize(density, *args, **kwds)
1171        return Z
1172
1173    def _magma_init_(self, magma):
1174        r"""
1175        EXAMPLES: We first coerce a square matrix.
1176       
1177        ::
1178       
1179            sage: magma(MatrixSpace(QQ,3))                      # optional - magma
1180            Full Matrix Algebra of degree 3 over Rational Field
1181       
1182        ::
1183       
1184            sage: magma(MatrixSpace(Integers(8),2,3))           # optional - magma
1185            Full RMatrixSpace of 2 by 3 matrices over IntegerRing(8)
1186        """
1187        K = self.base_ring()._magma_init_(magma)
1188        if self.__nrows == self.__ncols:
1189            s = 'MatrixAlgebra(%s,%s)'%(K, self.__nrows)
1190        else:
1191            s = 'RMatrixSpace(%s,%s,%s)'%(K, self.__nrows, self.__ncols)
1192        return s
1193
1194def dict_to_list(entries, nrows, ncols):
1195    """
1196    Given a dictionary of coordinate tuples, return the list given by
1197    reading off the nrows\*ncols matrix in row order.
1198   
1199    EXAMLES::
1200   
1201        sage: from sage.matrix.matrix_space import dict_to_list
1202        sage: d = {}
1203        sage: d[(0,0)] = 1
1204        sage: d[(1,1)] = 2
1205        sage: dict_to_list(d, 2, 2)
1206        [1, 0, 0, 2]
1207        sage: dict_to_list(d, 2, 3)
1208        [1, 0, 0, 0, 2, 0]
1209    """
1210    v = [0]*(nrows*ncols)
1211    for ij, y in entries.iteritems():
1212        i,j = ij
1213        v[i*ncols + j] = y
1214    return v
1215
1216def list_to_dict(entries, nrows, ncols, rows=True):
1217    """
1218    Given a list of entries, create a dictionary whose keys are
1219    coordinate tuples and values are the entries.
1220   
1221    EXAMPLES::
1222   
1223        sage: from sage.matrix.matrix_space import list_to_dict
1224        sage: d = list_to_dict([1,2,3,4],2,2)
1225        sage: d[(0,1)]
1226        2
1227        sage: d = list_to_dict([1,2,3,4],2,2,rows=False)
1228        sage: d[(0,1)]
1229        3
1230    """
1231    d = {}
1232    if ncols == 0 or nrows == 0:
1233        return d
1234    for i in range(len(entries)):
1235        x = entries[i]
1236        if x != 0:
1237            col = i % ncols
1238            row = i // ncols
1239            if rows:
1240                d[(row,col)] = x
1241            else:
1242                d[(col,row)] = x
1243    return d
1244
1245
1246
1247# sage: m = matrix(F, 0,0, sparse=False)
1248# sage: m.determinant()
1249# 0
1250
1251
1252def test_trivial_matrices_inverse(ring, sparse=True, checkrank=True):
1253    """
1254    Tests inversion, determinant and is_inverstible for trivial matrices.
1255   
1256    This function is a helper to check that the inversion of trivial matrices
1257    (of size 0x0, nx0, 0xn or 1x1) is handled consistently by the various
1258    implementation of matrices. The coherency is checked through a bunch of
1259    assertions. If an inconsistency is found, an AssertionError is raised
1260    which should make clear what is the problem.
1261
1262    INPUT:
1263
1264    - ``ring`` - a ring
1265
1266    - ``sparse`` - a boolean
1267
1268    - ``checkrank`` - a boolean
1269
1270    OUTPUT:
1271
1272    - nothing if everything is correct otherwise raise an AssertionError
1273   
1274    The methods determinant, is_invertible, rank and inverse are checked for
1275     - the 0x0 empty identity matrix
1276     - the 0x3 and 3x0 matrices
1277     - the 1x1 null matrix [0]
1278     - the 1x1 identity matrix [1]
1279
1280    If ``checkrank`` is ``False`` then the rank is not checked. This is used
1281    the check matrix over ring where echelon form is not implemented.
1282     
1283    TODO: must be adapted to category check framework when ready (see trac \#5274).
1284
1285    TESTS::
1286   
1287      sage: from sage.matrix.matrix_space import test_trivial_matrices_inverse as tinv
1288      sage: tinv(ZZ, sparse=True)
1289      sage: tinv(ZZ, sparse=False)
1290      sage: tinv(QQ, sparse=True)
1291      sage: tinv(QQ, sparse=False)
1292      sage: tinv(GF(11), sparse=True)
1293      sage: tinv(GF(11), sparse=False)
1294      sage: tinv(GF(2), sparse=True)
1295      sage: tinv(GF(2), sparse=False)
1296      sage: tinv(SR, sparse=True)
1297      sage: tinv(SR, sparse=False)
1298      sage: tinv(RDF, sparse=True)
1299      sage: tinv(RDF, sparse=False)
1300      sage: tinv(CDF, sparse=True)
1301      sage: tinv(CDF, sparse=False)
1302      sage: tinv(CyclotomicField(7), sparse=True)
1303      sage: tinv(CyclotomicField(7), sparse=False)
1304
1305   TODO: As soon as rank of sparse matrix over QQ['x,y'] is implemented,
1306   please remove the following test and the ``checkrank=False`` in the next one:
1307   
1308      sage: MatrixSpace(QQ['x,y'], 3, 3, sparse=True)(1).rank()
1309      Traceback (most recent call last):
1310      ...
1311      NotImplementedError: echelon form over Multivariate Polynomial Ring in x, y over Rational Field not yet implemented
1312      sage: tinv(QQ['x,y'], sparse=True, checkrank=False)
1313
1314     
1315   TODO: As soon as rank of dense matrix over QQ['x,y'] is implemented,
1316   please remove the following test and the ``checkrank=False`` in the next one:
1317     
1318      sage: MatrixSpace(QQ['x,y'], 3, 3, sparse=False)(1).rank()
1319      Traceback (most recent call last):
1320      ...
1321      RuntimeError: BUG: matrix pivots should have been set but weren't, matrix parent = 'Full MatrixSpace of 3 by 3 dense matrices over Multivariate Polynomial Ring in x, y over Rational Field'
1322     
1323      sage: tinv(QQ['x,y'], sparse=False, checkrank=False)
1324   
1325    """
1326    # Check that the empty 0x0 matrix is it's own inverse with det=1.
1327    ms00 = MatrixSpace(ring, 0, 0, sparse=sparse)
1328    m00  = ms00(0)
1329    assert(m00.determinant() == ring(1))
1330    assert(m00.is_invertible())
1331    assert(m00.inverse() == m00)
1332    if checkrank:
1333        assert(m00.rank() == 0)
1334
1335    # Check that the empty 0x3 and 3x0 matrices are not invertible and that
1336    # computing the detemininant raise the proper exception.
1337    for ms0 in [MatrixSpace(ring, 0, 3, sparse=sparse),
1338                MatrixSpace(ring, 3, 0, sparse=sparse)]:
1339        mn0  = ms0(0)
1340        assert(not mn0.is_invertible())
1341        try:
1342            d = mn0.determinant()
1343            print d
1344            res = False
1345        except ArithmeticError:
1346            res = True
1347        assert(res)
1348        try:
1349            mn0.inverse()
1350            res =False
1351        except ArithmeticError:
1352            res = True
1353        assert(res)
1354        if checkrank:
1355            assert(mn0.rank() == 0)
1356
1357    # Check that the null 1x1 matrix is not invertible and that det=0
1358    ms1 = MatrixSpace(ring, 1, 1, sparse=sparse)
1359    m0  = ms1(0)
1360    assert(not m0.is_invertible())
1361    assert(m0.determinant() == ring(0))
1362    try:
1363        m0.inverse()
1364        res = False
1365    except (ZeroDivisionError, RuntimeError):
1366        #FIXME: Make pynac throw a ZeroDivisionError on division by
1367        #zero instead of a runtime Error
1368        res = True
1369    assert(res)
1370    if checkrank:
1371        assert(m0.rank() == 0)
1372
1373    # Check that the identity 1x1 matrix is its own inverse with det=1
1374    m1  = ms1(1)
1375    assert(m1.is_invertible())
1376    assert(m1.determinant() == ring(1))
1377    inv = m1.inverse()
1378    assert(inv == m1)
1379    if checkrank:
1380        assert(m1.rank() == 1)
Note: See TracBrowser for help on using the repository browser.