Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/matrix/action.pyx
4056 views
1
"""
2
These are the actions used by the coercion model for matrix and vector multiplications.
3
4
AUTHORS:
5
-- Robert Bradshaw (2007-09): Initial version.
6
"""
7
8
#*****************************************************************************
9
# Copyright (C) 2007 Robert Bradshaw <[email protected]>
10
#
11
# Distributed under the terms of the GNU General Public License (GPL)
12
#
13
# http://www.gnu.org/licenses/
14
#*****************************************************************************
15
16
17
import operator
18
19
from matrix_space import MatrixSpace, is_MatrixSpace
20
from sage.modules.free_module import FreeModule, is_FreeModule
21
22
23
cdef class MatrixMulAction(Action):
24
def __init__(self, G, S, is_left):
25
if not is_MatrixSpace(G):
26
raise TypeError, "Not a matrix space: %s" % G
27
if G.base_ring() is not S.base_ring():
28
from sage.categories.pushout import pushout
29
base = pushout(G.base_ring(), S.base_ring())
30
else:
31
base = G.base_ring()
32
Action.__init__(self, G, S, is_left, operator.mul)
33
self._codomain = self._create_codomain(base)
34
self.fix_sparseness = G.is_sparse() != S.is_sparse()
35
36
def codomain(self):
37
return self._codomain
38
39
def domain(self):
40
"""
41
EXAMPLES:
42
sage: A = MatrixSpace(QQ, 2).get_action(MatrixSpace(ZZ['x'], 2)); A
43
Left action by Full MatrixSpace of 2 by 2 dense matrices over Rational Field on Full MatrixSpace of 2 by 2 dense matrices over Univariate Polynomial Ring in x over Integer Ring
44
sage: A.actor()
45
Full MatrixSpace of 2 by 2 dense matrices over Rational Field
46
sage: A.domain()
47
Full MatrixSpace of 2 by 2 dense matrices over Univariate Polynomial Ring in x over Integer Ring
48
sage: A.codomain()
49
Full MatrixSpace of 2 by 2 dense matrices over Univariate Polynomial Ring in x over Rational Field
50
"""
51
return self.S
52
53
54
cdef class MatrixMatrixAction(MatrixMulAction):
55
def __init__(self, G, S):
56
"""
57
EXAMPLES:
58
sage: R.<x> = ZZ[]
59
sage: from sage.matrix.action import MatrixMatrixAction
60
sage: A = MatrixMatrixAction(MatrixSpace(R, 3, 3), MatrixSpace(QQ, 3, 2)); A
61
Left action by Full MatrixSpace of 3 by 3 dense matrices over Univariate Polynomial Ring in x over Integer Ring on Full MatrixSpace of 3 by 2 dense matrices over Rational Field
62
sage: A.codomain()
63
Full MatrixSpace of 3 by 2 dense matrices over Univariate Polynomial Ring in x over Rational Field
64
sage: A(matrix(R, 3, 3, x), matrix(QQ, 3, 2, range(6)))
65
[ 0 x]
66
[2*x 3*x]
67
[4*x 5*x]
68
"""
69
if not is_MatrixSpace(S):
70
raise TypeError, "Not a matrix space: %s" % S
71
MatrixMulAction.__init__(self, G, S, True)
72
73
def _create_codomain(self, base):
74
"""
75
EXAMPLES:
76
sage: from sage.matrix.action import MatrixMatrixAction
77
sage: R.<x> = ZZ[]
78
sage: A = MatrixMatrixAction(MatrixSpace(R, 3, 3), MatrixSpace(QQ, 3, 2)); A
79
Left action by Full MatrixSpace of 3 by 3 dense matrices over Univariate Polynomial Ring in x over Integer Ring on Full MatrixSpace of 3 by 2 dense matrices over Rational Field
80
sage: A.codomain()
81
Full MatrixSpace of 3 by 2 dense matrices over Univariate Polynomial Ring in x over Rational Field
82
"""
83
if self.G.ncols() != self.S.nrows():
84
raise TypeError, "incompatible dimensions %s, %s" % (self.G.ncols(), self.S.nrows())
85
return MatrixSpace(base, self.G.nrows(), self.S.ncols(), sparse = self.G.is_sparse() and self.S.is_sparse())
86
87
cpdef _call_(self, g, s):
88
"""
89
EXAMPLES:
90
Respects compatible subdivisions:
91
sage: M = matrix(5, 5, prime_range(100))
92
sage: M.subdivide(2,3); M
93
[ 2 3 5| 7 11]
94
[13 17 19|23 29]
95
[--------+-----]
96
[31 37 41|43 47]
97
[53 59 61|67 71]
98
[73 79 83|89 97]
99
sage: N = matrix(5,2,[n^2 for n in range(10)])
100
sage: N.subdivide(3,1); N
101
[ 0| 1]
102
[ 4| 9]
103
[16|25]
104
[--+--]
105
[36|49]
106
[64|81]
107
sage: M*N
108
[ 1048| 1388]
109
[ 3056| 4117]
110
[-----+-----]
111
[ 5360| 7303]
112
[ 8168|11143]
113
[11056|15077]
114
115
Note that this is just like block matrix multiplication:
116
sage: M.subdivision(0,0) * N.subdivision(0,0) + M.subdivision(0,1) * N.subdivision(1,0)
117
[1048]
118
[3056]
119
120
If the subdivisions aren't compatible, ignore them.
121
sage: N.subdivide(1,1); N
122
[ 0| 1]
123
[--+--]
124
[ 4| 9]
125
[16|25]
126
[36|49]
127
[64|81]
128
sage: M*N
129
[ 1048 1388]
130
[ 3056 4117]
131
[ 5360 7303]
132
[ 8168 11143]
133
[11056 15077]
134
135
"""
136
cdef Matrix A = g #<Matrix>g
137
cdef Matrix B = s #<Matrix>s
138
if A._parent._base is not self._codomain._base:
139
A = A.change_ring(self._codomain._base)
140
if B._parent._base is not self._codomain._base:
141
B = B.change_ring(self._codomain._base)
142
if self.fix_sparseness:
143
if B.is_sparse_c():
144
B = B.dense_matrix()
145
else:
146
A = A.dense_matrix()
147
prod = A._matrix_times_matrix_(B)
148
if A._subdivisions is not None or B._subdivisions is not None:
149
Asubs = A.subdivisions()
150
Bsubs = B.subdivisions()
151
if Asubs[1] == Bsubs[0]:
152
prod.subdivide(Asubs[0], Bsubs[1])
153
return prod
154
155
156
cdef class MatrixVectorAction(MatrixMulAction):
157
def __init__(self, G, S):
158
"""
159
EXAMPLES:
160
sage: from sage.matrix.action import MatrixVectorAction
161
sage: A = MatrixVectorAction(MatrixSpace(QQ, 3, 3), VectorSpace(CDF, 4)); A
162
Traceback (most recent call last):
163
...
164
TypeError: incompatible dimensions 3, 4
165
"""
166
if not is_FreeModule(S):
167
raise TypeError, "Not a free module: %s" % S
168
MatrixMulAction.__init__(self, G, S, True)
169
170
def _create_codomain(self, base):
171
"""
172
EXAMPLES:
173
sage: from sage.matrix.action import MatrixVectorAction
174
sage: A = MatrixVectorAction(MatrixSpace(QQ, 5, 3), VectorSpace(CDF, 3)); A
175
Left action by Full MatrixSpace of 5 by 3 dense matrices over Rational Field on Vector space of dimension 3 over Complex Double Field
176
sage: A.codomain()
177
Vector space of dimension 5 over Complex Double Field
178
"""
179
if self.G.ncols() != self.S.degree():
180
raise TypeError, "incompatible dimensions %s, %s" % (self.G.ncols(), self.S.degree())
181
return FreeModule(base, self.G.nrows(), sparse = self.G.is_sparse())
182
183
cpdef _call_(self, g, s):
184
cdef Matrix A = g #<Matrix>g
185
cdef Vector v = s #<Vector>s
186
if A._parent._base is not self._codomain._base:
187
A = A.change_ring(self._codomain._base)
188
if v._parent._base is not self._codomain._base:
189
v = v.change_ring(self._codomain._base)
190
if self.fix_sparseness:
191
if A.is_sparse_c():
192
v = v.sparse_vector()
193
else:
194
v = v.dense_vector()
195
return A._matrix_times_vector_(v)
196
197
198
cdef class VectorMatrixAction(MatrixMulAction):
199
def __init__(self, G, S):
200
"""
201
EXAMPLES:
202
sage: from sage.matrix.action import VectorMatrixAction
203
sage: A = VectorMatrixAction(MatrixSpace(QQ, 5, 3), VectorSpace(CDF, 3)); A
204
Traceback (most recent call last):
205
...
206
TypeError: incompatible dimensions 5, 3
207
"""
208
if not is_FreeModule(S):
209
raise TypeError, "Not a free module: %s" % S
210
MatrixMulAction.__init__(self, G, S, False)
211
212
def _create_codomain(self, base):
213
"""
214
EXAMPLES:
215
sage: from sage.matrix.action import VectorMatrixAction
216
sage: A = VectorMatrixAction(MatrixSpace(QQ, 3, 5), VectorSpace(CDF, 3)); A
217
Right action by Full MatrixSpace of 3 by 5 dense matrices over Rational Field on Vector space of dimension 3 over Complex Double Field
218
sage: A.codomain()
219
Vector space of dimension 5 over Complex Double Field
220
"""
221
if self.G.nrows() != self.S.degree():
222
raise TypeError, "incompatible dimensions %s, %s" % (self.G.nrows(), self.S.degree())
223
return FreeModule(base, self.G.ncols(), sparse = self.G.is_sparse())
224
225
cpdef _call_(self, s, g):
226
cdef Matrix A = g #<Matrix>g
227
cdef Vector v = s #<Vector>s
228
if A._parent._base is not self._codomain._base:
229
A = A.change_ring(self._codomain._base)
230
if v._parent._base is not self._codomain._base:
231
v = v.change_ring(self._codomain._base)
232
if self.fix_sparseness:
233
if A.is_sparse_c():
234
v = v.sparse_vector()
235
else:
236
v = v.dense_vector()
237
return (<Matrix>A)._vector_times_matrix_(v) # v * A
238
239
240