Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/modules/vector_mod2_dense.pyx
4036 views
1
"""
2
Vectors with elements in GF(2).
3
4
AUTHOR:
5
6
- Martin Albrecht (2009-12): initial implementation
7
8
EXAMPLES::
9
10
sage: VS = GF(2)^3
11
sage: e = VS.random_element(); e
12
(1, 0, 0)
13
sage: f = VS.random_element(); f
14
(0, 1, 1)
15
sage: e + f
16
(1, 1, 1)
17
"""
18
19
##############################################################################
20
# Copyright (C) 2009 Martin Albrecht <[email protected]>
21
# Distributed under the terms of the GNU General Public License (GPL)
22
# The full text of the GPL is available at:
23
# http://www.gnu.org/licenses/
24
##############################################################################
25
26
include '../ext/interrupt.pxi'
27
include '../ext/stdsage.pxi'
28
29
from sage.rings.finite_rings.integer_mod cimport IntegerMod_int, IntegerMod_abstract
30
from sage.rings.integer cimport Integer
31
32
from sage.structure.element cimport Element, ModuleElement, RingElement, Vector
33
34
cimport free_module_element
35
from free_module_element import vector
36
37
from sage.libs.m4ri cimport *
38
39
cdef class Vector_mod2_dense(free_module_element.FreeModuleElement):
40
cdef _new_c(self):
41
"""
42
EXAMPLE::
43
44
sage: VS = VectorSpace(GF(2),3)
45
sage: VS([0,0,1])
46
(0, 0, 1)
47
sage: type(_)
48
<type 'sage.modules.vector_mod2_dense.Vector_mod2_dense'>
49
"""
50
cdef Vector_mod2_dense y
51
y = PY_NEW(Vector_mod2_dense)
52
y._init(self._degree, self._parent)
53
return y
54
55
cdef bint is_dense_c(self):
56
"""
57
EXAMPLE::
58
59
sage: VS = VectorSpace(GF(2),3)
60
sage: VS([0,0,1]).is_dense()
61
True
62
"""
63
return 1
64
65
cdef bint is_sparse_c(self):
66
"""
67
EXAMPLE::
68
69
sage: VS = VectorSpace(GF(2),3)
70
sage: VS([0,0,1]).is_sparse()
71
False
72
"""
73
return 0
74
75
def __copy__(self):
76
"""
77
EXAMPLE::
78
79
sage: VS = VectorSpace(GF(2),10^4)
80
sage: v = VS.random_element()
81
sage: w = copy(v)
82
sage: w == v
83
True
84
sage: v[:10]
85
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
86
sage: w[:10]
87
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
88
"""
89
cdef Vector_mod2_dense y = self._new_c()
90
if self._degree:
91
mzd_copy(y._entries, self._entries)
92
return y
93
94
cdef _init(self, Py_ssize_t degree, parent):
95
"""
96
EXAMPLE::
97
98
sage: VS = VectorSpace(GF(2),3)
99
sage: VS([0,0,1])
100
(0, 0, 1)
101
sage: type(_)
102
<type 'sage.modules.vector_mod2_dense.Vector_mod2_dense'>
103
"""
104
self._degree = degree
105
self._parent = parent
106
self._base_ring = parent.base_ring()
107
self._entries = mzd_init(1, degree)
108
if self._entries == NULL:
109
raise MemoryError("Allocation of Vector_mod2_dense failed.")
110
111
def __cinit__(self, parent=None, x=None, coerce=True, copy=True):
112
"""
113
EXAMPLE::
114
115
sage: VS = VectorSpace(GF(2),3)
116
sage: VS((0,0,1/3))
117
(0, 0, 1)
118
sage: type(_)
119
<type 'sage.modules.vector_mod2_dense.Vector_mod2_dense'>
120
"""
121
self._entries = NULL
122
self._is_mutable = 1
123
if not parent is None:
124
self._init(parent.degree(), parent)
125
126
def __init__(self, parent, x, coerce=True, copy=True):
127
"""
128
EXAMPLE::
129
130
sage: VS = VectorSpace(GF(2),3)
131
sage: VS((0,0,1/3))
132
(0, 0, 1)
133
sage: type(_)
134
<type 'sage.modules.vector_mod2_dense.Vector_mod2_dense'>
135
sage: VS((0,0,int(3)))
136
(0, 0, 1)
137
sage: VS((0,0,3))
138
(0, 0, 1)
139
sage: VS((0,0,GF(2)(1)))
140
(0, 0, 1)
141
142
TESTS:
143
144
Check that ticket #8601 is fixed::
145
146
sage: VS = VectorSpace(GF(2), 3)
147
sage: VS((-1,-2,-3))
148
(1, 0, 1)
149
sage: V = VectorSpace(GF(2), 2)
150
sage: V([1,3])
151
(1, 1)
152
sage: V([1,-3])
153
(1, 1)
154
"""
155
cdef Py_ssize_t i
156
cdef int xi
157
if isinstance(x, (list, tuple)):
158
if len(x) != self._degree:
159
raise TypeError("x must be a list of the right length")
160
for i from 0 <= i < self._degree:
161
if PY_TYPE_CHECK(x[i],IntegerMod_int) or PY_TYPE_CHECK(x[i],int) or PY_TYPE_CHECK(x[i],Integer):
162
xi = x[i]
163
# the if/else statement is because in some compilers, (-1)%2 is -1
164
mzd_write_bit(self._entries, 0, i, 0 if xi%2==0 else 1)
165
else:
166
mzd_write_bit(self._entries, 0, i, x[i]%2)
167
return
168
if x != 0:
169
raise TypeError("can't initialize vector from nonzero non-list")
170
else:
171
for i from 0 <= i < self._degree:
172
mzd_set_ui(self._entries, 0)
173
174
def __dealloc__(self):
175
"""
176
EXAMPLE::
177
178
sage: VS = VectorSpace(GF(2),10^3)
179
sage: import gc
180
sage: for i in range(10):
181
... v = VS.random_element()
182
... del v
183
... _ = gc.collect()
184
"""
185
if self._entries:
186
mzd_free(self._entries)
187
188
cdef int _cmp_c_impl(left, Element right) except -2:
189
"""
190
EXAMPLES::
191
sage: v = vector(GF(2), [0,0,0,0])
192
sage: v == 0
193
True
194
sage: v == 1
195
False
196
sage: v == v
197
True
198
sage: w = vector(GF(2), [1,0,0,0])
199
sage: w < v
200
False
201
sage: w > v
202
True
203
"""
204
if left._degree == 0:
205
return 0
206
return mzd_cmp(left._entries, (<Vector_mod2_dense>right)._entries)
207
208
def __len__(self):
209
"""
210
EXAMPLES::
211
212
sage: len(vector(GF(2),[0,0,1,1,1]))
213
5
214
"""
215
return self._degree
216
217
def __setitem__(self, i, value):
218
"""
219
EXAMPLES::
220
221
sage: VS = VectorSpace(GF(2),4)
222
sage: v = VS.random_element(); v
223
(1, 0, 0, 0)
224
sage: v[0] = 0; v
225
(0, 0, 0, 0)
226
sage: v[1:3] = [1, 1]; v
227
(0, 1, 1, 0)
228
sage: v[4] = 0
229
Traceback (most recent call last):
230
...
231
IndexError: Index '4' out of bound.
232
"""
233
if not self._is_mutable:
234
raise ValueError("vector is immutable; please change a copy instead (use copy())")
235
cdef IntegerMod_int m
236
cdef Py_ssize_t k, d, n
237
if isinstance(i, slice):
238
start, stop = i.start, i.stop
239
d = self.degree()
240
R = self.base_ring()
241
n = 0
242
for k from start <= k < stop:
243
if k >= d:
244
return
245
if k >= 0:
246
self[k] = R(value[n])
247
n = n + 1
248
else:
249
m = self.base_ring()(value)
250
if i < 0 or i >= self._degree:
251
raise IndexError("Index '%s' out of bound."%(i))
252
else:
253
mzd_write_bit(self._entries, 0, i, m)
254
255
def __getitem__(self, i):
256
"""
257
Returns `i`-th entry or slice of self.
258
259
EXAMPLES::
260
261
sage: v = vector(GF(2), [1,2,3]); v
262
(1, 0, 1)
263
sage: v[0]
264
1
265
sage: v[2]
266
1
267
sage: v[-2]
268
0
269
sage: v[0:2]
270
(1, 0)
271
sage: v[5]
272
Traceback (most recent call last):
273
...
274
IndexError: index '5' out of range
275
276
sage: v[-5]
277
Traceback (most recent call last):
278
...
279
IndexError: index '-2' out of range
280
"""
281
if isinstance(i, slice):
282
start, stop, step = i.indices(len(self))
283
return vector(self.base_ring(), self.list()[start:stop])
284
else:
285
if i < 0:
286
i += self._degree
287
if i < 0 or i >= self._degree:
288
raise IndexError("index '%s' out of range"%(i,))
289
return self._base_ring(mzd_read_bit(self._entries, 0, i))
290
291
def __reduce__(self):
292
"""
293
EXAMPLE::
294
295
sage: VS = VectorSpace(GF(2),10^4)
296
sage: e = VS.random_element()
297
sage: loads(dumps(e)) == e
298
True
299
"""
300
return unpickle_v0, (self._parent, self.list(), self._degree, self._is_mutable)
301
302
cpdef ModuleElement _add_(self, ModuleElement right):
303
"""
304
EXAMPLE::
305
306
sage: VS = VectorSpace(GF(2),10)
307
sage: e = VS([0,0,1,1,0,0,1,1,0,0])
308
sage: f = VS([0,1,0,1,0,1,0,1,0,1])
309
sage: e + f #indirect doctest
310
(0, 1, 1, 0, 0, 1, 1, 0, 0, 1)
311
"""
312
cdef Vector_mod2_dense z = self._new_c()
313
if self._degree:
314
mzd_add(z._entries, self._entries, (<Vector_mod2_dense>right)._entries)
315
return z
316
317
cpdef ModuleElement _sub_(self, ModuleElement right):
318
"""
319
EXAMPLE::
320
321
sage: VS = VectorSpace(GF(2),10)
322
sage: e = VS([0,0,1,1,0,0,1,1,0,0])
323
sage: f = VS([0,1,0,1,0,1,0,1,0,1])
324
sage: e - f #indirect doctest
325
(0, 1, 1, 0, 0, 1, 1, 0, 0, 1)
326
"""
327
cdef Vector_mod2_dense z = self._new_c()
328
if self._degree:
329
mzd_add(z._entries, self._entries, (<Vector_mod2_dense>right)._entries)
330
return z
331
332
cpdef Element _dot_product_(self, Vector right):
333
"""
334
EXAMPLES::
335
sage: VS = VectorSpace(GF(2),3)
336
sage: v = VS([1,1,1]); w = VS([0,0,0])
337
sage: v * w, w * v #indirect doctest
338
(0, 0)
339
sage: v = VS([1,1,1]); w = VS([0,1,0])
340
sage: v * w, w * v
341
(1, 1)
342
sage: v = VS([1,1,1]); w = VS([0,1,1])
343
sage: v * w, w * v
344
(0, 0)
345
sage: v = VS([1,1,1]); w = VS([1,1,1])
346
sage: v * w, w * v
347
(1, 1)
348
349
sage: VS = VectorSpace(GF(2),10^4)
350
sage: v = VS(0); w = VS(0)
351
sage: v[1337] = 1; w[1337] = 1
352
sage: v * w, w * v
353
(1, 1)
354
sage: v[9881] = 1; w[9881] = 1
355
sage: v * w, w * v
356
(0, 0)
357
sage: v[5172] = 1; w[6178] = 1
358
sage: v * w, w * v
359
(0, 0)
360
"""
361
cdef Py_ssize_t i
362
cdef IntegerMod_int n
363
cdef Vector_mod2_dense r = right
364
cdef m4ri_word tmp = 0
365
n = IntegerMod_int.__new__(IntegerMod_int)
366
IntegerMod_abstract.__init__(n, self.base_ring())
367
n.ivalue = 0
368
369
for i from 0 <= i < self._entries.width:
370
tmp ^= self._entries.rows[0][i] & r._entries.rows[0][i]
371
372
for i in range(64):
373
n.ivalue ^= <int>(tmp & 1)
374
tmp = tmp >> 1
375
376
return n
377
378
cpdef Vector _pairwise_product_(self, Vector right):
379
"""
380
EXAMPLE::
381
382
sage: VS = VectorSpace(GF(2),10)
383
sage: e = VS.random_element(); e
384
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
385
sage: f = VS.random_element(); f
386
(1, 1, 0, 1, 1, 1, 0, 0, 0, 1)
387
sage: e.pairwise_product(f) #indirect doctest
388
(1, 0, 0, 0, 1, 1, 0, 0, 0, 1)
389
"""
390
cdef Vector_mod2_dense z, r
391
r = right
392
z = self._new_c()
393
cdef Py_ssize_t i
394
for i from 0 <= i < self._entries.width:
395
z._entries.rows[0][i] = (self._entries.rows[0][i] & r._entries.rows[0][i])
396
return z
397
398
cpdef ModuleElement _rmul_(self, RingElement left):
399
"""
400
EXAMPLE::
401
402
sage: VS = VectorSpace(GF(2),10)
403
sage: e = VS.random_element(); e
404
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
405
sage: 0 * e
406
(0, 0, 0, 0, 0, 0, 0, 0, 0, 0)
407
sage: 1 * e
408
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
409
sage: 2 * e #indirect doctest
410
(0, 0, 0, 0, 0, 0, 0, 0, 0, 0)
411
"""
412
cdef IntegerMod_int a
413
414
if left:
415
return self.__copy__()
416
else:
417
return self._new_c()
418
419
420
cpdef ModuleElement _lmul_(self, RingElement right):
421
"""
422
EXAMPLE::
423
424
sage: VS = VectorSpace(GF(2),10)
425
sage: e = VS.random_element(); e
426
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
427
sage: e * 0 #indirect doctest
428
(0, 0, 0, 0, 0, 0, 0, 0, 0, 0)
429
sage: e * 1
430
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
431
sage: e * 2
432
(0, 0, 0, 0, 0, 0, 0, 0, 0, 0)
433
"""
434
return self._rmul_(right)
435
436
cpdef ModuleElement _neg_(self):
437
"""
438
EXAMPLE::
439
440
sage: VS = VectorSpace(GF(2),10)
441
sage: e = VS.random_element()
442
sage: -e == e
443
True
444
"""
445
return self.__copy__()
446
447
def n(self, *args, **kwargs):
448
"""
449
Returns a numerical approximation of ``self`` by calling the
450
:meth:`n()` method on all of its entries.
451
452
EXAMPLES::
453
454
sage: v = vector(GF(2), [1,2,3])
455
sage: v.n()
456
(1.00000000000000, 0.000000000000000, 1.00000000000000)
457
sage: _.parent()
458
Vector space of dimension 3 over Real Field with 53 bits of precision
459
sage: v.n(prec=75)
460
(1.000000000000000000000, 0.0000000000000000000000, 1.000000000000000000000)
461
sage: _.parent()
462
Vector space of dimension 3 over Real Field with 75 bits of precision
463
"""
464
return vector( [e.n(*args, **kwargs) for e in self] )
465
466
def list(self, copy=True):
467
"""
468
Return a list of entries in ``self``.
469
470
INPUT:
471
472
- ``copy`` - always ``True
473
474
EXAMPLE::
475
476
sage: VS = VectorSpace(GF(2),10)
477
sage: e = VS.random_element(); e
478
(1, 0, 0, 0, 1, 1, 1, 0, 0, 1)
479
sage: e.list()
480
[1, 0, 0, 0, 1, 1, 1, 0, 0, 1]
481
"""
482
cdef Py_ssize_t d = self._degree
483
cdef Py_ssize_t i
484
cdef list v = [0]*d
485
K = self.base_ring()
486
z = K.zero_element()
487
o = K.one_element()
488
cdef list switch = [z,o]
489
for i in range(d):
490
v[i] = switch[mzd_read_bit(self._entries, 0, i)]
491
return v
492
493
def unpickle_v0(parent, entries, degree, is_mutable):
494
"""
495
EXAMPLE::
496
497
sage: from sage.modules.vector_mod2_dense import unpickle_v0
498
sage: VS = VectorSpace(GF(2),10)
499
sage: unpickle_v0(VS, [0,1,2,3,4,5,6,7,8,9], 10, 0)
500
(0, 1, 0, 1, 0, 1, 0, 1, 0, 1)
501
"""
502
# If you think you want to change this function, don't.
503
cdef Vector_mod2_dense v
504
v = PY_NEW(Vector_mod2_dense)
505
v._init(degree, parent)
506
cdef int xi
507
508
for i from 0 <= i < degree:
509
if PY_TYPE_CHECK(entries[i],IntegerMod_int) or PY_TYPE_CHECK(entries[i],int) or PY_TYPE_CHECK(entries[i],Integer):
510
xi = entries[i]
511
mzd_write_bit(v._entries, 0, i, xi%2)
512
else:
513
mzd_write_bit(v._entries, 0, i, entries[i]%2)
514
v._is_mutable = int(is_mutable)
515
return v
516
517
518