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

403

404

405

406

407

408

409

410

411

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

455

456

457

458

459

460

461

462

463

464

465

466

467

468

469

470

471

472

473

474

475

476

477

478

479

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

521

522

523

524

525

526

527

528

529

530

531

532

533

534

535

536

537

538

539

540

541

542

543

544

545

546

547

548

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

641

642

643

644

645

646

647

648

649

650

651

652

653

654

655

656

657

658

659

660

661

662

663

664

665

666

667

668

669

670

671

672

673

674

675

676

677

678

679

680

681

682

683

684

685

686

687

688

689

690

691

692

693

694

695

696

697

698

699

700

701

702

703

704

705

706

707

708

709

710

711

712

713

714

715

716

717

718

719

720

721

722

723

724

725

726

727

728

729

730

731

732

733

734

735

736

737

738

739

740

741

742

743

744

745

746

747

748

749

750

751

752

753

754

755

756

757

758

759

760

761

762

763

764

765

766

767

768

769

770

771

772

773

774

775

776

777

778

779

780

781

782

783

784

785

786

787

788

789

790

791

792

793

794

795

796

797

798

799

800

801

802

803

804

805

806

807

808

809

810

811

812

813

814

815

816

817

818

819

820

821

822

823

824

825

826

827

828

829

830

831

832

833

834

835

836

837

838

839

840

841

842

843

844

845

846

847

848

849

850

851

852

853

854

855

856

857

858

859

860

861

862

863

864

""" 

Sparse rational matrices. 

  

AUTHORS: 

  

- William Stein (2007-02-21) 

- Soroosh Yazdani (2007-02-21) 

  

TESTS:: 

  

sage: a = matrix(QQ,2,range(4), sparse=True) 

sage: TestSuite(a).run() 

sage: matrix(QQ,0,0,sparse=True).inverse() 

[] 

""" 

  

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

# Copyright (C) 2007 William Stein <wstein@gmail.com> 

# 

# 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 

  

from cysignals.signals cimport sig_on, sig_off 

from cysignals.memory cimport sig_malloc, sig_free 

  

from collections import Iterator, Sequence 

  

from sage.data_structures.binary_search cimport * 

from sage.modules.vector_integer_sparse cimport * 

from sage.modules.vector_rational_sparse cimport * 

  

from cpython.sequence cimport * 

  

from sage.rings.rational cimport Rational 

from sage.rings.integer cimport Integer 

from .matrix cimport Matrix 

  

from sage.libs.gmp.mpz cimport * 

from sage.libs.gmp.mpq cimport * 

  

from sage.libs.flint.fmpq cimport fmpq_set_mpq 

from sage.libs.flint.fmpq_mat cimport fmpq_mat_entry 

  

from sage.rings.integer_ring import ZZ 

from sage.rings.rational_field import QQ 

  

cimport sage.structure.element 

  

import sage.matrix.matrix_space 

  

from .matrix_integer_sparse cimport Matrix_integer_sparse 

from .matrix_rational_dense cimport Matrix_rational_dense 

  

from sage.misc.misc import verbose 

  

cdef class Matrix_rational_sparse(Matrix_sparse): 

def __cinit__(self, parent, entries, copy, coerce): 

# set the parent, nrows, ncols, etc. 

Matrix_sparse.__init__(self, parent) 

  

self._matrix = <mpq_vector*> sig_malloc(parent.nrows()*sizeof(mpq_vector)) 

if self._matrix == NULL: 

raise MemoryError("error allocating sparse matrix") 

# initialize the rows 

for i from 0 <= i < parent.nrows(): 

mpq_vector_init(&self._matrix[i], self._ncols, 0) 

  

# record that rows have been initialized 

self._initialized = True 

  

  

def __dealloc__(self): 

self._dealloc() 

  

cdef _dealloc(self): 

cdef Py_ssize_t i 

if self._initialized: 

for i from 0 <= i < self._nrows: 

mpq_vector_clear(&self._matrix[i]) 

if self._matrix != NULL: 

sig_free(self._matrix) 

  

def __init__(self, parent, entries, copy, coerce): 

