Groebner Fans¶
Sage provides much of the functionality of gfan
, which is a
software package whose main function is to enumerate all reduced
Groebner bases of a polynomial ideal. The reduced Groebner bases
yield the maximal cones in the Groebner fan of the ideal. Several
subcomputations can be issued and additional tools are included.
Among these the highlights are:
Commands for computing tropical varieties.
Interactive walks in the Groebner fan of an ideal.
Commands for graphical renderings of Groebner fans and monomial ideals.
AUTHORS:
Anders Nedergaard Jensen: Wrote the
gfan
C++ program, which implements algorithms many of which were invented by Jensen, Komei Fukuda, and Rekha Thomas. All the underlying hard work of the Groebner fans functionality of Sage depends on this C++ program.William Stein (2006-04-20): Wrote first version of the Sage code for working with Groebner fans.
Tristram Bogart: the design of the Sage interface to
gfan
is joint work with Tristram Bogart, who also supplied numerous examples.Marshall Hampton (2008-03-25): Rewrote various functions to use
gfan-0.3
. This is still a work in progress, comments are appreciated on sage-devel@googlegroups.com (or personally at hamptonio@gmail.com).
EXAMPLES:
sage: x,y = QQ['x,y'].gens()
sage: i = ideal(x^2 - y^2 + 1)
sage: g = i.groebner_fan()
sage: g.reduced_groebner_bases()
[[x^2 - y^2 + 1], [-x^2 + y^2 - 1]]
>>> from sage.all import *
>>> x,y = QQ['x,y'].gens()
>>> i = ideal(x**Integer(2) - y**Integer(2) + Integer(1))
>>> g = i.groebner_fan()
>>> g.reduced_groebner_bases()
[[x^2 - y^2 + 1], [-x^2 + y^2 - 1]]
x,y = QQ['x,y'].gens() i = ideal(x^2 - y^2 + 1) g = i.groebner_fan() g.reduced_groebner_bases()
REFERENCES:
Anders N. Jensen; Gfan, a software system for Groebner fans; http://home.math.au.dk/jensen/software/gfan/gfan.html
- class sage.rings.polynomial.groebner_fan.GroebnerFan(I, is_groebner_basis=False, symmetry=None, verbose=False)[source]¶
Bases:
SageObject
This class is used to access capabilities of the program
Gfan
.In addition to computing Groebner fans,
Gfan
can compute other things in tropical geometry such as tropical prevarieties.INPUT:
I
– ideal in a multivariate polynomial ringis_groebner_basis
– boolean (default:False
); ifTrue
, then I.gens() must be a Groebner basis with respect to the standard degree lexicographic term ordersymmetry
– (default:None
) if notNone
, describes symmetries of the idealverbose
– (default:False
) ifTrue
, printout useful info during computations
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: I = R.ideal([x^2*y - z, y^2*z - x, z^2*x - y]) sage: G = I.groebner_fan(); G Groebner fan of the ideal: Ideal (x^2*y - z, y^2*z - x, x*z^2 - y) of Multivariate Polynomial Ring in x, y, z over Rational Field
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> I = R.ideal([x**Integer(2)*y - z, y**Integer(2)*z - x, z**Integer(2)*x - y]) >>> G = I.groebner_fan(); G Groebner fan of the ideal: Ideal (x^2*y - z, y^2*z - x, x*z^2 - y) of Multivariate Polynomial Ring in x, y, z over Rational Field
R.<x,y,z> = QQ[] I = R.ideal([x^2*y - z, y^2*z - x, z^2*x - y]) G = I.groebner_fan(); G
Here is an example of the use of the tropical_intersection command, and then using the RationalPolyhedralFan class to compute the Stanley-Reisner ideal of the tropical prevariety:
sage: R.<x,y,z> = QQ[] sage: I = R.ideal([(x+y+z)^3-1,(x+y+z)^3-x,(x+y+z)-3]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: PF.rays() [[-1, 0, 0], [0, -1, 0], [0, 0, -1], [1, 1, 1]] sage: RPF = PF.to_RationalPolyhedralFan() sage: RPF.Stanley_Reisner_ideal(PolynomialRing(QQ,4,'A, B, C, D')) Ideal (A*B, A*C, B*C*D) of Multivariate Polynomial Ring in A, B, C, D over Rational Field
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> I = R.ideal([(x+y+z)**Integer(3)-Integer(1),(x+y+z)**Integer(3)-x,(x+y+z)-Integer(3)]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> PF.rays() [[-1, 0, 0], [0, -1, 0], [0, 0, -1], [1, 1, 1]] >>> RPF = PF.to_RationalPolyhedralFan() >>> RPF.Stanley_Reisner_ideal(PolynomialRing(QQ,Integer(4),'A, B, C, D')) Ideal (A*B, A*C, B*C*D) of Multivariate Polynomial Ring in A, B, C, D over Rational Field
R.<x,y,z> = QQ[] I = R.ideal([(x+y+z)^3-1,(x+y+z)^3-x,(x+y+z)-3]) GF = I.groebner_fan() PF = GF.tropical_intersection() PF.rays() RPF = PF.to_RationalPolyhedralFan() RPF.Stanley_Reisner_ideal(PolynomialRing(QQ,4,'A, B, C, D'))
- buchberger()[source]¶
Return a lexicographic reduced Groebner basis for the ideal.
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: G = R.ideal([x - z^3, y^2 - x + x^2 - z^3*x]).groebner_fan() sage: G.buchberger() [-z^3 + y^2, -z^3 + x]
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> G = R.ideal([x - z**Integer(3), y**Integer(2) - x + x**Integer(2) - z**Integer(3)*x]).groebner_fan() >>> G.buchberger() [-z^3 + y^2, -z^3 + x]
R.<x,y,z> = PolynomialRing(QQ,3) G = R.ideal([x - z^3, y^2 - x + x^2 - z^3*x]).groebner_fan() G.buchberger()
- characteristic()[source]¶
Return the characteristic of the base ring.
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: i1 = ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) sage: gf = i1.groebner_fan() sage: gf.characteristic() 0
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> i1 = ideal(x*z + Integer(6)*y*z - z**Integer(2), x*y + Integer(6)*x*z + y*z - z**Integer(2), y**Integer(2) + x*z + y*z) >>> gf = i1.groebner_fan() >>> gf.characteristic() 0
R.<x,y,z> = PolynomialRing(QQ,3) i1 = ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) gf = i1.groebner_fan() gf.characteristic()
- dimension_of_homogeneity_space()[source]¶
Return the dimension of the homogeneity space.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: G.dimension_of_homogeneity_space() 0
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> G.dimension_of_homogeneity_space() 0
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() G.dimension_of_homogeneity_space()
- gfan(cmd='bases', I=None, format=None)[source]¶
Return the
gfan
output as a string given an inputcmd
.The default is to produce the list of reduced Groebner bases in
gfan
format.INPUT:
cmd
– string (default:'bases'
); GFan commandI
– ideal (default:None
)format
– boolean (default:None
); deprecated
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: gf = R.ideal([x^3-y,y^3-x-1]).groebner_fan() sage: gf.gfan() 'Q[x,y]\n{{\ny^9-1-y+3*y^3-3*y^6,\nx+1-y^3}\n,\n{\nx^3-y,\ny^3-1-x}\n,\n{\nx^9-1-x,\ny-x^3}\n}\n'
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> gf = R.ideal([x**Integer(3)-y,y**Integer(3)-x-Integer(1)]).groebner_fan() >>> gf.gfan() 'Q[x,y]\n{{\ny^9-1-y+3*y^3-3*y^6,\nx+1-y^3}\n,\n{\nx^3-y,\ny^3-1-x}\n,\n{\nx^9-1-x,\ny-x^3}\n}\n'
R.<x,y> = PolynomialRing(QQ,2) gf = R.ideal([x^3-y,y^3-x-1]).groebner_fan() gf.gfan()
- homogeneity_space()[source]¶
Return the homogeneity space of a the list of polynomials that define this Groebner fan.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: H = G.homogeneity_space()
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> H = G.homogeneity_space()
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() H = G.homogeneity_space()
- ideal()[source]¶
Return the ideal the was used to define this Groebner fan.
EXAMPLES:
sage: R.<x1,x2> = PolynomialRing(QQ,2) sage: gf = R.ideal([x1^3-x2,x2^3-2*x1-2]).groebner_fan() sage: gf.ideal() Ideal (x1^3 - x2, x2^3 - 2*x1 - 2) of Multivariate Polynomial Ring in x1, x2 over Rational Field
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x1', 'x2',)); (x1, x2,) = R._first_ngens(2) >>> gf = R.ideal([x1**Integer(3)-x2,x2**Integer(3)-Integer(2)*x1-Integer(2)]).groebner_fan() >>> gf.ideal() Ideal (x1^3 - x2, x2^3 - 2*x1 - 2) of Multivariate Polynomial Ring in x1, x2 over Rational Field
R.<x1,x2> = PolynomialRing(QQ,2) gf = R.ideal([x1^3-x2,x2^3-2*x1-2]).groebner_fan() gf.ideal()
- interactive(*args, **kwds)[source]¶
See the documentation for self[0].interactive(). This does not work with the notebook.
EXAMPLES:
sage: print("This is not easily doc-testable; please write a good one!") This is not easily doc-testable; please write a good one!
>>> from sage.all import * >>> print("This is not easily doc-testable; please write a good one!") This is not easily doc-testable; please write a good one!
print("This is not easily doc-testable; please write a good one!")
- maximal_total_degree_of_a_groebner_basis()[source]¶
Return the maximal total degree of any Groebner basis.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: G.maximal_total_degree_of_a_groebner_basis() 4
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> G.maximal_total_degree_of_a_groebner_basis() 4
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() G.maximal_total_degree_of_a_groebner_basis()
- minimal_total_degree_of_a_groebner_basis()[source]¶
Return the minimal total degree of any Groebner basis.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: G.minimal_total_degree_of_a_groebner_basis() 2
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> G.minimal_total_degree_of_a_groebner_basis() 2
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() G.minimal_total_degree_of_a_groebner_basis()
- mixed_volume()[source]¶
Return the mixed volume of the generators of this ideal.
This is not really an ideal property, it can depend on the generators used.
The generators must give a square system (as many polynomials as variables).
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: example_ideal = R.ideal([x^2-y-1,y^2-z-1,z^2-x-1]) sage: gf = example_ideal.groebner_fan() sage: mv = gf.mixed_volume() sage: mv 8 sage: R2.<x,y> = QQ[] sage: g1 = 1 - x + x^7*y^3 + 2*x^8*y^4 sage: g2 = 2 + y + 3*x^7*y^3 + x^8*y^4 sage: example2 = R2.ideal([g1,g2]) sage: example2.groebner_fan().mixed_volume() 15
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> example_ideal = R.ideal([x**Integer(2)-y-Integer(1),y**Integer(2)-z-Integer(1),z**Integer(2)-x-Integer(1)]) >>> gf = example_ideal.groebner_fan() >>> mv = gf.mixed_volume() >>> mv 8 >>> R2 = QQ['x, y']; (x, y,) = R2._first_ngens(2) >>> g1 = Integer(1) - x + x**Integer(7)*y**Integer(3) + Integer(2)*x**Integer(8)*y**Integer(4) >>> g2 = Integer(2) + y + Integer(3)*x**Integer(7)*y**Integer(3) + x**Integer(8)*y**Integer(4) >>> example2 = R2.ideal([g1,g2]) >>> example2.groebner_fan().mixed_volume() 15
R.<x,y,z> = QQ[] example_ideal = R.ideal([x^2-y-1,y^2-z-1,z^2-x-1]) gf = example_ideal.groebner_fan() mv = gf.mixed_volume() mv R2.<x,y> = QQ[] g1 = 1 - x + x^7*y^3 + 2*x^8*y^4 g2 = 2 + y + 3*x^7*y^3 + x^8*y^4 example2 = R2.ideal([g1,g2]) example2.groebner_fan().mixed_volume()
- number_of_reduced_groebner_bases()[source]¶
Return the number of reduced Groebner bases.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: G.number_of_reduced_groebner_bases() 3
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> G.number_of_reduced_groebner_bases() 3
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() G.number_of_reduced_groebner_bases()
- number_of_variables()[source]¶
Return the number of variables.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: G.number_of_variables() 2
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> G.number_of_variables() 2
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() G.number_of_variables()
sage: R = PolynomialRing(QQ,'x',10) sage: R.inject_variables(globals()) Defining x0, x1, x2, x3, x4, x5, x6, x7, x8, x9 sage: G = ideal([x0 - x9, sum(R.gens())]).groebner_fan() sage: G.number_of_variables() 10
>>> from sage.all import * >>> R = PolynomialRing(QQ,'x',Integer(10)) >>> R.inject_variables(globals()) Defining x0, x1, x2, x3, x4, x5, x6, x7, x8, x9 >>> G = ideal([x0 - x9, sum(R.gens())]).groebner_fan() >>> G.number_of_variables() 10
R = PolynomialRing(QQ,'x',10) R.inject_variables(globals()) G = ideal([x0 - x9, sum(R.gens())]).groebner_fan() G.number_of_variables()
>>> from sage.all import * >>> R = PolynomialRing(QQ,'x',Integer(10)) >>> R.inject_variables(globals()) Defining x0, x1, x2, x3, x4, x5, x6, x7, x8, x9 >>> G = ideal([x0 - x9, sum(R.gens())]).groebner_fan() >>> G.number_of_variables() 10
R = PolynomialRing(QQ,'x',10) R.inject_variables(globals()) G = ideal([x0 - x9, sum(R.gens())]).groebner_fan() G.number_of_variables()
- polyhedralfan()[source]¶
Return a polyhedral fan object corresponding to the reduced Groebner bases.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-1]).groebner_fan() sage: pf = gf.polyhedralfan() sage: pf.rays() [[0, 0, 1], [0, 1, 0], [1, 0, 0]]
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(1)]).groebner_fan() >>> pf = gf.polyhedralfan() >>> pf.rays() [[0, 0, 1], [0, 1, 0], [1, 0, 0]]
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-1]).groebner_fan() pf = gf.polyhedralfan() pf.rays()
- reduced_groebner_bases()[source]¶
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ, 3, order='lex') sage: G = R.ideal([x^2*y - z, y^2*z - x, z^2*x - y]).groebner_fan() sage: X = G.reduced_groebner_bases() sage: len(X) 33 sage: X[0] [z^15 - z, x - z^9, y - z^11] sage: X[0].ideal() Ideal (z^15 - z, x - z^9, y - z^11) of Multivariate Polynomial Ring in x, y, z over Rational Field sage: X[:5] [[z^15 - z, x - z^9, y - z^11], [y^2 - z^8, x - z^9, y*z^4 - z, -y + z^11], [y^3 - z^5, x - y^2*z, y^2*z^3 - y, y*z^4 - z, -y^2 + z^8], [y^4 - z^2, x - y^2*z, y^2*z^3 - y, y*z^4 - z, -y^3 + z^5], [y^9 - z, y^6*z - y, x - y^2*z, -y^4 + z^2]] sage: R3.<x,y,z> = PolynomialRing(GF(2477),3) sage: gf = R3.ideal([300*x^3-y,y^2-z,z^2-12]).groebner_fan() sage: gf.reduced_groebner_bases() [[z^2 - 12, y^2 - z, x^3 + 933*y], [y^4 - 12, x^3 + 933*y, -y^2 + z], [x^6 - 1062*z, z^2 - 12, -300*x^3 + y], [x^12 + 200, -300*x^3 + y, -828*x^6 + z]]
>>> from sage.all import * >>> R = PolynomialRing(QQ, Integer(3), order='lex', names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> G = R.ideal([x**Integer(2)*y - z, y**Integer(2)*z - x, z**Integer(2)*x - y]).groebner_fan() >>> X = G.reduced_groebner_bases() >>> len(X) 33 >>> X[Integer(0)] [z^15 - z, x - z^9, y - z^11] >>> X[Integer(0)].ideal() Ideal (z^15 - z, x - z^9, y - z^11) of Multivariate Polynomial Ring in x, y, z over Rational Field >>> X[:Integer(5)] [[z^15 - z, x - z^9, y - z^11], [y^2 - z^8, x - z^9, y*z^4 - z, -y + z^11], [y^3 - z^5, x - y^2*z, y^2*z^3 - y, y*z^4 - z, -y^2 + z^8], [y^4 - z^2, x - y^2*z, y^2*z^3 - y, y*z^4 - z, -y^3 + z^5], [y^9 - z, y^6*z - y, x - y^2*z, -y^4 + z^2]] >>> R3 = PolynomialRing(GF(Integer(2477)),Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([Integer(300)*x**Integer(3)-y,y**Integer(2)-z,z**Integer(2)-Integer(12)]).groebner_fan() >>> gf.reduced_groebner_bases() [[z^2 - 12, y^2 - z, x^3 + 933*y], [y^4 - 12, x^3 + 933*y, -y^2 + z], [x^6 - 1062*z, z^2 - 12, -300*x^3 + y], [x^12 + 200, -300*x^3 + y, -828*x^6 + z]]
R.<x,y,z> = PolynomialRing(QQ, 3, order='lex') G = R.ideal([x^2*y - z, y^2*z - x, z^2*x - y]).groebner_fan() X = G.reduced_groebner_bases() len(X) X[0] X[0].ideal() X[:5] R3.<x,y,z> = PolynomialRing(GF(2477),3) gf = R3.ideal([300*x^3-y,y^2-z,z^2-12]).groebner_fan() gf.reduced_groebner_bases()
- render(file=None, larger=False, shift=0, rgbcolor=(0, 0, 0), polyfill=True, scale_colors=True)[source]¶
Render a Groebner fan as sage graphics or save as an xfig file.
More precisely, the output is a drawing of the Groebner fan intersected with a triangle. The corners of the triangle are (1,0,0) to the right, (0,1,0) to the left and (0,0,1) at the top. If there are more than three variables in the ring we extend these coordinates with zeros.
INPUT:
file
– a filename if you prefer the output saved to a file; this will be in xfig formatshift
– shift the positions of the variables in the drawing. For example, with shift=1, the corners will be b (right), c (left), and d (top). The shifting is done modulo the number of variables in the polynomial ring. The default is 0.larger
– boolean (default:False
); ifTrue
, make the triangle larger so that the shape of the Groebner region appears. Affects the xfig file but probably not the sage graphics (?).rgbcolor
– this will not affect the saved xfig file, only the sage graphics producedpolyfill
– whether or not to fill the cones with a color determined by the highest degree in each reduced Groebner basis for that conescale_colors
– ifTrue
, this will normalize color values to try to maximize the range
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x,z]).groebner_fan() sage: test_render = G.render() # needs sage.plot
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x,z]).groebner_fan() >>> test_render = G.render() # needs sage.plot
R.<x,y,z> = PolynomialRing(QQ,3) G = R.ideal([y^3 - x^2, y^2 - 13*x,z]).groebner_fan() test_render = G.render() # needs sage.plot
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: G = R.ideal([x^2*y - z, y^2*z - x, z^2*x - y]).groebner_fan() sage: test_render = G.render(larger=True) # needs sage.plot
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> G = R.ideal([x**Integer(2)*y - z, y**Integer(2)*z - x, z**Integer(2)*x - y]).groebner_fan() >>> test_render = G.render(larger=True) # needs sage.plot
R.<x,y,z> = PolynomialRing(QQ,3) G = R.ideal([x^2*y - z, y^2*z - x, z^2*x - y]).groebner_fan() test_render = G.render(larger=True) # needs sage.plot
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> G = R.ideal([x**Integer(2)*y - z, y**Integer(2)*z - x, z**Integer(2)*x - y]).groebner_fan() >>> test_render = G.render(larger=True) # needs sage.plot
R.<x,y,z> = PolynomialRing(QQ,3) G = R.ideal([x^2*y - z, y^2*z - x, z^2*x - y]).groebner_fan() test_render = G.render(larger=True) # needs sage.plot
- render3d(verbose=False)[source]¶
For a Groebner fan of an ideal in a ring with four variables, this function intersects the fan with the standard simplex perpendicular to (1,1,1,1), creating a 3d polytope, which is then projected into 3 dimensions. The edges of this projected polytope are returned as lines.
EXAMPLES:
sage: R4.<w,x,y,z> = PolynomialRing(QQ,4) sage: gf = R4.ideal([w^2-x,x^2-y,y^2-z,z^2-x]).groebner_fan() sage: three_d = gf.render3d() # needs sage.plot
>>> from sage.all import * >>> R4 = PolynomialRing(QQ,Integer(4), names=('w', 'x', 'y', 'z',)); (w, x, y, z,) = R4._first_ngens(4) >>> gf = R4.ideal([w**Integer(2)-x,x**Integer(2)-y,y**Integer(2)-z,z**Integer(2)-x]).groebner_fan() >>> three_d = gf.render3d() # needs sage.plot
R4.<w,x,y,z> = PolynomialRing(QQ,4) gf = R4.ideal([w^2-x,x^2-y,y^2-z,z^2-x]).groebner_fan() three_d = gf.render3d() # needs sage.plot
- ring()[source]¶
Return the multivariate polynomial ring.
EXAMPLES:
sage: R.<x1,x2> = PolynomialRing(QQ,2) sage: gf = R.ideal([x1^3-x2,x2^3-x1-2]).groebner_fan() sage: gf.ring() Multivariate Polynomial Ring in x1, x2 over Rational Field
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x1', 'x2',)); (x1, x2,) = R._first_ngens(2) >>> gf = R.ideal([x1**Integer(3)-x2,x2**Integer(3)-x1-Integer(2)]).groebner_fan() >>> gf.ring() Multivariate Polynomial Ring in x1, x2 over Rational Field
R.<x1,x2> = PolynomialRing(QQ,2) gf = R.ideal([x1^3-x2,x2^3-x1-2]).groebner_fan() gf.ring()
- tropical_basis(check=True, verbose=False)[source]¶
Return a tropical basis for the tropical curve associated to this ideal.
INPUT:
check
– boolean (default:True
); ifTrue
raises aValueError
exception if this ideal does not define a tropical curve (i.e., the condition that R/I has dimension equal to 1 + the dimension of the homogeneity space is not satisfied)
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3, order='lex') sage: G = R.ideal([y^3-3*x^2, z^3-x-y-2*y^3+2*x^2]).groebner_fan() sage: G Groebner fan of the ideal: Ideal (-3*x^2 + y^3, 2*x^2 - x - 2*y^3 - y + z^3) of Multivariate Polynomial Ring in x, y, z over Rational Field sage: G.tropical_basis() [-3*x^2 + y^3, 2*x^2 - x - 2*y^3 - y + z^3, 3/4*x + y^3 + 3/4*y - 3/4*z^3]
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), order='lex', names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> G = R.ideal([y**Integer(3)-Integer(3)*x**Integer(2), z**Integer(3)-x-y-Integer(2)*y**Integer(3)+Integer(2)*x**Integer(2)]).groebner_fan() >>> G Groebner fan of the ideal: Ideal (-3*x^2 + y^3, 2*x^2 - x - 2*y^3 - y + z^3) of Multivariate Polynomial Ring in x, y, z over Rational Field >>> G.tropical_basis() [-3*x^2 + y^3, 2*x^2 - x - 2*y^3 - y + z^3, 3/4*x + y^3 + 3/4*y - 3/4*z^3]
R.<x,y,z> = PolynomialRing(QQ,3, order='lex') G = R.ideal([y^3-3*x^2, z^3-x-y-2*y^3+2*x^2]).groebner_fan() G G.tropical_basis()
- tropical_intersection(parameters=[], symmetry_generators=[], *args, **kwds)[source]¶
Return information about the tropical intersection of the polynomials defining the ideal.
This is the common refinement of the outward-pointing normal fans of the Newton polytopes of the generators of the ideal. Note that some people use the inward-pointing normal fans.
INPUT:
parameters
– (optional) list of variables to be considered as parameterssymmetry_generators
– (optional) generators of the symmetry group
OUTPUT: a TropicalPrevariety object
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: I = R.ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) sage: gf = I.groebner_fan() sage: pf = gf.tropical_intersection() sage: pf.rays() [[-2, 1, 1]] sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: f1 = x*y*z - 1 sage: f2 = f1*(x^2 + y^2 + z^2) sage: f3 = f2*(x + y + z - 1) sage: I = R.ideal([f1,f2,f3]) sage: gf = I.groebner_fan() sage: pf = gf.tropical_intersection(symmetry_generators = '(1,2,0),(1,0,2)') sage: pf.rays() [[-2, 1, 1], [1, -2, 1], [1, 1, -2]] sage: R.<x,y,z> = QQ[] sage: I = R.ideal([(x+y+z)^2-1,(x+y+z)-x,(x+y+z)-3]) sage: GF = I.groebner_fan() sage: TI = GF.tropical_intersection() sage: TI.rays() [[-1, 0, 0], [0, -1, -1], [1, 1, 1]] sage: GF = I.groebner_fan() sage: TI = GF.tropical_intersection(parameters=[y]) sage: TI.rays() [[-1, 0, 0]]
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> I = R.ideal(x*z + Integer(6)*y*z - z**Integer(2), x*y + Integer(6)*x*z + y*z - z**Integer(2), y**Integer(2) + x*z + y*z) >>> gf = I.groebner_fan() >>> pf = gf.tropical_intersection() >>> pf.rays() [[-2, 1, 1]] >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> f1 = x*y*z - Integer(1) >>> f2 = f1*(x**Integer(2) + y**Integer(2) + z**Integer(2)) >>> f3 = f2*(x + y + z - Integer(1)) >>> I = R.ideal([f1,f2,f3]) >>> gf = I.groebner_fan() >>> pf = gf.tropical_intersection(symmetry_generators = '(1,2,0),(1,0,2)') >>> pf.rays() [[-2, 1, 1], [1, -2, 1], [1, 1, -2]] >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> I = R.ideal([(x+y+z)**Integer(2)-Integer(1),(x+y+z)-x,(x+y+z)-Integer(3)]) >>> GF = I.groebner_fan() >>> TI = GF.tropical_intersection() >>> TI.rays() [[-1, 0, 0], [0, -1, -1], [1, 1, 1]] >>> GF = I.groebner_fan() >>> TI = GF.tropical_intersection(parameters=[y]) >>> TI.rays() [[-1, 0, 0]]
R.<x,y,z> = PolynomialRing(QQ,3) I = R.ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) gf = I.groebner_fan() pf = gf.tropical_intersection() pf.rays() R.<x,y,z> = PolynomialRing(QQ,3) f1 = x*y*z - 1 f2 = f1*(x^2 + y^2 + z^2) f3 = f2*(x + y + z - 1) I = R.ideal([f1,f2,f3]) gf = I.groebner_fan() pf = gf.tropical_intersection(symmetry_generators = '(1,2,0),(1,0,2)') pf.rays() R.<x,y,z> = QQ[] I = R.ideal([(x+y+z)^2-1,(x+y+z)-x,(x+y+z)-3]) GF = I.groebner_fan() TI = GF.tropical_intersection() TI.rays() GF = I.groebner_fan() TI = GF.tropical_intersection(parameters=[y]) TI.rays()
- weight_vectors()[source]¶
Return the weight vectors corresponding to the reduced Groebner bases.
EXAMPLES:
sage: r3.<x,y,z> = PolynomialRing(QQ,3) sage: g = r3.ideal([x^3+y,y^3-z,z^2-x]).groebner_fan() sage: g.weight_vectors() [(3, 7, 1), (5, 1, 2), (7, 1, 4), (5, 1, 4), (1, 1, 1), (1, 4, 8), (1, 4, 10)] sage: r4.<x,y,z,w> = PolynomialRing(QQ,4) sage: g4 = r4.ideal([x^3+y,y^3-z,z^2-x,z^3 - w]).groebner_fan() sage: len(g4.weight_vectors()) 23
>>> from sage.all import * >>> r3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = r3._first_ngens(3) >>> g = r3.ideal([x**Integer(3)+y,y**Integer(3)-z,z**Integer(2)-x]).groebner_fan() >>> g.weight_vectors() [(3, 7, 1), (5, 1, 2), (7, 1, 4), (5, 1, 4), (1, 1, 1), (1, 4, 8), (1, 4, 10)] >>> r4 = PolynomialRing(QQ,Integer(4), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = r4._first_ngens(4) >>> g4 = r4.ideal([x**Integer(3)+y,y**Integer(3)-z,z**Integer(2)-x,z**Integer(3) - w]).groebner_fan() >>> len(g4.weight_vectors()) 23
r3.<x,y,z> = PolynomialRing(QQ,3) g = r3.ideal([x^3+y,y^3-z,z^2-x]).groebner_fan() g.weight_vectors() r4.<x,y,z,w> = PolynomialRing(QQ,4) g4 = r4.ideal([x^3+y,y^3-z,z^2-x,z^3 - w]).groebner_fan() len(g4.weight_vectors())
- class sage.rings.polynomial.groebner_fan.InitialForm(cone, rays, initial_forms)[source]¶
Bases:
SageObject
A system of initial forms from a polynomial system.
To each form is associated a cone and a list of polynomials (the initial form system itself).
This class is intended for internal use inside of the
TropicalPrevariety
class.EXAMPLES:
sage: from sage.rings.polynomial.groebner_fan import InitialForm sage: R.<x,y> = QQ[] sage: inform = InitialForm([0], [[-1, 0]], [y^2 - 1, y^2 - 2, y^2 - 3]) sage: inform._cone [0]
>>> from sage.all import * >>> from sage.rings.polynomial.groebner_fan import InitialForm >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> inform = InitialForm([Integer(0)], [[-Integer(1), Integer(0)]], [y**Integer(2) - Integer(1), y**Integer(2) - Integer(2), y**Integer(2) - Integer(3)]) >>> inform._cone [0]
from sage.rings.polynomial.groebner_fan import InitialForm R.<x,y> = QQ[] inform = InitialForm([0], [[-1, 0]], [y^2 - 1, y^2 - 2, y^2 - 3]) inform._cone
- cone()[source]¶
The cone associated with the initial form system.
EXAMPLES:
sage: R.<x,y> = QQ[] sage: I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: pfi0 = PF.initial_form_systems()[0] sage: pfi0.cone() [0]
>>> from sage.all import * >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> I = R.ideal([(x+y)**Integer(2)-Integer(1),(x+y)**Integer(2)-Integer(2),(x+y)**Integer(2)-Integer(3)]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> pfi0 = PF.initial_form_systems()[Integer(0)] >>> pfi0.cone() [0]
R.<x,y> = QQ[] I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) GF = I.groebner_fan() PF = GF.tropical_intersection() pfi0 = PF.initial_form_systems()[0] pfi0.cone()
- initial_forms()[source]¶
The initial forms (polynomials).
EXAMPLES:
sage: R.<x,y> = QQ[] sage: I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: pfi0 = PF.initial_form_systems()[0] sage: pfi0.initial_forms() [y^2 - 1, y^2 - 2, y^2 - 3]
>>> from sage.all import * >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> I = R.ideal([(x+y)**Integer(2)-Integer(1),(x+y)**Integer(2)-Integer(2),(x+y)**Integer(2)-Integer(3)]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> pfi0 = PF.initial_form_systems()[Integer(0)] >>> pfi0.initial_forms() [y^2 - 1, y^2 - 2, y^2 - 3]
R.<x,y> = QQ[] I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) GF = I.groebner_fan() PF = GF.tropical_intersection() pfi0 = PF.initial_form_systems()[0] pfi0.initial_forms()
- internal_ray()[source]¶
A ray internal to the cone associated with the initial form system.
EXAMPLES:
sage: R.<x,y> = QQ[] sage: I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: pfi0 = PF.initial_form_systems()[0] sage: pfi0.internal_ray() (-1, 0)
>>> from sage.all import * >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> I = R.ideal([(x+y)**Integer(2)-Integer(1),(x+y)**Integer(2)-Integer(2),(x+y)**Integer(2)-Integer(3)]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> pfi0 = PF.initial_form_systems()[Integer(0)] >>> pfi0.internal_ray() (-1, 0)
R.<x,y> = QQ[] I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) GF = I.groebner_fan() PF = GF.tropical_intersection() pfi0 = PF.initial_form_systems()[0] pfi0.internal_ray()
- rays()[source]¶
The rays of the cone associated with the initial form system.
EXAMPLES:
sage: R.<x,y> = QQ[] sage: I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: pfi0 = PF.initial_form_systems()[0] sage: pfi0.rays() [[-1, 0]]
>>> from sage.all import * >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> I = R.ideal([(x+y)**Integer(2)-Integer(1),(x+y)**Integer(2)-Integer(2),(x+y)**Integer(2)-Integer(3)]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> pfi0 = PF.initial_form_systems()[Integer(0)] >>> pfi0.rays() [[-1, 0]]
R.<x,y> = QQ[] I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) GF = I.groebner_fan() PF = GF.tropical_intersection() pfi0 = PF.initial_form_systems()[0] pfi0.rays()
- class sage.rings.polynomial.groebner_fan.PolyhedralCone(gfan_polyhedral_cone, ring=Rational Field)[source]¶
Bases:
SageObject
Convert polymake/gfan data on a polyhedral cone into a sage class.
Currently (18-03-2008) needs a lot of work.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf[0].groebner_cone() sage: a.facets() [[0, 0, 1], [0, 1, 0], [1, 0, 0]]
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf[Integer(0)].groebner_cone() >>> a.facets() [[0, 0, 1], [0, 1, 0], [1, 0, 0]]
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf[0].groebner_cone() a.facets()
- ambient_dim()[source]¶
Return the ambient dimension of the Groebner cone.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf[0].groebner_cone() sage: a.ambient_dim() 3
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf[Integer(0)].groebner_cone() >>> a.ambient_dim() 3
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf[0].groebner_cone() a.ambient_dim()
- dim()[source]¶
Return the dimension of the Groebner cone.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf[0].groebner_cone() sage: a.dim() 3
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf[Integer(0)].groebner_cone() >>> a.dim() 3
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf[0].groebner_cone() a.dim()
- facets()[source]¶
Return the inward facet normals of the Groebner cone.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf[0].groebner_cone() sage: a.facets() [[0, 0, 1], [0, 1, 0], [1, 0, 0]]
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf[Integer(0)].groebner_cone() >>> a.facets() [[0, 0, 1], [0, 1, 0], [1, 0, 0]]
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf[0].groebner_cone() a.facets()
- lineality_dim()[source]¶
Return the lineality dimension of the Groebner cone. This is just the difference between the ambient dimension and the dimension of the cone.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf[0].groebner_cone() sage: a.lineality_dim() 0
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf[Integer(0)].groebner_cone() >>> a.lineality_dim() 0
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf[0].groebner_cone() a.lineality_dim()
- relative_interior_point()[source]¶
Return a point in the relative interior of the Groebner cone.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf[0].groebner_cone() sage: a.relative_interior_point() [1, 1, 1]
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf[Integer(0)].groebner_cone() >>> a.relative_interior_point() [1, 1, 1]
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf[0].groebner_cone() a.relative_interior_point()
- class sage.rings.polynomial.groebner_fan.PolyhedralFan(gfan_polyhedral_fan, parameter_indices=None)[source]¶
Bases:
SageObject
Convert polymake/gfan data on a polyhedral fan into a sage class.
INPUT:
gfan_polyhedral_fan
– output from gfan of a polyhedral fan
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: i2 = ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) sage: gf2 = i2.groebner_fan(verbose=False) sage: pf = gf2.polyhedralfan() sage: pf.rays() [[-1, 0, 1], [-1, 1, 0], [1, -2, 1], [1, 1, -2], [2, -1, -1]]
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> i2 = ideal(x*z + Integer(6)*y*z - z**Integer(2), x*y + Integer(6)*x*z + y*z - z**Integer(2), y**Integer(2) + x*z + y*z) >>> gf2 = i2.groebner_fan(verbose=False) >>> pf = gf2.polyhedralfan() >>> pf.rays() [[-1, 0, 1], [-1, 1, 0], [1, -2, 1], [1, 1, -2], [2, -1, -1]]
R.<x,y,z> = PolynomialRing(QQ,3) i2 = ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) gf2 = i2.groebner_fan(verbose=False) pf = gf2.polyhedralfan() pf.rays()
- ambient_dim()[source]¶
Return the ambient dimension of the Groebner fan.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf.polyhedralfan() sage: a.ambient_dim() 3
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf.polyhedralfan() >>> a.ambient_dim() 3
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf.polyhedralfan() a.ambient_dim()
- cones()[source]¶
A dictionary of cones in which the keys are the cone dimensions. For each dimension, the value is a list of the cones, where each element consists of a list of ray indices.
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: f = 1+x+y+x*y sage: I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: PF.cones() {1: [[0], [1], [2], [3], [4], [5]], 2: [[0, 1], [0, 2], [0, 3], [0, 4], [1, 2], [1, 3], [2, 4], [3, 4], [1, 5], [2, 5], [3, 5], [4, 5]]}
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> f = Integer(1)+x+y+x*y >>> I = R.ideal([f+z*f, Integer(2)*f+z*f, Integer(3)*f+z**Integer(2)*f]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> PF.cones() {1: [[0], [1], [2], [3], [4], [5]], 2: [[0, 1], [0, 2], [0, 3], [0, 4], [1, 2], [1, 3], [2, 4], [3, 4], [1, 5], [2, 5], [3, 5], [4, 5]]}
R.<x,y,z> = QQ[] f = 1+x+y+x*y I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) GF = I.groebner_fan() PF = GF.tropical_intersection() PF.cones()
- dim()[source]¶
Return the dimension of the Groebner fan.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf.polyhedralfan() sage: a.dim() 3
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf.polyhedralfan() >>> a.dim() 3
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf.polyhedralfan() a.dim()
- f_vector()[source]¶
The f-vector of the fan.
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: f = 1+x+y+x*y sage: I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: PF.f_vector() [1, 6, 12]
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> f = Integer(1)+x+y+x*y >>> I = R.ideal([f+z*f, Integer(2)*f+z*f, Integer(3)*f+z**Integer(2)*f]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> PF.f_vector() [1, 6, 12]
R.<x,y,z> = QQ[] f = 1+x+y+x*y I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) GF = I.groebner_fan() PF = GF.tropical_intersection() PF.f_vector()
- is_simplicial()[source]¶
Whether the fan is simplicial or not.
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: f = 1+x+y+x*y sage: I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: PF.is_simplicial() True
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> f = Integer(1)+x+y+x*y >>> I = R.ideal([f+z*f, Integer(2)*f+z*f, Integer(3)*f+z**Integer(2)*f]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> PF.is_simplicial() True
R.<x,y,z> = QQ[] f = 1+x+y+x*y I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) GF = I.groebner_fan() PF = GF.tropical_intersection() PF.is_simplicial()
- lineality_dim()[source]¶
Return the lineality dimension of the fan. This is the dimension of the largest subspace contained in the fan.
EXAMPLES:
sage: R3.<x,y,z> = PolynomialRing(QQ,3) sage: gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() sage: a = gf.polyhedralfan() sage: a.lineality_dim() 0
>>> from sage.all import * >>> R3 = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R3._first_ngens(3) >>> gf = R3.ideal([x**Integer(8)-y**Integer(4),y**Integer(4)-z**Integer(2),z**Integer(2)-Integer(2)]).groebner_fan() >>> a = gf.polyhedralfan() >>> a.lineality_dim() 0
R3.<x,y,z> = PolynomialRing(QQ,3) gf = R3.ideal([x^8-y^4,y^4-z^2,z^2-2]).groebner_fan() a = gf.polyhedralfan() a.lineality_dim()
- maximal_cones()[source]¶
A dictionary of the maximal cones in which the keys are the cone dimensions. For each dimension, the value is a list of the maximal cones, where each element consists of a list of ray indices.
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: f = 1+x+y+x*y sage: I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: PF.maximal_cones() {2: [[0, 1], [0, 2], [0, 3], [0, 4], [1, 2], [1, 3], [2, 4], [3, 4], [1, 5], [2, 5], [3, 5], [4, 5]]}
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> f = Integer(1)+x+y+x*y >>> I = R.ideal([f+z*f, Integer(2)*f+z*f, Integer(3)*f+z**Integer(2)*f]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> PF.maximal_cones() {2: [[0, 1], [0, 2], [0, 3], [0, 4], [1, 2], [1, 3], [2, 4], [3, 4], [1, 5], [2, 5], [3, 5], [4, 5]]}
R.<x,y,z> = QQ[] f = 1+x+y+x*y I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) GF = I.groebner_fan() PF = GF.tropical_intersection() PF.maximal_cones()
- rays()[source]¶
A list of rays of the polyhedral fan.
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: i2 = ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) sage: gf2 = i2.groebner_fan(verbose=False) sage: pf = gf2.polyhedralfan() sage: pf.rays() [[-1, 0, 1], [-1, 1, 0], [1, -2, 1], [1, 1, -2], [2, -1, -1]]
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> i2 = ideal(x*z + Integer(6)*y*z - z**Integer(2), x*y + Integer(6)*x*z + y*z - z**Integer(2), y**Integer(2) + x*z + y*z) >>> gf2 = i2.groebner_fan(verbose=False) >>> pf = gf2.polyhedralfan() >>> pf.rays() [[-1, 0, 1], [-1, 1, 0], [1, -2, 1], [1, 1, -2], [2, -1, -1]]
R.<x,y,z> = PolynomialRing(QQ,3) i2 = ideal(x*z + 6*y*z - z^2, x*y + 6*x*z + y*z - z^2, y^2 + x*z + y*z) gf2 = i2.groebner_fan(verbose=False) pf = gf2.polyhedralfan() pf.rays()
- to_RationalPolyhedralFan()[source]¶
Convert to the RationalPolyhedralFan class, which is more actively maintained. While the information in each class is essentially the same, the methods and implementation are different.
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: f = 1+x+y+x*y sage: I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: fan = PF.to_RationalPolyhedralFan() sage: [tuple(q.facet_normals()) for q in fan] [(M(0, -1, 0), M(-1, 0, 0)), (M(0, 0, -1), M(-1, 0, 0)), (M(0, 0, 1), M(-1, 0, 0)), (M(0, 1, 0), M(-1, 0, 0)), (M(0, 0, -1), M(0, -1, 0)), (M(0, 0, 1), M(0, -1, 0)), (M(0, 1, 0), M(0, 0, -1)), (M(0, 1, 0), M(0, 0, 1)), (M(1, 0, 0), M(0, -1, 0)), (M(1, 0, 0), M(0, 0, -1)), (M(1, 0, 0), M(0, 0, 1)), (M(1, 0, 0), M(0, 1, 0))]
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> f = Integer(1)+x+y+x*y >>> I = R.ideal([f+z*f, Integer(2)*f+z*f, Integer(3)*f+z**Integer(2)*f]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> fan = PF.to_RationalPolyhedralFan() >>> [tuple(q.facet_normals()) for q in fan] [(M(0, -1, 0), M(-1, 0, 0)), (M(0, 0, -1), M(-1, 0, 0)), (M(0, 0, 1), M(-1, 0, 0)), (M(0, 1, 0), M(-1, 0, 0)), (M(0, 0, -1), M(0, -1, 0)), (M(0, 0, 1), M(0, -1, 0)), (M(0, 1, 0), M(0, 0, -1)), (M(0, 1, 0), M(0, 0, 1)), (M(1, 0, 0), M(0, -1, 0)), (M(1, 0, 0), M(0, 0, -1)), (M(1, 0, 0), M(0, 0, 1)), (M(1, 0, 0), M(0, 1, 0))]
R.<x,y,z> = QQ[] f = 1+x+y+x*y I = R.ideal([f+z*f, 2*f+z*f, 3*f+z^2*f]) GF = I.groebner_fan() PF = GF.tropical_intersection() fan = PF.to_RationalPolyhedralFan() [tuple(q.facet_normals()) for q in fan]
Here we use the RationalPolyhedralFan’s Gale_transform method on a tropical prevariety.
sage: fan.Gale_transform() [ 1 0 0 0 0 1 -2] [ 0 1 0 0 1 0 -2] [ 0 0 1 1 0 0 -2]
>>> from sage.all import * >>> fan.Gale_transform() [ 1 0 0 0 0 1 -2] [ 0 1 0 0 1 0 -2] [ 0 0 1 1 0 0 -2]
fan.Gale_transform()
- class sage.rings.polynomial.groebner_fan.ReducedGroebnerBasis(groebner_fan, gens, gfan_gens)[source]¶
Bases:
SageObject
,list
A class for representing reduced Groebner bases as produced by
gfan
.INPUT:
groebner_fan
– a GroebnerFan object from an idealgens
– the generators of the idealgfan_gens
– the generators as a gfan string
EXAMPLES:
sage: R.<a,b> = PolynomialRing(QQ,2) sage: gf = R.ideal([a^2-b^2,b-a-1]).groebner_fan() sage: from sage.rings.polynomial.groebner_fan import ReducedGroebnerBasis sage: ReducedGroebnerBasis(gf,gf[0],gf[0]._gfan_gens()) [b - 1/2, a + 1/2]
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('a', 'b',)); (a, b,) = R._first_ngens(2) >>> gf = R.ideal([a**Integer(2)-b**Integer(2),b-a-Integer(1)]).groebner_fan() >>> from sage.rings.polynomial.groebner_fan import ReducedGroebnerBasis >>> ReducedGroebnerBasis(gf,gf[Integer(0)],gf[Integer(0)]._gfan_gens()) [b - 1/2, a + 1/2]
R.<a,b> = PolynomialRing(QQ,2) gf = R.ideal([a^2-b^2,b-a-1]).groebner_fan() from sage.rings.polynomial.groebner_fan import ReducedGroebnerBasis ReducedGroebnerBasis(gf,gf[0],gf[0]._gfan_gens())
- groebner_cone(restrict=False)[source]¶
Return defining inequalities for the full-dimensional Groebner cone associated to this marked minimal reduced Groebner basis.
INPUT:
restrict
– boolean (default:False
); ifTrue
, add an inequality for each coordinate, so that the cone is restricted to the positive orthant
OUTPUT: tuple of integer vectors
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: poly_cone = G[1].groebner_cone() sage: poly_cone.facets() [[-1, 2], [1, -1]] sage: [g.groebner_cone().facets() for g in G] [[[0, 1], [1, -2]], [[-1, 2], [1, -1]], [[-1, 1], [1, 0]]] sage: G[1].groebner_cone(restrict=True).facets() [[-1, 2], [1, -1]]
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> poly_cone = G[Integer(1)].groebner_cone() >>> poly_cone.facets() [[-1, 2], [1, -1]] >>> [g.groebner_cone().facets() for g in G] [[[0, 1], [1, -2]], [[-1, 2], [1, -1]], [[-1, 1], [1, 0]]] >>> G[Integer(1)].groebner_cone(restrict=True).facets() [[-1, 2], [1, -1]]
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() poly_cone = G[1].groebner_cone() poly_cone.facets() [g.groebner_cone().facets() for g in G] G[1].groebner_cone(restrict=True).facets()
- ideal()[source]¶
Return the ideal generated by this basis.
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: G = R.ideal([x - z^3, y^2 - 13*x]).groebner_fan() sage: G[0].ideal() Ideal (-13*z^3 + y^2, -z^3 + x) of Multivariate Polynomial Ring in x, y, z over Rational Field
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> G = R.ideal([x - z**Integer(3), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> G[Integer(0)].ideal() Ideal (-13*z^3 + y^2, -z^3 + x) of Multivariate Polynomial Ring in x, y, z over Rational Field
R.<x,y,z> = PolynomialRing(QQ,3) G = R.ideal([x - z^3, y^2 - 13*x]).groebner_fan() G[0].ideal()
- interactive(latex=False, flippable=False, wall=False, inequalities=False, weight=False)[source]¶
Do an interactive walk of the Groebner fan starting at this reduced Groebner basis.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ,2) sage: G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() sage: G[0].interactive() # not tested Initializing gfan interactive mode ********************************************* * Press control-C to return to Sage * ********************************************* ....
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> G = R.ideal([y**Integer(3) - x**Integer(2), y**Integer(2) - Integer(13)*x]).groebner_fan() >>> G[Integer(0)].interactive() # not tested Initializing gfan interactive mode ********************************************* * Press control-C to return to Sage * ********************************************* ....
R.<x,y> = PolynomialRing(QQ,2) G = R.ideal([y^3 - x^2, y^2 - 13*x]).groebner_fan() G[0].interactive() # not tested
- class sage.rings.polynomial.groebner_fan.TropicalPrevariety(gfan_polyhedral_fan, polynomial_system, poly_ring, parameters=None)[source]¶
Bases:
PolyhedralFan
This class is a subclass of the PolyhedralFan class, with some additional methods for tropical prevarieties.
INPUT:
gfan_polyhedral_fan
– output fromgfan
of a polyhedral fanpolynomial_system
– list of polynomialspoly_ring
– the polynomial ring of the list of polynomialsparameters
– (optional) list of variables to be considered as parameters
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: I = R.ideal([(x+y+z)^2-1,(x+y+z)-x,(x+y+z)-3]) sage: GF = I.groebner_fan() sage: TI = GF.tropical_intersection() sage: TI._polynomial_system [x^2 + 2*x*y + y^2 + 2*x*z + 2*y*z + z^2 - 1, y + z, x + y + z - 3]
>>> from sage.all import * >>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3) >>> I = R.ideal([(x+y+z)**Integer(2)-Integer(1),(x+y+z)-x,(x+y+z)-Integer(3)]) >>> GF = I.groebner_fan() >>> TI = GF.tropical_intersection() >>> TI._polynomial_system [x^2 + 2*x*y + y^2 + 2*x*z + 2*y*z + z^2 - 1, y + z, x + y + z - 3]
R.<x,y,z> = QQ[] I = R.ideal([(x+y+z)^2-1,(x+y+z)-x,(x+y+z)-3]) GF = I.groebner_fan() TI = GF.tropical_intersection() TI._polynomial_system
- initial_form_systems()[source]¶
Return a list of systems of initial forms for each cone in the tropical prevariety.
EXAMPLES:
sage: R.<x,y> = QQ[] sage: I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) sage: GF = I.groebner_fan() sage: PF = GF.tropical_intersection() sage: pfi = PF.initial_form_systems() sage: for q in pfi: ....: print(q.initial_forms()) [y^2 - 1, y^2 - 2, y^2 - 3] [x^2 - 1, x^2 - 2, x^2 - 3] [x^2 + 2*x*y + y^2, x^2 + 2*x*y + y^2, x^2 + 2*x*y + y^2]
>>> from sage.all import * >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> I = R.ideal([(x+y)**Integer(2)-Integer(1),(x+y)**Integer(2)-Integer(2),(x+y)**Integer(2)-Integer(3)]) >>> GF = I.groebner_fan() >>> PF = GF.tropical_intersection() >>> pfi = PF.initial_form_systems() >>> for q in pfi: ... print(q.initial_forms()) [y^2 - 1, y^2 - 2, y^2 - 3] [x^2 - 1, x^2 - 2, x^2 - 3] [x^2 + 2*x*y + y^2, x^2 + 2*x*y + y^2, x^2 + 2*x*y + y^2]
R.<x,y> = QQ[] I = R.ideal([(x+y)^2-1,(x+y)^2-2,(x+y)^2-3]) GF = I.groebner_fan() PF = GF.tropical_intersection() pfi = PF.initial_form_systems() for q in pfi: print(q.initial_forms())
- sage.rings.polynomial.groebner_fan.ideal_to_gfan_format(input_ring, polys)[source]¶
Return the ideal in gfan’s notation.
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: polys = [x^2*y - z, y^2*z - x, z^2*x - y] sage: from sage.rings.polynomial.groebner_fan import ideal_to_gfan_format sage: ideal_to_gfan_format(R, polys) 'Q[x, y, z]{x^2*y-z,y^2*z-x,x*z^2-y}'
>>> from sage.all import * >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> polys = [x**Integer(2)*y - z, y**Integer(2)*z - x, z**Integer(2)*x - y] >>> from sage.rings.polynomial.groebner_fan import ideal_to_gfan_format >>> ideal_to_gfan_format(R, polys) 'Q[x, y, z]{x^2*y-z,y^2*z-x,x*z^2-y}'
R.<x,y,z> = PolynomialRing(QQ,3) polys = [x^2*y - z, y^2*z - x, z^2*x - y] from sage.rings.polynomial.groebner_fan import ideal_to_gfan_format ideal_to_gfan_format(R, polys)
- sage.rings.polynomial.groebner_fan.max_degree(list_of_polys)[source]¶
Compute the maximum degree of a list of polynomials.
EXAMPLES:
sage: from sage.rings.polynomial.groebner_fan import max_degree sage: R.<x,y> = PolynomialRing(QQ,2) sage: p_list = [x^2-y,x*y^10-x] sage: max_degree(p_list) 11.0
>>> from sage.all import * >>> from sage.rings.polynomial.groebner_fan import max_degree >>> R = PolynomialRing(QQ,Integer(2), names=('x', 'y',)); (x, y,) = R._first_ngens(2) >>> p_list = [x**Integer(2)-y,x*y**Integer(10)-x] >>> max_degree(p_list) 11.0
from sage.rings.polynomial.groebner_fan import max_degree R.<x,y> = PolynomialRing(QQ,2) p_list = [x^2-y,x*y^10-x] max_degree(p_list)
- sage.rings.polynomial.groebner_fan.prefix_check(str_list)[source]¶
Check if any strings in a list are prefixes of another string in the list.
EXAMPLES:
sage: from sage.rings.polynomial.groebner_fan import prefix_check sage: prefix_check(['z1','z1z1']) False sage: prefix_check(['z1','zz1']) True
>>> from sage.all import * >>> from sage.rings.polynomial.groebner_fan import prefix_check >>> prefix_check(['z1','z1z1']) False >>> prefix_check(['z1','zz1']) True
from sage.rings.polynomial.groebner_fan import prefix_check prefix_check(['z1','z1z1']) prefix_check(['z1','zz1'])
- sage.rings.polynomial.groebner_fan.ring_to_gfan_format(input_ring)[source]¶
Convert a ring to gfan’s format.
EXAMPLES:
sage: R.<w,x,y,z> = QQ[] sage: from sage.rings.polynomial.groebner_fan import ring_to_gfan_format sage: ring_to_gfan_format(R) 'Q[w, x, y, z]' sage: R2.<x,y> = GF(2)[] sage: ring_to_gfan_format(R2) 'Z/2Z[x, y]'
>>> from sage.all import * >>> R = QQ['w, x, y, z']; (w, x, y, z,) = R._first_ngens(4) >>> from sage.rings.polynomial.groebner_fan import ring_to_gfan_format >>> ring_to_gfan_format(R) 'Q[w, x, y, z]' >>> R2 = GF(Integer(2))['x, y']; (x, y,) = R2._first_ngens(2) >>> ring_to_gfan_format(R2) 'Z/2Z[x, y]'
R.<w,x,y,z> = QQ[] from sage.rings.polynomial.groebner_fan import ring_to_gfan_format ring_to_gfan_format(R) R2.<x,y> = GF(2)[] ring_to_gfan_format(R2)
- sage.rings.polynomial.groebner_fan.verts_for_normal(normal, poly)[source]¶
Return the exponents of the vertices of a Newton polytope that make up the supporting hyperplane for the given outward normal.
EXAMPLES:
sage: from sage.rings.polynomial.groebner_fan import verts_for_normal sage: R.<x,y,z> = PolynomialRing(QQ,3) sage: f1 = x*y*z - 1 sage: f2 = f1*(x^2 + y^2 + 1) sage: verts_for_normal([1,1,1],f2) [(3, 1, 1), (1, 3, 1)]
>>> from sage.all import * >>> from sage.rings.polynomial.groebner_fan import verts_for_normal >>> R = PolynomialRing(QQ,Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3) >>> f1 = x*y*z - Integer(1) >>> f2 = f1*(x**Integer(2) + y**Integer(2) + Integer(1)) >>> verts_for_normal([Integer(1),Integer(1),Integer(1)],f2) [(3, 1, 1), (1, 3, 1)]
from sage.rings.polynomial.groebner_fan import verts_for_normal R.<x,y,z> = PolynomialRing(QQ,3) f1 = x*y*z - 1 f2 = f1*(x^2 + y^2 + 1) verts_for_normal([1,1,1],f2)