Zariski-Van Kampen method implementation¶
This file contains functions to compute the fundamental group of
the complement of a curve in the complex affine or projective plane,
using Zariski-Van Kampen approach. It depends on the package sirocco
.
The current implementation allows to compute a presentation of the fundamental group of curves over the rationals or number fields with a fixed embedding on \(\QQbar\).
Instead of computing a representation of the braid monodromy, we choose several base points and a system of paths joining them that generate all the necessary loops around the points of the discriminant. The group is generated by the free groups over these points, and braids over these paths give relations between these generators. This big group presentation is simplified at the end.
AUTHORS:
Miguel Marco (2015-09-30): Initial version
EXAMPLES:
sage: # needs sirocco
sage: from sage.schemes.curves.zariski_vankampen import fundamental_group, braid_monodromy
sage: R.<x, y> = QQ[]
sage: f = y^3 + x^3 - 1
sage: braid_monodromy(f)
([s1*s0, s1*s0, s1*s0], {0: 0, 1: 0, 2: 0}, {}, 3)
sage: fundamental_group(f)
Finitely presented group < x0 | >
>>> from sage.all import *
>>> # needs sirocco
>>> from sage.schemes.curves.zariski_vankampen import fundamental_group, braid_monodromy
>>> R = QQ['x, y']; (x, y,) = R._first_ngens(2)
>>> f = y**Integer(3) + x**Integer(3) - Integer(1)
>>> braid_monodromy(f)
([s1*s0, s1*s0, s1*s0], {0: 0, 1: 0, 2: 0}, {}, 3)
>>> fundamental_group(f)
Finitely presented group < x0 | >
# needs sirocco from sage.schemes.curves.zariski_vankampen import fundamental_group, braid_monodromy R.<x, y> = QQ[] f = y^3 + x^3 - 1 braid_monodromy(f) fundamental_group(f)
- sage.schemes.curves.zariski_vankampen.braid2rels(L)[source]¶
Return a minimal set of relations of the group
F / [(b * F([j])) / F([j]) for j in (1..d)]
whereF = FreeGroup(d)
andb
is a conjugate of a positive braid . One starts from the non-trivial relations determined by the positive braid and transform them in relations determined byb
.INPUT:
L
– tuple whose first element is a positive braid and the second element is a list of permutation braids
OUTPUT:
A list of Tietze words for a minimal set of relations of
F / [(g * b) / g for g in F.gens()]
.EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import braid2rels sage: B.<s0, s1, s2> = BraidGroup(4) sage: L = ((s1*s0)^2, [s2]) sage: braid2rels(L) [(4, 1, -2, -1), (2, -4, -2, 1)]
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import braid2rels >>> B = BraidGroup(Integer(4), names=('s0', 's1', 's2',)); (s0, s1, s2,) = B._first_ngens(3) >>> L = ((s1*s0)**Integer(2), [s2]) >>> braid2rels(L) [(4, 1, -2, -1), (2, -4, -2, 1)]
from sage.schemes.curves.zariski_vankampen import braid2rels B.<s0, s1, s2> = BraidGroup(4) L = ((s1*s0)^2, [s2]) braid2rels(L)
- sage.schemes.curves.zariski_vankampen.braid_from_piecewise(strands)[source]¶
Compute the braid corresponding to the piecewise linear curves strands.
INPUT:
strands
– list of lists of tuples(t, c1, c2)
, wheret
is a number between 0 and 1, andc1
andc2
are rationals or algebraic reals
OUTPUT: the braid formed by the piecewise linear strands
EXAMPLES:
sage: # needs sirocco sage: from sage.schemes.curves.zariski_vankampen import braid_from_piecewise sage: paths = [[(0, 0, 1), (0.2, -1, -0.5), (0.8, -1, 0), (1, 0, -1)], ....: [(0, -1, 0), (0.5, 0, -1), (1, 1, 0)], ....: [(0, 1, 0), (0.5, 1, 1), (1, 0, 1)]] sage: braid_from_piecewise(paths) s0*s1
>>> from sage.all import * >>> # needs sirocco >>> from sage.schemes.curves.zariski_vankampen import braid_from_piecewise >>> paths = [[(Integer(0), Integer(0), Integer(1)), (RealNumber('0.2'), -Integer(1), -RealNumber('0.5')), (RealNumber('0.8'), -Integer(1), Integer(0)), (Integer(1), Integer(0), -Integer(1))], ... [(Integer(0), -Integer(1), Integer(0)), (RealNumber('0.5'), Integer(0), -Integer(1)), (Integer(1), Integer(1), Integer(0))], ... [(Integer(0), Integer(1), Integer(0)), (RealNumber('0.5'), Integer(1), Integer(1)), (Integer(1), Integer(0), Integer(1))]] >>> braid_from_piecewise(paths) s0*s1
# needs sirocco from sage.schemes.curves.zariski_vankampen import braid_from_piecewise paths = [[(0, 0, 1), (0.2, -1, -0.5), (0.8, -1, 0), (1, 0, -1)], [(0, -1, 0), (0.5, 0, -1), (1, 1, 0)], [(0, 1, 0), (0.5, 1, 1), (1, 0, 1)]] braid_from_piecewise(paths)
- sage.schemes.curves.zariski_vankampen.braid_monodromy(f, arrangement=(), vertical=False)[source]¶
Compute the braid monodromy of a projection of the curve defined by a polynomial.
INPUT:
f
– a polynomial with two variables, over a number field with an embedding in the complex numbersarrangement
– tuple (default:()
); an optional tuple of polynomials whose product equalsf
vertical
– boolean (default:False
); if set toTrue
,arrangements
contains more than one polynomial, some of them are of degree \(1\) in \(x\) and degree \(0\) in \(y\), and none of the other components have vertical asymptotes, then these components are marked as vertical and not used for the computation of the braid monodromy. The other ones are marked as horizontal. If a vertical component does not pass through a singular points of the projection of the horizontal components a trivial braid is added to the list.
OUTPUT:
A list of braids, images by the braid monodromy of a geometric basis of the complement of the discriminant of \(f\) in \(\CC\).
A dictionary:
i
, index of a strand is sent to the index of the corresponding factor inarrangement
.Another dictionary
dv
, only relevant ifvertical
isTrue
. Ifj
is the index of a braid corresponding to a vertical line with indexi
inarrangement
, thendv[j] = i
.A nonnegative integer: the number of strands of the braids, only necessary if the list of braids is empty.
Note
The projection over the \(x\) axis is used if there are no vertical asymptotes. Otherwise, a linear change of variables is done to fall into the previous case except if the only vertical asymptotes are lines and
vertical=True
.EXAMPLES:
sage: # needs sirocco sage: from sage.schemes.curves.zariski_vankampen import braid_monodromy sage: R.<x, y> = QQ[] sage: f = (x^2 - y^3) * (x + 3*y - 5) sage: bm = braid_monodromy(f); bm ([s1*s0*(s1*s2)^2*s0*s2^2*s0^-1*(s2^-1*s1^-1)^2*s0^-1*s1^-1, s1*s0*(s1*s2)^2*(s0*s2^-1*s1*s2*s1*s2^-1)^2*(s2^-1*s1^-1)^2*s0^-1*s1^-1, s1*s0*(s1*s2)^2*s2*s1^-1*s2^-1*s1^-1*s0^-1*s1^-1, s1*s0*s2*s0^-1*s2*s1^-1], {0: 0, 1: 0, 2: 0, 3: 0}, {}, 4) sage: flist = (x^2 - y^3, x + 3*y - 5) sage: bm1 = braid_monodromy(f, arrangement=flist) sage: bm1[0] == bm[0] True sage: bm1[1] {0: 0, 1: 1, 2: 0, 3: 0} sage: braid_monodromy(R(1)) ([], {}, {}, 0) sage: braid_monodromy(x*y^2 - 1) ([s0*s1*s0^-1*s1*s0*s1^-1*s0^-1, s0*s1*s0^-1, s0], {0: 0, 1: 0, 2: 0}, {}, 3) sage: L = [x, y, x - 1, x -y] sage: braid_monodromy(prod(L), arrangement=L, vertical=True) ([s^2, 1], {0: 1, 1: 3}, {0: 0, 1: 2}, 2)
>>> from sage.all import * >>> # needs sirocco >>> from sage.schemes.curves.zariski_vankampen import braid_monodromy >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> f = (x**Integer(2) - y**Integer(3)) * (x + Integer(3)*y - Integer(5)) >>> bm = braid_monodromy(f); bm ([s1*s0*(s1*s2)^2*s0*s2^2*s0^-1*(s2^-1*s1^-1)^2*s0^-1*s1^-1, s1*s0*(s1*s2)^2*(s0*s2^-1*s1*s2*s1*s2^-1)^2*(s2^-1*s1^-1)^2*s0^-1*s1^-1, s1*s0*(s1*s2)^2*s2*s1^-1*s2^-1*s1^-1*s0^-1*s1^-1, s1*s0*s2*s0^-1*s2*s1^-1], {0: 0, 1: 0, 2: 0, 3: 0}, {}, 4) >>> flist = (x**Integer(2) - y**Integer(3), x + Integer(3)*y - Integer(5)) >>> bm1 = braid_monodromy(f, arrangement=flist) >>> bm1[Integer(0)] == bm[Integer(0)] True >>> bm1[Integer(1)] {0: 0, 1: 1, 2: 0, 3: 0} >>> braid_monodromy(R(Integer(1))) ([], {}, {}, 0) >>> braid_monodromy(x*y**Integer(2) - Integer(1)) ([s0*s1*s0^-1*s1*s0*s1^-1*s0^-1, s0*s1*s0^-1, s0], {0: 0, 1: 0, 2: 0}, {}, 3) >>> L = [x, y, x - Integer(1), x -y] >>> braid_monodromy(prod(L), arrangement=L, vertical=True) ([s^2, 1], {0: 1, 1: 3}, {0: 0, 1: 2}, 2)
# needs sirocco from sage.schemes.curves.zariski_vankampen import braid_monodromy R.<x, y> = QQ[] f = (x^2 - y^3) * (x + 3*y - 5) bm = braid_monodromy(f); bm flist = (x^2 - y^3, x + 3*y - 5) bm1 = braid_monodromy(f, arrangement=flist) bm1[0] == bm[0] bm1[1] braid_monodromy(R(1)) braid_monodromy(x*y^2 - 1) L = [x, y, x - 1, x -y] braid_monodromy(prod(L), arrangement=L, vertical=True)
- sage.schemes.curves.zariski_vankampen.conjugate_positive_form(braid)[source]¶
For a
braid
which is conjugate to a product of disjoint positive braids a list of such decompositions is given.INPUT:
braid
– a braid \(\sigma\)
OUTPUT:
A list of \(r\) lists. Each such list is another list with two elements, a positive braid \(\alpha_i\) and a list of permutation braids \(\gamma_{1}^{i},\dots,\gamma_{r}^{n_i}\) such that if \(\gamma_i=\prod_{j=1}^{n_i} \gamma_j^i\) then the braids \(\tau_i=\gamma_i\alpha_i\gamma_i^{-1}\) pairwise commute and \(\alpha=\prod_{i=1}^{r} \tau_i\).
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import conjugate_positive_form sage: B = BraidGroup(4) sage: t = B((1, 3, 2, -3, 1, 1)) sage: conjugate_positive_form(t) [[(s1*s0)^2, [s2]]] sage: B = BraidGroup(5) sage: t = B((1, 2, 3, 4, -1, -2, 3, 3, 2, -4)) sage: L = conjugate_positive_form(t); L [[s1^2, [s3*s2]], [s1*s2, [s0]]] sage: s = B.one() sage: for a, l in L: ....: b = prod(l) ....: s *= b * a / b sage: s == t True sage: s1 = B.gen(1)^3 sage: conjugate_positive_form(s1) [[s1^3, []]]
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import conjugate_positive_form >>> B = BraidGroup(Integer(4)) >>> t = B((Integer(1), Integer(3), Integer(2), -Integer(3), Integer(1), Integer(1))) >>> conjugate_positive_form(t) [[(s1*s0)^2, [s2]]] >>> B = BraidGroup(Integer(5)) >>> t = B((Integer(1), Integer(2), Integer(3), Integer(4), -Integer(1), -Integer(2), Integer(3), Integer(3), Integer(2), -Integer(4))) >>> L = conjugate_positive_form(t); L [[s1^2, [s3*s2]], [s1*s2, [s0]]] >>> s = B.one() >>> for a, l in L: ... b = prod(l) ... s *= b * a / b >>> s == t True >>> s1 = B.gen(Integer(1))**Integer(3) >>> conjugate_positive_form(s1) [[s1^3, []]]
from sage.schemes.curves.zariski_vankampen import conjugate_positive_form B = BraidGroup(4) t = B((1, 3, 2, -3, 1, 1)) conjugate_positive_form(t) B = BraidGroup(5) t = B((1, 2, 3, 4, -1, -2, 3, 3, 2, -4)) L = conjugate_positive_form(t); L s = B.one() for a, l in L: b = prod(l) s *= b * a / b s == t s1 = B.gen(1)^3 conjugate_positive_form(s1)
- sage.schemes.curves.zariski_vankampen.corrected_voronoi_diagram()[source]¶
Compute a Voronoi diagram of a set of points with rational coordinates. The given points are granted to lie one in each bounded region.
INPUT:
points
– tuple of complex numbers
OUTPUT:
A Voronoi diagram constructed from rational approximations of the points, with the guarantee that each bounded region contains exactly one of the input points.
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import corrected_voronoi_diagram sage: points = (2, I, 0.000001, 0, 0.000001*I) sage: V = corrected_voronoi_diagram(points) sage: V The Voronoi diagram of 9 points of dimension 2 in the Rational Field sage: V.regions() {P(-7, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices and 2 rays, P(0, -7): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices and 2 rays, P(0, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices, P(0, 1): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 5 vertices, P(0, 1/1000000): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices, P(0, 7): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 3 vertices and 2 rays, P(1/1000000, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 5 vertices, P(2, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 5 vertices, P(7, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 2 vertices and 2 rays}
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import corrected_voronoi_diagram >>> points = (Integer(2), I, RealNumber('0.000001'), Integer(0), RealNumber('0.000001')*I) >>> V = corrected_voronoi_diagram(points) >>> V The Voronoi diagram of 9 points of dimension 2 in the Rational Field >>> V.regions() {P(-7, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices and 2 rays, P(0, -7): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices and 2 rays, P(0, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices, P(0, 1): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 5 vertices, P(0, 1/1000000): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices, P(0, 7): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 3 vertices and 2 rays, P(1/1000000, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 5 vertices, P(2, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 5 vertices, P(7, 0): A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 2 vertices and 2 rays}
from sage.schemes.curves.zariski_vankampen import corrected_voronoi_diagram points = (2, I, 0.000001, 0, 0.000001*I) V = corrected_voronoi_diagram(points) V V.regions()
- sage.schemes.curves.zariski_vankampen.discrim(pols)[source]¶
Return the points in the discriminant of the product of the polynomials of a list or tuple
pols
.The result is the set of values of the first variable for which two roots in the second variable coincide.
INPUT:
pols
– list or tuple of polynomials in two variables with coefficients in a number field with a fixed embedding in \(\QQbar\)
OUTPUT: a tuple with the roots of the discriminant in \(\QQbar\)
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import discrim sage: R.<x, y> = QQ[] sage: flist = (y^3 + x^3 - 1, 2 * x + y) sage: sorted((discrim(flist))) [-0.522757958574711?, -0.500000000000000? - 0.866025403784439?*I, -0.500000000000000? + 0.866025403784439?*I, 0.2613789792873551? - 0.4527216721561923?*I, 0.2613789792873551? + 0.4527216721561923?*I, 1]
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import discrim >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> flist = (y**Integer(3) + x**Integer(3) - Integer(1), Integer(2) * x + y) >>> sorted((discrim(flist))) [-0.522757958574711?, -0.500000000000000? - 0.866025403784439?*I, -0.500000000000000? + 0.866025403784439?*I, 0.2613789792873551? - 0.4527216721561923?*I, 0.2613789792873551? + 0.4527216721561923?*I, 1]
from sage.schemes.curves.zariski_vankampen import discrim R.<x, y> = QQ[] flist = (y^3 + x^3 - 1, 2 * x + y) sorted((discrim(flist)))
- sage.schemes.curves.zariski_vankampen.fieldI(field)[source]¶
Return the (either double or trivial) extension of a number field which contains
I
.INPUT:
field
– a number field with an embedding in \(\QQbar\)
OUTPUT: the extension
F
offield
containingI
with an embedding in \(\QQbar\)EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import fieldI sage: p = QQ[x](x^5 + 2 * x + 1) sage: a0 = p.roots(QQbar, multiplicities=False)[0] sage: F0.<a> = NumberField(p, embedding=a0) sage: fieldI(F0) Number Field in prim with defining polynomial x^10 + 5*x^8 + 14*x^6 - 2*x^5 - 10*x^4 + 20*x^3 - 11*x^2 - 14*x + 10 with prim = 0.4863890359345430? + 1.000000000000000?*I sage: F0 = CyclotomicField(5) sage: fieldI(F0) Number Field in prim with defining polynomial x^8 - 2*x^7 + 7*x^6 - 10*x^5 + 16*x^4 - 10*x^3 - 2*x^2 + 4*x + 1 with prim = -0.3090169943749474? + 0.04894348370484643?*I sage: fieldI(QuadraticField(3)) Number Field in prim with defining polynomial x^4 - 4*x^2 + 16 with prim = -1.732050807568878? + 1.000000000000000?*I sage: fieldI(QuadraticField(-3)) Number Field in prim with defining polynomial x^4 + 8*x^2 + 4 with prim = 0.?e-18 - 0.732050807568878?*I
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import fieldI >>> p = QQ[x](x**Integer(5) + Integer(2) * x + Integer(1)) >>> a0 = p.roots(QQbar, multiplicities=False)[Integer(0)] >>> F0 = NumberField(p, embedding=a0, names=('a',)); (a,) = F0._first_ngens(1) >>> fieldI(F0) Number Field in prim with defining polynomial x^10 + 5*x^8 + 14*x^6 - 2*x^5 - 10*x^4 + 20*x^3 - 11*x^2 - 14*x + 10 with prim = 0.4863890359345430? + 1.000000000000000?*I >>> F0 = CyclotomicField(Integer(5)) >>> fieldI(F0) Number Field in prim with defining polynomial x^8 - 2*x^7 + 7*x^6 - 10*x^5 + 16*x^4 - 10*x^3 - 2*x^2 + 4*x + 1 with prim = -0.3090169943749474? + 0.04894348370484643?*I >>> fieldI(QuadraticField(Integer(3))) Number Field in prim with defining polynomial x^4 - 4*x^2 + 16 with prim = -1.732050807568878? + 1.000000000000000?*I >>> fieldI(QuadraticField(-Integer(3))) Number Field in prim with defining polynomial x^4 + 8*x^2 + 4 with prim = 0.?e-18 - 0.732050807568878?*I
from sage.schemes.curves.zariski_vankampen import fieldI p = QQ[x](x^5 + 2 * x + 1) a0 = p.roots(QQbar, multiplicities=False)[0] F0.<a> = NumberField(p, embedding=a0) fieldI(F0) F0 = CyclotomicField(5) fieldI(F0) fieldI(QuadraticField(3)) fieldI(QuadraticField(-3))
If
I
is already in the field, the result is the field itself:sage: from sage.schemes.curves.zariski_vankampen import fieldI sage: p = QQ[x](x^4 + 1) sage: a0 = p.roots(QQbar, multiplicities=False)[0] sage: F0.<a> = NumberField(p, embedding=a0) sage: F1 = fieldI(F0) sage: F0 == F1 True sage: QuadraticField(-1) == fieldI(QuadraticField(-1)) True
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import fieldI >>> p = QQ[x](x**Integer(4) + Integer(1)) >>> a0 = p.roots(QQbar, multiplicities=False)[Integer(0)] >>> F0 = NumberField(p, embedding=a0, names=('a',)); (a,) = F0._first_ngens(1) >>> F1 = fieldI(F0) >>> F0 == F1 True >>> QuadraticField(-Integer(1)) == fieldI(QuadraticField(-Integer(1))) True
from sage.schemes.curves.zariski_vankampen import fieldI p = QQ[x](x^4 + 1) a0 = p.roots(QQbar, multiplicities=False)[0] F0.<a> = NumberField(p, embedding=a0) F1 = fieldI(F0) F0 == F1 QuadraticField(-1) == fieldI(QuadraticField(-1))
- sage.schemes.curves.zariski_vankampen.followstrand(f, factors, x0, x1, y0a, prec=53)[source]¶
Return a piecewise linear approximation of the homotopy continuation of the root
y0a
fromx0
tox1
.INPUT:
f
– an irreducible polynomial in two variablesfactors
– list of irreducible polynomials in two variablesx0
– a complex value, where the homotopy startsx1
– a complex value, where the homotopy endsy0a
– an approximate solution of the polynomial \(F(y) = f(x_0, y)\)prec
– the precision to use
OUTPUT:
A list of values \((t, y_{tr}, y_{ti})\) such that:
t
is a real number between zero and one\(f(t \cdot x_1 + (1-t) \cdot x_0, y_{tr} + I \cdot y_{ti})\) is zero (or a good enough approximation)
the piecewise linear path determined by the points has a tubular neighborhood where the actual homotopy continuation path lies, and no other root of
f
, nor any root of the polynomials infactors
, intersects it.
EXAMPLES:
sage: # needs sirocco sage: from sage.schemes.curves.zariski_vankampen import followstrand sage: R.<x, y> = QQ[] sage: f = x^2 + y^3 sage: x0 = CC(1, 0) sage: x1 = CC(1, 0.5) sage: followstrand(f, [], x0, x1, -1.0) # abs tol 1e-15 [(0.0, -1.0, 0.0), (0.7500000000000001, -1.015090921153253, -0.24752813818386948), (1.0, -1.026166099551513, -0.32768940253604323)] sage: fup = f.subs({y: y - 1/10}) sage: fdown = f.subs({y: y + 1/10}) sage: followstrand(f, [fup, fdown], x0, x1, -1.0) # abs tol 1e-15 [(0.0, -1.0, 0.0), (0.5303300858899107, -1.0076747107983448, -0.17588022709184917), (0.7651655429449553, -1.015686131039112, -0.25243563967299404), (1.0, -1.026166099551513, -0.3276894025360433)]
>>> from sage.all import * >>> # needs sirocco >>> from sage.schemes.curves.zariski_vankampen import followstrand >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> f = x**Integer(2) + y**Integer(3) >>> x0 = CC(Integer(1), Integer(0)) >>> x1 = CC(Integer(1), RealNumber('0.5')) >>> followstrand(f, [], x0, x1, -RealNumber('1.0')) # abs tol 1e-15 [(0.0, -1.0, 0.0), (0.7500000000000001, -1.015090921153253, -0.24752813818386948), (1.0, -1.026166099551513, -0.32768940253604323)] >>> fup = f.subs({y: y - Integer(1)/Integer(10)}) >>> fdown = f.subs({y: y + Integer(1)/Integer(10)}) >>> followstrand(f, [fup, fdown], x0, x1, -RealNumber('1.0')) # abs tol 1e-15 [(0.0, -1.0, 0.0), (0.5303300858899107, -1.0076747107983448, -0.17588022709184917), (0.7651655429449553, -1.015686131039112, -0.25243563967299404), (1.0, -1.026166099551513, -0.3276894025360433)]
# needs sirocco from sage.schemes.curves.zariski_vankampen import followstrand R.<x, y> = QQ[] f = x^2 + y^3 x0 = CC(1, 0) x1 = CC(1, 0.5) followstrand(f, [], x0, x1, -1.0) # abs tol 1e-15 fup = f.subs({y: y - 1/10}) fdown = f.subs({y: y + 1/10}) followstrand(f, [fup, fdown], x0, x1, -1.0) # abs tol 1e-15
- sage.schemes.curves.zariski_vankampen.fundamental_group(f, simplified=True, projective=False, puiseux=True)[source]¶
Return a presentation of the fundamental group of the complement of the algebraic set defined by the polynomial
f
.INPUT:
f
– a polynomial in two variables, with coefficients in either the rationals or a number field with a fixed embedding in \(\QQbar\)simplified
– boolean (default:True
); if set toTrue
the presentation will be simplified (see below)projective
– boolean (default:False
); if set toTrue
, the fundamental group of the complement of the projective completion of the curve will be computed, otherwise, the fundamental group of the complement in the affine plane will be computedpuiseux
– boolean (default:True
); if set toTrue
, a presentation of the fundamental group with the homotopy type of the complement of the affine curve is computed. If the Euler characteristic does not match, the homotopy type is obtained with a wedge of 2-spheres. One relation is added ifprojective
is set toTrue
.
If
projective
isFalse
andpuiseux
isTrue
, a Zariski-VanKampen presentation is returned.OUTPUT:
A presentation of the fundamental group of the complement of the curve defined by
f
.EXAMPLES:
sage: # needs sirocco sage: from sage.schemes.curves.zariski_vankampen import fundamental_group, braid_monodromy sage: R.<x, y> = QQ[] sage: f = x^2 + y^3 sage: fundamental_group(f) Finitely presented group < x0, x1 | x0*x1^-1*x0^-1*x1^-1*x0*x1 > sage: fundamental_group(f, simplified=False, puiseux=False).sorted_presentation() Finitely presented group < x0, x1, x2 | x2^-1*x1^-1*x0*x1, x2^-1*x0*x1*x0^-1, x1^-1*x0^-1*x1^-1*x0*x1*x0 > sage: fundamental_group(f, projective=True) Finitely presented group < x0 | x0^3 > sage: fundamental_group(f).sorted_presentation() Finitely presented group < x0, x1 | x1^-1*x0^-1*x1^-1*x0*x1*x0 >
>>> from sage.all import * >>> # needs sirocco >>> from sage.schemes.curves.zariski_vankampen import fundamental_group, braid_monodromy >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> f = x**Integer(2) + y**Integer(3) >>> fundamental_group(f) Finitely presented group < x0, x1 | x0*x1^-1*x0^-1*x1^-1*x0*x1 > >>> fundamental_group(f, simplified=False, puiseux=False).sorted_presentation() Finitely presented group < x0, x1, x2 | x2^-1*x1^-1*x0*x1, x2^-1*x0*x1*x0^-1, x1^-1*x0^-1*x1^-1*x0*x1*x0 > >>> fundamental_group(f, projective=True) Finitely presented group < x0 | x0^3 > >>> fundamental_group(f).sorted_presentation() Finitely presented group < x0, x1 | x1^-1*x0^-1*x1^-1*x0*x1*x0 >
# needs sirocco from sage.schemes.curves.zariski_vankampen import fundamental_group, braid_monodromy R.<x, y> = QQ[] f = x^2 + y^3 fundamental_group(f) fundamental_group(f, simplified=False, puiseux=False).sorted_presentation() fundamental_group(f, projective=True) fundamental_group(f).sorted_presentation()
sage: # needs sirocco sage: from sage.schemes.curves.zariski_vankampen import fundamental_group sage: R.<x, y> = QQ[] sage: f = y^3 + x^3 sage: fundamental_group(f).sorted_presentation() Finitely presented group < x0, x1, x2 | x2^-1*x1^-1*x0^-1*x2*x0*x1, x2^-1*x1^-1*x2*x0*x1*x0^-1 >
>>> from sage.all import * >>> # needs sirocco >>> from sage.schemes.curves.zariski_vankampen import fundamental_group >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> f = y**Integer(3) + x**Integer(3) >>> fundamental_group(f).sorted_presentation() Finitely presented group < x0, x1, x2 | x2^-1*x1^-1*x0^-1*x2*x0*x1, x2^-1*x1^-1*x2*x0*x1*x0^-1 >
# needs sirocco from sage.schemes.curves.zariski_vankampen import fundamental_group R.<x, y> = QQ[] f = y^3 + x^3 fundamental_group(f).sorted_presentation()
It is also possible to have coefficients in a number field with a fixed embedding in \(\QQbar\):
sage: from sage.schemes.curves.zariski_vankampen import fundamental_group sage: zeta = QQbar['x']('x^2 + x+ 1').roots(multiplicities=False)[0] sage: zeta -0.50000000000000000? - 0.866025403784439?*I sage: F = NumberField(zeta.minpoly(), 'zeta', embedding=zeta) sage: F.inject_variables() Defining zeta sage: R.<x, y> = F[] sage: f = y^3 + x^3 + zeta * x + 1 sage: fundamental_group(f) # needs sirocco Finitely presented group < x0 | >
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import fundamental_group >>> zeta = QQbar['x']('x^2 + x+ 1').roots(multiplicities=False)[Integer(0)] >>> zeta -0.50000000000000000? - 0.866025403784439?*I >>> F = NumberField(zeta.minpoly(), 'zeta', embedding=zeta) >>> F.inject_variables() Defining zeta >>> R = F['x, y']; (x, y,) = R._first_ngens(2) >>> f = y**Integer(3) + x**Integer(3) + zeta * x + Integer(1) >>> fundamental_group(f) # needs sirocco Finitely presented group < x0 | >
from sage.schemes.curves.zariski_vankampen import fundamental_group zeta = QQbar['x']('x^2 + x+ 1').roots(multiplicities=False)[0] zeta F = NumberField(zeta.minpoly(), 'zeta', embedding=zeta) F.inject_variables() R.<x, y> = F[] f = y^3 + x^3 + zeta * x + 1 fundamental_group(f) # needs sirocco
We compute the fundamental group of the complement of a quartic using the
puiseux
option:sage: # optional - sirocco sage: from sage.schemes.curves.zariski_vankampen import fundamental_group sage: R.<x, y> = QQ[] sage: f = x^2 * y^2 + x^2 + y^2 - 2 * x * y * (x + y + 1) sage: g = fundamental_group(f); g.sorted_presentation() Finitely presented group < x0, x1 | x1^-2*x0^2, (x1^-1*x0)^3 > sage: g = fundamental_group(f, projective=True) sage: g.order(), g.abelian_invariants() (12, (4,)) sage: fundamental_group(y * (y - 1)) Finitely presented group < x0, x1 | >
>>> from sage.all import * >>> # optional - sirocco >>> from sage.schemes.curves.zariski_vankampen import fundamental_group >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> f = x**Integer(2) * y**Integer(2) + x**Integer(2) + y**Integer(2) - Integer(2) * x * y * (x + y + Integer(1)) >>> g = fundamental_group(f); g.sorted_presentation() Finitely presented group < x0, x1 | x1^-2*x0^2, (x1^-1*x0)^3 > >>> g = fundamental_group(f, projective=True) >>> g.order(), g.abelian_invariants() (12, (4,)) >>> fundamental_group(y * (y - Integer(1))) Finitely presented group < x0, x1 | >
# optional - sirocco from sage.schemes.curves.zariski_vankampen import fundamental_group R.<x, y> = QQ[] f = x^2 * y^2 + x^2 + y^2 - 2 * x * y * (x + y + 1) g = fundamental_group(f); g.sorted_presentation() g = fundamental_group(f, projective=True) g.order(), g.abelian_invariants() fundamental_group(y * (y - 1))
- sage.schemes.curves.zariski_vankampen.fundamental_group_arrangement(flist, simplified=True, projective=False, puiseux=True, vertical=False, braid_data=None)[source]¶
Compute the fundamental group of the complement of a curve defined by a list of polynomials with the extra information about the correspondence of the generators and meridians of the elements of the list.
INPUT:
flist
– a tuple of polynomial with two variables, over a number field with an embedding in the complex numberssimplified
– boolean (default:True
); if set toTrue
the presentation will be simplified (see below)projective
– boolean (default:False
); if set toTrue
, the fundamental group of the complement of the projective completion of the curve will be computed, otherwise, the fundamental group of the complement in the affine plane will be computedpuiseux
– boolean (default:True
); if set toTrue
a presentation of the fundamental group with the homotopy type of the complement of the affine curve will be computed, adding one relation ifprojective
is set toTrue
.vertical
– boolean (default:False
); if set toTrue
, whenever no curve has vertical asymptotes the computation of braid monodromy is simpler if some lines are verticalbraid_data
– tuple (default:None
); if it is not the default it is the output offundamental_group_from_braid_mon
previously computed
OUTPUT:
A list of braids. The braids correspond to paths based in the same point; each of this paths is the conjugated of a loop around one of the points in the discriminant of the projection of
f
.A dictionary attaching to
j
a tuple a list of elements of the group which are meridians of the curve in positionj
. Ifprojective
isFalse
and the \(y\)-degree of the horizontal components coincide with the total degree, another key is added to give a meridian of the line at infinity.
EXAMPLES:
sage: # needs sirocco sage: from sage.schemes.curves.zariski_vankampen import braid_monodromy sage: from sage.schemes.curves.zariski_vankampen import fundamental_group_arrangement sage: R.<x, y> = QQ[] sage: flist = [x^2 - y^3, x + 3 * y - 5] sage: g, dic = fundamental_group_arrangement(flist) sage: g.sorted_presentation() Finitely presented group < x0, x1, x2 | x2^-1*x1^-1*x2*x1, x2^-1*x0^-1*x2^-1*x0*x2*x0, x1^-1*x0^-1*x1*x0 > sage: dic {0: [x0, x2], 1: [x1], 2: [x0^-1*x2^-1*x1^-1*x0^-1]} sage: g, dic = fundamental_group_arrangement(flist, simplified=False, puiseux=False) sage: g.sorted_presentation(), dic (Finitely presented group < x0, x1, x2, x3 | 1, 1, 1, 1, 1, 1, 1, x3^-1*x2^-1*x1^-1*x2*x3*x2^-1*x1*x2, x3^-1*x2^-1*x1^-1*x0^-1*x1*x2*x3*x2, x3^-1*x2^-1*x1^-1*x0^-1*x1*x2*x1^-1*x0*x1*x2, x3^-1*x2^-1*x1^-1*x2*x3*x2^-1*x1*x2, x3^-1*x1^-1*x0*x1, x1^-1*x0^-1*x1*x0, x1^-1*x0^-1*x1*x0, x1^-1*x0^-1*x1*x0, x1^-1*x0^-1*x1*x0 >, {0: [x0, x2, x3], 1: [x1], 2: [x3^-1*x2^-1*x1^-1*x0^-1]}) sage: fundamental_group_arrangement(flist, projective=True) (Finitely presented group < x | >, {0: [x], 1: [x^-3]}) sage: fundamental_group_arrangement([]) (Finitely presented group < | >, {}) sage: g, dic = fundamental_group_arrangement([x * y]) sage: g.sorted_presentation(), dic (Finitely presented group < x0, x1 | x1^-1*x0^-1*x1*x0 >, {0: [x0, x1], 1: [x1^-1*x0^-1]}) sage: fundamental_group_arrangement([y + x^2]) (Finitely presented group < x | >, {0: [x]}) sage: fundamental_group_arrangement([y^2 + x], projective=True) (Finitely presented group < x | x^2 >, {0: [x]}) sage: L = [x, y, x - 1, x -y] sage: G, dic =fundamental_group_arrangement(L) sage: G.sorted_presentation() Finitely presented group < x0, x1, x2, x3 | x3^-1*x2^-1*x3*x2, x3^-1*x1^-1*x0^-1*x1*x3*x0, x3^-1*x1^-1*x3*x0*x1*x0^-1, x2^-1*x0^-1*x2*x0 > sage: dic {0: [x1], 1: [x3], 2: [x2], 3: [x0], 4: [x3^-1*x2^-1*x1^-1*x0^-1]} sage: fundamental_group_arrangement(L, vertical=True) (Finitely presented group < x0, x1, x2, x3 | x3*x0*x3^-1*x0^-1, x3*x1*x3^-1*x1^-1, x1*x2*x0*x2^-1*x1^-1*x0^-1, x1*x2*x0*x1^-1*x0^-1*x2^-1 >, {0: [x2], 1: [x0], 2: [x3], 3: [x1], 4: [x3^-1*x2^-1*x1^-1*x0^-1]})
>>> from sage.all import * >>> # needs sirocco >>> from sage.schemes.curves.zariski_vankampen import braid_monodromy >>> from sage.schemes.curves.zariski_vankampen import fundamental_group_arrangement >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> flist = [x**Integer(2) - y**Integer(3), x + Integer(3) * y - Integer(5)] >>> g, dic = fundamental_group_arrangement(flist) >>> g.sorted_presentation() Finitely presented group < x0, x1, x2 | x2^-1*x1^-1*x2*x1, x2^-1*x0^-1*x2^-1*x0*x2*x0, x1^-1*x0^-1*x1*x0 > >>> dic {0: [x0, x2], 1: [x1], 2: [x0^-1*x2^-1*x1^-1*x0^-1]} >>> g, dic = fundamental_group_arrangement(flist, simplified=False, puiseux=False) >>> g.sorted_presentation(), dic (Finitely presented group < x0, x1, x2, x3 | 1, 1, 1, 1, 1, 1, 1, x3^-1*x2^-1*x1^-1*x2*x3*x2^-1*x1*x2, x3^-1*x2^-1*x1^-1*x0^-1*x1*x2*x3*x2, x3^-1*x2^-1*x1^-1*x0^-1*x1*x2*x1^-1*x0*x1*x2, x3^-1*x2^-1*x1^-1*x2*x3*x2^-1*x1*x2, x3^-1*x1^-1*x0*x1, x1^-1*x0^-1*x1*x0, x1^-1*x0^-1*x1*x0, x1^-1*x0^-1*x1*x0, x1^-1*x0^-1*x1*x0 >, {0: [x0, x2, x3], 1: [x1], 2: [x3^-1*x2^-1*x1^-1*x0^-1]}) >>> fundamental_group_arrangement(flist, projective=True) (Finitely presented group < x | >, {0: [x], 1: [x^-3]}) >>> fundamental_group_arrangement([]) (Finitely presented group < | >, {}) >>> g, dic = fundamental_group_arrangement([x * y]) >>> g.sorted_presentation(), dic (Finitely presented group < x0, x1 | x1^-1*x0^-1*x1*x0 >, {0: [x0, x1], 1: [x1^-1*x0^-1]}) >>> fundamental_group_arrangement([y + x**Integer(2)]) (Finitely presented group < x | >, {0: [x]}) >>> fundamental_group_arrangement([y**Integer(2) + x], projective=True) (Finitely presented group < x | x^2 >, {0: [x]}) >>> L = [x, y, x - Integer(1), x -y] >>> G, dic =fundamental_group_arrangement(L) >>> G.sorted_presentation() Finitely presented group < x0, x1, x2, x3 | x3^-1*x2^-1*x3*x2, x3^-1*x1^-1*x0^-1*x1*x3*x0, x3^-1*x1^-1*x3*x0*x1*x0^-1, x2^-1*x0^-1*x2*x0 > >>> dic {0: [x1], 1: [x3], 2: [x2], 3: [x0], 4: [x3^-1*x2^-1*x1^-1*x0^-1]} >>> fundamental_group_arrangement(L, vertical=True) (Finitely presented group < x0, x1, x2, x3 | x3*x0*x3^-1*x0^-1, x3*x1*x3^-1*x1^-1, x1*x2*x0*x2^-1*x1^-1*x0^-1, x1*x2*x0*x1^-1*x0^-1*x2^-1 >, {0: [x2], 1: [x0], 2: [x3], 3: [x1], 4: [x3^-1*x2^-1*x1^-1*x0^-1]})
# needs sirocco from sage.schemes.curves.zariski_vankampen import braid_monodromy from sage.schemes.curves.zariski_vankampen import fundamental_group_arrangement R.<x, y> = QQ[] flist = [x^2 - y^3, x + 3 * y - 5] g, dic = fundamental_group_arrangement(flist) g.sorted_presentation() dic g, dic = fundamental_group_arrangement(flist, simplified=False, puiseux=False) g.sorted_presentation(), dic fundamental_group_arrangement(flist, projective=True) fundamental_group_arrangement([]) g, dic = fundamental_group_arrangement([x * y]) g.sorted_presentation(), dic fundamental_group_arrangement([y + x^2]) fundamental_group_arrangement([y^2 + x], projective=True) L = [x, y, x - 1, x -y] G, dic =fundamental_group_arrangement(L) G.sorted_presentation() dic fundamental_group_arrangement(L, vertical=True)
- sage.schemes.curves.zariski_vankampen.fundamental_group_from_braid_mon(bm, degree=None, simplified=True, projective=False, puiseux=True, vertical=[])[source]¶
Return a presentation of the fundamental group computed from a braid monodromy.
INPUT:
bm
– list of braidsdegree
– integer (default:None
); only needed if the braid monodromy is an empty listsimplified
– boolean (default:True
); if set toTrue
the presentation will be simplified (see below)projective
– boolean (default:False
); if set toTrue
, the fundamental group of the complement of the projective completion of the curve will be computed, otherwise, the fundamental group of the complement in the affine plane will be computedpuiseux
– boolean (default:True
); if set toTrue
a presentation of the fundamental group with the homotopy type of the complement of the affine curve will be computed, adding one relation ifprojective
is set toTrue
.vertical
– list of integers (default:[]
); the indices in[0 .. r - 1]
of the braids that surround a vertical line
If
projective
isFalse
andpuiseux
isTrue
, a Zariski-VanKampen presentation is returned.OUTPUT:
A presentation of the fundamental group of the complement of the union of the curve with some vertical lines from its braid monodromy.
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import fundamental_group_from_braid_mon sage: B.<s0, s1, s2> = BraidGroup(4) sage: bm = [s1*s2*s0*s1*s0^-1*s1^-1*s0^-1, ....: s0*s1^2*s0*s2*s1*(s0^-1*s1^-1)^2*s0^-1, ....: (s0*s1)^2] sage: g = fundamental_group_from_braid_mon(bm, projective=True); g # needs sirocco Finitely presented group < x1, x3 | x3^2*x1^2, x1^-1*x3^-1*x1*x3^-1*x1^-1*x3^-1 > sage: print(g.order(), g.abelian_invariants()) # needs sirocco 12 (4,) sage: B2 = BraidGroup(2) sage: bm = [B2(3 * [1])] sage: g = fundamental_group_from_braid_mon(bm, vertical=[0]); g # needs sirocco Finitely presented group < x0, x1, x2 | x2*x0*x1*x2^-1*x1^-1*x0^-1, x2*x0*x1*x0*x1^-1*x0^-1*x2^-1*x1^-1 > sage: fundamental_group_from_braid_mon([]) is None # needs sirocco True sage: fundamental_group_from_braid_mon([], degree=2) # needs sirocco Finitely presented group < x0, x1 | > sage: fundamental_group_from_braid_mon([SymmetricGroup(1).one()]) # needs sirocco Finitely presented group < x | >
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import fundamental_group_from_braid_mon >>> B = BraidGroup(Integer(4), names=('s0', 's1', 's2',)); (s0, s1, s2,) = B._first_ngens(3) >>> bm = [s1*s2*s0*s1*s0**-Integer(1)*s1**-Integer(1)*s0**-Integer(1), ... s0*s1**Integer(2)*s0*s2*s1*(s0**-Integer(1)*s1**-Integer(1))**Integer(2)*s0**-Integer(1), ... (s0*s1)**Integer(2)] >>> g = fundamental_group_from_braid_mon(bm, projective=True); g # needs sirocco Finitely presented group < x1, x3 | x3^2*x1^2, x1^-1*x3^-1*x1*x3^-1*x1^-1*x3^-1 > >>> print(g.order(), g.abelian_invariants()) # needs sirocco 12 (4,) >>> B2 = BraidGroup(Integer(2)) >>> bm = [B2(Integer(3) * [Integer(1)])] >>> g = fundamental_group_from_braid_mon(bm, vertical=[Integer(0)]); g # needs sirocco Finitely presented group < x0, x1, x2 | x2*x0*x1*x2^-1*x1^-1*x0^-1, x2*x0*x1*x0*x1^-1*x0^-1*x2^-1*x1^-1 > >>> fundamental_group_from_braid_mon([]) is None # needs sirocco True >>> fundamental_group_from_braid_mon([], degree=Integer(2)) # needs sirocco Finitely presented group < x0, x1 | > >>> fundamental_group_from_braid_mon([SymmetricGroup(Integer(1)).one()]) # needs sirocco Finitely presented group < x | >
from sage.schemes.curves.zariski_vankampen import fundamental_group_from_braid_mon B.<s0, s1, s2> = BraidGroup(4) bm = [s1*s2*s0*s1*s0^-1*s1^-1*s0^-1, s0*s1^2*s0*s2*s1*(s0^-1*s1^-1)^2*s0^-1, (s0*s1)^2] g = fundamental_group_from_braid_mon(bm, projective=True); g # needs sirocco print(g.order(), g.abelian_invariants()) # needs sirocco B2 = BraidGroup(2) bm = [B2(3 * [1])] g = fundamental_group_from_braid_mon(bm, vertical=[0]); g # needs sirocco fundamental_group_from_braid_mon([]) is None # needs sirocco fundamental_group_from_braid_mon([], degree=2) # needs sirocco fundamental_group_from_braid_mon([SymmetricGroup(1).one()]) # needs sirocco
- sage.schemes.curves.zariski_vankampen.geometric_basis(G, E, EC0, p, dual_graph, vertical_regions={})[source]¶
Return a geometric basis, based on a vertex.
INPUT:
G
– a graph with the bounded edges of a Voronoi DiagramE
– a subgraph ofG
which is a cycle containing the bounded edges touching an unbounded region of a Voronoi DiagramEC0
– a counterclockwise orientation of the vertices ofE
p
– a vertex ofE
dual_graph
– a dual graph for a plane embedding ofG
such thatE
is the boundary of the non-bounded component of the complement. The edges are labelled as the dual edges and the vertices are labelled by a tuple whose first element is the an integer for the position and the second one is the cyclic ordered list of vertices in the regionvertical_regions
– dictionary (default: \({}\)); its keys are the vertices ofdual_graph
to fix regions associated with vertical lines
OUTPUT: a geometric basis and a dictionary
The geometric basis is formed by a list of sequences of paths. Each path is a ist of vertices, that form a closed path in
G
, based atp
, that goes to a region, surrounds it, and comes back by the same path it came. The concatenation of all these paths is equivalent toE
.The dictionary associates to each vertical line the index of the generator of the geometric basis associated to it.
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import geometric_basis, corrected_voronoi_diagram, voronoi_cells sage: points = (0, -1, I, 1, -I) sage: V = corrected_voronoi_diagram(points) sage: G, E, p, EC, DG, VR = voronoi_cells(V, vertical_lines=frozenset((0 .. 4))) sage: gb, vd = geometric_basis(G, E, EC, p, DG, vertical_regions=VR) sage: gb [[A vertex at (5/2, -5/2), A vertex at (5/2, 5/2), A vertex at (-5/2, 5/2), A vertex at (-1/2, 1/2), A vertex at (-1/2, -1/2), A vertex at (1/2, -1/2), A vertex at (1/2, 1/2), A vertex at (-1/2, 1/2), A vertex at (-5/2, 5/2), A vertex at (5/2, 5/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (5/2, 5/2), A vertex at (-5/2, 5/2), A vertex at (-1/2, 1/2), A vertex at (1/2, 1/2), A vertex at (5/2, 5/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (5/2, 5/2), A vertex at (1/2, 1/2), A vertex at (1/2, -1/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (1/2, -1/2), A vertex at (-1/2, -1/2), A vertex at (-1/2, 1/2), A vertex at (-5/2, 5/2), A vertex at (-5/2, -5/2), A vertex at (-1/2, -1/2), A vertex at (1/2, -1/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (1/2, -1/2), A vertex at (-1/2, -1/2), A vertex at (-5/2, -5/2), A vertex at (5/2, -5/2)]] sage: vd {0: 0, 1: 3, 2: 1, 3: 2, 4: 4}
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import geometric_basis, corrected_voronoi_diagram, voronoi_cells >>> points = (Integer(0), -Integer(1), I, Integer(1), -I) >>> V = corrected_voronoi_diagram(points) >>> G, E, p, EC, DG, VR = voronoi_cells(V, vertical_lines=frozenset((ellipsis_iter(Integer(0) ,Ellipsis, Integer(4))))) >>> gb, vd = geometric_basis(G, E, EC, p, DG, vertical_regions=VR) >>> gb [[A vertex at (5/2, -5/2), A vertex at (5/2, 5/2), A vertex at (-5/2, 5/2), A vertex at (-1/2, 1/2), A vertex at (-1/2, -1/2), A vertex at (1/2, -1/2), A vertex at (1/2, 1/2), A vertex at (-1/2, 1/2), A vertex at (-5/2, 5/2), A vertex at (5/2, 5/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (5/2, 5/2), A vertex at (-5/2, 5/2), A vertex at (-1/2, 1/2), A vertex at (1/2, 1/2), A vertex at (5/2, 5/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (5/2, 5/2), A vertex at (1/2, 1/2), A vertex at (1/2, -1/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (1/2, -1/2), A vertex at (-1/2, -1/2), A vertex at (-1/2, 1/2), A vertex at (-5/2, 5/2), A vertex at (-5/2, -5/2), A vertex at (-1/2, -1/2), A vertex at (1/2, -1/2), A vertex at (5/2, -5/2)], [A vertex at (5/2, -5/2), A vertex at (1/2, -1/2), A vertex at (-1/2, -1/2), A vertex at (-5/2, -5/2), A vertex at (5/2, -5/2)]] >>> vd {0: 0, 1: 3, 2: 1, 3: 2, 4: 4}
from sage.schemes.curves.zariski_vankampen import geometric_basis, corrected_voronoi_diagram, voronoi_cells points = (0, -1, I, 1, -I) V = corrected_voronoi_diagram(points) G, E, p, EC, DG, VR = voronoi_cells(V, vertical_lines=frozenset((0 .. 4))) gb, vd = geometric_basis(G, E, EC, p, DG, vertical_regions=VR) gb vd
- sage.schemes.curves.zariski_vankampen.newton(f, x0, i0)[source]¶
Return the interval Newton operator.
INPUT:
f
– a univariate polynomialx0
– a numberI0
– an interval
OUTPUT:
The interval \(x_0-\frac{f(x_0)}{f'(I_0)}\)
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import newton sage: R.<x> = QQbar[] sage: f = x^3 + x sage: x0 = 1/10 sage: I0 = RIF((-1/5,1/5)) sage: n = newton(f, x0, I0) sage: n 0.0? sage: n.real().endpoints() (-0.0147727272727274, 0.00982142857142862) sage: n.imag().endpoints() (0.000000000000000, -0.000000000000000)
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import newton >>> R = QQbar['x']; (x,) = R._first_ngens(1) >>> f = x**Integer(3) + x >>> x0 = Integer(1)/Integer(10) >>> I0 = RIF((-Integer(1)/Integer(5),Integer(1)/Integer(5))) >>> n = newton(f, x0, I0) >>> n 0.0? >>> n.real().endpoints() (-0.0147727272727274, 0.00982142857142862) >>> n.imag().endpoints() (0.000000000000000, -0.000000000000000)
from sage.schemes.curves.zariski_vankampen import newton R.<x> = QQbar[] f = x^3 + x x0 = 1/10 I0 = RIF((-1/5,1/5)) n = newton(f, x0, I0) n n.real().endpoints() n.imag().endpoints()
- sage.schemes.curves.zariski_vankampen.orient_circuit(circuit, convex=False, precision=53, verbose=False)[source]¶
Reverse a circuit if it goes clockwise; otherwise leave it unchanged.
INPUT:
circuit
– a circuit in the graph of a Voronoi Diagram, given by a list of edgesconvex
– boolean (default:False
); if set toTrue
a simpler computation is madeprecision
– bits of precision (default: 53)verbose
– boolean (default:False
); for testing purposes
OUTPUT:
The same circuit if it goes counterclockwise, and its reversed otherwise, given as the ordered list of vertices with identic extremities.
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import orient_circuit sage: points = [(-4, 0), (4, 0), (0, 4), (0, -4), (0, 0)] sage: V = VoronoiDiagram(points) sage: E = Graph() sage: for reg in V.regions().values(): ....: if reg.rays() or reg.lines(): ....: E = E.union(reg.vertex_graph()) sage: E.vertices(sort=True) [A vertex at (-2, -2), A vertex at (-2, 2), A vertex at (2, -2), A vertex at (2, 2)] sage: cir = E.eulerian_circuit() sage: cir [(A vertex at (-2, -2), A vertex at (2, -2), None), (A vertex at (2, -2), A vertex at (2, 2), None), (A vertex at (2, 2), A vertex at (-2, 2), None), (A vertex at (-2, 2), A vertex at (-2, -2), None)] sage: cir_oriented = orient_circuit(cir); cir_oriented (A vertex at (-2, -2), A vertex at (2, -2), A vertex at (2, 2), A vertex at (-2, 2), A vertex at (-2, -2)) sage: cirinv = list(reversed([(c[1],c[0],c[2]) for c in cir])) sage: cirinv [(A vertex at (-2, -2), A vertex at (-2, 2), None), (A vertex at (-2, 2), A vertex at (2, 2), None), (A vertex at (2, 2), A vertex at (2, -2), None), (A vertex at (2, -2), A vertex at (-2, -2), None)] sage: orient_circuit(cirinv) == cir_oriented True sage: cir_oriented == orient_circuit(cir, convex=True) True sage: P0=[(1,1/2),(0,1),(1,1)]; P1=[(0,3/2),(-1,0)] sage: Q=Polyhedron(P0).vertices() sage: Q = [Q[2], Q[0], Q[1]] + [_ for _ in reversed(Polyhedron(P1).vertices())] sage: Q [A vertex at (1, 1/2), A vertex at (0, 1), A vertex at (1, 1), A vertex at (0, 3/2), A vertex at (-1, 0)] sage: E = Graph() sage: for v, w in zip(Q, Q[1:] + [Q[0]]): ....: E.add_edge((v, w)) sage: cir = orient_circuit(E.eulerian_circuit(), precision=1, verbose=True) 2 sage: cir (A vertex at (1, 1/2), A vertex at (0, 1), A vertex at (1, 1), A vertex at (0, 3/2), A vertex at (-1, 0), A vertex at (1, 1/2))
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import orient_circuit >>> points = [(-Integer(4), Integer(0)), (Integer(4), Integer(0)), (Integer(0), Integer(4)), (Integer(0), -Integer(4)), (Integer(0), Integer(0))] >>> V = VoronoiDiagram(points) >>> E = Graph() >>> for reg in V.regions().values(): ... if reg.rays() or reg.lines(): ... E = E.union(reg.vertex_graph()) >>> E.vertices(sort=True) [A vertex at (-2, -2), A vertex at (-2, 2), A vertex at (2, -2), A vertex at (2, 2)] >>> cir = E.eulerian_circuit() >>> cir [(A vertex at (-2, -2), A vertex at (2, -2), None), (A vertex at (2, -2), A vertex at (2, 2), None), (A vertex at (2, 2), A vertex at (-2, 2), None), (A vertex at (-2, 2), A vertex at (-2, -2), None)] >>> cir_oriented = orient_circuit(cir); cir_oriented (A vertex at (-2, -2), A vertex at (2, -2), A vertex at (2, 2), A vertex at (-2, 2), A vertex at (-2, -2)) >>> cirinv = list(reversed([(c[Integer(1)],c[Integer(0)],c[Integer(2)]) for c in cir])) >>> cirinv [(A vertex at (-2, -2), A vertex at (-2, 2), None), (A vertex at (-2, 2), A vertex at (2, 2), None), (A vertex at (2, 2), A vertex at (2, -2), None), (A vertex at (2, -2), A vertex at (-2, -2), None)] >>> orient_circuit(cirinv) == cir_oriented True >>> cir_oriented == orient_circuit(cir, convex=True) True >>> P0=[(Integer(1),Integer(1)/Integer(2)),(Integer(0),Integer(1)),(Integer(1),Integer(1))]; P1=[(Integer(0),Integer(3)/Integer(2)),(-Integer(1),Integer(0))] >>> Q=Polyhedron(P0).vertices() >>> Q = [Q[Integer(2)], Q[Integer(0)], Q[Integer(1)]] + [_ for _ in reversed(Polyhedron(P1).vertices())] >>> Q [A vertex at (1, 1/2), A vertex at (0, 1), A vertex at (1, 1), A vertex at (0, 3/2), A vertex at (-1, 0)] >>> E = Graph() >>> for v, w in zip(Q, Q[Integer(1):] + [Q[Integer(0)]]): ... E.add_edge((v, w)) >>> cir = orient_circuit(E.eulerian_circuit(), precision=Integer(1), verbose=True) 2 >>> cir (A vertex at (1, 1/2), A vertex at (0, 1), A vertex at (1, 1), A vertex at (0, 3/2), A vertex at (-1, 0), A vertex at (1, 1/2))
from sage.schemes.curves.zariski_vankampen import orient_circuit points = [(-4, 0), (4, 0), (0, 4), (0, -4), (0, 0)] V = VoronoiDiagram(points) E = Graph() for reg in V.regions().values(): if reg.rays() or reg.lines(): E = E.union(reg.vertex_graph()) E.vertices(sort=True) cir = E.eulerian_circuit() cir cir_oriented = orient_circuit(cir); cir_oriented cirinv = list(reversed([(c[1],c[0],c[2]) for c in cir])) cirinv orient_circuit(cirinv) == cir_oriented cir_oriented == orient_circuit(cir, convex=True) P0=[(1,1/2),(0,1),(1,1)]; P1=[(0,3/2),(-1,0)] Q=Polyhedron(P0).vertices() Q = [Q[2], Q[0], Q[1]] + [_ for _ in reversed(Polyhedron(P1).vertices())] Q E = Graph() for v, w in zip(Q, Q[1:] + [Q[0]]): E.add_edge((v, w)) cir = orient_circuit(E.eulerian_circuit(), precision=1, verbose=True) cir
- sage.schemes.curves.zariski_vankampen.populate_roots_interval_cache(inputs)[source]¶
Call
roots_interval()
to the inputs that have not been computed previously, and cache them.INPUT:
inputs
– list of tuples(f, x0)
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import populate_roots_interval_cache, roots_interval_cache, fieldI sage: R.<x,y> = QQ[] sage: K=fieldI(QQ) sage: f = y^5 - x^2 sage: f = f.change_ring(K) sage: (f, 3) in roots_interval_cache False sage: populate_roots_interval_cache([(f, 3)]) sage: (f, 3) in roots_interval_cache True sage: roots_interval_cache[(f, 3)] {-1.255469441943070? - 0.9121519421827974?*I: -2.? - 1.?*I, -1.255469441943070? + 0.9121519421827974?*I: -2.? + 1.?*I, 0.4795466549853897? - 1.475892845355996?*I: 1.? - 2.?*I, 0.4795466549853897? + 1.475892845355996?*I: 1.? + 2.?*I, 14421467174121563/9293107134194871: 2.? + 0.?*I}
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import populate_roots_interval_cache, roots_interval_cache, fieldI >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> K=fieldI(QQ) >>> f = y**Integer(5) - x**Integer(2) >>> f = f.change_ring(K) >>> (f, Integer(3)) in roots_interval_cache False >>> populate_roots_interval_cache([(f, Integer(3))]) >>> (f, Integer(3)) in roots_interval_cache True >>> roots_interval_cache[(f, Integer(3))] {-1.255469441943070? - 0.9121519421827974?*I: -2.? - 1.?*I, -1.255469441943070? + 0.9121519421827974?*I: -2.? + 1.?*I, 0.4795466549853897? - 1.475892845355996?*I: 1.? - 2.?*I, 0.4795466549853897? + 1.475892845355996?*I: 1.? + 2.?*I, 14421467174121563/9293107134194871: 2.? + 0.?*I}
from sage.schemes.curves.zariski_vankampen import populate_roots_interval_cache, roots_interval_cache, fieldI R.<x,y> = QQ[] K=fieldI(QQ) f = y^5 - x^2 f = f.change_ring(K) (f, 3) in roots_interval_cache populate_roots_interval_cache([(f, 3)]) (f, 3) in roots_interval_cache roots_interval_cache[(f, 3)]
- sage.schemes.curves.zariski_vankampen.roots_interval_cached(f, x0)[source]¶
Cached version of
roots_interval()
.
- sage.schemes.curves.zariski_vankampen.strand_components(f, pols, p1)[source]¶
Compute only the assignment from strands to elements of
flist
.INPUT:
f
– a reduced polynomial with two variables, over a number field with an embedding in the complex numberspols
– a list of polynomials with two variables whose product equalsf
p1
– a Gauss rational
OUTPUT:
A list and a dictionary. The first one is an ordered list of pairs consisting of
(z,i)
wherez
is a root off(p_1,y)
and \(i\) is the position of the polynomial in the list whose root isz
. The second one attaches a number \(i\) (strand) to a number \(j\) (a polynomial in the list).
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import strand_components sage: R.<x, y> = QQ[] sage: flist = [x^2 - y^3, x + 3 * y - 5] sage: strand_components(prod(flist), flist, 1) ([(-0.500000000000000? - 0.866025403784439?*I, 0), (-0.500000000000000? + 0.866025403784439?*I, 0), (1, 0), (1.333333333333334?, 1)], {0: 0, 1: 0, 2: 0, 3: 1})
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import strand_components >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> flist = [x**Integer(2) - y**Integer(3), x + Integer(3) * y - Integer(5)] >>> strand_components(prod(flist), flist, Integer(1)) ([(-0.500000000000000? - 0.866025403784439?*I, 0), (-0.500000000000000? + 0.866025403784439?*I, 0), (1, 0), (1.333333333333334?, 1)], {0: 0, 1: 0, 2: 0, 3: 1})
from sage.schemes.curves.zariski_vankampen import strand_components R.<x, y> = QQ[] flist = [x^2 - y^3, x + 3 * y - 5] strand_components(prod(flist), flist, 1)
- sage.schemes.curves.zariski_vankampen.vertical_lines_in_braidmon(pols)[source]¶
Return the vertical lines in
pols
, unless one of the other components has a vertical asymptote.INPUT:
pols
– a list of polynomials with two variables whose product equalsf
OUTPUT:
A list with the indices of the vertical lines in
flist
if there is no other component with vertical asymptote; otherwise it returns an empty list.EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import vertical_lines_in_braidmon sage: R.<x, y> = QQ[] sage: flist = [x^2 - y^3, x, x + 3 * y - 5, 1 - x] sage: vertical_lines_in_braidmon(flist) [1, 3] sage: flist += [x * y - 1] sage: vertical_lines_in_braidmon(flist) [] sage: vertical_lines_in_braidmon([]) []
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import vertical_lines_in_braidmon >>> R = QQ['x, y']; (x, y,) = R._first_ngens(2) >>> flist = [x**Integer(2) - y**Integer(3), x, x + Integer(3) * y - Integer(5), Integer(1) - x] >>> vertical_lines_in_braidmon(flist) [1, 3] >>> flist += [x * y - Integer(1)] >>> vertical_lines_in_braidmon(flist) [] >>> vertical_lines_in_braidmon([]) []
from sage.schemes.curves.zariski_vankampen import vertical_lines_in_braidmon R.<x, y> = QQ[] flist = [x^2 - y^3, x, x + 3 * y - 5, 1 - x] vertical_lines_in_braidmon(flist) flist += [x * y - 1] vertical_lines_in_braidmon(flist) vertical_lines_in_braidmon([])
- sage.schemes.curves.zariski_vankampen.voronoi_cells(V, vertical_lines=frozenset())[source]¶
Compute the graph, the boundary graph, a base point, a positive orientation of the boundary graph, and the dual graph of a corrected Voronoi diagram.
INPUT:
V
– a corrected Voronoi diagramvertical_lines
– frozenset (default:frozenset()
); indices of the vertical lines
OUTPUT:
G
– the graph of the 1-skeleton ofV
E
– the subgraph of the boundaryp
– a vertex inE
EC
– list of vertices (representing a counterclockwise orientation ofE
) with identical first and last elements)DG
– the dual graph ofV
, where the vertices are labelled by the compact regions ofV
and the edges by their dual edgesvertical_regions
– dictionary for the regions associated with vertical lines
EXAMPLES:
sage: from sage.schemes.curves.zariski_vankampen import corrected_voronoi_diagram, voronoi_cells sage: points = (2, I, 0.000001, 0, 0.000001*I) sage: V = corrected_voronoi_diagram(points) sage: G, E, p, EC, DG, VR = voronoi_cells(V, vertical_lines=frozenset((1,))) sage: Gv = G.vertices(sort=True) sage: Ge = G.edges(sort=True) sage: len(Gv), len(Ge) (12, 16) sage: Ev = E.vertices(sort=True); Ev [A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (-7/2, -7/2), A vertex at (-7/2, 1/2000000), A vertex at (1/2000000, -7/2), A vertex at (2000001/2000000, -24500001/7000000), A vertex at (11/4, 4), A vertex at (9/2, -9/2), A vertex at (9/2, 9/2)] sage: Ev.index(p) 7 sage: EC (A vertex at (9/2, -9/2), A vertex at (9/2, 9/2), A vertex at (11/4, 4), A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (-7/2, 1/2000000), A vertex at (-7/2, -7/2), A vertex at (1/2000000, -7/2), A vertex at (2000001/2000000, -24500001/7000000), A vertex at (9/2, -9/2)) sage: len(DG.vertices(sort=True)), len(DG.edges(sort=True)) (5, 7) sage: edg = DG.edges(sort=True)[0]; edg ((0, (A vertex at (9/2, -9/2), A vertex at (9/2, 9/2), A vertex at (11/4, 4), A vertex at (2000001/2000000, 500001/1000000), A vertex at (2000001/2000000, -24500001/7000000), A vertex at (9/2, -9/2))), (1, (A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (1000001/2000000, 1000001/2000000), A vertex at (2000001/2000000, 500001/1000000), A vertex at (11/4, 4), A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000))), (A vertex at (2000001/2000000, 500001/1000000), A vertex at (11/4, 4), None)) sage: edg[-1] in Ge True sage: VR {1: (A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (1000001/2000000, 1000001/2000000), A vertex at (2000001/2000000, 500001/1000000), A vertex at (11/4, 4), A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000))}
>>> from sage.all import * >>> from sage.schemes.curves.zariski_vankampen import corrected_voronoi_diagram, voronoi_cells >>> points = (Integer(2), I, RealNumber('0.000001'), Integer(0), RealNumber('0.000001')*I) >>> V = corrected_voronoi_diagram(points) >>> G, E, p, EC, DG, VR = voronoi_cells(V, vertical_lines=frozenset((Integer(1),))) >>> Gv = G.vertices(sort=True) >>> Ge = G.edges(sort=True) >>> len(Gv), len(Ge) (12, 16) >>> Ev = E.vertices(sort=True); Ev [A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (-7/2, -7/2), A vertex at (-7/2, 1/2000000), A vertex at (1/2000000, -7/2), A vertex at (2000001/2000000, -24500001/7000000), A vertex at (11/4, 4), A vertex at (9/2, -9/2), A vertex at (9/2, 9/2)] >>> Ev.index(p) 7 >>> EC (A vertex at (9/2, -9/2), A vertex at (9/2, 9/2), A vertex at (11/4, 4), A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (-7/2, 1/2000000), A vertex at (-7/2, -7/2), A vertex at (1/2000000, -7/2), A vertex at (2000001/2000000, -24500001/7000000), A vertex at (9/2, -9/2)) >>> len(DG.vertices(sort=True)), len(DG.edges(sort=True)) (5, 7) >>> edg = DG.edges(sort=True)[Integer(0)]; edg ((0, (A vertex at (9/2, -9/2), A vertex at (9/2, 9/2), A vertex at (11/4, 4), A vertex at (2000001/2000000, 500001/1000000), A vertex at (2000001/2000000, -24500001/7000000), A vertex at (9/2, -9/2))), (1, (A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (1000001/2000000, 1000001/2000000), A vertex at (2000001/2000000, 500001/1000000), A vertex at (11/4, 4), A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000))), (A vertex at (2000001/2000000, 500001/1000000), A vertex at (11/4, 4), None)) >>> edg[-Integer(1)] in Ge True >>> VR {1: (A vertex at (-49000001/14000000, 1000001/2000000), A vertex at (1000001/2000000, 1000001/2000000), A vertex at (2000001/2000000, 500001/1000000), A vertex at (11/4, 4), A vertex at (-4, 4), A vertex at (-49000001/14000000, 1000001/2000000))}
from sage.schemes.curves.zariski_vankampen import corrected_voronoi_diagram, voronoi_cells points = (2, I, 0.000001, 0, 0.000001*I) V = corrected_voronoi_diagram(points) G, E, p, EC, DG, VR = voronoi_cells(V, vertical_lines=frozenset((1,))) Gv = G.vertices(sort=True) Ge = G.edges(sort=True) len(Gv), len(Ge) Ev = E.vertices(sort=True); Ev Ev.index(p) EC len(DG.vertices(sort=True)), len(DG.edges(sort=True)) edg = DG.edges(sort=True)[0]; edg edg[-1] in Ge VR