""" 

Create a sparse matrix over the rational numbers. 

  

INPUT: 

  

- ``parent`` -- a matrix space 

  

- ``entries`` -- can be one of the following: 

  

* a Python dictionary whose items have the 

form ``(i, j): x``, where ``0 <= i < nrows``, 

``0 <= j < ncols``, and ``x`` is coercible to 

a rational. The ``i,j`` entry of ``self`` is 

set to ``x``. The ``x``'s can be ``0``. 

* Alternatively, entries can be a list of *all* 

the entries of the sparse matrix, read 

row-by-row from top to bottom (so they would 

be mostly 0). 

  

- ``copy`` -- ignored 

  

- ``coerce`` -- ignored 

""" 

cdef Py_ssize_t i, j, k 

cdef Rational z 

cdef PyObject** X 

  

if entries is None: 

return 

# fill in entries in the dict case 

if isinstance(entries, dict): 

R = self.base_ring() 

for ij, x in entries.iteritems(): 

z = R(x) 

if z != 0: 

i, j = ij # nothing better to do since this is user input, which may be bogus. 

if i < 0 or j < 0 or i >= self._nrows or j >= self._ncols: 

raise IndexError("invalid entries list") 

mpq_vector_set_entry(&self._matrix[i], j, z.value) 

elif isinstance(entries, (Iterator, Sequence)): 

if not isinstance(entries, (list, tuple)): 

entries = list(entries) 

# Dense input format -- fill in entries 

if len(entries) != self._nrows * self._ncols: 

raise TypeError("list of entries must be a dictionary of (i,j):x or a dense list of n * m elements") 

seq = PySequence_Fast(entries,"expected a list") 

X = PySequence_Fast_ITEMS(seq) 

k = 0 

R = self.base_ring() 

# Get fast access to the entries list. 

for i from 0 <= i < self._nrows: 

for j from 0 <= j < self._ncols: 

z = R(<object>X[k]) 

if z != 0: 

mpq_vector_set_entry(&self._matrix[i], j, z.value) 

k = k + 1 

else: 

# fill in entries in the scalar case 

z = Rational(entries) 

if z == 0: 

return 

if self._nrows != self._ncols: 

raise TypeError("matrix must be square to initialize with a scalar.") 

for i from 0 <= i < self._nrows: 

mpq_vector_set_entry(&self._matrix[i], i, z.value) 

  

  

cdef set_unsafe(self, Py_ssize_t i, Py_ssize_t j, x): 

mpq_vector_set_entry(&self._matrix[i], j, (<Rational> x).value) 

  

cdef get_unsafe(self, Py_ssize_t i, Py_ssize_t j): 

cdef Rational x 

x = Rational() 

mpq_vector_get_entry(x.value, &self._matrix[i], j) 

return x 

  

def add_to_entry(self, Py_ssize_t i, Py_ssize_t j, elt): 

r""" 

Add ``elt`` to the entry at position ``(i, j)``. 

  

EXAMPLES:: 

  

sage: m = matrix(QQ, 2, 2, sparse=True) 

sage: m.add_to_entry(0, 0, -1/3) 

sage: m 

[-1/3 0] 

[ 0 0] 

sage: m.add_to_entry(0, 0, 1/3) 

sage: m 

[0 0] 

[0 0] 

sage: m.nonzero_positions() 

[] 

""" 

if not isinstance(elt, Rational): 

elt = Rational(elt) 

if i < 0: 

i += self._nrows 

if i < 0 or i >= self._nrows: 

raise IndexError("row index out of range") 

if j < 0: 

j += self._ncols 

if j < 0 or j >= self._ncols: 

raise IndexError("column index out of range") 

  

cdef mpq_t z 

mpq_init(z) 

mpq_vector_get_entry(z, &self._matrix[i], j) 

mpq_add(z, z, (<Rational>elt).value) 

mpq_vector_set_entry(&self._matrix[i], j, z) 

mpq_clear(z) 

  

  

######################################################################## 

# LEVEL 2 functionality 

# * def _pickle 

# * def _unpickle 

# * cdef _add_ 

# * cdef _sub_ 

# * cdef _mul_ 

# * cpdef _cmp_ 

# * __neg__ 

# * __invert__ 

# * __copy__ 

# * _multiply_classical 

# * _matrix_times_matrix_ 

# * _list -- list of underlying elements (need not be a copy) 

# * x _dict -- sparse dictionary of underlying elements (need not be a copy) 

  

cdef sage.structure.element.Matrix _matrix_times_matrix_(self, sage.structure.element.Matrix _right): 

