Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

381

382

383

384

385

386

387

388

389

390

391

392

393

394

395

396

397

398

399

400

401

402

""" 

Actions used by the coercion model for matrix and vector multiplications 

  

.. WARNING:: 

  

The class :class:`MatrixMulAction` and its descendants extends the class 

:class:`Action`. As a consequence objects from these classes only keep weak 

references to the underlying sets which are acted upon. This decision was 

made in :trac:`715` in order to allow garbage collection within the coercion 

framework, where actions are mainly used, and avoid memory leaks. 

  

To ensure that the underlying set of such an object does not get garbage 

collected, it is sufficient to explicitly create a strong reference to it 

before creating the action. 

  

:: 

  

sage: MSQ = MatrixSpace(QQ, 2) 

sage: MSZ = MatrixSpace(ZZ['x'], 2) 

sage: A = MSQ.get_action(MSZ) 

sage: A 

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 

sage: import gc 

sage: _ = gc.collect() 

sage: A 

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 

  

.. NOTE:: 

  

The :func:`MatrixSpace` function caches the objects it creates. Therefore, 

the underlying set ``MSZ`` in the above example will not be garbage 

collected, even if it is not strongly ref'ed. Nonetheless, there is no 

guarantee that the set that is acted upon will always be cached in such a 

way, so that following the above example is good practice. 

  

EXAMPLES: 

  

An action requires a common parent for the base rings, so the following 

doesn't work (see :trac:`17859`):: 

  

sage: vector(QQ, [1]) * matrix(Zmod(2), [[1]]) 

Traceback (most recent call last): 

... 

TypeError: unsupported operand parent(s) for *: 'Vector space of 

dimension 1 over Rational Field' and 'Full MatrixSpace of 1 by 1 

dense matrices over Ring of integers modulo 2' 

  

AUTHOR: 

  

- Robert Bradshaw (2007-09): Initial version. 

""" 

  

#***************************************************************************** 

# Copyright (C) 2007 Robert Bradshaw <robertwb@math.washington.edu> 

# 

# This program is free software: you can redistribute it and/or modify 

# it under the terms of the GNU General Public License as published by 

# the Free Software Foundation, either version 2 of the License, or 

# (at your option) any later version. 

# http://www.gnu.org/licenses/ 

#***************************************************************************** 

from __future__ import absolute_import 

  

import operator 

  

from .matrix_space import MatrixSpace, is_MatrixSpace 

from sage.modules.free_module import FreeModule, is_FreeModule 

from sage.structure.element cimport coercion_model 

  

  

cdef class MatrixMulAction(Action): 

def __init__(self, G, S, is_left): 

if not is_MatrixSpace(G): 

raise TypeError("Not a matrix space: %s" % G) 

if G.base_ring() is not S.base_ring(): 

base = coercion_model.common_parent(G.base_ring(), S.base_ring()) 

else: 

base = G.base_ring() 

Action.__init__(self, G, S, is_left, operator.mul) 

self._codomain = self._create_codomain(base) 

self.fix_sparseness = G.is_sparse() != S.is_sparse() 

  

def codomain(self): 

return self._codomain 

  

def domain(self): 

""" 

EXAMPLES: 

  

By :trac:`715`, there only is a weak reference on the underlying set, 

so that it can be garbage collected if only the action itself is 

explicitly referred to. Hence, we first assign the involved matrix 

spaces to a variable:: 

  

sage: MSQ = MatrixSpace(QQ, 2) 

sage: MSZ = MatrixSpace(ZZ['x'], 2) 

sage: A = MSQ.get_action(MSZ); A 

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 

sage: A.actor() 

Full MatrixSpace of 2 by 2 dense matrices over Rational Field 

sage: A.domain() 

Full MatrixSpace of 2 by 2 dense matrices over Univariate Polynomial Ring in x over Integer Ring 

sage: A.codomain() 

Full MatrixSpace of 2 by 2 dense matrices over Univariate Polynomial Ring in x over Rational Field 

  

.. NOTE:: 

  

The :func:`MatrixSpace` function caches the object it creates. 

Therefore, the underlying set ``MSZ`` in the above example will not 

be garbage collected, even if it is not strongly ref'ed. 

Nonetheless, there is no guarantee that the set that is acted upon 

will always be cached in such a way, so that following the above 

example is good practice. 

  

""" 

return self.underlying_set() 

  

  

cdef class MatrixMatrixAction(MatrixMulAction): 

