Path: blob/master/sage/modules/vector_real_double_dense.pyx
4045 views
r"""1Dense real double vectors using a NumPy backend.23EXAMPLES:4sage: v = vector(RDF,[1, pi, sqrt(2)])5sage: v6(1.0, 3.14159265359, 1.41421356237)7sage: type(v)8<type 'sage.modules.vector_real_double_dense.Vector_real_double_dense'>9sage: parent(v)10Vector space of dimension 3 over Real Double Field11sage: v[0] = 512sage: v13(5.0, 3.14159265359, 1.41421356237)14sage: loads(dumps(v)) == v15True1617AUTHORS:18-- Jason Grout, Oct 2008: switch to numpy backend, factored out19Vector_double_dense class20"""2122##############################################################################23# Copyright (C) 2008 Jason Grout <[email protected]>24# Distributed under the terms of the GNU General Public License (GPL)25# The full text of the GPL is available at:26# http://www.gnu.org/licenses/27##############################################################################28from sage.rings.real_double import RDF2930cimport numpy as cnumpy3132numpy=None33scipy=None3435cdef class Vector_real_double_dense(vector_double_dense.Vector_double_dense):36"""37Vectors over the Real Double Field. These are supposed to be fast38vector operations using C doubles. Most operations are implemented39using numpy which will call the underlying BLAS, if needed, on the40system.4142EXAMPLES:43sage: v = vector(RDF, [1,2,3,4]); v44(1.0, 2.0, 3.0, 4.0)45sage: v*v4630.047"""48def __cinit__(self, parent, entries, coerce=True, copy=True):49global numpy50if numpy is None:51import numpy52self._numpy_dtype = numpy.dtype('float64')53self._numpy_dtypeint = cnumpy.NPY_DOUBLE54self._python_dtype = float55# TODO: Make RealDoubleElement instead of RDF for speed56self._sage_dtype = RDF57self.__create_vector__()58return5960def stats_skew(self):61"""62Computes the skewness of a data set.6364For normally distributed data, the skewness should be about650. A skewness value > 0 means that there is more weight in the66left tail of the distribution. (Paragraph from the scipy.stats67docstring.)6869EXAMPLE:70sage: v = vector(RDF, range(9))71sage: v.stats_skew()720.073"""74global scipy75if scipy is None:76import scipy77import scipy.stats78return self._sage_dtype(scipy.stats.skew(self._vector_numpy))798081def __reduce__(self):82"""83Pickling8485EXAMPLE:86sage: a = vector(RDF, range(9))87sage: loads(dumps(a)) == a88True89"""90return (unpickle_v1, (self._parent, self.list(), self._degree, self._is_mutable))919293# For backwards compatibility, we must keep the function unpickle_v094def unpickle_v0(parent, entries, degree):95"""96Create a real double vector containing the entries.9798EXAMPLE:99sage: v = vector(RDF, [1,2,3])100sage: w = sage.modules.vector_real_double_dense.unpickle_v0(v.parent(), list(v), v.degree())101sage: v == w102True103"""104return unpickle_v1(parent, entries, degree)105106def unpickle_v1(parent, entries, degree, is_mutable=None):107"""108Create a real double vector with the given parent, entries,109degree, and mutability.110111EXAMPLE:112sage: v = vector(RDF, [1,2,3])113sage: w = sage.modules.vector_real_double_dense.unpickle_v1(v.parent(), list(v), v.degree(), v.is_mutable())114sage: v == w115True116"""117cdef Vector_real_double_dense v = Vector_real_double_dense(parent, entries)118if is_mutable is not None:119v._is_mutable = is_mutable120return v121122123