Direct low-level access to SINGULAR’s Groebner basis engine via libSINGULAR

AUTHOR:

  • Martin Albrecht (2007-08-08): initial version

EXAMPLES:

sage: x,y,z = QQ['x,y,z'].gens()
sage: I = ideal(x^5 + y^4 + z^3 - 1,  x^3 + y^3 + z^2 - 1)
sage: I.groebner_basis('libsingular:std')
[y^6 + x*y^4 + 2*y^3*z^2 + x*z^3 + z^4 - 2*y^3 - 2*z^2 - x + 1,
 x^2*y^3 - y^4 + x^2*z^2 - z^3 - x^2 + 1, x^3 + y^3 + z^2 - 1]
>>> from sage.all import *
>>> x,y,z = QQ['x,y,z'].gens()
>>> I = ideal(x**Integer(5) + y**Integer(4) + z**Integer(3) - Integer(1),  x**Integer(3) + y**Integer(3) + z**Integer(2) - Integer(1))
>>> I.groebner_basis('libsingular:std')
[y^6 + x*y^4 + 2*y^3*z^2 + x*z^3 + z^4 - 2*y^3 - 2*z^2 - x + 1,
 x^2*y^3 - y^4 + x^2*z^2 - z^3 - x^2 + 1, x^3 + y^3 + z^2 - 1]

We compute a Groebner basis for cyclic 6, which is a standard benchmark and test ideal:

sage: R.<x,y,z,t,u,v> = QQ['x,y,z,t,u,v']
sage: I = sage.rings.ideal.Cyclic(R,6)
sage: B = I.groebner_basis('libsingular:std')
sage: len(B)
45
>>> from sage.all import *
>>> R = QQ['x,y,z,t,u,v']; (x, y, z, t, u, v,) = R._first_ngens(6)
>>> I = sage.rings.ideal.Cyclic(R,Integer(6))
>>> B = I.groebner_basis('libsingular:std')
>>> len(B)
45

Two examples from the Mathematica documentation (done in Sage):

  • We compute a Groebner basis:

    sage: R.<x,y> = PolynomialRing(QQ, order='lex')
    sage: ideal(x^2 - 2*y^2, x*y - 3).groebner_basis('libsingular:slimgb')
    [x - 2/3*y^3, y^4 - 9/2]
    
    >>> from sage.all import *
    >>> R = PolynomialRing(QQ, order='lex', names=('x', 'y',)); (x, y,) = R._first_ngens(2)
    >>> ideal(x**Integer(2) - Integer(2)*y**Integer(2), x*y - Integer(3)).groebner_basis('libsingular:slimgb')
    [x - 2/3*y^3, y^4 - 9/2]
    
  • We show that three polynomials have no common root:

    sage: R.<x,y> = QQ[]
    sage: ideal(x+y, x^2 - 1, y^2 - 2*x).groebner_basis('libsingular:slimgb')
    [1]
    
    >>> from sage.all import *
    >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2)
    >>> ideal(x+y, x**Integer(2) - Integer(1), y**Integer(2) - Integer(2)*x).groebner_basis('libsingular:slimgb')
    [1]
    
sage.rings.polynomial.multi_polynomial_ideal_libsingular.interred_libsingular(I)[source]

SINGULAR’s interred() command.

INPUT:

  • I – a Sage ideal

EXAMPLES:

sage: P.<x,y,z> = PolynomialRing(ZZ)
sage: I = ideal( x^2 - 3*y, y^3 - x*y, z^3 - x, x^4 - y*z + 1 )
sage: I.interreduced_basis()
[y*z^2 - 81*x*y - 9*y - z, z^3 - x, x^2 - 3*y, 9*y^2 - y*z + 1]

sage: P.<x,y,z> = PolynomialRing(QQ)
sage: I = ideal( x^2 - 3*y, y^3 - x*y, z^3 - x, x^4 - y*z + 1 )
sage: I.interreduced_basis()
[y*z^2 - 81*x*y - 9*y - z, z^3 - x, x^2 - 3*y, y^2 - 1/9*y*z + 1/9]
>>> from sage.all import *
>>> P = PolynomialRing(ZZ, names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> I = ideal( x**Integer(2) - Integer(3)*y, y**Integer(3) - x*y, z**Integer(3) - x, x**Integer(4) - y*z + Integer(1) )
>>> I.interreduced_basis()
[y*z^2 - 81*x*y - 9*y - z, z^3 - x, x^2 - 3*y, 9*y^2 - y*z + 1]

>>> P = PolynomialRing(QQ, names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> I = ideal( x**Integer(2) - Integer(3)*y, y**Integer(3) - x*y, z**Integer(3) - x, x**Integer(4) - y*z + Integer(1) )
>>> I.interreduced_basis()
[y*z^2 - 81*x*y - 9*y - z, z^3 - x, x^2 - 3*y, y^2 - 1/9*y*z + 1/9]
sage.rings.polynomial.multi_polynomial_ideal_libsingular.kbase_libsingular(I, degree=None)[source]

SINGULAR’s kbase() algorithm.

INPUT:

  • I – a groebner basis of an ideal

  • degree – integer (default: None); if not None, return only the monomials of the given degree

OUTPUT:

Computes a vector space basis (consisting of monomials) of the quotient ring by the ideal, resp. of a free module by the module, in case it is finite dimensional and if the input is a standard basis with respect to the ring ordering. If the input is not a standard basis, the leading terms of the input are used and the result may have no meaning.

With two arguments: computes the part of a vector space basis of the respective quotient with degree of the monomials equal to the second argument. Here, the quotient does not need to be finite dimensional.

EXAMPLES:

sage: R.<x,y> = PolynomialRing(QQ, order='lex')
sage: I = R.ideal(x^2-2*y^2, x*y-3)
sage: I.normal_basis()  # indirect doctest
[y^3, y^2, y, 1]
sage: J = R.ideal(x^2-2*y^2)
sage: [J.normal_basis(d) for d in (0..4)]  # indirect doctest
[[1], [y, x], [y^2, x*y], [y^3, x*y^2], [y^4, x*y^3]]
>>> from sage.all import *
>>> R = PolynomialRing(QQ, order='lex', names=('x', 'y',)); (x, y,) = R._first_ngens(2)
>>> I = R.ideal(x**Integer(2)-Integer(2)*y**Integer(2), x*y-Integer(3))
>>> I.normal_basis()  # indirect doctest
[y^3, y^2, y, 1]
>>> J = R.ideal(x**Integer(2)-Integer(2)*y**Integer(2))
>>> [J.normal_basis(d) for d in (ellipsis_iter(Integer(0),Ellipsis,Integer(4)))]  # indirect doctest
[[1], [y, x], [y^2, x*y], [y^3, x*y^2], [y^4, x*y^3]]
sage.rings.polynomial.multi_polynomial_ideal_libsingular.slimgb_libsingular(I)[source]

SINGULAR’s slimgb() algorithm.

INPUT:

  • I – a Sage ideal

sage.rings.polynomial.multi_polynomial_ideal_libsingular.std_libsingular(I)[source]

SINGULAR’s std() algorithm.

INPUT:

  • I – a Sage ideal