""" 

Action of a matrix on another matrix. 

  

EXAMPLES: 

  

By :trac:`715`, there only is a weak reference on the underlying set, 

so that it can be garbage collected if only the action itself is 

explicitly referred to. Hence, we first assign the involved matrix 

spaces to a variable:: 

  

sage: R.<x> = ZZ[] 

sage: MSR = MatrixSpace(R, 3, 3) 

sage: MSQ = MatrixSpace(QQ, 3, 2) 

sage: from sage.matrix.action import MatrixMatrixAction 

sage: A = MatrixMatrixAction(MSR, MSQ); A 

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 

sage: A.codomain() 

Full MatrixSpace of 3 by 2 dense matrices over Univariate Polynomial Ring in x over Rational Field 

sage: A(matrix(R, 3, 3, x), matrix(QQ, 3, 2, range(6))) 

[ 0 x] 

[2*x 3*x] 

[4*x 5*x] 

  

.. NOTE:: 

  

The :func:`MatrixSpace` function caches the object it creates. 

Therefore, the underlying set ``MSZ`` in the above example will not 

be garbage collected, even if it is not strongly ref'ed. 

Nonetheless, there is no guarantee that the set that is acted upon 

will always be cached in such a way, so that following the above 

example is good practice. 

""" 

def __init__(self, G, S): 

""" 

TESTS: 

  

Check that multiplication for matrices with different backends are not allowed:: 

  

sage: M1 = MatrixSpace(ZZ, 2, implementation='flint') 

sage: M2 = MatrixSpace(ZZ, 2, implementation='generic') 

sage: M3 = MatrixSpace(ZZ, 2, implementation='gap') 

sage: M4 = MatrixSpace(ZZ, 2, sparse=True) 

sage: M = [M1, M2, M3, M4] 

  

sage: coercions = '' 

sage: for M1 in M: 

....: for M2 in M: 

....: try: 

....: s = M1.an_element() * M2.an_element() 

....: coercions += 'X' 

....: except TypeError: 

....: coercions += ' ' 

....: coercions += '\n' 

sage: print(coercions) 

X X 

X 

X 

X X 

""" 

if not is_MatrixSpace(S): 

raise TypeError("Not a matrix space: %s" % S) 

  

MatrixMulAction.__init__(self, G, S, True) 

  

# disallow multiplication on different backends (same size and rings) 

if G.base_ring() is S.base_ring() and \ 

G.is_sparse() == S.is_sparse() and \ 

G._matrix_class != S._matrix_class: 

raise TypeError("no matrix multiplication between different implementations") 

  

# disallow multiplication (sparse) x (dense) when the densification is not the default 

# implementation 

if self.fix_sparseness: 

if G.is_sparse(): 

if not S._has_default_implementation(): 

raise TypeError("matrix multiplication not allowed") 

else: 

if not G._has_default_implementation(): 

raise TypeError("matrix multiplication not allowed") 

  

def _create_codomain(self, base): 

""" 

EXAMPLES: 

  

By :trac:`715`, there only is a weak reference on the underlying set, 

so that it can be garbage collected if only the action itself is 

explicitly referred to. Hence, we first assign the involved matrix 

spaces to a variable:: 

  

sage: from sage.matrix.action import MatrixMatrixAction 

sage: R.<x> = ZZ[] 

sage: MSR = MatrixSpace(R, 3, 3) 

sage: MSQ = MatrixSpace(QQ, 3, 2) 

sage: A = MatrixMatrixAction(MSR, MSQ); A 

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 

sage: A.codomain() 

Full MatrixSpace of 3 by 2 dense matrices over Univariate Polynomial Ring in x over Rational Field 

  

.. NOTE:: 

  

The :func:`MatrixSpace` function caches the object it creates. 

Therefore, the underlying set ``MSZ`` in the above example will not 

be garbage collected, even if it is not strongly ref'ed. 

Nonetheless, there is no guarantee that the set that is acted upon 

will always be cached in such a way, so that following the above 

example is good practice. 

  

""" 

if self.G.ncols() != self.underlying_set().nrows(): 

raise TypeError("incompatible dimensions %s, %s" % 

(self.G.ncols(), self.underlying_set().nrows())) 

return MatrixSpace(base, self.G.nrows(), self.underlying_set().ncols(), 

sparse = self.G.is_sparse() and self.underlying_set().is_sparse()) 

  

cpdef _call_(self, g, s): 

""" 

EXAMPLES: 

  

Respects compatible subdivisions:: 

  

sage: M = matrix(5, 5, prime_range(100)) 

sage: M.subdivide(2,3); M 

[ 2 3 5| 7 11] 

[13 17 19|23 29] 

[--------+-----] 

[31 37 41|43 47] 

[53 59 61|67 71] 

[73 79 83|89 97] 

sage: N = matrix(5,2,[n^2 for n in range(10)]) 

sage: N.subdivide(3,1); N 

[ 0| 1] 

[ 4| 9] 

[16|25] 

[--+--] 

[36|49] 

[64|81] 

sage: M*N 

[ 1048| 1388] 

[ 3056| 4117] 

[-----+-----] 

[ 5360| 7303] 

[ 8168|11143] 

[11056|15077] 

  

Note that this is just like block matrix multiplication:: 

  

sage: M.subdivision(0,0) * N.subdivision(0,0) + M.subdivision(0,1) * N.subdivision(1,0) 

[1048] 

[3056] 

  

If the subdivisions aren't compatible, ignore them. 

:: 

  

sage: N.subdivide(1,1); N 

[ 0| 1] 

[--+--] 

[ 4| 9] 

[16|25] 

[36|49] 

[64|81] 

sage: M*N 

[ 1048 1388] 

[ 3056 4117] 

[ 5360 7303] 

[ 8168 11143] 

[11056 15077] 

  

""" 

