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

""" 

Vectors with rational entries. 

  

AUTHOR: 

  

- William Stein (2007) 

- Soroosh Yazdani (2007) 

  

EXAMPLES:: 

  

sage: v = vector(QQ,[1,2,3,4,5]) 

sage: v 

(1, 2, 3, 4, 5) 

sage: 3*v 

(3, 6, 9, 12, 15) 

sage: v/2 

(1/2, 1, 3/2, 2, 5/2) 

sage: -v 

(-1, -2, -3, -4, -5) 

sage: v - v 

(0, 0, 0, 0, 0) 

sage: v + v 

(2, 4, 6, 8, 10) 

sage: v * v 

55 

  

We make a large zero vector:: 

  

sage: k = QQ^100000; k 

Vector space of dimension 100000 over Rational Field 

sage: v = k(0) 

sage: v[:10] 

(0, 0, 0, 0, 0, 0, 0, 0, 0, 0) 

  

TESTS:: 

  

sage: v = vector(QQ, [1,2/5,-3/8,4]) 

sage: loads(dumps(v)) == v 

True 

  

sage: w = vector(QQ, [-1,0,0,0]) 

sage: w.set_immutable() 

sage: isinstance(hash(w), int) 

True 

""" 

  

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

# 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, print_function 

  

from cysignals.memory cimport check_allocarray, sig_free 

from cysignals.signals cimport sig_on, sig_off 

  

from sage.structure.element cimport Element, ModuleElement, RingElement, Vector 

  

from sage.rings.integer cimport Integer 

from sage.rings.rational cimport Rational 

  

cimport sage.modules.free_module_element as free_module_element 

from .free_module_element import vector 

  

from sage.libs.gmp.mpq cimport * 

  

  

cdef inline _Rational_from_mpq(mpq_t e): 

cdef Rational z = Rational.__new__(Rational) 

mpq_set(z.value, e) 

return z 

  

cdef class Vector_rational_dense(free_module_element.FreeModuleElement): 

cdef bint is_dense_c(self): 

return 1 

cdef bint is_sparse_c(self): 

return 0 

  

def __copy__(self): 

cdef Vector_rational_dense y 

y = self._new_c() 

cdef Py_ssize_t i 

for i in range(self._degree): 

mpq_set(y._entries[i], self._entries[i]) 

return y 

  

cdef int _init(self, Py_ssize_t degree, Parent parent) except -1: 

""" 

Initialize the C data structures in this vector. After calling 

this, ``self`` is a zero vector of degree ``degree`` with 

parent ``parent``. 

  

Only if you call ``__new__`` without a ``parent`` argument, it 

is needed to call this function manually. The only reason to do 

that is for efficiency: calling ``__new__`` without any 

additional arguments besides the type and then calling ``_init`` 

is (slightly) more efficient than calling ``__new__`` with a 

``parent`` argument. 

  

TESTS: 

  

Check implicitly that :trac:`10257` works:: 

  

sage: from sage.modules.vector_rational_dense import Vector_rational_dense 

sage: Vector_rational_dense(QQ^(sys.maxsize)) 

Traceback (most recent call last): 

... 

MemoryError: failed to allocate ... bytes 

sage: try: 

....: # Note: some malloc() implementations (on OS X 

....: # for example) print stuff when an allocation 

....: # fails. # We catch this with the ... in the 

....: # doctest result. The * is needed because a 

....: # result cannot start with ... 

....: print("*") 

....: Vector_rational_dense(QQ^(2^56)) 

....: except (MemoryError, OverflowError): 

....: print("allocation failed") 

*... 

allocation failed 

""" 

# Assign variables only when the array is fully initialized 

cdef mpq_t* entries = <mpq_t*>check_allocarray(degree, sizeof(mpq_t)) 

cdef Py_ssize_t i 

sig_on() 

for i in range(degree): 

mpq_init(entries[i]) 

sig_off() 

self._entries = entries 

self._degree = degree 

self._parent = parent 

  

def __cinit__(self, parent=None, x=None, coerce=True, copy=True): 

self._entries = NULL 

self._is_mutable = 1 

if parent is None: 

self._degree = 0 

return 

self._init(parent.degree(), <Parent?>parent) 

  

def __init__(self, parent, x, coerce=True, copy=True): 

cdef Py_ssize_t i 

cdef Rational z 

if isinstance(x, (list, tuple)): 

if len(x) != self._degree: 

raise TypeError("entries must be a list of length %s"%self._degree) 

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

z = Rational(x[i]) 

mpq_set(self._entries[i], z.value) 

return 

if x != 0: 

raise TypeError("can't initialize vector from nonzero non-list") 

  

def __dealloc__(self): 

cdef Py_ssize_t i 

if self._entries: 

# Do *not* use sig_on() here, since __dealloc__ 

# cannot raise exceptions! 

for i in range(self._degree): 

mpq_clear(self._entries[i]) 

sig_free(self._entries) 

  

cpdef int _cmp_(left, right) except -2: 

""" 

EXAMPLES:: 

  

sage: v = vector(QQ, [0,0,0,0]) 

sage: v == 0 

True 

sage: v == 1 

False 

sage: v == v 

True 

sage: w = vector(QQ, [-1,3/2,0,0]) 

sage: w < v 

True 

sage: w > v 

False 

sage: w = vector(QQ, [-1,0,0,0]) 

sage: w == w 

True 

""" 

cdef Py_ssize_t i 

cdef int c 

for i from 0 <= i < left.degree(): 

c = mpq_cmp(left._entries[i], (<Vector_rational_dense>right)._entries[i]) 

if c < 0: 

return -1 

elif c > 0: 

return 1 

return 0 

  

cdef get_unsafe(self, Py_ssize_t i): 

