Wrappers for Giac functions#
We provide a python function to compute and convert to sage a Groebner
basis using the giacpy_sage
module.
AUTHORS:
Martin Albrecht (2015-07-01): initial version
Han Frederic (2015-07-01): initial version
EXAMPLES:
sage: from sage.libs.giac import groebner_basis as gb_giac # random
sage: P = PolynomialRing(QQ, 6, 'x')
sage: I = sage.rings.ideal.Cyclic(P)
sage: B = gb_giac(I.gens()) # random
sage: B
Polynomial Sequence with 45 Polynomials in 6 Variables
- class sage.libs.giac.GiacSettingsDefaultContext#
Bases:
object
Context preserve libgiac settings.
- sage.libs.giac.groebner_basis(gens, proba_epsilon=None, threads=None, prot=False, elim_variables=None, *args, **kwds)#
Compute a Groebner Basis of an ideal using
giacpy_sage
. The result is automatically converted to sage.Supported term orders of the underlying polynomial ring are
lex
,deglex
,degrevlex
and block orders with 2degrevlex
blocks.INPUT:
gens
- an ideal (or a list) of polynomials over a prime field of characteristic 0 or p<2^31proba_epsilon
- (default: None) majoration of the probabilityof a wrong answer when probabilistic algorithms are allowed.
if
proba_epsilon
is None, the value ofsage.structure.proof.all.polynomial()
is taken. If it is false then the globalgiacpy_sage.giacsettings.proba_epsilon
is used.if
proba_epsilon
is 0, probabilistic algorithms are disabled.
threads
- (default: None) Maximal number of threads allowed for giac. If None, the globalgiacpy_sage.giacsettings.threads
is considered.prot
- (default: False) if True print detailled informationselim_variables
- (default: None) a list of variables to eliminate from the ideal.if
elim_variables
is None, a Groebner basis with respect to the term ordering of the parent polynomial ring of the polynomialsgens
is computed.if
elim_variables
is a list of variables, a Groebner basis of the elimination ideal with respect to adegrevlex
term order is computed, regardless of the term order of the polynomial ring.
OUTPUT:
Polynomial sequence of the reduced Groebner basis.
EXAMPLES:
sage: from sage.libs.giac import groebner_basis as gb_giac sage: P = PolynomialRing(GF(previous_prime(2**31)), 6, 'x') sage: I = sage.rings.ideal.Cyclic(P) sage: B=gb_giac(I.gens());B // Groebner basis computation time ... Polynomial Sequence with 45 Polynomials in 6 Variables sage: B.is_groebner() True
Elimination ideals can be computed by passing
elim_variables
:sage: P = PolynomialRing(GF(previous_prime(2**31)), 5, 'x') sage: I = sage.rings.ideal.Cyclic(P) sage: B = gb_giac(I.gens(), elim_variables=[P.gen(0), P.gen(2)]) // Groebner basis computation time ... sage: B.is_groebner() True sage: B.ideal() == I.elimination_ideal([P.gen(0), P.gen(2)]) True
Computations over QQ can benefit from
a probabilistic lifting:
sage: P = PolynomialRing(QQ,5, 'x') sage: I = ideal([P.random_element(3,7) for j in range(5)]) sage: B1 = gb_giac(I.gens(),1e-16) # long time (1s) ... sage: sage.structure.proof.all.polynomial(True) sage: B2 = gb_giac(I.gens()) # long time (4s) // Groebner basis computation time... sage: B1 == B2 # long time True sage: B1.is_groebner() # long time (20s) True
multi threaded operations:
sage: P = PolynomialRing(QQ, 8, 'x') sage: I = sage.rings.ideal.Cyclic(P) sage: time B = gb_giac(I.gens(),1e-6,threads=2) # doctest: +SKIP ... Time: CPU 168.98 s, Wall: 94.13 s
You can get detailled information by setting
prot=True
sage: I = sage.rings.ideal.Katsura(P) sage: gb_giac(I,prot=True) # random, long time (3s) 9381383 begin computing basis modulo 535718473 9381501 begin new iteration zmod, number of pairs: 8, base size: 8 ...end, basis size 74 prime number 1 G=Vector [0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18,... ...creating reconstruction #0 ... ++++++++basis size 74 checking pairs for i=0, j= checking pairs for i=1, j=2,6,12,17,19,24,29,34,39,42,43,48,56,61,64,69, ... checking pairs for i=72, j=73, checking pairs for i=73, j= Number of critical pairs to check 373 +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++... Successful... check of 373 critical pairs 12380865 end final check Polynomial Sequence with 74 Polynomials in 8 Variables
- sage.libs.giac.local_giacsettings(func)#
Decorator to preserve Giac’s proba_epsilon and threads settings.
EXAMPLES:
sage: def testf(a,b): ....: giacsettings.proba_epsilon = a/100 ....: giacsettings.threads = b+2 ....: return (giacsettings.proba_epsilon, giacsettings.threads) sage: from sage.libs.giac.giac import giacsettings sage: from sage.libs.giac import local_giacsettings sage: gporig, gtorig = (giacsettings.proba_epsilon,giacsettings.threads) sage: gp, gt = local_giacsettings(testf)(giacsettings.proba_epsilon,giacsettings.threads) sage: gporig == giacsettings.proba_epsilon True sage: gtorig == giacsettings.threads True sage: gp<gporig, gt-gtorig (True, 2)