cdef Matrix_rational_sparse right, ans 

right = _right 

  

cdef mpq_vector* v 

  

# Build a table that gives the nonzero positions in each column of right 

nonzero_positions_in_columns = [set([]) for _ in range(right._ncols)] 

cdef Py_ssize_t i, j, k 

for i from 0 <= i < right._nrows: 

v = &(right._matrix[i]) 

for j from 0 <= j < right._matrix[i].num_nonzero: 

nonzero_positions_in_columns[v.positions[j]].add(i) 

  

ans = self.new_matrix(self._nrows, right._ncols) 

  

# Now do the multiplication, getting each row completely before filling it in. 

cdef mpq_t x, y, s 

mpq_init(x) 

mpq_init(y) 

mpq_init(s) 

for i from 0 <= i < self._nrows: 

v = &self._matrix[i] 

for j from 0 <= j < right._ncols: 

mpq_set_si(s, 0, 1) 

c = nonzero_positions_in_columns[j] 

for k from 0 <= k < v.num_nonzero: 

if v.positions[k] in c: 

mpq_vector_get_entry(y, &right._matrix[v.positions[k]], j) 

mpq_mul(x, v.entries[k], y) 

mpq_add(s, s, x) 

mpq_vector_set_entry(&ans._matrix[i], j, s) 

  

mpq_clear(x) 

mpq_clear(y) 

mpq_clear(s) 

return ans 

  

def _matrix_times_matrix_dense(self, sage.structure.element.Matrix _right): 

""" 

Do the sparse matrix multiply, but return a dense matrix as the result. 

  

EXAMPLES:: 

  

sage: a = matrix(QQ, 2, [1,2,3,4], sparse=True) 

sage: b = matrix(QQ, 2, 3, [1..6], sparse=True) 

sage: a * b 

[ 9 12 15] 

[19 26 33] 

sage: c = a._matrix_times_matrix_dense(b); c 

[ 9 12 15] 

[19 26 33] 

sage: type(c) 

<type 'sage.matrix.matrix_rational_dense.Matrix_rational_dense'> 

""" 

cdef Matrix_rational_sparse right 

cdef Matrix_rational_dense ans 

right = _right 

  

cdef mpq_vector* v 

  

# Build a table that gives the nonzero positions in each column of right 

nonzero_positions_in_columns = [set([]) for _ in range(right._ncols)] 

cdef Py_ssize_t i, j, k 

for i in range(right._nrows): 

v = &(right._matrix[i]) 

for j in range(right._matrix[i].num_nonzero): 

nonzero_positions_in_columns[v.positions[j]].add(i) 

  

ans = self.new_matrix(self._nrows, right._ncols, sparse=False) 

  

# Now do the multiplication, getting each row completely before filling it in. 

cdef mpq_t x, y, s 

mpq_init(x) 

mpq_init(y) 

mpq_init(s) 

for i in range(self._nrows): 

v = &self._matrix[i] 

for j in range(right._ncols): 

mpq_set_si(s, 0, 1) 

c = nonzero_positions_in_columns[j] 

for k in range(v.num_nonzero): 

if v.positions[k] in c: 

mpq_vector_get_entry(y, &right._matrix[v.positions[k]], j) 

mpq_mul(x, v.entries[k], y) 

mpq_add(s, s, x) 

fmpq_set_mpq(fmpq_mat_entry(ans._matrix, i, j), s) 

  

mpq_clear(x) 

mpq_clear(y) 

mpq_clear(s) 

return ans 

  

  

######################################################################## 

# def _pickle(self): 

# def _unpickle(self, data, int version): # use version >= 0 

# cpdef _add_(self, right): 

# cdef _mul_(self, Matrix right): 

# cpdef int _cmp_(self, Matrix right) except -2: 

# def __neg__(self): 

# def __invert__(self): 

# def __copy__(self): 

# def _multiply_classical(left, matrix.Matrix _right): 

# def _list(self): 

  

# TODO 

## cpdef _lmul_(self, Element right): 

## """ 

## EXAMPLES: 

## sage: a = matrix(QQ,2,range(6)) 

## sage: (3/4) * a 

## [ 0 3/4 3/2] 

## [ 9/4 3 15/4] 

## """ 

  

def _dict(self): 

