Package sage :: Package rings :: Package polynomial :: Module laurent_polynomial_ring
[hide private]
[frames] | no frames]

Module laurent_polynomial_ring

source code


Ring of Laurent Polynomials.

If R is a commutative ring, then the ring of laurent polynomials in n variables over R is $R[x_1^{\pm 1}, x_2^{\pm 1}, \ldots, x_n^{\pm 1}].$
We implement it as a quotient ring $R[x_1, xx_1, x_2, xx_2, \ldots, x_n, xx_n] / (x_1*xx_1 - 1, x_2*xx_2 - 1, \ldots, x_n, xx_n - 1)$.

AUTHORS:
  -- David Roe (2008-2-23)



Classes [hide private]
  LaurentPolynomialRing_generic
  LaurentPolynomialRing_mpair
Functions [hide private]
 
is_LaurentPolynomialRing(R)
Returns True if and only if R is a Laurent polynomial ring.
source code
 
LaurentPolynomialRing(base_ring, arg1=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a..., arg2=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a..., sparse=False, order='degrevlex', names=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a..., name=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a...)
Return the globally unique univariate or multivariate laurent polynomial ring with given properties and variable name or names.
source code
 
_get_from_cache(key)
EXAMPLES:...
source code
 
_save_in_cache(key, R)
EXAMPLES:...
source code
 
_single_variate(base_ring, names, sparse)
EXAMPLES:...
source code
 
_multi_variate(base_ring, names, n, sparse, order)
EXAMPLES:...
source code
Variables [hide private]
  _cache = {}
Function Details [hide private]

is_LaurentPolynomialRing(R)

source code 

Returns True if and only if R is a Laurent polynomial ring.

EXAMPLES:
    sage: P = PolynomialRing(QQ,2,'x')
    sage: is_LaurentPolynomialRing(P)
    False

    sage: R = LaurentPolynomialRing(QQ,3,'x')
    sage: is_LaurentPolynomialRing(R)
    True

