Path: blob/master/sage/modules/vector_rational_dense.pyx
4045 views
"""1Vectors with rational entries.23AUTHOR:45- William Stein (2007)6- Soroosh Yazdani (2007)78EXAMPLES::910sage: v = vector(QQ,[1,2,3,4,5])11sage: v12(1, 2, 3, 4, 5)13sage: 3*v14(3, 6, 9, 12, 15)15sage: v/216(1/2, 1, 3/2, 2, 5/2)17sage: -v18(-1, -2, -3, -4, -5)19sage: v - v20(0, 0, 0, 0, 0)21sage: v + v22(2, 4, 6, 8, 10)23sage: v * v24552526We make a large zero vector::2728sage: k = QQ^100000; k29Vector space of dimension 100000 over Rational Field30sage: v = k(0)31sage: v[:10]32(0, 0, 0, 0, 0, 0, 0, 0, 0, 0)3334TESTS::3536sage: v = vector(QQ, [1,2/5,-3/8,4])37sage: loads(dumps(v)) == v38True39"""4041###############################################################################42# Sage: System for Algebra and Geometry Experimentation43# Copyright (C) 2007 William Stein <[email protected]>44# Distributed under the terms of the GNU General Public License (GPL)45# http://www.gnu.org/licenses/46###############################################################################4748include '../ext/interrupt.pxi'49include '../ext/stdsage.pxi'5051from sage.structure.element cimport Element, ModuleElement, RingElement, Vector5253from sage.rings.integer cimport Integer54from sage.rings.rational cimport Rational5556cimport free_module_element57from free_module_element import vector5859cdef inline _Rational_from_mpq(mpq_t e):60cdef Rational z = PY_NEW(Rational)61mpq_set(z.value, e)62return z6364cdef class Vector_rational_dense(free_module_element.FreeModuleElement):65cdef bint is_dense_c(self):66return 167cdef bint is_sparse_c(self):68return 06970cdef _new_c(self):71cdef Vector_rational_dense y72y = PY_NEW(Vector_rational_dense)73y._init(self._degree, self._parent)74return y7576def __copy__(self):77cdef Vector_rational_dense y78y = self._new_c()79cdef Py_ssize_t i80for i from 0 <= i < self._degree:81mpq_init(y._entries[i])82mpq_set(y._entries[i], self._entries[i])83return y8485cdef _init(self, Py_ssize_t degree, parent):86self._degree = degree87self._parent = parent88self._entries = <mpq_t *> sage_malloc(sizeof(mpq_t) * degree)89if self._entries == NULL:90raise MemoryError9192def __cinit__(self, parent=None, x=None, coerce=True,copy=True):93self._entries = NULL94self._is_mutable = 195if not parent is None:96self._init(parent.degree(), parent)9798def __init__(self, parent, x, coerce=True, copy=True):99cdef Py_ssize_t i100cdef Rational z101# we have to do this to avoid a garbage collection error in dealloc102for i from 0 <= i < self._degree:103mpq_init(self._entries[i])104if isinstance(x, (list, tuple)):105if len(x) != self._degree:106raise TypeError("entries must be a list of length %s"%self._degree)107for i from 0 <= i < self._degree:108z = Rational(x[i])109mpq_set(self._entries[i], z.value)110return111if x != 0:112raise TypeError("can't initialize vector from nonzero non-list")113114def __dealloc__(self):115cdef Py_ssize_t i116if self._entries:117sig_on()118for i from 0 <= i < self._degree:119#print "clearing gmp's entry %s"%i120mpq_clear(self._entries[i])121sig_off()122#print "clearing python entries"123sage_free(self._entries)124125cdef int _cmp_c_impl(left, Element right) except -2:126"""127EXAMPLES::128129sage: v = vector(QQ, [0,0,0,0])130sage: v == 0131True132sage: v == 1133False134sage: v == v135True136sage: w = vector(QQ, [-1,3/2,0,0])137sage: w < v138True139sage: w > v140False141"""142cdef Py_ssize_t i143cdef int c144for i from 0 <= i < left.degree():145c = mpq_cmp(left._entries[i], (<Vector_rational_dense>right)._entries[i])146if c < 0:147return -1148elif c > 0:149return 1150return 0151152def __len__(self):153return self._degree154155def __setitem__(self, i, value):156if not self._is_mutable:157raise ValueError("vector is immutable; please change a copy instead (use copy())")158cdef Rational z159cdef Py_ssize_t k, d, n160if isinstance(i, slice):161start, stop = i.start, i.stop162d = self.degree()163R = self.base_ring()164n = 0165for k from start <= k < stop:166if k >= d:167return168if k >= 0:169self[k] = R(value[n])170n = n + 1171else:172if i < 0 or i >= self._degree:173raise IndexError174else:175z = Rational(value)176mpq_set(self._entries[i], z.value)177178def __getitem__(self, i):179"""180Returns `i`-th entry or slice of self.181182EXAMPLES::183184sage: v = vector([1/2,2/3,3/4]); v185(1/2, 2/3, 3/4)186sage: v[0]1871/2188sage: v[2]1893/4190sage: v[-2]1912/3192sage: v[5]193Traceback (most recent call last):194...195IndexError: index out of range196sage: v[-5]197Traceback (most recent call last):198...199IndexError: index out of range200"""201cdef Rational z = PY_NEW(Rational)202if isinstance(i, slice):203start, stop, step = i.indices(len(self))204return vector(self.base_ring(), self.list()[start:stop])205else:206if i < 0:207i += self._degree208if i < 0 or i >= self._degree:209raise IndexError('index out of range')210else:211mpq_set(z.value, self._entries[i])212return z213214def list(self,copy=True):215"""216The list of entries of the vector.217218INPUT:219220- ``copy``, ignored optional argument.221222EXAMPLES::223224sage: v = vector(QQ,[1,2,3,4])225sage: a = v.list(copy=False); a226[1, 2, 3, 4]227sage: a[0] = 0228sage: v229(1, 2, 3, 4)230"""231cdef int i232return [_Rational_from_mpq(self._entries[i]) for i in233xrange(self._degree)]234235def __reduce__(self):236return (unpickle_v1, (self._parent, self.list(), self._degree, self._is_mutable))237238cpdef ModuleElement _add_(self, ModuleElement right):239cdef Vector_rational_dense z, r240r = right241z = self._new_c()242cdef Py_ssize_t i243for i from 0 <= i < self._degree:244mpq_init(z._entries[i])245mpq_add(z._entries[i], self._entries[i], r._entries[i])246return z247248249cpdef ModuleElement _sub_(self, ModuleElement right):250cdef Vector_rational_dense z, r251r = right252z = self._new_c()253cdef Py_ssize_t i254for i from 0 <= i < self._degree:255mpq_init(z._entries[i])256mpq_sub(z._entries[i], self._entries[i], r._entries[i])257return z258259cpdef Element _dot_product_(self, Vector right):260"""261Dot product of dense vectors over the rationals.262263EXAMPLES::264265sage: v = vector(QQ, [1,2,-3]); w = vector(QQ,[4,3,2])266sage: v*w2674268sage: w*v2694270"""271cdef Vector_rational_dense r = right272cdef Rational z273z = PY_NEW(Rational)274cdef mpq_t t275mpq_init(t)276mpq_set_si(z.value, 0, 1)277cdef Py_ssize_t i278for i from 0 <= i < self._degree:279mpq_mul(t, self._entries[i], r._entries[i])280mpq_add(z.value, z.value, t)281mpq_clear(t)282return z283284285cpdef Vector _pairwise_product_(self, Vector right):286"""287EXAMPLES::288289sage: v = vector(QQ, [1,2,-3]); w = vector(QQ,[4,3,2])290sage: v.pairwise_product(w)291(4, 6, -6)292"""293cdef Vector_rational_dense z, r294r = right295z = self._new_c()296cdef Py_ssize_t i297for i from 0 <= i < self._degree:298mpq_init(z._entries[i])299mpq_mul(z._entries[i], self._entries[i], r._entries[i])300return z301302cpdef ModuleElement _rmul_(self, RingElement left):303cdef Vector_rational_dense z304cdef Rational a305if PY_TYPE_CHECK(left, Rational):306a = <Rational>left307elif PY_TYPE_CHECK(left, Integer):308a = <Rational>PY_NEW(Rational)309mpq_set_z(a.value, (<Integer>left).value)310else:311# should not happen312raise TypeError("Cannot convert %s to %s" % (type(left).__name__, Rational.__name__))313z = self._new_c()314cdef Py_ssize_t i315for i from 0 <= i < self._degree:316mpq_init(z._entries[i])317mpq_mul(z._entries[i], self._entries[i], a.value)318return z319320321cpdef ModuleElement _lmul_(self, RingElement right):322cdef Vector_rational_dense z323cdef Rational a324if PY_TYPE_CHECK(right, Rational):325a = <Rational>right326elif PY_TYPE_CHECK(right, Integer):327a = <Rational>PY_NEW(Rational)328mpq_set_z(a.value, (<Integer>right).value)329else:330# should not happen331raise TypeError("Cannot convert %s to %s" % (type(right).__name__, Rational.__name__))332z = self._new_c()333cdef Py_ssize_t i334for i from 0 <= i < self._degree:335mpq_init(z._entries[i])336mpq_mul(z._entries[i], self._entries[i], a.value)337return z338339cpdef ModuleElement _neg_(self):340cdef Vector_rational_dense z341z = self._new_c()342cdef Py_ssize_t i343for i from 0 <= i < self._degree:344mpq_init(z._entries[i])345mpq_neg(z._entries[i], self._entries[i])346return z347348349def unpickle_v0(parent, entries, degree):350# If you think you want to change this function, don't.351# Instead make a new version with a name like352# make_FreeModuleElement_generic_dense_v1353# and changed the reduce method below.354cdef Vector_rational_dense v355v = PY_NEW(Vector_rational_dense)356v._init(degree, parent)357cdef Rational z358for i from 0 <= i < degree:359z = Rational(entries[i])360mpq_init(v._entries[i])361mpq_set(v._entries[i], z.value)362return v363364def unpickle_v1(parent, entries, degree, is_mutable):365cdef Vector_rational_dense v366v = PY_NEW(Vector_rational_dense)367v._init(degree, parent)368cdef Rational z369for i from 0 <= i < degree:370z = Rational(entries[i])371mpq_init(v._entries[i])372mpq_set(v._entries[i], z.value)373v._is_mutable = is_mutable374return v375376377