""" 

Unsafe version of the dict method, mainly for internal use. 

This may return the dict of elements, but as an *unsafe* 

reference to the underlying dict of the object. It might 

be dangerous if you change entries of the returned dict. 

""" 

d = self.fetch('dict') 

if not d is None: 

return d 

  

cdef Py_ssize_t i, j, k 

d = {} 

for i from 0 <= i < self._nrows: 

for j from 0 <= j < self._matrix[i].num_nonzero: 

x = Rational() 

mpq_set((<Rational>x).value, self._matrix[i].entries[j]) 

d[(int(i),int(self._matrix[i].positions[j]))] = x 

self.cache('dict', d) 

return d 

  

  

######################################################################## 

# LEVEL 3 functionality (Optional) 

# * cdef _sub_ 

# * __deepcopy__ 

# * __invert__ 

# * Matrix windows -- only if you need strassen for that base 

# * Other functions (list them here): 

######################################################################## 

  

def _nonzero_positions_by_row(self, copy=True): 

""" 

Returns the list of pairs (i,j) such that self[i,j] != 0. 

  

It is safe to change the resulting list (unless you give the option copy=False). 

  

EXAMPLES:: 

  

sage: M = Matrix(QQ, [[0,0,0,1,0,0,0,0],[0,1,0,0,0,0,1,0]], sparse=True); M 

[0 0 0 1 0 0 0 0] 

[0 1 0 0 0 0 1 0] 

sage: M.nonzero_positions() 

[(0, 3), (1, 1), (1, 6)] 

  

""" 

x = self.fetch('nonzero_positions') 

if not x is None: 

if copy: 

return list(x) 

return x 

nzp = [] 

cdef Py_ssize_t i, j 

for i from 0 <= i < self._nrows: 

for j from 0 <= j < self._matrix[i].num_nonzero: 

nzp.append((i,self._matrix[i].positions[j])) 

self.cache('nonzero_positions', nzp) 

if copy: 

return list(nzp) 

return nzp 

  

def height(self): 

""" 

Return the height of this matrix, which is the least common 

multiple of all numerators and denominators of elements of 

this matrix. 

  

OUTPUT: 

  

-- Integer 

  

EXAMPLES:: 

  

sage: b = matrix(QQ,2,range(6), sparse=True); b[0,0]=-5007/293; b 

[-5007/293 1 2] 

[ 3 4 5] 

sage: b.height() 

5007 

""" 

cdef Integer z = Integer.__new__(Integer) 

self.mpz_height(z.value) 

return z 

  

cdef int mpz_height(self, mpz_t height) except -1: 

cdef mpz_t x, h 

mpz_init(x) 

mpz_init_set_si(h, 0) 

cdef int i, j 

sig_on() 

for i from 0 <= i < self._nrows: 

for j from 0 <= j < self._matrix[i].num_nonzero: 

mpq_get_num(x, self._matrix[i].entries[j]) 

mpz_abs(x, x) 

if mpz_cmp(h,x) < 0: 

mpz_set(h,x) 

mpq_get_den(x, self._matrix[i].entries[j]) 

mpz_abs(x, x) 

if mpz_cmp(h,x) < 0: 

mpz_set(h,x) 

sig_off() 

mpz_set(height, h) 

mpz_clear(h) 

mpz_clear(x) 

return 0 

  

cdef int mpz_denom(self, mpz_t d) except -1: 

mpz_set_si(d,1) 

cdef Py_ssize_t i, j 

cdef mpq_vector *v 

  

sig_on() 

for i from 0 <= i < self._nrows: 

for j from 0 <= j < self._matrix[i].num_nonzero: 

mpz_lcm(d, d, mpq_denref(self._matrix[i].entries[j])) 

sig_off() 

return 0 

  

  

def denominator(self): 

""" 

Return the denominator of this matrix. 

  

OUTPUT: 

  

-- Sage Integer 

  

EXAMPLES:: 

  

sage: b = matrix(QQ,2,range(6)); b[0,0]=-5007/293; b 

[-5007/293 1 2] 

[ 3 4 5] 

sage: b.denominator() 

293 

""" 

cdef Integer z = Integer.__new__(Integer) 

self.mpz_denom(z.value) 

return z 

  

def _clear_denom(self): 