LaurentPolynomialRing(base_ring, arg1=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a..., arg2=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a..., sparse=False, order='degrevlex', names=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a..., name=['4ti2-20061025', 'R-2.6.0', 'atlas-3.7.37', 'atlas-3.8.1', 'a...)

source code 

Return the globally unique univariate or multivariate laurent polynomial
ring with given properties and variable name or names.

There are four ways to call the polynomial ring constructor:
      1. LaurentPolynomialRing(base_ring, name,    sparse=False)
      2. LaurentPolynomialRing(base_ring, names,   order='degrevlex')
      3. LaurentPolynomialRing(base_ring, name, n, order='degrevlex')
      4. LaurentPolynomialRing(base_ring, n, name, order='degrevlex')

The optional arguments sparse and order *must* be explicitly
named, and the other arguments must be given positionally.

INPUT:
     base_ring -- a commutative ring
     name -- a string
     names -- a list or tuple of names, or a comma separated string
     n -- a positive integer
     sparse -- bool (default: False), whether or not elements are sparse
     order -- string or TermOrder, e.g.,
             'degrevlex' (default) -- degree reverse lexicographic
             'lex'  -- lexicographic
             'deglex' -- degree lexicographic
             TermOrder('deglex',3) + TermOrder('deglex',3) -- block ordering
             
OUTPUT:
    LaurentPolynomialRing(base_ring, name, sparse=False) returns a univariate
    laurent polynomial ring; all other input formats return a multivariate
    laurent polynomial ring.

UNIQUENESS and IMMUTABILITY: In SAGE there is exactly one
single-variate laurent polynomial ring over each base ring in each choice
of variable and sparseness.  There is also exactly one multivariate
laurent polynomial ring over each base ring for each choice of names of
variables and term order.  

        sage: R.<x,y> = LaurentPolynomialRing(QQ,2); R
        Multivariate Laurent Polynomial Ring in x, y over Rational Field
        sage: f = x^2 - 2*y^-2

    You can't just globally change the names of those variables.
    This is because objects all over SAGE could have pointers to
    that polynomial ring.
        sage: R._assign_names(['z','w'])
        Traceback (most recent call last):
        ...
        ValueError: variable names cannot be changed after object creation.

        
EXAMPLES:
1. LaurentPolynomialRing(base_ring, name,    sparse=False):
    sage: LaurentPolynomialRing(QQ, 'w')
    Univariate Laurent Polynomial Ring in w over Rational Field

Use the diamond brackets notation to make the variable
ready for use after you define the ring:
    sage: R.<w> = LaurentPolynomialRing(QQ)
    sage: (1 + w)^3
    w^3 + 3*w^2 + 3*w + 1
    
You must specify a name:
    sage: LaurentPolynomialRing(QQ)
    Traceback (most recent call last):
    ...
    TypeError: You must specify the names of the variables.

    sage: R.<abc> = LaurentPolynomialRing(QQ, sparse=True); R
    Univariate Laurent Polynomial Ring in abc over Rational Field

    sage: R.<w> = LaurentPolynomialRing(PolynomialRing(GF(7),'k')); R
    Univariate Laurent Polynomial Ring in w over Univariate Polynomial Ring in k over Finite Field of size 7

Rings with different variables are different:
    sage: LaurentPolynomialRing(QQ, 'x') == LaurentPolynomialRing(QQ, 'y')
    False
    
2. LaurentPolynomialRing(base_ring, names,   order='degrevlex'):
    sage: R = LaurentPolynomialRing(QQ, 'a,b,c'); R
    Multivariate Laurent Polynomial Ring in a, b, c over Rational Field

    sage: S = LaurentPolynomialRing(QQ, ['a','b','c']); S
    Multivariate Laurent Polynomial Ring in a, b, c over Rational Field

    sage: T = LaurentPolynomialRing(QQ, ('a','b','c')); T
    Multivariate Laurent Polynomial Ring in a, b, c over Rational Field

All three rings are identical.
    sage: (R is S) and  (S is T)
    True

There is a unique laurent polynomial ring with each term order:
    sage: R = LaurentPolynomialRing(QQ, 'x,y,z', order='degrevlex'); R
    Multivariate Laurent Polynomial Ring in x, y, z over Rational Field
    sage: S = LaurentPolynomialRing(QQ, 'x,y,z', order='invlex'); S
    Multivariate Laurent Polynomial Ring in x, y, z over Rational Field
    sage: S is LaurentPolynomialRing(QQ, 'x,y,z', order='invlex')
    True
    sage: R == S
    False


3. LaurentPolynomialRing(base_ring, name, n, order='degrevlex'):

If you specify a single name as a string and a number of
variables, then variables labeled with numbers are created.
    sage: LaurentPolynomialRing(QQ, 'x', 10)
    Multivariate Laurent Polynomial Ring in x0, x1, x2, x3, x4, x5, x6, x7, x8, x9 over Rational Field
    
    sage: LaurentPolynomialRing(GF(7), 'y', 5)
    Multivariate Laurent Polynomial Ring in y0, y1, y2, y3, y4 over Finite Field of size 7

    sage: LaurentPolynomialRing(QQ, 'y', 3, sparse=True)
    Multivariate Laurent Polynomial Ring in y0, y1, y2 over Rational Field

You can call \code{injvar}, which is a convenient shortcut for \code{inject_variables()}.
    sage: R = LaurentPolynomialRing(GF(7),15,'w'); R
    Multivariate Laurent Polynomial Ring in w0, w1, w2, w3, w4, w5, w6, w7, w8, w9, w10, w11, w12, w13, w14 over Finite Field of size 7        
    sage: R.injvar()
    Defining w0, w1, w2, w3, w4, w5, w6, w7, w8, w9, w10, w11, w12, w13, w14
    sage: (w0 + 2*w8 + w13)^2
    w0^2 + 4*w0*w8 + 4*w8^2 + 2*w0*w13 + 4*w8*w13 + w13^2

_get_from_cache(key)

source code 

EXAMPLES:
    sage: from sage.rings.polynomial.laurent_polynomial_ring import _get_from_cache
    sage: L = LaurentPolynomialRing(QQ,2,'x')
    sage: L2 = _get_from_cache( (QQ,('x0','x1'),2,False,TermOrder('degrevlex')) ); L2
    Multivariate Laurent Polynomial Ring in x0, x1 over Rational Field
    sage: L is L2
    True

_save_in_cache(key, R)

source code 

EXAMPLES:
    sage: from sage.rings.polynomial.laurent_polynomial_ring import _save_in_cache, _get_from_cache
    sage: L = LaurentPolynomialRing(QQ,2,'x')
    sage: _save_in_cache('testkey', L)
    sage: _get_from_cache('testkey')
    Multivariate Laurent Polynomial Ring in x0, x1 over Rational Field
    sage: _ is L
    True

_single_variate(base_ring, names, sparse)

source code 

EXAMPLES:
    sage: from sage.rings.polynomial.laurent_polynomial_ring import _single_variate
    sage: _single_variate(QQ, ('x',), False)
    Univariate Laurent Polynomial Ring in x over Rational Field

_multi_variate(base_ring, names, n, sparse, order)

source code 

EXAMPLES:
    sage: from sage.rings.polynomial.laurent_polynomial_ring import _multi_variate
    sage: _multi_variate(QQ, ('x','y'), 2, False, 'degrevlex')
    Multivariate Laurent Polynomial Ring in x, y over Rational Field