cdef Matrix A = g #<Matrix>g 

cdef Matrix B = s #<Matrix>s 

if A._parent._base is not self._codomain._base: 

A = A.change_ring(self._codomain._base) 

if B._parent._base is not self._codomain._base: 

B = B.change_ring(self._codomain._base) 

if self.fix_sparseness: 

if B.is_sparse_c(): 

B = B.dense_matrix() 

else: 

A = A.dense_matrix() 

assert type(A) == type(B), (type(A), type(B)) 

prod = A._matrix_times_matrix_(B) 

if A._subdivisions is not None or B._subdivisions is not None: 

Asubs = A.subdivisions() 

Bsubs = B.subdivisions() 

if Asubs[1] == Bsubs[0]: 

prod.subdivide(Asubs[0], Bsubs[1]) 

return prod 

  

  

cdef class MatrixVectorAction(MatrixMulAction): 

def __init__(self, G, S): 

""" 

EXAMPLES:: 

  

sage: from sage.matrix.action import MatrixVectorAction 

sage: A = MatrixVectorAction(MatrixSpace(QQ, 3, 3), VectorSpace(CDF, 4)); A 

Traceback (most recent call last): 

... 

TypeError: incompatible dimensions 3, 4 

""" 

if not is_FreeModule(S): 

raise TypeError("Not a free module: %s" % S) 

MatrixMulAction.__init__(self, G, S, True) 

  

def _create_codomain(self, base): 

""" 

EXAMPLES:: 

  

sage: from sage.matrix.action import MatrixVectorAction 

sage: M = MatrixSpace(QQ, 5, 3) 

sage: V = VectorSpace(CDF, 3) # strong reference prevents garbage collection 

sage: A = MatrixVectorAction(M, V); A 

Left action by Full MatrixSpace of 5 by 3 dense matrices over Rational Field on Vector space of dimension 3 over Complex Double Field 

sage: A.codomain() 

Vector space of dimension 5 over Complex Double Field 

""" 

if self.G.ncols() != self.underlying_set().degree(): 

raise TypeError("incompatible dimensions %s, %s" % (self.G.ncols(), 

self.underlying_set().degree())) 

return FreeModule(base, self.G.nrows(), sparse = self.G.is_sparse()) 

  

cpdef _call_(self, g, s): 

cdef Matrix A = g #<Matrix>g 

cdef Vector v = s #<Vector>s 

if A._parent._base is not self._codomain._base: 

A = A.change_ring(self._codomain._base) 

if v._parent._base is not self._codomain._base: 

v = v.change_ring(self._codomain._base) 

if self.fix_sparseness: 

if A.is_sparse_c(): 

v = v.sparse_vector() 

else: 

v = v.dense_vector() 

return A._matrix_times_vector_(v) 

  

  

cdef class VectorMatrixAction(MatrixMulAction): 

def __init__(self, G, S): 

""" 

EXAMPLES:: 

  

sage: from sage.matrix.action import VectorMatrixAction 

sage: A = VectorMatrixAction(MatrixSpace(QQ, 5, 3), VectorSpace(CDF, 3)); A 

Traceback (most recent call last): 

... 

TypeError: incompatible dimensions 5, 3 

""" 

if not is_FreeModule(S): 

raise TypeError("Not a free module: %s" % S) 

MatrixMulAction.__init__(self, G, S, False) 

  

def _create_codomain(self, base): 

""" 

EXAMPLES:: 

  

sage: from sage.matrix.action import VectorMatrixAction 

sage: M = MatrixSpace(QQ, 3, 5) 

sage: V = VectorSpace(CDF, 3) 

sage: A = VectorMatrixAction(M, V) 

sage: A 

Right action by Full MatrixSpace of 3 by 5 dense matrices over Rational Field on Vector space of dimension 3 over Complex Double Field 

sage: A.codomain() 

Vector space of dimension 5 over Complex Double Field 

""" 

if self.G.nrows() != self.underlying_set().degree(): 

raise TypeError("incompatible dimensions %s, %s" % (self.G.nrows(), 

self.underlying_set().degree())) 

return FreeModule(base, self.G.ncols(), sparse = self.G.is_sparse()) 

  

cpdef _call_(self, s, g): 

cdef Matrix A = g #<Matrix>g 

cdef Vector v = s #<Vector>s 

if A._parent._base is not self._codomain._base: 

A = A.change_ring(self._codomain._base) 

if v._parent._base is not self._codomain._base: 

v = v.change_ring(self._codomain._base) 

if self.fix_sparseness: 

if A.is_sparse_c(): 

v = v.sparse_vector() 

else: 

v = v.dense_vector() 

return (<Matrix>A)._vector_times_matrix_(v) # v * A