""" 

INPUT: 

  

self -- a matrix 

  

OUTPUT: 

  

D*self, D 

  

The product D*self is a matrix over ZZ 

  

EXAMPLES:: 

  

sage: a = matrix(QQ,3,[-2/7, -1/4, -2, 0, 1/7, 1, 0, 1/2, 1/5],sparse=True) 

sage: a.denominator() 

140 

sage: a._clear_denom() 

( 

[ -40 -35 -280] 

[ 0 20 140] 

[ 0 70 28], 140 

) 

""" 

cdef Integer D 

cdef Py_ssize_t i, j 

cdef Matrix_integer_sparse A 

cdef mpz_t t 

cdef mpq_vector* v 

  

D = Integer() 

self.mpz_denom(D.value) 

  

MZ = sage.matrix.matrix_space.MatrixSpace(ZZ, self._nrows, self._ncols, sparse=True) 

A = MZ.zero_matrix().__copy__() 

  

mpz_init(t) 

sig_on() 

for i from 0 <= i < self._nrows: 

v = &(self._matrix[i]) 

for j from 0 <= j < v.num_nonzero: 

mpz_divexact(t, D.value, mpq_denref(v.entries[j])) 

mpz_mul(t, t, mpq_numref(v.entries[j])) 

mpz_vector_set_entry(&(A._matrix[i]), v.positions[j], t) 

sig_off() 

mpz_clear(t) 

A._initialized = 1 

return A, D 

  

  

################################################ 

# Echelon form 

################################################ 

def echelonize(self, height_guess=None, proof=True, **kwds): 

""" 

Transform the matrix ``self`` into reduced row echelon form 

in place. 

  

INPUT: 

  

``height_guess``, ``proof``, ``**kwds`` -- all passed to the multimodular 

algorithm; ignored by the p-adic algorithm. 

  

OUTPUT: 

  

Nothing. The matrix ``self`` is transformed into reduced row 

echelon form in place. 

  

ALGORITHM: a multimodular algorithm. 

  

EXAMPLES:: 

  

sage: a = matrix(QQ, 4, range(16), sparse=True); a[0,0] = 1/19; a[0,1] = 1/5; a 

[1/19 1/5 2 3] 

[ 4 5 6 7] 

[ 8 9 10 11] 

[ 12 13 14 15] 

sage: a.echelonize(); a 

[ 1 0 0 -76/157] 

[ 0 1 0 -5/157] 

[ 0 0 1 238/157] 

[ 0 0 0 0] 

  

:trac:`10319` has been fixed:: 

  

sage: m = Matrix(QQ, [1], sparse=True); m.echelonize() 

sage: m = Matrix(QQ, [1], sparse=True); m.echelonize(); m 

[1] 

""" 

  

x = self.fetch('in_echelon_form') 

if not x is None: return # already known to be in echelon form 

self.check_mutability() 

  

pivots = self._echelonize_multimodular(height_guess, proof, **kwds) 

  

self.cache('in_echelon_form', True) 

self.cache('echelon_form', self) 

self.cache('pivots', pivots) 

self.cache('rank', len(pivots)) 

  

def echelon_form(self, algorithm='default', 

height_guess=None, proof=True, **kwds): 

""" 

INPUT: 

  

``height_guess``, ``proof``, ``**kwds`` -- all passed to the multimodular 

algorithm; ignored by the p-adic algorithm. 

  

OUTPUT: 

  

self is no in reduced row echelon form. 

  

EXAMPLES:: 

  

sage: a = matrix(QQ, 4, range(16), sparse=True); a[0,0] = 1/19; a[0,1] = 1/5; a 

[1/19 1/5 2 3] 

[ 4 5 6 7] 

[ 8 9 10 11] 

[ 12 13 14 15] 

sage: a.echelon_form() 

[ 1 0 0 -76/157] 

[ 0 1 0 -5/157] 

[ 0 0 1 238/157] 

[ 0 0 0 0] 

""" 

label = 'echelon_form_%s'%algorithm 

x = self.fetch(label) 

if not x is None: 

return x 

if self.fetch('in_echelon_form'): return self 

  

E, pivots = self._echelon_form_multimodular(height_guess, proof=proof) 

  

self.cache(label, E) 

self.cache('pivots', pivots) 

return E 

  

# Multimodular echelonization algorithms 

def _echelonize_multimodular(self, height_guess=None, proof=True, **kwds): 

cdef Py_ssize_t i, j 

cdef Matrix_rational_sparse E 

