Elements of free modules#
AUTHORS:
William Stein
Josh Kantor
Thomas Feulner (2012-11): Added
FreeModuleElement.hamming_weight()
andFreeModuleElement_generic_sparse.hamming_weight()
Jeroen Demeyer (2015-02-24): Implement fast Cython methods
get_unsafe
andset_unsafe
similar to other places in Sage (trac ticket #17562)
EXAMPLES: We create a vector space over \(\QQ\) and a subspace of this space.
sage: V = QQ^5
sage: W = V.span([V.1, V.2])
Arithmetic operations always return something in the ambient space, since there is a canonical map from \(W\) to \(V\) but not from \(V\) to \(W\).
sage: parent(W.0 + V.1)
Vector space of dimension 5 over Rational Field
sage: parent(V.1 + W.0)
Vector space of dimension 5 over Rational Field
sage: W.0 + V.1
(0, 2, 0, 0, 0)
sage: W.0 - V.0
(-1, 1, 0, 0, 0)
Next we define modules over \(\ZZ\) and a finite field.
sage: K = ZZ^5
sage: M = GF(7)^5
Arithmetic between the \(\QQ\) and \(\ZZ\) modules is defined, and the result is always over \(\QQ\), since there is a canonical coercion map to \(\QQ\).
sage: K.0 + V.1
(1, 1, 0, 0, 0)
sage: parent(K.0 + V.1)
Vector space of dimension 5 over Rational Field
Since there is no canonical coercion map to the finite field from \(\QQ\) the following arithmetic is not defined:
sage: V.0 + M.0
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for +: 'Vector space of dimension 5 over Rational Field' and 'Vector space of dimension 5 over Finite Field of size 7'
However, there is a map from \(\ZZ\) to the finite field, so the following is defined, and the result is in the finite field.
sage: w = K.0 + M.0; w
(2, 0, 0, 0, 0)
sage: parent(w)
Vector space of dimension 5 over Finite Field of size 7
sage: parent(M.0 + K.0)
Vector space of dimension 5 over Finite Field of size 7
Matrix vector multiply:
sage: MS = MatrixSpace(QQ,3)
sage: A = MS([0,1,0,1,0,0,0,0,1])
sage: V = QQ^3
sage: v = V([1,2,3])
sage: v * A
(2, 1, 3)
- class sage.modules.free_module_element.FreeModuleElement#
Bases:
sage.structure.element.Vector
An element of a generic free module.
- Mod(p)#
EXAMPLES:
sage: V = vector(ZZ, [5, 9, 13, 15]) sage: V.Mod(7) (5, 2, 6, 1) sage: parent(V.Mod(7)) Vector space of dimension 4 over Ring of integers modulo 7
- additive_order()#
Return the additive order of self.
EXAMPLES:
sage: v = vector(Integers(4), [1,2]) sage: v.additive_order() 4
sage: v = vector([1,2,3]) sage: v.additive_order() +Infinity
sage: v = vector(Integers(30), [6, 15]); v (6, 15) sage: v.additive_order() 10 sage: 10*v (0, 0)
- apply_map(phi, R=None, sparse=None)#
Apply the given map phi (an arbitrary Python function or callable object) to this free module element. If R is not given, automatically determine the base ring of the resulting element.
- INPUT:
- sparse – True or False will control whether the result
is sparse. By default, the result is sparse iff self is sparse.
phi
- arbitrary Python function or callable objectR
- (optional) ring
OUTPUT: a free module element over R
EXAMPLES:
sage: m = vector([1,x,sin(x+1)]) sage: m.apply_map(lambda x: x^2) (1, x^2, sin(x + 1)^2) sage: m.apply_map(sin) (sin(1), sin(x), sin(sin(x + 1)))
sage: m = vector(ZZ, 9, range(9)) sage: k.<a> = GF(9) sage: m.apply_map(k) (0, 1, 2, 0, 1, 2, 0, 1, 2)
In this example, we explicitly specify the codomain.
sage: s = GF(3) sage: f = lambda x: s(x) sage: n = m.apply_map(f, k); n (0, 1, 2, 0, 1, 2, 0, 1, 2) sage: n.parent() Vector space of dimension 9 over Finite Field in a of size 3^2
If your map sends 0 to a non-zero value, then your resulting vector is not mathematically sparse:
sage: v = vector([0] * 6 + [1], sparse=True); v (0, 0, 0, 0, 0, 0, 1) sage: v2 = v.apply_map(lambda x: x+1); v2 (1, 1, 1, 1, 1, 1, 2)
but it’s still represented with a sparse data type:
sage: parent(v2) Ambient sparse free module of rank 7 over the principal ideal domain Integer Ring
This data type is inefficient for dense vectors, so you may want to specify sparse=False:
sage: v2 = v.apply_map(lambda x: x+1, sparse=False); v2 (1, 1, 1, 1, 1, 1, 2) sage: parent(v2) Ambient free module of rank 7 over the principal ideal domain Integer Ring
Or if you have a map that will result in mostly zeroes, you may want to specify sparse=True:
sage: v = vector(srange(10)) sage: v2 = v.apply_map(lambda x: 0 if x else 1, sparse=True); v2 (1, 0, 0, 0, 0, 0, 0, 0, 0, 0) sage: parent(v2) Ambient sparse free module of rank 10 over the principal ideal domain Integer Ring
- change_ring(R)#
Change the base ring of this vector.
EXAMPLES:
sage: v = vector(QQ['x,y'], [1..5]); v.change_ring(GF(3)) (1, 2, 0, 1, 2)
- column()#
Return a matrix with a single column and the same entries as the vector
self
.OUTPUT:
A matrix over the same ring as the vector (or free module element), with a single column. The entries of the column are identical to those of the vector, and in the same order.
EXAMPLES:
sage: v = vector(ZZ, [1,2,3]) sage: w = v.column(); w [1] [2] [3] sage: w.parent() Full MatrixSpace of 3 by 1 dense matrices over Integer Ring sage: x = vector(FiniteField(13), [2,4,8,16]) sage: x.column() [2] [4] [8] [3]
There is more than one way to get one-column matrix from a vector. The
column
method is about equally efficient to making a row and then taking a transpose. Notice that supplying a vector to the matrix constructor demonstrates Sage’s preference for rows.sage: x = vector(RDF, [sin(i*pi/20) for i in range(10)]) sage: x.column() == matrix(x).transpose() True sage: x.column() == x.row().transpose() True
Sparse or dense implementations are preserved.
sage: d = vector(RR, [1.0, 2.0, 3.0]) sage: s = vector(CDF, {2:5.0+6.0*I}) sage: dm = d.column() sage: sm = s.column() sage: all([d.is_dense(), dm.is_dense(), s.is_sparse(), sm.is_sparse()]) True
- conjugate()#
Returns a vector where every entry has been replaced by its complex conjugate.
OUTPUT:
A vector of the same length, over the same ring, but with each entry replaced by the complex conjugate, as implemented by the
conjugate()
method for elements of the base ring, which is presently always complex conjugation.EXAMPLES:
sage: v = vector(CDF, [2.3 - 5.4*I, -1.7 + 3.6*I]) sage: w = v.conjugate(); w (2.3 + 5.4*I, -1.7 - 3.6*I) sage: w.parent() Vector space of dimension 2 over Complex Double Field
Even if conjugation seems nonsensical over a certain ring, this method for vectors cooperates silently.
sage: u = vector(ZZ, range(6)) sage: u.conjugate() (0, 1, 2, 3, 4, 5)
Sage implements a few specialized subfields of the complex numbers, such as the cyclotomic fields. This example uses such a field containing a primitive 7-th root of unity named
a
.sage: F.<a> = CyclotomicField(7) sage: v = vector(F, [a^i for i in range(7)]) sage: v (1, a, a^2, a^3, a^4, a^5, -a^5 - a^4 - a^3 - a^2 - a - 1) sage: v.conjugate() (1, -a^5 - a^4 - a^3 - a^2 - a - 1, a^5, a^4, a^3, a^2, a)
Sparse vectors are returned as such.
sage: v = vector(CC, {1: 5 - 6*I, 3: -7*I}); v (0.000000000000000, 5.00000000000000 - 6.00000000000000*I, 0.000000000000000, -7.00000000000000*I) sage: v.is_sparse() True sage: vc = v.conjugate(); vc (0.000000000000000, 5.00000000000000 + 6.00000000000000*I, 0.000000000000000, 7.00000000000000*I) sage: vc.conjugate() (0.000000000000000, 5.00000000000000 - 6.00000000000000*I, 0.000000000000000, -7.00000000000000*I)
- coordinate_ring()#
Return the ring from which the coefficients of this vector come.
This is different from
base_ring()
, which returns the ring of scalars.EXAMPLES:
sage: M = (ZZ^2) * (1/2) sage: v = M([0,1/2]) sage: v.base_ring() Integer Ring sage: v.coordinate_ring() Rational Field
- cross_product(right)#
Return the cross product of self and right, which is only defined for vectors of length 3 or 7.
INPUT:
right
- A vector of the same size asself
, either degree three or degree seven.
OUTPUT:
The cross product (vector product) of
self
andright
, a vector of the same size ofself
andright
.This product is performed under the assumption that the basis vectors are orthonormal. See the method
cross_product()
of vector fields for more general cases.EXAMPLES:
sage: v = vector([1,2,3]); w = vector([0,5,-9]) sage: v.cross_product(v) (0, 0, 0) sage: u = v.cross_product(w); u (-33, 9, 5) sage: u.dot_product(v) 0 sage: u.dot_product(w) 0
The cross product is defined for degree seven vectors as well: see Wikipedia article Cross_product. The 3-D cross product is achieved using the quaternions, whereas the 7-D cross product is achieved using the octonions.
sage: u = vector(QQ, [1, -1/3, 57, -9, 56/4, -4,1]) sage: v = vector(QQ, [37, 55, -99/57, 9, -12, 11/3, 4/98]) sage: u.cross_product(v) (1394815/2793, -2808401/2793, 39492/49, -48737/399, -9151880/2793, 62513/2793, -326603/171)
The degree seven cross product is anticommutative.
sage: u.cross_product(v) + v.cross_product(u) (0, 0, 0, 0, 0, 0, 0)
The degree seven cross product is distributive across addition.
sage: v = vector([-12, -8/9, 42, 89, -37, 60/99, 73]) sage: u = vector([31, -42/7, 97, 80, 30/55, -32, 64]) sage: w = vector([-25/4, 40, -89, -91, -72/7, 79, 58]) sage: v.cross_product(u + w) - (v.cross_product(u) + v.cross_product(w)) (0, 0, 0, 0, 0, 0, 0)
The degree seven cross product respects scalar multiplication.
sage: v = vector([2, 17, -11/5, 21, -6, 2/17, 16]) sage: u = vector([-8, 9, -21, -6, -5/3, 12, 99]) sage: (5*v).cross_product(u) - 5*(v.cross_product(u)) (0, 0, 0, 0, 0, 0, 0) sage: v.cross_product(5*u) - 5*(v.cross_product(u)) (0, 0, 0, 0, 0, 0, 0) sage: (5*v).cross_product(u) - (v.cross_product(5*u)) (0, 0, 0, 0, 0, 0, 0)
The degree seven cross product respects the scalar triple product.
sage: v = vector([2,6,-7/4,-9/12,-7,12,9]) sage: u = vector([22,-7,-9/11,12,15,15/7,11]) sage: w = vector([-11,17,19,-12/5,44,21/56,-8]) sage: v.dot_product(u.cross_product(w)) - w.dot_product(v.cross_product(u)) 0
AUTHOR:
Billy Wonderly (2010-05-11), Added 7-D Cross Product
- cross_product_matrix()#
Return the matrix which describes a cross product between
self
and some other vector.This operation is sometimes written using the hat operator: see Wikipedia article Hat_operator#Cross_product. It is only defined for vectors of length 3 or 7. For a vector \(v\) the cross product matrix \(\hat{v}\) is a matrix which satisfies \(\hat{v} \cdot w = v \times w\) and also \(w \cdot \hat{v} = w \times v\) for all vectors \(w\). The basis vectors are assumed to be orthonormal.
OUTPUT:
The cross product matrix of this vector.
EXAMPLES:
sage: v = vector([1, 2, 3]) sage: vh = v.cross_product_matrix() sage: vh [ 0 -3 2] [ 3 0 -1] [-2 1 0] sage: w = random_vector(3, x=1, y=100) sage: vh*w == v.cross_product(w) True sage: w*vh == w.cross_product(v) True sage: vh.is_alternating() True
- curl(variables=None)#
Return the curl of this two-dimensional or three-dimensional vector function.
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: vector([-y, x, 0]).curl() (0, 0, 2) sage: vector([y, -x, x*y*z]).curl() (x*z, -y*z, -2) sage: vector([y^2, 0, 0]).curl() (0, 0, -2*y) sage: (R^3).random_element().curl().div() 0
For rings where the variable order is not well defined, it must be defined explicitly:
sage: v = vector(SR, [-y, x, 0]) sage: v.curl() Traceback (most recent call last): ... ValueError: Unable to determine ordered variable names for Symbolic Ring sage: v.curl([x, y, z]) (0, 0, 2)
Note that callable vectors have well defined variable orderings:
sage: v(x, y, z) = (-y, x, 0) sage: v.curl() (x, y, z) |--> (0, 0, 2)
In two-dimensions, this returns a scalar value:
sage: R.<x,y> = QQ[] sage: vector([-y, x]).curl() 2
See also
curl()
of vector fields on Euclidean spaces (and more generally pseudo-Riemannian manifolds), in particular for computing the curl in curvilinear coordinates.
- degree()#
Return the degree of this vector, which is simply the number of entries.
EXAMPLES:
sage: sage.modules.free_module_element.FreeModuleElement(QQ^389).degree() 389 sage: vector([1,2/3,8]).degree() 3
- denominator()#
Return the least common multiple of the denominators of the entries of self.
EXAMPLES:
sage: v = vector([1/2,2/5,3/14]) sage: v.denominator() 70 sage: 2*5*7 70
sage: M = (ZZ^2)*(1/2) sage: M.basis()[0].denominator() 2
- dense_vector()#
Return dense version of self. If self is dense, just return self; otherwise, create and return correspond dense vector.
EXAMPLES:
sage: vector([-1,0,3,0,0,0]).dense_vector().is_dense() True sage: vector([-1,0,3,0,0,0],sparse=True).dense_vector().is_dense() True sage: vector([-1,0,3,0,0,0],sparse=True).dense_vector() (-1, 0, 3, 0, 0, 0)
- derivative(*args)#
Derivative with respect to variables supplied in args.
Multiple variables and iteration counts may be supplied; see documentation for the global derivative() function for more details.
diff()
is an alias of this function.EXAMPLES:
sage: v = vector([1,x,x^2]) sage: v.derivative(x) (0, 1, 2*x) sage: type(v.derivative(x)) == type(v) True sage: v = vector([1,x,x^2], sparse=True) sage: v.derivative(x) (0, 1, 2*x) sage: type(v.derivative(x)) == type(v) True sage: v.derivative(x,x) (0, 0, 2)
- dict(copy=True)#
Return dictionary of nonzero entries of
self
.More precisely, this returns a dictionary whose keys are indices of basis elements in the support of
self
and whose values are the corresponding coefficients.INPUT:
copy
– (default:True
) ifself
is internally represented by a dictionaryd
, then make a copy ofd
; ifFalse
, then this can cause undesired behavior by mutatingd
OUTPUT:
Python dictionary
EXAMPLES:
sage: v = vector([0,0,0,0,1/2,0,3/14]) sage: v.dict() {4: 1/2, 6: 3/14} sage: sorted(v.support()) [4, 6]
In some cases, when
copy=False
, we get back a dangerous reference:sage: v = vector({0:5, 2:3/7}, sparse=True) sage: v.dict(copy=False) {0: 5, 2: 3/7} sage: v.dict(copy=False)[0] = 18 sage: v (18, 0, 3/7)
- diff(*args)#
Derivative with respect to variables supplied in args.
Multiple variables and iteration counts may be supplied; see documentation for the global derivative() function for more details.
diff()
is an alias of this function.EXAMPLES:
sage: v = vector([1,x,x^2]) sage: v.derivative(x) (0, 1, 2*x) sage: type(v.derivative(x)) == type(v) True sage: v = vector([1,x,x^2], sparse=True) sage: v.derivative(x) (0, 1, 2*x) sage: type(v.derivative(x)) == type(v) True sage: v.derivative(x,x) (0, 0, 2)
- div(variables=None)#
Return the divergence of this vector function.
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: vector([x, y, z]).div() 3 sage: vector([x*y, y*z, z*x]).div() x + y + z sage: R.<x,y,z,w> = QQ[] sage: vector([x*y, y*z, z*x]).div([x, y, z]) x + y + z sage: vector([x*y, y*z, z*x]).div([z, x, y]) 0 sage: vector([x*y, y*z, z*x]).div([x, y, w]) y + z sage: vector(SR, [x*y, y*z, z*x]).div() Traceback (most recent call last): ... ValueError: Unable to determine ordered variable names for Symbolic Ring sage: vector(SR, [x*y, y*z, z*x]).div([x, y, z]) x + y + z
See also
divergence()
of vector fields on Euclidean spaces (and more generally pseudo-Riemannian manifolds), in particular for computing the divergence in curvilinear coordinates.
- dot_product(right)#
Return the dot product of
self
andright
, which is the sum of the product of the corresponding entries.INPUT:
right
– a vector of the same degree asself
. It does not need to belong to the same parent asself
, so long as the necessary products and sums are defined.
OUTPUT:
If
self
andright
are the vectors \(\vec{x}\) and \(\vec{y}\), of degree \(n\), then this method returns\[\sum_{i=1}^{n}x_iy_i\]Note
The
inner_product()
is a more general version of this method, and thehermitian_inner_product()
method may be more appropriate if your vectors have complex entries.EXAMPLES:
sage: V = FreeModule(ZZ, 3) sage: v = V([1,2,3]) sage: w = V([4,5,6]) sage: v.dot_product(w) 32
sage: R.<x> = QQ[] sage: v = vector([x,x^2,3*x]); w = vector([2*x,x,3+x]) sage: v*w x^3 + 5*x^2 + 9*x sage: (x*2*x) + (x^2*x) + (3*x*(3+x)) x^3 + 5*x^2 + 9*x sage: w*v x^3 + 5*x^2 + 9*x
The vectors may be from different vector spaces, provided the necessary operations make sense. Notice that coercion will generate a result of the same type, even if the order of the arguments is reversed.:
sage: v = vector(ZZ, [1,2,3]) sage: w = vector(FiniteField(3), [0,1,2]) sage: ip = w.dot_product(v); ip 2 sage: ip.parent() Finite Field of size 3 sage: ip = v.dot_product(w); ip 2 sage: ip.parent() Finite Field of size 3
The dot product of a vector with itself is the 2-norm, squared.
sage: v = vector(QQ, [3, 4, 7]) sage: v.dot_product(v) - v.norm()^2 0
- element()#
Simply returns self. This is useful, since for many objects, self.element() returns a vector corresponding to self.
EXAMPLES:
sage: v = vector([1/2,2/5,0]); v (1/2, 2/5, 0) sage: v.element() (1/2, 2/5, 0)
- get(i)#
Like
__getitem__
but without bounds checking: \(i\) must satisfy0 <= i < self.degree
.EXAMPLES:
sage: vector(SR, [1/2,2/5,0]).get(0) 1/2
- hamming_weight()#
Return the number of positions
i
such thatself[i] != 0
.EXAMPLES:
sage: vector([-1,0,3,0,0,0,0.01]).hamming_weight() 3
- hermitian_inner_product(right)#
Returns the dot product, but with the entries of the first vector conjugated beforehand.
INPUT:
right
- a vector of the same degree asself
OUTPUT:
If
self
andright
are the vectors \(\vec{x}\) and \(\vec{y}\) of degree \(n\) then this routine computes\[\sum_{i=1}^{n}\overline{x}_i{y}_i\]where the bar indicates complex conjugation.
Note
If your vectors do not contain complex entries, then
dot_product()
will return the same result without the overhead of conjugating elements ofself
.If you are not computing a weighted inner product, and your vectors do not have complex entries, then the
dot_product()
will return the same result.EXAMPLES:
sage: v = vector(CDF, [2+3*I, 5-4*I]) sage: w = vector(CDF, [6-4*I, 2+3*I]) sage: v.hermitian_inner_product(w) -2.0 - 3.0*I
Sage implements a few specialized fields over the complex numbers, such as cyclotomic fields and quadratic number fields. So long as the base rings have a conjugate method, then the Hermitian inner product will be available.
sage: Q.<a> = QuadraticField(-7) sage: a^2 -7 sage: v = vector(Q, [3+a, 5-2*a]) sage: w = vector(Q, [6, 4+3*a]) sage: v.hermitian_inner_product(w) 17*a - 4
The Hermitian inner product should be additive in each argument (we only need to test one), linear in each argument (with conjugation on the first scalar), and anti-commutative.
sage: alpha = CDF(5.0 + 3.0*I) sage: u = vector(CDF, [2+4*I, -3+5*I, 2-7*I]) sage: v = vector(CDF, [-1+3*I, 5+4*I, 9-2*I]) sage: w = vector(CDF, [8+3*I, -4+7*I, 3-6*I]) sage: (u+v).hermitian_inner_product(w) == u.hermitian_inner_product(w) + v.hermitian_inner_product(w) True sage: (alpha*u).hermitian_inner_product(w) == alpha.conjugate()*u.hermitian_inner_product(w) True sage: u.hermitian_inner_product(alpha*w) == alpha*u.hermitian_inner_product(w) True sage: u.hermitian_inner_product(v) == v.hermitian_inner_product(u).conjugate() True
For vectors with complex entries, the Hermitian inner product has a more natural relationship with the 2-norm (which is the default for the
norm()
method). The norm squared equals the Hermitian inner product of the vector with itself.sage: v = vector(CDF, [-0.66+0.47*I, -0.60+0.91*I, -0.62-0.87*I, 0.53+0.32*I]) sage: abs(v.norm()^2 - v.hermitian_inner_product(v)) < 1.0e-10 True
- inner_product(right)#
Returns the inner product of
self
andright
, possibly using an inner product matrix from the parent ofself
.INPUT:
right
- a vector of the same degree asself
OUTPUT:
If the parent vector space does not have an inner product matrix defined, then this is the usual dot product (
dot_product()
). Ifself
andright
are considered as single column matrices, \(\vec{x}\) and \(\vec{y}\), and \(A\) is the inner product matrix, then this method computes\[\left(\vec{x}\right)^tA\vec{y}\]where \(t\) indicates the transpose.
Note
If your vectors have complex entries, the
hermitian_inner_product()
may be more appropriate for your purposes.EXAMPLES:
sage: v = vector(QQ, [1,2,3]) sage: w = vector(QQ, [-1,2,-3]) sage: v.inner_product(w) -6 sage: v.inner_product(w) == v.dot_product(w) True
The vector space or free module that is the parent to
self
can have an inner product matrix defined, which will be used by this method. This matrix will be passed through to subspaces.sage: ipm = matrix(ZZ,[[2,0,-1], [0,2,0], [-1,0,6]]) sage: M = FreeModule(ZZ, 3, inner_product_matrix = ipm) sage: v = M([1,0,0]) sage: v.inner_product(v) 2 sage: K = M.span_of_basis([[0/2,-1/2,-1/2], [0,1/2,-1/2],[2,0,0]]) sage: (K.0).inner_product(K.0) 2 sage: w = M([1,3,-1]) sage: v = M([2,-4,5]) sage: w.row()*ipm*v.column() == w.inner_product(v) True
Note that the inner product matrix comes from the parent of
self
. So if a vector is not an element of the correct parent, the result could be a source of confusion.sage: V = VectorSpace(QQ, 2, inner_product_matrix=[[1,2],[2,1]]) sage: v = V([12, -10]) sage: w = vector(QQ, [10,12]) sage: v.inner_product(w) 88 sage: w.inner_product(v) 0 sage: w = V(w) sage: w.inner_product(v) 88
Note
The use of an inner product matrix makes no restrictions on the nature of the matrix. In particular, in this context it need not be Hermitian and positive-definite (as it is in the example above).
- integral(*args, **kwds)#
Returns a symbolic integral of the vector, component-wise.
integrate()
is an alias of the function.EXAMPLES:
sage: t=var('t') sage: r=vector([t,t^2,sin(t)]) sage: r.integral(t) (1/2*t^2, 1/3*t^3, -cos(t)) sage: integrate(r,t) (1/2*t^2, 1/3*t^3, -cos(t)) sage: r.integrate(t,0,1) (1/2, 1/3, -cos(1) + 1)
- integrate(*args, **kwds)#
Returns a symbolic integral of the vector, component-wise.
integrate()
is an alias of the function.EXAMPLES:
sage: t=var('t') sage: r=vector([t,t^2,sin(t)]) sage: r.integral(t) (1/2*t^2, 1/3*t^3, -cos(t)) sage: integrate(r,t) (1/2*t^2, 1/3*t^3, -cos(t)) sage: r.integrate(t,0,1) (1/2, 1/3, -cos(1) + 1)
- is_dense()#
Return True if this is a dense vector, which is just a statement about the data structure, not the number of nonzero entries.
EXAMPLES:
sage: vector([1/2,2/5,0]).is_dense() True sage: vector([1/2,2/5,0],sparse=True).is_dense() False
- is_sparse()#
Return True if this is a sparse vector, which is just a statement about the data structure, not the number of nonzero entries.
EXAMPLES:
sage: vector([1/2,2/5,0]).is_sparse() False sage: vector([1/2,2/5,0],sparse=True).is_sparse() True
- is_vector()#
Return True, since this is a vector.
EXAMPLES:
sage: vector([1/2,2/5,0]).is_vector() True
- items()#
Return an iterator over
self
.EXAMPLES:
sage: v = vector([1,2/3,pi]) sage: v.items() <generator object at ...> sage: list(v.items()) [(0, 1), (1, 2/3), (2, pi)]
- iteritems()#
Return an iterator over
self
.EXAMPLES:
sage: v = vector([1,2/3,pi]) sage: v.items() <generator object at ...> sage: list(v.items()) [(0, 1), (1, 2/3), (2, pi)]
- lift()#
Lift
self
to the cover ring.OUTPUT:
Return a lift of self to the covering ring of the base ring \(R\), which is by definition the ring returned by calling
cover_ring()
on \(R\), or just \(R\) itself if thecover_ring()
method is not defined.EXAMPLES:
sage: V = vector(Integers(7), [5, 9, 13, 15]) ; V (5, 2, 6, 1) sage: V.lift() (5, 2, 6, 1) sage: parent(V.lift()) Ambient free module of rank 4 over the principal ideal domain Integer Ring
If the base ring does not have a cover method, return a copy of the vector:
sage: W = vector(QQ, [1, 2, 3]) sage: W1 = W.lift() sage: W is W1 False sage: parent(W1) Vector space of dimension 3 over Rational Field
- lift_centered()#
Lift to a congruent, centered vector.
INPUT:
self
A vector with coefficients in \(Integers(n)\).
OUTPUT:
The unique integer vector \(v\) such that foreach \(i\), \(Mod(v[i],n) = Mod(self[i],n)\) and \(-n/2 < v[i] \leq n/2\).
EXAMPLES:
sage: V = vector(Integers(7), [5, 9, 13, 15]) ; V (5, 2, 6, 1) sage: V.lift_centered() (-2, 2, -1, 1) sage: parent(V.lift_centered()) Ambient free module of rank 4 over the principal ideal domain Integer Ring
- list(copy=True)#
Return list of elements of self.
INPUT:
copy – bool, whether returned list is a copy that is safe to change, is ignored.
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: v = vector([x,y,z], sparse=True) sage: type(v) <class 'sage.modules.free_module_element.FreeModuleElement_generic_sparse'> sage: a = v.list(); a [x, y, z] sage: a[0] = x*y; v (x, y, z)
The optional argument
copy
is ignored:sage: a = v.list(copy=False); a [x, y, z] sage: a[0] = x*y; v (x, y, z)
- list_from_positions(positions)#
Return list of elements chosen from this vector using the given positions of this vector.
INPUT:
positions – iterable of ints
EXAMPLES:
sage: v = vector([1,2/3,pi]) sage: v.list_from_positions([0,0,0,2,1]) [1, 1, 1, pi, 2/3]
- monic()#
Return this vector divided through by the first nonzero entry of this vector.
EXAMPLES:
sage: v = vector(QQ, [0, 4/3, 5, 1, 2]) sage: v.monic() (0, 1, 15/4, 3/4, 3/2) sage: v = vector(QQ, []) sage: v.monic() ()
- monomial_coefficients(copy=True)#
Return dictionary of nonzero entries of
self
.More precisely, this returns a dictionary whose keys are indices of basis elements in the support of
self
and whose values are the corresponding coefficients.INPUT:
copy
– (default:True
) ifself
is internally represented by a dictionaryd
, then make a copy ofd
; ifFalse
, then this can cause undesired behavior by mutatingd
OUTPUT:
Python dictionary
EXAMPLES:
sage: v = vector([0,0,0,0,1/2,0,3/14]) sage: v.dict() {4: 1/2, 6: 3/14} sage: sorted(v.support()) [4, 6]
In some cases, when
copy=False
, we get back a dangerous reference:sage: v = vector({0:5, 2:3/7}, sparse=True) sage: v.dict(copy=False) {0: 5, 2: 3/7} sage: v.dict(copy=False)[0] = 18 sage: v (18, 0, 3/7)
- nintegral(*args, **kwds)#
Returns a numeric integral of the vector, component-wise, and the result of the nintegral command on each component of the input.
nintegrate()
is an alias of the function.EXAMPLES:
sage: t=var('t') sage: r=vector([t,t^2,sin(t)]) sage: vec,answers=r.nintegral(t,0,1) sage: vec (0.5, 0.3333333333333334, 0.4596976941318602) sage: type(vec) <class 'sage.modules.vector_real_double_dense.Vector_real_double_dense'> sage: answers [(0.5, 5.55111512312578...e-15, 21, 0), (0.3333333333333..., 3.70074341541719...e-15, 21, 0), (0.45969769413186..., 5.10366964392284...e-15, 21, 0)] sage: r=vector([t,0,1], sparse=True) sage: r.nintegral(t,0,1) ((0.5, 0.0, 1.0), {0: (0.5, 5.55111512312578...e-15, 21, 0), 2: (1.0, 1.11022302462515...e-14, 21, 0)})
- nintegrate(*args, **kwds)#
Returns a numeric integral of the vector, component-wise, and the result of the nintegral command on each component of the input.
nintegrate()
is an alias of the function.EXAMPLES:
sage: t=var('t') sage: r=vector([t,t^2,sin(t)]) sage: vec,answers=r.nintegral(t,0,1) sage: vec (0.5, 0.3333333333333334, 0.4596976941318602) sage: type(vec) <class 'sage.modules.vector_real_double_dense.Vector_real_double_dense'> sage: answers [(0.5, 5.55111512312578...e-15, 21, 0), (0.3333333333333..., 3.70074341541719...e-15, 21, 0), (0.45969769413186..., 5.10366964392284...e-15, 21, 0)] sage: r=vector([t,0,1], sparse=True) sage: r.nintegral(t,0,1) ((0.5, 0.0, 1.0), {0: (0.5, 5.55111512312578...e-15, 21, 0), 2: (1.0, 1.11022302462515...e-14, 21, 0)})
- nonzero_positions()#
Return the sorted list of integers
i
such thatself[i] != 0
.EXAMPLES:
sage: vector([-1,0,3,0,0,0,0.01]).nonzero_positions() [0, 2, 6]
- norm(p='__two__')#
Return the \(p\)-norm of
self
.INPUT:
p
- default: 2 -p
can be a real number greater than 1, infinity (oo
orInfinity
), or a symbolic expression.\(p=1\): the taxicab (Manhattan) norm
\(p=2\): the usual Euclidean norm (the default)
\(p=\infty\): the maximum entry (in absolute value)
Note
See also
sage.misc.functional.norm()
EXAMPLES:
sage: v = vector([1,2,-3]) sage: v.norm(5) 276^(1/5)
The default is the usual Euclidean norm.
sage: v.norm() sqrt(14) sage: v.norm(2) sqrt(14)
The infinity norm is the maximum size (in absolute value) of the entries.
sage: v.norm(Infinity) 3 sage: v.norm(oo) 3
Real or symbolic values may be used for
p
.sage: v=vector(RDF,[1,2,3]) sage: v.norm(5) 3.077384885394063 sage: v.norm(pi/2) #abs tol 1e-15 4.216595864704748 sage: _=var('a b c d p'); v=vector([a, b, c, d]) sage: v.norm(p) (abs(a)^p + abs(b)^p + abs(c)^p + abs(d)^p)^(1/p)
Notice that the result may be a symbolic expression, owing to the necessity of taking a square root (in the default case). These results can be converted to numerical values if needed.
sage: v = vector(ZZ, [3,4]) sage: nrm = v.norm(); nrm 5 sage: nrm.parent() Rational Field sage: v = vector(QQ, [3, 5]) sage: nrm = v.norm(); nrm sqrt(34) sage: nrm.parent() Symbolic Ring sage: numeric = N(nrm); numeric 5.83095189484... sage: numeric.parent() Real Field with 53 bits of precision
- normalized(p='__two__')#
Return the input vector divided by the p-norm.
INPUT:
“p” - default: 2 - p value for the norm
EXAMPLES:
sage: v = vector(QQ, [4, 1, 3, 2]) sage: v.normalized() (2/15*sqrt(30), 1/30*sqrt(30), 1/10*sqrt(30), 1/15*sqrt(30)) sage: sum(v.normalized(1)) 1
Note that normalizing the vector may change the base ring:
sage: v.base_ring() == v.normalized().base_ring() False sage: u = vector(RDF, [-3, 4, 6, 9]) sage: u.base_ring() == u.normalized().base_ring() True
- numerical_approx(prec=None, digits=None, algorithm=None)#
Return a numerical approximation of
self
withprec
bits (or decimaldigits
) of precision, by approximating all entries.INPUT:
prec
– precision in bitsdigits
– precision in decimal digits (only used ifprec
is not given)algorithm
– which algorithm to use to compute the approximation of the entries (the accepted algorithms depend on the object)
If neither
prec
nordigits
is given, the default precision is 53 bits (roughly 16 digits).EXAMPLES:
sage: v = vector(RealField(212), [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 53 bits of precision sage: numerical_approx(v) (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 53 bits of precision sage: v.n(prec=75) (1.000000000000000000000, 2.000000000000000000000, 3.000000000000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 75 bits of precision sage: numerical_approx(v, digits=3) (1.00, 2.00, 3.00) sage: _.parent() Vector space of dimension 3 over Real Field with 14 bits of precision
Both functional and object-oriented usage is possible.
sage: u = vector(QQ, [1/2, 1/3, 1/4]) sage: u.n() (0.500000000000000, 0.333333333333333, 0.250000000000000) sage: u.numerical_approx() (0.500000000000000, 0.333333333333333, 0.250000000000000) sage: n(u) (0.500000000000000, 0.333333333333333, 0.250000000000000) sage: N(u) (0.500000000000000, 0.333333333333333, 0.250000000000000) sage: numerical_approx(u) (0.500000000000000, 0.333333333333333, 0.250000000000000)
Precision (bits) and digits (decimal) may be specified. When both are given,
prec
wins.sage: u = vector(QQ, [1/2, 1/3, 1/4]) sage: n(u, prec=15) (0.5000, 0.3333, 0.2500) sage: n(u, digits=5) (0.50000, 0.33333, 0.25000) sage: n(u, prec=30, digits=100) (0.50000000, 0.33333333, 0.25000000)
These are some legacy doctests that were part of various specialized versions of the numerical approximation routine that were removed as part of trac ticket #12195.
sage: v = vector(ZZ, [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 53 bits of precision sage: v.n(prec=75) (1.000000000000000000000, 2.000000000000000000000, 3.000000000000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 75 bits of precision sage: v = vector(RDF, [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 53 bits of precision sage: v = vector(CDF, [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Complex Field with 53 bits of precision sage: v = vector(Integers(8), [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 53 bits of precision sage: v.n(prec=75) (1.000000000000000000000, 2.000000000000000000000, 3.000000000000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 75 bits of precision sage: v = vector(QQ, [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 53 bits of precision sage: v.n(prec=75) (1.000000000000000000000, 2.000000000000000000000, 3.000000000000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 75 bits of precision
sage: v = vector(GF(2), [1,2,3]) sage: v.n() (1.00000000000000, 0.000000000000000, 1.00000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 53 bits of precision sage: v.n(prec=75) (1.000000000000000000000, 0.0000000000000000000000, 1.000000000000000000000) sage: _.parent() Vector space of dimension 3 over Real Field with 75 bits of precision
- numpy(dtype='object')#
Convert self to a numpy array.
INPUT:
dtype
– the numpy dtypeof the returned array
EXAMPLES:
sage: v = vector([1,2,3]) sage: v.numpy() array([1, 2, 3], dtype=object) sage: v.numpy() * v.numpy() array([1, 4, 9], dtype=object) sage: vector(QQ, [1, 2, 5/6]).numpy() array([1, 2, 5/6], dtype=object)
By default the
object
dtype is used. Alternatively, the desired dtype can be passed in as a parameter:sage: v = vector(QQ, [1, 2, 5/6]) sage: v.numpy() array([1, 2, 5/6], dtype=object) sage: v.numpy(dtype=float) array([1. , 2. , 0.83333333]) sage: v.numpy(dtype=int) array([1, 2, 0]) sage: import numpy sage: v.numpy(dtype=numpy.uint8) array([1, 2, 0], dtype=uint8)
Passing a dtype of None will let numpy choose a native type, which can be more efficient but may have unintended consequences:
sage: v.numpy(dtype=None) array([1. , 2. , 0.83333333]) sage: w = vector(ZZ, [0, 1, 2^63 -1]); w (0, 1, 9223372036854775807) sage: wn = w.numpy(dtype=None); wn array([ 0, 1, 9223372036854775807]...) sage: wn.dtype dtype('int64') sage: w.dot_product(w) 85070591730234615847396907784232501250 sage: wn.dot(wn) # overflow 2
Numpy can give rather obscure errors; we wrap these to give a bit of context:
sage: vector([1, 1/2, QQ['x'].0]).numpy(dtype=float) Traceback (most recent call last): ... ValueError: Could not convert vector over Univariate Polynomial Ring in x over Rational Field to numpy array of type <... 'float'>: setting an array element with a sequence.
- outer_product(right)#
Returns a matrix, the outer product of two vectors
self
andright
.INPUT:
right
- a vector (or free module element) of any size, whose elements are compatible (with regard to multiplication) with the elements ofself
.
OUTPUT:
The outer product of two vectors \(x\) and \(y\) (respectively
self
andright
) can be described several ways. If we interpret \(x\) as a \(m\times 1\) matrix and interpret \(y\) as a \(1\times n\) matrix, then the outer product is the \(m\times n\) matrix from the usual matrix product \(xy\). Notice how this is the “opposite” in some ways from an inner product (which would require \(m=n\)).If we just consider vectors, use each entry of \(x\) to create a scalar multiples of the vector \(y\) and use these vectors as the rows of a matrix. Or use each entry of \(y\) to create a scalar multiples of \(x\) and use these vectors as the columns of a matrix.
EXAMPLES:
sage: u = vector(QQ, [1/2, 1/3, 1/4, 1/5]) sage: v = vector(ZZ, [60, 180, 600]) sage: u.outer_product(v) [ 30 90 300] [ 20 60 200] [ 15 45 150] [ 12 36 120] sage: M = v.outer_product(u); M [ 30 20 15 12] [ 90 60 45 36] [300 200 150 120] sage: M.parent() Full MatrixSpace of 3 by 4 dense matrices over Rational Field
The more general
sage.matrix.matrix2.tensor_product()
is an operation on a pair of matrices. If we construct a pair of vectors as a column vector and a row vector, then an outer product and a tensor product are identical. Thus \(tensor_product\) is a synonym for this method.sage: u = vector(QQ, [1/2, 1/3, 1/4, 1/5]) sage: v = vector(ZZ, [60, 180, 600]) sage: u.tensor_product(v) == (u.column()).tensor_product(v.row()) True
The result is always a dense matrix, no matter if the two vectors are, or are not, dense.
sage: d = vector(ZZ,[4,5], sparse=False) sage: s = vector(ZZ, [1,2,3], sparse=True) sage: dd = d.outer_product(d) sage: ds = d.outer_product(s) sage: sd = s.outer_product(d) sage: ss = s.outer_product(s) sage: all([dd.is_dense(), ds.is_dense(), sd.is_dense(), dd.is_dense()]) True
Vectors with no entries do the right thing.
sage: v = vector(ZZ, []) sage: z = v.outer_product(v) sage: z.parent() Full MatrixSpace of 0 by 0 dense matrices over Integer Ring
There is a fair amount of latitude in the value of the
right
vector, and the matrix that results can have entries from a new ring large enough to contain the result. If you know better, you can sometimes bring the result down to a less general ring.sage: R.<t> = ZZ[] sage: v = vector(R, [12, 24*t]) sage: w = vector(QQ, [1/2, 1/3, 1/4]) sage: op = v.outer_product(w) sage: op [ 6 4 3] [12*t 8*t 6*t] sage: op.base_ring() Univariate Polynomial Ring in t over Rational Field sage: m = op.change_ring(R); m [ 6 4 3] [12*t 8*t 6*t] sage: m.base_ring() Univariate Polynomial Ring in t over Integer Ring
But some inputs are not compatible, even if vectors.
sage: w = vector(GF(5), [1,2]) sage: v = vector(GF(7), [1,2,3,4]) sage: z = w.outer_product(v) Traceback (most recent call last): ... TypeError: unsupported operand parent(s) for *: 'Full MatrixSpace of 2 by 1 dense matrices over Finite Field of size 5' and 'Full MatrixSpace of 1 by 4 dense matrices over Finite Field of size 7'
And some inputs don’t make any sense at all.
sage: w=vector(QQ, [5,10]) sage: z=w.outer_product(6) Traceback (most recent call last): ... TypeError: right operand in an outer product must be a vector, not an element of Integer Ring
- pairwise_product(right)#
Return the pairwise product of self and right, which is a vector of the products of the corresponding entries.
INPUT:
right
- vector of the same degree as self. It need not be in the same vector space as self, as long as the coefficients can be multiplied.
EXAMPLES:
sage: V = FreeModule(ZZ, 3) sage: v = V([1,2,3]) sage: w = V([4,5,6]) sage: v.pairwise_product(w) (4, 10, 18) sage: sum(v.pairwise_product(w)) == v.dot_product(w) True
sage: W = VectorSpace(GF(3),3) sage: w = W([0,1,2]) sage: w.pairwise_product(v) (0, 2, 0) sage: w.pairwise_product(v).parent() Vector space of dimension 3 over Finite Field of size 3
Implicit coercion is well defined (regardless of order), so we get 2 even if we do the dot product in the other order.
sage: v.pairwise_product(w).parent() Vector space of dimension 3 over Finite Field of size 3
- plot(plot_type=None, start=None, **kwds)#
INPUT:
plot_type
- (default: ‘arrow’ if v has 3 or fewer components,otherwise ‘step’) type of plot. Options are:
‘arrow’ to draw an arrow
‘point’ to draw a point at the coordinates specified by the vector
‘step’ to draw a step function representing the coordinates of the vector.
Both ‘arrow’ and ‘point’ raise exceptions if the vector has more than 3 dimensions.
start
- (default: origin in correct dimension) may be a tuple, list, or vector.
EXAMPLES:
The following both plot the given vector:
sage: v = vector(RDF, (1,2)) sage: A = plot(v) sage: B = v.plot() sage: A+B # should just show one vector Graphics object consisting of 2 graphics primitives
Examples of the plot types:
sage: A = plot(v, plot_type='arrow') sage: B = plot(v, plot_type='point', color='green', size=20) sage: C = plot(v, plot_type='step') # calls v.plot_step() sage: A+B+C Graphics object consisting of 3 graphics primitives
You can use the optional arguments for
plot_step()
:sage: eps = 0.1 sage: plot(v, plot_type='step', eps=eps, xmax=5, hue=0) Graphics object consisting of 1 graphics primitive
Three-dimensional examples:
sage: v = vector(RDF, (1,2,1)) sage: plot(v) # defaults to an arrow plot Graphics3d Object
sage: plot(v, plot_type='arrow') Graphics3d Object
sage: from sage.plot.plot3d.shapes2 import frame3d sage: plot(v, plot_type='point')+frame3d((0,0,0), v.list()) Graphics3d Object
sage: plot(v, plot_type='step') # calls v.plot_step() Graphics object consisting of 1 graphics primitive
sage: plot(v, plot_type='step', eps=eps, xmax=5, hue=0) Graphics object consisting of 1 graphics primitive
With greater than three coordinates, it defaults to a step plot:
sage: v = vector(RDF, (1,2,3,4)) sage: plot(v) Graphics object consisting of 1 graphics primitive
One dimensional vectors are plotted along the horizontal axis of the coordinate plane:
sage: plot(vector([1])) Graphics object consisting of 1 graphics primitive
An optional start argument may also be specified by a tuple, list, or vector:
sage: u = vector([1,2]); v = vector([2,5]) sage: plot(u, start=v) Graphics object consisting of 1 graphics primitive
- plot_step(xmin=0, xmax=1, eps=None, res=None, connect=True, **kwds)#
INPUT:
xmin
- (default: 0) start x position to start plottingxmax
- (default: 1) stop x position to stop plottingeps
- (default: determined by xmax) we view this vector as defining a function at the points xmin, xmin + eps, xmin + 2*eps, …,res
- (default: all points) total number of points to include in the graphconnect
- (default: True) if True draws a line; otherwise draw a list of points.
EXAMPLES:
sage: eps=0.1 sage: v = vector(RDF, [sin(n*eps) for n in range(100)]) sage: v.plot_step(eps=eps, xmax=5, hue=0) Graphics object consisting of 1 graphics primitive
- row()#
Return a matrix with a single row and the same entries as the vector
self
.OUTPUT:
A matrix over the same ring as the vector (or free module element), with a single row. The entries of the row are identical to those of the vector, and in the same order.
EXAMPLES:
sage: v = vector(ZZ, [1,2,3]) sage: w = v.row(); w [1 2 3] sage: w.parent() Full MatrixSpace of 1 by 3 dense matrices over Integer Ring sage: x = vector(FiniteField(13), [2,4,8,16]) sage: x.row() [2 4 8 3]
There is more than one way to get one-row matrix from a vector, but the
row
method is more efficient than making a column and then taking a transpose. Notice that supplying a vector to the matrix constructor demonstrates Sage’s preference for rows.sage: x = vector(RDF, [sin(i*pi/20) for i in range(10)]) sage: x.row() == matrix(x) True sage: x.row() == x.column().transpose() True
Sparse or dense implementations are preserved.
sage: d = vector(RR, [1.0, 2.0, 3.0]) sage: s = vector(CDF, {2:5.0+6.0*I}) sage: dm = d.row() sage: sm = s.row() sage: all([d.is_dense(), dm.is_dense(), s.is_sparse(), sm.is_sparse()]) True
- set(i, value)#
Like
__setitem__
but without type or bounds checking: \(i\) must satisfy0 <= i < self.degree
andvalue
must be an element of the coordinate ring.EXAMPLES:
sage: v = vector(SR, [1/2,2/5,0]); v (1/2, 2/5, 0) sage: v.set(2, pi); v (1/2, 2/5, pi)
- sparse_vector()#
Return sparse version of self. If self is sparse, just return self; otherwise, create and return correspond sparse vector.
EXAMPLES:
sage: vector([-1,0,3,0,0,0]).sparse_vector().is_sparse() True sage: vector([-1,0,3,0,0,0]).sparse_vector().is_sparse() True sage: vector([-1,0,3,0,0,0]).sparse_vector() (-1, 0, 3, 0, 0, 0)
- subs(in_dict=None, **kwds)#
EXAMPLES:
sage: var('a,b,d,e') (a, b, d, e) sage: v = vector([a, b, d, e]) sage: v.substitute(a=1) (1, b, d, e) sage: v.subs(a=b, b=d) (b, d, d, e)
- support()#
Return the integers
i
such thatself[i] != 0
. This is the same as thenonzero_positions
function.EXAMPLES:
sage: vector([-1,0,3,0,0,0,0.01]).support() [0, 2, 6]
- tensor_product(right)#
Returns a matrix, the outer product of two vectors
self
andright
.INPUT:
right
- a vector (or free module element) of any size, whose elements are compatible (with regard to multiplication) with the elements ofself
.
OUTPUT:
The outer product of two vectors \(x\) and \(y\) (respectively
self
andright
) can be described several ways. If we interpret \(x\) as a \(m\times 1\) matrix and interpret \(y\) as a \(1\times n\) matrix, then the outer product is the \(m\times n\) matrix from the usual matrix product \(xy\). Notice how this is the “opposite” in some ways from an inner product (which would require \(m=n\)).If we just consider vectors, use each entry of \(x\) to create a scalar multiples of the vector \(y\) and use these vectors as the rows of a matrix. Or use each entry of \(y\) to create a scalar multiples of \(x\) and use these vectors as the columns of a matrix.
EXAMPLES:
sage: u = vector(QQ, [1/2, 1/3, 1/4, 1/5]) sage: v = vector(ZZ, [60, 180, 600]) sage: u.outer_product(v) [ 30 90 300] [ 20 60 200] [ 15 45 150] [ 12 36 120] sage: M = v.outer_product(u); M [ 30 20 15 12] [ 90 60 45 36] [300 200 150 120] sage: M.parent() Full MatrixSpace of 3 by 4 dense matrices over Rational Field
The more general
sage.matrix.matrix2.tensor_product()
is an operation on a pair of matrices. If we construct a pair of vectors as a column vector and a row vector, then an outer product and a tensor product are identical. Thus \(tensor_product\) is a synonym for this method.sage: u = vector(QQ, [1/2, 1/3, 1/4, 1/5]) sage: v = vector(ZZ, [60, 180, 600]) sage: u.tensor_product(v) == (u.column()).tensor_product(v.row()) True
The result is always a dense matrix, no matter if the two vectors are, or are not, dense.
sage: d = vector(ZZ,[4,5], sparse=False) sage: s = vector(ZZ, [1,2,3], sparse=True) sage: dd = d.outer_product(d) sage: ds = d.outer_product(s) sage: sd = s.outer_product(d) sage: ss = s.outer_product(s) sage: all([dd.is_dense(), ds.is_dense(), sd.is_dense(), dd.is_dense()]) True
Vectors with no entries do the right thing.
sage: v = vector(ZZ, []) sage: z = v.outer_product(v) sage: z.parent() Full MatrixSpace of 0 by 0 dense matrices over Integer Ring
There is a fair amount of latitude in the value of the
right
vector, and the matrix that results can have entries from a new ring large enough to contain the result. If you know better, you can sometimes bring the result down to a less general ring.sage: R.<t> = ZZ[] sage: v = vector(R, [12, 24*t]) sage: w = vector(QQ, [1/2, 1/3, 1/4]) sage: op = v.outer_product(w) sage: op [ 6 4 3] [12*t 8*t 6*t] sage: op.base_ring() Univariate Polynomial Ring in t over Rational Field sage: m = op.change_ring(R); m [ 6 4 3] [12*t 8*t 6*t] sage: m.base_ring() Univariate Polynomial Ring in t over Integer Ring
But some inputs are not compatible, even if vectors.
sage: w = vector(GF(5), [1,2]) sage: v = vector(GF(7), [1,2,3,4]) sage: z = w.outer_product(v) Traceback (most recent call last): ... TypeError: unsupported operand parent(s) for *: 'Full MatrixSpace of 2 by 1 dense matrices over Finite Field of size 5' and 'Full MatrixSpace of 1 by 4 dense matrices over Finite Field of size 7'
And some inputs don’t make any sense at all.
sage: w=vector(QQ, [5,10]) sage: z=w.outer_product(6) Traceback (most recent call last): ... TypeError: right operand in an outer product must be a vector, not an element of Integer Ring
- class sage.modules.free_module_element.FreeModuleElement_generic_dense#
Bases:
sage.modules.free_module_element.FreeModuleElement
A generic dense element of a free module.
- function(*args)#
Return a vector over a callable symbolic expression ring.
EXAMPLES:
sage: x,y=var('x,y') sage: v=vector([x,y,x*sin(y)]) sage: w=v.function([x,y]); w (x, y) |--> (x, y, x*sin(y)) sage: w.coordinate_ring() Callable function ring with arguments (x, y) sage: w(1,2) (1, 2, sin(2)) sage: w(2,1) (2, 1, 2*sin(1)) sage: w(y=1,x=2) (2, 1, 2*sin(1))
sage: x,y=var('x,y') sage: v=vector([x,y,x*sin(y)]) sage: w=v.function([x]); w x |--> (x, y, x*sin(y)) sage: w.coordinate_ring() Callable function ring with argument x sage: w(4) (4, y, 4*sin(y))
- list(copy=True)#
Return list of elements of self.
INPUT:
copy – bool, return list of underlying entries
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: v = vector([x,y,z]) sage: type(v) <class 'sage.modules.free_module_element.FreeModuleElement_generic_dense'> sage: a = v.list(); a [x, y, z] sage: a[0] = x*y; v (x, y, z) sage: a = v.list(copy=False); a [x, y, z] sage: a[0] = x*y; v (x*y, y, z)
- class sage.modules.free_module_element.FreeModuleElement_generic_sparse#
Bases:
sage.modules.free_module_element.FreeModuleElement
A generic sparse free module element is a dictionary with keys ints i and entries in the base ring.
- denominator()#
Return the least common multiple of the denominators of the entries of self.
EXAMPLES:
sage: v = vector([1/2,2/5,3/14], sparse=True) sage: v.denominator() 70
- dict(copy=True)#
Return dictionary of nonzero entries of
self
.More precisely, this returns a dictionary whose keys are indices of basis elements in the support of
self
and whose values are the corresponding coefficients.INPUT:
copy
– (default:True
) ifself
is internally represented by a dictionaryd
, then make a copy ofd
; ifFalse
, then this can cause undesired behavior by mutatingd
OUTPUT:
Python dictionary
EXAMPLES:
sage: v = vector([0,0,0,0,1/2,0,3/14], sparse=True) sage: v.dict() {4: 1/2, 6: 3/14} sage: sorted(v.support()) [4, 6]
- hamming_weight()#
Returns the number of positions
i
such thatself[i] != 0
.EXAMPLES:
sage: v = vector({1: 1, 3: -2}) sage: w = vector({1: 4, 3: 2}) sage: v+w (0, 5, 0, 0) sage: (v+w).hamming_weight() 1
- items()#
Return an iterator over the entries of
self
.EXAMPLES:
sage: v = vector([1,2/3,pi], sparse=True) sage: next(v.items()) (0, 1) sage: list(v.items()) [(0, 1), (1, 2/3), (2, pi)]
- iteritems()#
Return an iterator over the entries of
self
.EXAMPLES:
sage: v = vector([1,2/3,pi], sparse=True) sage: next(v.items()) (0, 1) sage: list(v.items()) [(0, 1), (1, 2/3), (2, pi)]
- list(copy=True)#
Return list of elements of
self
.INPUT:
copy
– ignored for sparse vectors
EXAMPLES:
sage: R.<x> = QQ[] sage: M = FreeModule(R, 3, sparse=True) * (1/x) sage: v = M([-x^2, 3/x, 0]) sage: type(v) <class 'sage.modules.free_module_element.FreeModuleElement_generic_sparse'> sage: a = v.list() sage: a [-x^2, 3/x, 0] sage: [parent(c) for c in a] [Fraction Field of Univariate Polynomial Ring in x over Rational Field, Fraction Field of Univariate Polynomial Ring in x over Rational Field, Fraction Field of Univariate Polynomial Ring in x over Rational Field]
- monomial_coefficients(copy=True)#
Return dictionary of nonzero entries of
self
.More precisely, this returns a dictionary whose keys are indices of basis elements in the support of
self
and whose values are the corresponding coefficients.INPUT:
copy
– (default:True
) ifself
is internally represented by a dictionaryd
, then make a copy ofd
; ifFalse
, then this can cause undesired behavior by mutatingd
OUTPUT:
Python dictionary
EXAMPLES:
sage: v = vector([0,0,0,0,1/2,0,3/14], sparse=True) sage: v.dict() {4: 1/2, 6: 3/14} sage: sorted(v.support()) [4, 6]
- nonzero_positions()#
Returns the list of numbers
i
such thatself[i] != 0
.EXAMPLES:
sage: v = vector({1: 1, 3: -2}) sage: w = vector({1: 4, 3: 2}) sage: v+w (0, 5, 0, 0) sage: (v+w).nonzero_positions() [1]
- numerical_approx(prec=None, digits=None, algorithm=None)#
Return a numerical approximation of
self
withprec
bits (or decimaldigits
) of precision, by approximating all entries.INPUT:
prec
– precision in bitsdigits
– precision in decimal digits (only used ifprec
is not given)algorithm
– which algorithm to use to compute the approximation of the entries (the accepted algorithms depend on the object)
If neither
prec
nordigits
is given, the default precision is 53 bits (roughly 16 digits).EXAMPLES:
sage: v = vector(RealField(200), [1,2,3], sparse=True) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Sparse vector space of dimension 3 over Real Field with 53 bits of precision sage: v.n(prec=75) (1.000000000000000000000, 2.000000000000000000000, 3.000000000000000000000) sage: _.parent() Sparse vector space of dimension 3 over Real Field with 75 bits of precision
- sage.modules.free_module_element.free_module_element(arg0, arg1=None, arg2=None, sparse=None, immutable=False)#
Return a vector or free module element with specified entries.
CALL FORMATS:
This constructor can be called in several different ways. In each case,
sparse=True
orsparse=False
as well asimmutable=True
orimmutable=False
can be supplied as an option.free_module_element()
is an alias forvector()
.vector(object)
vector(ring, object)
vector(object, ring)
vector(ring, degree, object)
vector(ring, degree)
INPUT:
object
– a list, dictionary, or other iterable containing the entries of the vector, including any object that is palatable to theSequence
constructorring
– a base ring (or field) for the vector space or free module, which contains all of the elementsdegree
– an integer specifying the number of entries in the vector or free module elementsparse
– boolean, whether the result should be a sparse vectorimmutable
– boolean (default:False
); whether the result should be an immutable vector
In call format 4, an error is raised if the
degree
does not match the length ofobject
so this call can provide some safeguards. Note however that using this format whenobject
is a dictionary is unlikely to work properly.OUTPUT:
An element of the ambient vector space or free module with the given base ring and implied or specified dimension or rank, containing the specified entries and with correct degree.
In call format 5, no entries are specified, so the element is populated with all zeros.
If the
sparse
option is not supplied, the output will generally have a dense representation. The exception is ifobject
is a dictionary, then the representation will be sparse.EXAMPLES:
sage: v = vector([1,2,3]); v (1, 2, 3) sage: v.parent() Ambient free module of rank 3 over the principal ideal domain Integer Ring sage: v = vector([1,2,3/5]); v (1, 2, 3/5) sage: v.parent() Vector space of dimension 3 over Rational Field
All entries must canonically coerce to some common ring:
sage: v = vector([17, GF(11)(5), 19/3]); v Traceback (most recent call last): ... TypeError: unable to find a common ring for all elements
sage: v = vector([17, GF(11)(5), 19]); v (6, 5, 8) sage: v.parent() Vector space of dimension 3 over Finite Field of size 11 sage: v = vector([17, GF(11)(5), 19], QQ); v (17, 5, 19) sage: v.parent() Vector space of dimension 3 over Rational Field sage: v = vector((1,2,3), QQ); v (1, 2, 3) sage: v.parent() Vector space of dimension 3 over Rational Field sage: v = vector(QQ, (1,2,3)); v (1, 2, 3) sage: v.parent() Vector space of dimension 3 over Rational Field sage: v = vector(vector([1,2,3])); v (1, 2, 3) sage: v.parent() Ambient free module of rank 3 over the principal ideal domain Integer Ring
You can also use
free_module_element
, which is the same asvector
.sage: free_module_element([1/3, -4/5]) (1/3, -4/5)
We make a vector mod 3 out of a vector over \(\ZZ\).
sage: vector(vector([1,2,3]), GF(3)) (1, 2, 0)
The degree of a vector may be specified:
sage: vector(QQ, 4, [1,1/2,1/3,1/4]) (1, 1/2, 1/3, 1/4)
But it is an error if the degree and size of the list of entries are mismatched:
sage: vector(QQ, 5, [1,1/2,1/3,1/4]) Traceback (most recent call last): ... ValueError: incompatible degrees in vector constructor
Providing no entries populates the vector with zeros, but of course, you must specify the degree since it is not implied. Here we use a finite field as the base ring.
sage: w = vector(FiniteField(7), 4); w (0, 0, 0, 0) sage: w.parent() Vector space of dimension 4 over Finite Field of size 7
The fastest method to construct a zero vector is to call the
zero_vector()
method directly on a free module or vector space, since vector(…) must do a small amount of type checking. Almost as fast as thezero_vector()
method is thezero_vector()
constructor, which defaults to the integers.sage: vector(ZZ, 5) # works fine (0, 0, 0, 0, 0) sage: (ZZ^5).zero_vector() # very tiny bit faster (0, 0, 0, 0, 0) sage: zero_vector(ZZ, 5) # similar speed to vector(...) (0, 0, 0, 0, 0) sage: z = zero_vector(5); z (0, 0, 0, 0, 0) sage: z.parent() Ambient free module of rank 5 over the principal ideal domain Integer Ring
Here we illustrate the creation of sparse vectors by using a dictionary:
sage: vector({1:1.1, 3:3.14}) (0.000000000000000, 1.10000000000000, 0.000000000000000, 3.14000000000000)
With no degree given, a dictionary of entries implicitly declares a degree by the largest index (key) present. So you can provide a terminal element (perhaps a zero?) to set the degree. But it is probably safer to just include a degree in your construction.
sage: v = vector(QQ, {0:1/2, 4:-6, 7:0}); v (1/2, 0, 0, 0, -6, 0, 0, 0) sage: v.degree() 8 sage: v.is_sparse() True sage: w = vector(QQ, 8, {0:1/2, 4:-6}) sage: w == v True
It is an error to specify a negative degree.
sage: vector(RR, -4, [1.0, 2.0, 3.0, 4.0]) Traceback (most recent call last): ... ValueError: cannot specify the degree of a vector as a negative integer (-4)
It is an error to create a zero vector but not provide a ring as the first argument.
sage: vector('junk', 20) Traceback (most recent call last): ... TypeError: first argument must be base ring of zero vector, not junk
And it is an error to specify an index in a dictionary that is greater than or equal to a requested degree.
sage: vector(ZZ, 10, {3:4, 7:-2, 10:637}) Traceback (most recent call last): ... ValueError: dictionary of entries has a key (index) exceeding the requested degree
A 1-dimensional numpy array of type float or complex may be passed to vector. Unless an explicit ring is given, the result will be a vector in the appropriate dimensional vector space over the real double field or the complex double field. The data in the array must be contiguous, so column-wise slices of numpy matrices will raise an exception.
sage: import numpy sage: x = numpy.random.randn(10) sage: y = vector(x) sage: parent(y) Vector space of dimension 10 over Real Double Field sage: parent(vector(RDF, x)) Vector space of dimension 10 over Real Double Field sage: parent(vector(CDF, x)) Vector space of dimension 10 over Complex Double Field sage: parent(vector(RR, x)) Vector space of dimension 10 over Real Field with 53 bits of precision sage: v = numpy.random.randn(10) * complex(0,1) sage: w = vector(v) sage: parent(w) Vector space of dimension 10 over Complex Double Field
Multi-dimensional arrays are not supported:
sage: import numpy as np sage: a = np.array([[1, 2, 3], [4, 5, 6]], np.float64) sage: vector(a) Traceback (most recent call last): ... TypeError: cannot convert 2-dimensional array to a vector
If any of the arguments to vector have Python type int, long, real, or complex, they will first be coerced to the appropriate Sage objects. This fixes trac ticket #3847.
sage: v = vector([int(0)]); v (0) sage: v[0].parent() Integer Ring sage: v = vector(range(10)); v (0, 1, 2, 3, 4, 5, 6, 7, 8, 9) sage: v[3].parent() Integer Ring sage: v = vector([float(23.4), int(2), complex(2+7*I), 1]); v (23.4, 2.0, 2.0 + 7.0*I, 1.0) sage: v[1].parent() Complex Double Field
If the argument is a vector, it doesn’t change the base ring. This fixes trac ticket #6643:
sage: K.<sqrt3> = QuadraticField(3) sage: u = vector(K, (1/2, sqrt3/2) ) sage: vector(u).base_ring() Number Field in sqrt3 with defining polynomial x^2 - 3 with sqrt3 = 1.732050807568878? sage: v = vector(K, (0, 1) ) sage: vector(v).base_ring() Number Field in sqrt3 with defining polynomial x^2 - 3 with sqrt3 = 1.732050807568878?
Constructing a vector from a numpy array behaves as expected:
sage: import numpy sage: a=numpy.array([1,2,3]) sage: v=vector(a); v (1, 2, 3) sage: parent(v) Ambient free module of rank 3 over the principal ideal domain Integer Ring
Complex numbers can be converted naturally to a sequence of length 2. And then to a vector.
sage: c = CDF(2 + 3*I) sage: v = vector(c); v (2.0, 3.0)
A generator, or other iterable, may also be supplied as input. Anything that can be converted to a
Sequence
is a possible input.sage: type(i^2 for i in range(3)) <... 'generator'> sage: v = vector(i^2 for i in range(3)); v (0, 1, 4)
An empty list, without a ring given, will default to the integers.
sage: x = vector([]); x () sage: x.parent() Ambient free module of rank 0 over the principal ideal domain Integer Ring
The
immutable
switch allows to create an immutable vector.sage: v = vector(QQ, {0:1/2, 4:-6, 7:0}, immutable=True); v (1/2, 0, 0, 0, -6, 0, 0, 0) sage: v.is_immutable() True
The
immutable
switch works regardless of the type of valid input to the constructor.sage: v = vector(ZZ, 4, immutable=True) sage: v.is_immutable() True sage: w = vector(ZZ, [1,2,3]) sage: v = vector(w, ZZ, immutable=True) sage: v.is_immutable() True sage: v = vector(QQ, w, immutable=True) sage: v.is_immutable() True sage: import numpy as np sage: w = np.array([1, 2, pi], float) sage: v = vector(w, immutable=True) sage: v.is_immutable() True sage: w = np.array([i, 2, 3], complex) sage: v = vector(w, immutable=True) sage: v.is_immutable() True
- sage.modules.free_module_element.is_FreeModuleElement(x)#
EXAMPLES:
sage: sage.modules.free_module_element.is_FreeModuleElement(0) False sage: sage.modules.free_module_element.is_FreeModuleElement(vector([1,2,3])) True
- sage.modules.free_module_element.make_FreeModuleElement_generic_dense(parent, entries, degree)#
EXAMPLES:
sage: sage.modules.free_module_element.make_FreeModuleElement_generic_dense(QQ^3, [1,2,-3/7], 3) (1, 2, -3/7)
- sage.modules.free_module_element.make_FreeModuleElement_generic_dense_v1(parent, entries, degree, is_mutable)#
EXAMPLES:
sage: v = sage.modules.free_module_element.make_FreeModuleElement_generic_dense_v1(QQ^3, [1,2,-3/7], 3, True); v (1, 2, -3/7) sage: v[0] = 10; v (10, 2, -3/7) sage: v = sage.modules.free_module_element.make_FreeModuleElement_generic_dense_v1(QQ^3, [1,2,-3/7], 3, False); v (1, 2, -3/7) sage: v[0] = 10 Traceback (most recent call last): ... ValueError: vector is immutable; please change a copy instead (use copy())
- sage.modules.free_module_element.make_FreeModuleElement_generic_sparse(parent, entries, degree)#
EXAMPLES:
sage: v = sage.modules.free_module_element.make_FreeModuleElement_generic_sparse(QQ^3, {2:5/2}, 3); v (0, 0, 5/2)
- sage.modules.free_module_element.make_FreeModuleElement_generic_sparse_v1(parent, entries, degree, is_mutable)#
EXAMPLES:
sage: v = sage.modules.free_module_element.make_FreeModuleElement_generic_sparse_v1(QQ^3, {2:5/2}, 3, False); v (0, 0, 5/2) sage: v.is_mutable() False
- sage.modules.free_module_element.prepare(v, R, degree=None)#
Converts an object describing elements of a vector into a list of entries in a common ring.
INPUT:
v
- a dictionary with non-negative integers as keys, or a list or other object that can be converted by theSequence
constructorR
- a ring containing all the entries, possibly given asNone
degree
- a requested size for the list when the input is a dictionary, otherwise ignored
OUTPUT:
A pair.
The first item is a list of the values specified in the object
v
. If the object is a dictionary , entries are placed in the list according to the indices that were their keys in the dictionary, and the remainder of the entries are zero. The value ofdegree
is assumed to be larger than any index provided in the dictionary and will be used as the number of entries in the returned list.The second item returned is a ring that contains all of the entries in the list. If
R
is given, the entries are coerced in. Otherwise a common ring is found. For more details, see theSequence
object. Whenv
has no elements andR
isNone
, the ring returned is the integers.EXAMPLES:
sage: from sage.modules.free_module_element import prepare sage: prepare([1,2/3,5],None) ([1, 2/3, 5], Rational Field) sage: prepare([1,2/3,5],RR) ([1.00000000000000, 0.666666666666667, 5.00000000000000], Real Field with 53 bits of precision) sage: prepare({1:4, 3:-2}, ZZ, 6) ([0, 4, 0, -2, 0, 0], Integer Ring) sage: prepare({3:1, 5:3}, QQ, 6) ([0, 0, 0, 1, 0, 3], Rational Field) sage: prepare([1,2/3,'10',5],RR) ([1.00000000000000, 0.666666666666667, 10.0000000000000, 5.00000000000000], Real Field with 53 bits of precision) sage: prepare({},QQ, 0) ([], Rational Field) sage: prepare([1,2/3,'10',5],None) Traceback (most recent call last): ... TypeError: unable to find a common ring for all elements
Some objects can be converted to sequences even if they are not always thought of as sequences.
sage: c = CDF(2+3*I) sage: prepare(c, None) ([2.0, 3.0], Real Double Field)
This checks a bug listed at trac ticket #10595. Without good evidence for a ring, the default is the integers.
sage: prepare([], None) ([], Integer Ring)
- sage.modules.free_module_element.random_vector(ring, degree=None, *args, **kwds)#
Returns a vector (or module element) with random entries.
INPUT:
ring – default:
ZZ
- the base ring for the entriesdegree – a non-negative integer for the number of entries in the vector
sparse – default:
False
- whether to use a sparse implementationargs, kwds - additional arguments and keywords are passed to the
random_element()
method of the ring
OUTPUT:
A vector, or free module element, with
degree
elements fromring
, chosen randomly from the ring according to the ring’srandom_element()
method.Note
See below for examples of how random elements are generated by some common base rings.
EXAMPLES:
First, module elements over the integers. The default distribution is tightly clustered around -1, 0, 1. Uniform distributions can be specified by giving bounds, though the upper bound is never met. See
sage.rings.integer_ring.IntegerRing_class.random_element()
for several other variants.sage: random_vector(10).parent() Ambient free module of rank 10 over the principal ideal domain Integer Ring sage: random_vector(20).parent() Ambient free module of rank 20 over the principal ideal domain Integer Ring sage: v = random_vector(ZZ, 20, x=4) sage: all(i in range(4) for i in v) True sage: v = random_vector(ZZ, 20, x=-20, y=100) sage: all(i in range(-20, 100) for i in v) True
If the ring is not specified, the default is the integers, and parameters for the random distribution may be passed without using keywords. This is a random vector with 20 entries uniformly distributed between -20 and 100.
sage: random_vector(20, -20, 100).parent() Ambient free module of rank 20 over the principal ideal domain Integer Ring
Now over the rationals. Note that bounds on the numerator and denominator may be specified. See
sage.rings.rational_field.RationalField.random_element()
for documentation.sage: random_vector(QQ, 10).parent() Vector space of dimension 10 over Rational Field sage: v = random_vector(QQ, 10, num_bound=15, den_bound=5) sage: v.parent() Vector space of dimension 10 over Rational Field sage: all(q.numerator() <= 15 and q.denominator() <= 5 for q in v) True
Inexact rings may be used as well. The reals have uniform distributions, with the range \((-1,1)\) as the default. More at:
sage.rings.real_mpfr.RealField_class.random_element()
sage: v = random_vector(RR, 5) sage: v.parent() Vector space of dimension 5 over Real Field with 53 bits of precision sage: all(-1 <= r <= 1 for r in v) True sage: v = random_vector(RR, 5, min=8, max=14) sage: v.parent() Vector space of dimension 5 over Real Field with 53 bits of precision sage: all(8 <= r <= 14 for r in v) True
Any ring with a
random_element()
method may be used.sage: F = FiniteField(23) sage: hasattr(F, 'random_element') True sage: v = random_vector(F, 10) sage: v.parent() Vector space of dimension 10 over Finite Field of size 23
The default implementation is a dense representation, equivalent to setting
sparse=False
.sage: v = random_vector(10) sage: v.is_sparse() False sage: w = random_vector(ZZ, 20, sparse=True) sage: w.is_sparse() True
The elements are chosen using the ring’s
random_element
method:sage: from sage.misc.randstate import current_randstate sage: seed = current_randstate().seed() sage: set_random_seed(seed) sage: v1 = random_vector(ZZ, 20, distribution="1/n") sage: v2 = random_vector(ZZ, 15, x=-1000, y=1000) sage: v3 = random_vector(QQ, 10) sage: v4 = random_vector(FiniteField(17), 10) sage: v5 = random_vector(RR, 10) sage: set_random_seed(seed) sage: w1 = vector(ZZ.random_element(distribution="1/n") for _ in range(20)) sage: w2 = vector(ZZ.random_element(x=-1000, y=1000) for _ in range(15)) sage: w3 = vector(QQ.random_element() for _ in range(10)) sage: w4 = vector(FiniteField(17).random_element() for _ in range(10)) sage: w5 = vector(RR.random_element() for _ in range(10)) sage: [v1, v2, v3, v4, v5] == [w1, w2, w3, w4, w5] True
Inputs get checked before constructing the vector.
sage: random_vector('junk') Traceback (most recent call last): ... TypeError: degree of a random vector must be an integer, not None sage: random_vector('stuff', 5) Traceback (most recent call last): ... TypeError: elements of a vector, or module element, must come from a ring, not stuff sage: random_vector(ZZ, -9) Traceback (most recent call last): ... ValueError: degree of a random vector must be non-negative, not -9
- sage.modules.free_module_element.vector(arg0, arg1=None, arg2=None, sparse=None, immutable=False)#
Return a vector or free module element with specified entries.
CALL FORMATS:
This constructor can be called in several different ways. In each case,
sparse=True
orsparse=False
as well asimmutable=True
orimmutable=False
can be supplied as an option.free_module_element()
is an alias forvector()
.vector(object)
vector(ring, object)
vector(object, ring)
vector(ring, degree, object)
vector(ring, degree)
INPUT:
object
– a list, dictionary, or other iterable containing the entries of the vector, including any object that is palatable to theSequence
constructorring
– a base ring (or field) for the vector space or free module, which contains all of the elementsdegree
– an integer specifying the number of entries in the vector or free module elementsparse
– boolean, whether the result should be a sparse vectorimmutable
– boolean (default:False
); whether the result should be an immutable vector
In call format 4, an error is raised if the
degree
does not match the length ofobject
so this call can provide some safeguards. Note however that using this format whenobject
is a dictionary is unlikely to work properly.OUTPUT:
An element of the ambient vector space or free module with the given base ring and implied or specified dimension or rank, containing the specified entries and with correct degree.
In call format 5, no entries are specified, so the element is populated with all zeros.
If the
sparse
option is not supplied, the output will generally have a dense representation. The exception is ifobject
is a dictionary, then the representation will be sparse.EXAMPLES:
sage: v = vector([1,2,3]); v (1, 2, 3) sage: v.parent() Ambient free module of rank 3 over the principal ideal domain Integer Ring sage: v = vector([1,2,3/5]); v (1, 2, 3/5) sage: v.parent() Vector space of dimension 3 over Rational Field
All entries must canonically coerce to some common ring:
sage: v = vector([17, GF(11)(5), 19/3]); v Traceback (most recent call last): ... TypeError: unable to find a common ring for all elements
sage: v = vector([17, GF(11)(5), 19]); v (6, 5, 8) sage: v.parent() Vector space of dimension 3 over Finite Field of size 11 sage: v = vector([17, GF(11)(5), 19], QQ); v (17, 5, 19) sage: v.parent() Vector space of dimension 3 over Rational Field sage: v = vector((1,2,3), QQ); v (1, 2, 3) sage: v.parent() Vector space of dimension 3 over Rational Field sage: v = vector(QQ, (1,2,3)); v (1, 2, 3) sage: v.parent() Vector space of dimension 3 over Rational Field sage: v = vector(vector([1,2,3])); v (1, 2, 3) sage: v.parent() Ambient free module of rank 3 over the principal ideal domain Integer Ring
You can also use
free_module_element
, which is the same asvector
.sage: free_module_element([1/3, -4/5]) (1/3, -4/5)
We make a vector mod 3 out of a vector over \(\ZZ\).
sage: vector(vector([1,2,3]), GF(3)) (1, 2, 0)
The degree of a vector may be specified:
sage: vector(QQ, 4, [1,1/2,1/3,1/4]) (1, 1/2, 1/3, 1/4)
But it is an error if the degree and size of the list of entries are mismatched:
sage: vector(QQ, 5, [1,1/2,1/3,1/4]) Traceback (most recent call last): ... ValueError: incompatible degrees in vector constructor
Providing no entries populates the vector with zeros, but of course, you must specify the degree since it is not implied. Here we use a finite field as the base ring.
sage: w = vector(FiniteField(7), 4); w (0, 0, 0, 0) sage: w.parent() Vector space of dimension 4 over Finite Field of size 7
The fastest method to construct a zero vector is to call the
zero_vector()
method directly on a free module or vector space, since vector(…) must do a small amount of type checking. Almost as fast as thezero_vector()
method is thezero_vector()
constructor, which defaults to the integers.sage: vector(ZZ, 5) # works fine (0, 0, 0, 0, 0) sage: (ZZ^5).zero_vector() # very tiny bit faster (0, 0, 0, 0, 0) sage: zero_vector(ZZ, 5) # similar speed to vector(...) (0, 0, 0, 0, 0) sage: z = zero_vector(5); z (0, 0, 0, 0, 0) sage: z.parent() Ambient free module of rank 5 over the principal ideal domain Integer Ring
Here we illustrate the creation of sparse vectors by using a dictionary:
sage: vector({1:1.1, 3:3.14}) (0.000000000000000, 1.10000000000000, 0.000000000000000, 3.14000000000000)
With no degree given, a dictionary of entries implicitly declares a degree by the largest index (key) present. So you can provide a terminal element (perhaps a zero?) to set the degree. But it is probably safer to just include a degree in your construction.
sage: v = vector(QQ, {0:1/2, 4:-6, 7:0}); v (1/2, 0, 0, 0, -6, 0, 0, 0) sage: v.degree() 8 sage: v.is_sparse() True sage: w = vector(QQ, 8, {0:1/2, 4:-6}) sage: w == v True
It is an error to specify a negative degree.
sage: vector(RR, -4, [1.0, 2.0, 3.0, 4.0]) Traceback (most recent call last): ... ValueError: cannot specify the degree of a vector as a negative integer (-4)
It is an error to create a zero vector but not provide a ring as the first argument.
sage: vector('junk', 20) Traceback (most recent call last): ... TypeError: first argument must be base ring of zero vector, not junk
And it is an error to specify an index in a dictionary that is greater than or equal to a requested degree.
sage: vector(ZZ, 10, {3:4, 7:-2, 10:637}) Traceback (most recent call last): ... ValueError: dictionary of entries has a key (index) exceeding the requested degree
A 1-dimensional numpy array of type float or complex may be passed to vector. Unless an explicit ring is given, the result will be a vector in the appropriate dimensional vector space over the real double field or the complex double field. The data in the array must be contiguous, so column-wise slices of numpy matrices will raise an exception.
sage: import numpy sage: x = numpy.random.randn(10) sage: y = vector(x) sage: parent(y) Vector space of dimension 10 over Real Double Field sage: parent(vector(RDF, x)) Vector space of dimension 10 over Real Double Field sage: parent(vector(CDF, x)) Vector space of dimension 10 over Complex Double Field sage: parent(vector(RR, x)) Vector space of dimension 10 over Real Field with 53 bits of precision sage: v = numpy.random.randn(10) * complex(0,1) sage: w = vector(v) sage: parent(w) Vector space of dimension 10 over Complex Double Field
Multi-dimensional arrays are not supported:
sage: import numpy as np sage: a = np.array([[1, 2, 3], [4, 5, 6]], np.float64) sage: vector(a) Traceback (most recent call last): ... TypeError: cannot convert 2-dimensional array to a vector
If any of the arguments to vector have Python type int, long, real, or complex, they will first be coerced to the appropriate Sage objects. This fixes trac ticket #3847.
sage: v = vector([int(0)]); v (0) sage: v[0].parent() Integer Ring sage: v = vector(range(10)); v (0, 1, 2, 3, 4, 5, 6, 7, 8, 9) sage: v[3].parent() Integer Ring sage: v = vector([float(23.4), int(2), complex(2+7*I), 1]); v (23.4, 2.0, 2.0 + 7.0*I, 1.0) sage: v[1].parent() Complex Double Field
If the argument is a vector, it doesn’t change the base ring. This fixes trac ticket #6643:
sage: K.<sqrt3> = QuadraticField(3) sage: u = vector(K, (1/2, sqrt3/2) ) sage: vector(u).base_ring() Number Field in sqrt3 with defining polynomial x^2 - 3 with sqrt3 = 1.732050807568878? sage: v = vector(K, (0, 1) ) sage: vector(v).base_ring() Number Field in sqrt3 with defining polynomial x^2 - 3 with sqrt3 = 1.732050807568878?
Constructing a vector from a numpy array behaves as expected:
sage: import numpy sage: a=numpy.array([1,2,3]) sage: v=vector(a); v (1, 2, 3) sage: parent(v) Ambient free module of rank 3 over the principal ideal domain Integer Ring
Complex numbers can be converted naturally to a sequence of length 2. And then to a vector.
sage: c = CDF(2 + 3*I) sage: v = vector(c); v (2.0, 3.0)
A generator, or other iterable, may also be supplied as input. Anything that can be converted to a
Sequence
is a possible input.sage: type(i^2 for i in range(3)) <... 'generator'> sage: v = vector(i^2 for i in range(3)); v (0, 1, 4)
An empty list, without a ring given, will default to the integers.
sage: x = vector([]); x () sage: x.parent() Ambient free module of rank 0 over the principal ideal domain Integer Ring
The
immutable
switch allows to create an immutable vector.sage: v = vector(QQ, {0:1/2, 4:-6, 7:0}, immutable=True); v (1/2, 0, 0, 0, -6, 0, 0, 0) sage: v.is_immutable() True
The
immutable
switch works regardless of the type of valid input to the constructor.sage: v = vector(ZZ, 4, immutable=True) sage: v.is_immutable() True sage: w = vector(ZZ, [1,2,3]) sage: v = vector(w, ZZ, immutable=True) sage: v.is_immutable() True sage: v = vector(QQ, w, immutable=True) sage: v.is_immutable() True sage: import numpy as np sage: w = np.array([1, 2, pi], float) sage: v = vector(w, immutable=True) sage: v.is_immutable() True sage: w = np.array([i, 2, 3], complex) sage: v = vector(w, immutable=True) sage: v.is_immutable() True
- sage.modules.free_module_element.zero_vector(arg0, arg1=None)#
Returns a vector or free module element with a specified number of zeros.
CALL FORMATS:
zero_vector(degree)
zero_vector(ring, degree)
INPUT:
degree
- the number of zero entries in the vector or free module elementring
- defaultZZ
- the base ring of the vector space or module containing the constructed zero vector
OUTPUT:
A vector or free module element with
degree
entries, all equal to zero and belonging to the ring if specified. If no ring is given, a free module element overZZ
is returned.EXAMPLES:
A zero vector over the field of rationals.
sage: v = zero_vector(QQ, 5); v (0, 0, 0, 0, 0) sage: v.parent() Vector space of dimension 5 over Rational Field
A free module zero element.
sage: w = zero_vector(Integers(6), 3); w (0, 0, 0) sage: w.parent() Ambient free module of rank 3 over Ring of integers modulo 6
If no ring is given, the integers are used.
sage: u = zero_vector(9); u (0, 0, 0, 0, 0, 0, 0, 0, 0) sage: u.parent() Ambient free module of rank 9 over the principal ideal domain Integer Ring
Non-integer degrees produce an error.
sage: zero_vector(5.6) Traceback (most recent call last): ... TypeError: Attempt to coerce non-integral RealNumber to Integer
Negative degrees also give an error.
sage: zero_vector(-3) Traceback (most recent call last): ... ValueError: rank (=-3) must be nonnegative
Garbage instead of a ring will be recognized as such.
sage: zero_vector(x^2, 5) Traceback (most recent call last): ... TypeError: first argument must be a ring