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
""" 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 class Vector_rational_dense(free_module_element.FreeModuleElement): cdef bint is_dense_c(self): cdef bint is_sparse_c(self):
def __copy__(self): cdef Vector_rational_dense y cdef Py_ssize_t i
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 Py_ssize_t i
def __cinit__(self, parent=None, x=None, coerce=True, copy=True):
def __init__(self, parent, x, coerce=True, copy=True): cdef Py_ssize_t i cdef Rational z
def __dealloc__(self): cdef Py_ssize_t i # Do *not* use sig_on() here, since __dealloc__ # cannot raise exceptions!
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 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) """
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
def __reduce__(self):
cpdef _add_(self, right): cdef Vector_rational_dense z, r cdef Py_ssize_t i
cpdef _sub_(self, right): cdef Vector_rational_dense z, r cdef Py_ssize_t i
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 Rational z cdef mpq_t t cdef Py_ssize_t i
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 cdef Py_ssize_t i
cpdef _rmul_(self, Element left): cdef Vector_rational_dense z cdef Rational a else: # should not happen raise TypeError("Cannot convert %s to %s" % (type(left).__name__, Rational.__name__)) cdef Py_ssize_t i
cpdef _lmul_(self, Element right): cdef Vector_rational_dense z cdef Rational a else: # should not happen raise TypeError("Cannot convert %s to %s" % (type(right).__name__, Rational.__name__)) cdef Py_ssize_t i
cpdef _neg_(self): cdef Vector_rational_dense z cdef Py_ssize_t i
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 cdef Rational z cdef Py_ssize_t i |