cdef mpq_vector* v 

cdef mpq_vector* w 

E, pivots = self._echelon_form_multimodular(height_guess, proof=proof, **kwds) 

self.clear_cache() 

# Get rid of self's data 

self._dealloc() 

# Copy E's data to self's data. 

self._matrix = <mpq_vector*> sig_malloc(E._nrows * sizeof(mpq_vector)) 

if self._matrix == NULL: 

raise MemoryError("error allocating sparse matrix") 

for i in range(E._nrows): 

v = &self._matrix[i] 

w = &E._matrix[i] 

mpq_vector_init(v, E._ncols, w.num_nonzero) 

for j in range(w.num_nonzero): 

mpq_set(v.entries[j], w.entries[j]) 

v.positions[j] = w.positions[j] 

return pivots 

  

  

def _echelon_form_multimodular(self, height_guess=None, proof=True): 

""" 

Returns reduced row-echelon form using a multi-modular 

algorithm. Does not change self. 

  

INPUT: 

  

- height_guess -- integer or None 

- proof -- boolean (default: True) 

""" 

from .misc import matrix_rational_echelon_form_multimodular 

cdef Matrix E 

E, pivots = matrix_rational_echelon_form_multimodular(self, 

height_guess=height_guess, proof=proof) 

E._parent = self._parent 

return E, pivots 

  

  

def set_row_to_multiple_of_row(self, i, j, s): 

""" 

Set row i equal to s times row j. 

  

EXAMPLES:: 

  

sage: a = matrix(QQ,2,3,range(6), sparse=True); a 

[0 1 2] 

[3 4 5] 

sage: a.set_row_to_multiple_of_row(1,0,-3) 

sage: a 

[ 0 1 2] 

[ 0 -3 -6] 

""" 

self.check_row_bounds_and_mutability(i, j) 

cdef Rational _s 

_s = Rational(s) 

mpq_vector_scalar_multiply(&self._matrix[i], &self._matrix[j], _s.value) 

  

def dense_matrix(self): 

""" 

Return dense version of this matrix. 

  

EXAMPLES:: 

  

sage: a = matrix(QQ,2,[1..4],sparse=True); type(a) 

<type 'sage.matrix.matrix_rational_sparse.Matrix_rational_sparse'> 

sage: type(a.dense_matrix()) 

<type 'sage.matrix.matrix_rational_dense.Matrix_rational_dense'> 

sage: a.dense_matrix() 

[1 2] 

[3 4] 

  

Check that subdivisions are preserved when converting between 

dense and sparse matrices:: 

  

sage: a.subdivide([1,1], [2]) 

sage: b = a.dense_matrix().sparse_matrix().dense_matrix() 

sage: b.subdivisions() == a.subdivisions() 

True 

""" 

cdef Matrix_rational_dense B 

cdef mpq_vector* v 

  

B = self.matrix_space(sparse=False).zero_matrix().__copy__() 

for i from 0 <= i < self._nrows: 

v = &(self._matrix[i]) 

for j from 0 <= j < v.num_nonzero: 

fmpq_set_mpq(fmpq_mat_entry(B._matrix, i, v.positions[j]), v.entries[j]) 

B.subdivide(self.subdivisions()) 

return B 

  

## def _set_row_to_negative_of_row_of_A_using_subset_of_columns(self, Py_ssize_t i, Matrix A, 

## Py_ssize_t r, cols): 

## B = self.__copy__() 

## B.x_set_row_to_negative_of_row_of_A_using_subset_of_columns(i, A, r, cols) 

## cdef Py_ssize_t l 

## l = 0 

## for z in range(self.ncols()): 

## self[i,z] = 0 

## for k in cols: 

## self.set_unsafe(i,l,-A.get_unsafe(r,k)) #self[i,l] = -A[r,k] 

## l += 1 

## if self != B: 

## print("correct =\n", self.str()) 

## print("wrong = \n", B.str()) 

## print("diff = \n", (self-B).str()) 

  

def _set_row_to_negative_of_row_of_A_using_subset_of_columns(self, Py_ssize_t i, Matrix A, 

Py_ssize_t r, cols, 

cols_index=None): 