""" 

EXAMPLES:: 

  

sage: v = vector([1/2,2/3,3/4]); v 

(1/2, 2/3, 3/4) 

sage: v[0] 

1/2 

sage: v[2] 

3/4 

sage: v[-2] 

2/3 

sage: v[5] 

Traceback (most recent call last): 

... 

IndexError: vector index out of range 

sage: v[-5] 

Traceback (most recent call last): 

... 

IndexError: vector index out of range 

""" 

cdef Rational z = Rational.__new__(Rational) 

mpq_set(z.value, self._entries[i]) 

return z 

  

cdef int set_unsafe(self, Py_ssize_t i, value) except -1: 

""" 

EXAMPLES:: 

  

sage: v = vector(QQ, [1/2,2/5,0]); v 

(1/2, 2/5, 0) 

sage: v.set(2, -15/17); v 

(1/2, 2/5, -15/17) 

""" 

mpq_set(self._entries[i], (<Rational>value).value) 

  

  

def list(self,copy=True): 

""" 

The list of entries of the vector. 

  

INPUT: 

  

- ``copy``, ignored optional argument. 

  

EXAMPLES:: 

  

sage: v = vector(QQ,[1,2,3,4]) 

sage: a = v.list(copy=False); a 

[1, 2, 3, 4] 

sage: a[0] = 0 

sage: v 

(1, 2, 3, 4) 

""" 

cdef int i 

return [_Rational_from_mpq(self._entries[i]) for i in 

xrange(self._degree)] 

  

def __reduce__(self): 

return (unpickle_v1, (self._parent, self.list(), self._degree, self._is_mutable)) 

  

cpdef _add_(self, right): 

cdef Vector_rational_dense z, r 

r = right 

z = self._new_c() 

cdef Py_ssize_t i 

for i in range(self._degree): 

mpq_add(z._entries[i], self._entries[i], r._entries[i]) 

return z 

  

  

cpdef _sub_(self, right): 

cdef Vector_rational_dense z, r 

r = right 

z = self._new_c() 

cdef Py_ssize_t i 

for i in range(self._degree): 

mpq_sub(z._entries[i], self._entries[i], r._entries[i]) 

return z 

  

cpdef _dot_product_(self, Vector right): 

""" 

Dot product of dense vectors over the rationals. 

  

EXAMPLES:: 

  

sage: v = vector(QQ, [1,2,-3]); w = vector(QQ,[4,3,2]) 

sage: v*w 

4 

sage: w*v 

4 

""" 

cdef Vector_rational_dense r = right 

cdef Rational z 

z = Rational.__new__(Rational) 

cdef mpq_t t 

mpq_init(t) 

mpq_set_si(z.value, 0, 1) 

cdef Py_ssize_t i 

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

mpq_mul(t, self._entries[i], r._entries[i]) 

mpq_add(z.value, z.value, t) 

mpq_clear(t) 

return z 

  

  

cpdef _pairwise_product_(self, Vector right): 

""" 

EXAMPLES:: 

  

sage: v = vector(QQ, [1,2,-3]); w = vector(QQ,[4,3,2]) 

sage: v.pairwise_product(w) 

(4, 6, -6) 

""" 

cdef Vector_rational_dense z, r 

r = right 

z = self._new_c() 

cdef Py_ssize_t i 

for i in range(self._degree): 

mpq_mul(z._entries[i], self._entries[i], r._entries[i]) 

return z 

  

cpdef _rmul_(self, Element left): 

cdef Vector_rational_dense z 

cdef Rational a 

if isinstance(left, Rational): 

a = <Rational>left 

elif isinstance(left, Integer): 

a = <Rational>Rational.__new__(Rational) 

mpq_set_z(a.value, (<Integer>left).value) 

else: 

# should not happen 

raise TypeError("Cannot convert %s to %s" % (type(left).__name__, Rational.__name__)) 

z = self._new_c() 

cdef Py_ssize_t i 

for i in range(self._degree): 

mpq_mul(z._entries[i], self._entries[i], a.value) 

return z 

  

cpdef _lmul_(self, Element right): 

cdef Vector_rational_dense z 

cdef Rational a 

if isinstance(right, Rational): 

a = <Rational>right 

elif isinstance(right, Integer): 

a = <Rational>Rational.__new__(Rational) 

mpq_set_z(a.value, (<Integer>right).value) 

else: 

# should not happen 

raise TypeError("Cannot convert %s to %s" % (type(right).__name__, Rational.__name__)) 

z = self._new_c() 

cdef Py_ssize_t i 

for i in range(self._degree): 

mpq_mul(z._entries[i], self._entries[i], a.value) 

return z 

  

cpdef _neg_(self): 

cdef Vector_rational_dense z 

z = self._new_c() 

cdef Py_ssize_t i 

for i in range(self._degree): 

mpq_neg(z._entries[i], self._entries[i]) 

return z 

  

  

def unpickle_v0(parent, entries, degree): 

# If you think you want to change this function, don't. 

# Instead make a new version with a name like 

# make_FreeModuleElement_generic_dense_v1 

# and changed the reduce method below. 

cdef Vector_rational_dense v 

v = Vector_rational_dense.__new__(Vector_rational_dense) 

v._init(degree, parent) 

cdef Rational z 

cdef Py_ssize_t i 

for i in range(degree): 

z = Rational(entries[i]) 

mpq_set(v._entries[i], z.value) 

return v 

  

def unpickle_v1(parent, entries, degree, is_mutable): 

cdef Vector_rational_dense v 

v = Vector_rational_dense.__new__(Vector_rational_dense) 

v._init(degree, parent) 

cdef Rational z 

cdef Py_ssize_t i 

for i in range(degree): 

z = Rational(entries[i]) 

mpq_set(v._entries[i], z.value) 

v._is_mutable = is_mutable 

return v