Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
""" 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]
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
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]
- 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] """
#***************************************************************************** # Copyright (C) 2007 Martin Albrecht <malb@informatik.uni-bremen.de> # Copyright (C) 2007 William Stein <wstein@gmail.com> # # This program is free software: you can redistribute it and/or modify # it under the terms of the GNU General Public License as published by # the Free Software Foundation, either version 2 of the License, or # (at your option) any later version. # http://www.gnu.org/licenses/ #*****************************************************************************
from cysignals.signals cimport sig_on, sig_off
from sage.libs.singular.decl cimport tHomog, number, IDELEMS, p_Copy, rChangeCurrRing from sage.libs.singular.decl cimport idInit, id_Delete, currRing, Sy_bit, OPT_REDSB from sage.libs.singular.decl cimport scKBase, poly, testHomog, idSkipZeroes, id_RankFreeModule, kStd from sage.libs.singular.decl cimport OPT_REDTAIL, singular_options, kInterRed, t_rep_gb, p_GetCoeff from sage.libs.singular.decl cimport pp_Mult_nn, p_Delete, n_Delete from sage.libs.singular.decl cimport rIsPluralRing from sage.libs.singular.decl cimport n_unknown, n_Zp, n_Q, n_R, n_GF, n_long_R, n_algExt,n_transExt,n_long_C, n_Z, n_Zn, n_Znm, n_Z2m, n_CF
from sage.rings.polynomial.multi_polynomial_libsingular cimport new_MP from sage.rings.polynomial.plural cimport new_NCP
from sage.rings.polynomial.multi_polynomial_ideal import MPolynomialIdeal from sage.rings.polynomial.multi_polynomial_libsingular cimport MPolynomial_libsingular from sage.rings.polynomial.multi_polynomial_libsingular cimport MPolynomialRing_libsingular from sage.structure.sequence import Sequence
from sage.rings.polynomial.plural cimport NCPolynomialRing_plural, NCPolynomial_plural
cdef object singular_ideal_to_sage_sequence(ideal *i, ring *r, object parent): """ convert a SINGULAR ideal to a Sage Sequence (the format Sage stores a Groebner basis in)
INPUT:
- ``i`` -- a SINGULAR ideal - ``r`` -- a SINGULAR ring - ``parent`` -- a Sage ring matching r """ cdef int j cdef MPolynomial_libsingular p cdef NCPolynomial_plural p_nc
else:
cdef ideal *sage_ideal_to_singular_ideal(I) except NULL: """ convert a Sage ideal to a SINGULAR ideal
INPUT:
- ``I`` -- a Sage ideal in a ring of type :class:`~sage.rings.polynomial.multi_polynomial_libsingular.MPolynomialRing_libsingular` or a list of generators.
TESTS:
We test conversion::
sage: P.<x,y,z> = QQ[] sage: sage.libs.singular.function_factory.ff.std(Sequence([x,y,z])) [z, y, x] sage: sage.libs.singular.function_factory.ff.std(Ideal([x,y,z])) [z, y, x] """ cdef ideal *result cdef ring *r cdef ideal *i
else: raise TypeError("Ring must be of type 'MPolynomialRing_libsingular'")
else: id_Delete(&i, r) raise TypeError("All generators must be of type MPolynomial_libsingular.")
def kbase_libsingular(I): """ SINGULAR's ``kbase()`` algorithm.
INPUT:
- ``I`` -- a groebner basis of an ideal
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.
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() [y^3, y^2, y, 1] """
global singular_options
cdef ideal *result
def std_libsingular(I): """ SINGULAR's ``std()`` algorithm.
INPUT:
- ``I`` -- a Sage ideal """ global singular_options
cdef ideal *result
def slimgb_libsingular(I): """ SINGULAR's ``slimgb()`` algorithm.
INPUT:
- ``I`` -- a Sage ideal """ global singular_options
cdef ideal *i cdef ring *r cdef ideal *result
id_Delete(&i, r) raise TypeError("ordering must be global for slimgb")
id_Delete(&i, r) raise TypeError
def interred_libsingular(I): """ 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] """ global singular_options
cdef ideal *i cdef ideal *result cdef ring *r cdef number *n cdef poly *p cdef int j cdef int bck
return Sequence([], check=False, immutable=True) pass
# divide head by coefficients
|