""" 

Set row i of self to -(row r of A), but where we only take the 

given column positions in that row of A. Note that we *DO* 

zero out the other entries of self's row i. 

  

INPUT: 

  

- i -- integer, index into the rows of self 

- A -- a sparse matrix 

- r -- integer, index into rows of A 

- cols -- a *sorted* list of integers. 

- cols_index -- (optional). But set it to this to vastly speed up 

calls to this function:: 

  

dict([(cols[i], i) for i in range(len(cols))]) 

  

EXAMPLES:: 

  

sage: a = matrix(QQ,2,3,range(6), sparse=True); a 

[0 1 2] 

[3 4 5] 

  

Note that the row is zeroed out before being set in the sparse case. :: 

  

sage: a._set_row_to_negative_of_row_of_A_using_subset_of_columns(0,a,1,[1,2]) 

sage: a 

[-4 -5 0] 

[ 3 4 5] 

""" 

# this function exists just because it is useful for modular symbols presentations. 

self.check_row_bounds_and_mutability(i,i) 

if r < 0 or r >= A.nrows(): 

raise IndexError("invalid row") 

  

if not A.is_sparse(): 

A = A.sparse_matrix() 

  

if A.base_ring() != QQ: 

A = A.change_ring(QQ) 

  

cdef Matrix_rational_sparse _A 

_A = A 

  

cdef Py_ssize_t l, n 

  

cdef mpq_vector *v 

cdef mpq_vector *w 

v = &self._matrix[i] 

w = &_A._matrix[r] 

  

  

if cols_index is None: 

cols_index = dict([(cols[i], i) for i in range(len(cols))]) 

  

_cols = set(cols) 

pos = [i for i from 0 <= i < w.num_nonzero if w.positions[i] in _cols] 

n = len(pos) 

  

mpq_vector_clear(v) 

allocate_mpq_vector(v, n) 

v.num_nonzero = n 

v.degree = self._ncols 

  

  

for l from 0 <= l < n: 

v.positions[l] = cols_index[w.positions[pos[l]]] 

mpq_mul(v.entries[l], w.entries[pos[l]], minus_one) 

  

def _right_kernel_matrix(self, **kwds): 

r""" 

Returns a pair that includes a matrix of basis vectors 

for the right kernel of ``self``. 

  

INPUT: 

  

- ``kwds`` - these are provided for consistency with other versions 

of this method. Here they are ignored as there is no optional 

behavior available. 

  

OUTPUT: 

  

Returns a pair. First item is the string 'computed-iml-rational' 

that identifies the nature of the basis vectors. 

  

Second item is a matrix whose rows are a basis for the right kernel, 

over the rationals, as computed by the IML library. Notice that the 

IML library returns a matrix that is in the 'pivot' format, once the 

whole matrix is multiplied by -1. So the 'computed' format is very 

close to the 'pivot' format. 

  

EXAMPLES:: 

  

sage: A = matrix(QQ, [ 

....: [1, 0, 1, -3, 1], 

....: [-5, 1, 0, 7, -3], 

....: [0, -1, -4, 6, -2], 

....: [4, -1, 0, -6, 2]], 

....: sparse=True) 

sage: result = A._right_kernel_matrix() 

sage: result[0] 

'computed-iml-rational' 

sage: result[1] 

[-1 2 -2 -1 0] 

[ 1 2 0 0 -1] 

sage: X = result[1].transpose() 

sage: A*X == zero_matrix(QQ, 4, 2) 

True 

  

Computed result is the negative of the pivot basis, which 

is just slightly more efficient to compute. :: 

  

sage: A.right_kernel_matrix(basis='pivot') == -A.right_kernel_matrix(basis='computed') 

True 

  

TESTS: 

  

We test three trivial cases. :: 

  

sage: A = matrix(QQ, 0, 2, sparse=True) 

sage: A._right_kernel_matrix()[1] 

[1 0] 

[0 1] 

sage: A = matrix(QQ, 2, 0, sparse=True) 

sage: A._right_kernel_matrix()[1].parent() 

Full MatrixSpace of 0 by 0 dense matrices over Rational Field 

sage: A = zero_matrix(QQ, 4, 3, sparse=True) 

sage: A._right_kernel_matrix()[1] 

[1 0 0] 

[0 1 0] 

[0 0 1] 

""" 

return self.dense_matrix()._right_kernel_matrix() 

  

  

######################### 

  

# used for a function above 

cdef mpq_t minus_one 

mpq_init(minus_one) 

mpq_set_si(minus_one, -1,1)