def _test_result(number_field, prec=53, epsilon=1e-10): """ sage: from sage.rings.complex_field import ComplexField sage: CF = ComplexField() sage: x = var('x') sage: nf = NumberField(x**2 + 1, 'x', embedding = CF(1.0j)) sage: _test_result(nf) sage: nf = NumberField(x**2 + 7, 'x', embedding = CF(2.64575j)) sage: _test_result(nf) sage: nf = NumberField(x**8 + 6 * x ** 4 + x + 23, 'x', embedding = CF(0.7747 + 1.25937j)) sage: _test_result(nf, 212, epsilon = 1e-30) """ CIF = ComplexIntervalField(prec) RIF = RealIntervalField(prec) real_number_field, x_expression, y_expression = ( field_containing_real_and_imaginary_part_of_number_field(number_field)) x_val = x_expression.lift()(RIF(real_number_field.gen_embedding())) y_val = y_expression.lift()(RIF(real_number_field.gen_embedding())) z_val = CIF(x_val, y_val) diff = z_val - CIF(number_field.gen_embedding()) if not abs(diff) < epsilon: raise Exception("Test failed")
def _singularities(self, dom, include_apparent=True, multiplicities=False): dop = self if include_apparent else self.desingularize() # TBI if isinstance(dom, ComplexBallField): # TBI dom1 = ComplexIntervalField(dom.precision()) else: dom1 = dom lc = dop.leading_coefficient() sing = lc.roots(dom1, multiplicities=multiplicities) if dom1 is not dom: sing = [dom(s) for s in sing] return sing
def _singularities(self, dom, include_apparent=True, multiplicities=False): if not multiplicities: rts = self._singularities(dom, include_apparent, multiplicities=True) return [s for s, _ in rts] if isinstance(dom, ComplexBallField): # TBI dom1 = ComplexIntervalField(dom.precision()) rts = self._singularities(dom1, include_apparent, multiplicities) return [(dom(s), m) for s, m in rts] dop = self if include_apparent else self.desingularize() # TBI lc = dop.leading_coefficient() try: return lc.roots(dom) except NotImplementedError: return lc.change_ring(QQbar).roots(dom)
def interval_roots(p, rts, prec): """ We are given a squarefree polynomial p, a list of estimated roots, and a precision. We attempt to verify that the estimated roots are in fact distinct roots of the polynomial, using interval arithmetic of precision prec. If we succeed, we return a list of intervals bounding the roots; if we fail, we return None. EXAMPLES:: sage: x = polygen(ZZ) sage: p = x^3 - 1 sage: rts = [CC.zeta(3)^i for i in range(0, 3)] sage: from sage.rings.polynomial.complex_roots import interval_roots sage: interval_roots(p, rts, 53) [1, -0.500000000000000? + 0.866025403784439?*I, -0.500000000000000? - 0.866025403784439?*I] sage: interval_roots(p, rts, 200) [1, -0.500000000000000000000000000000000000000000000000000000000000? + 0.866025403784438646763723170752936183471402626905190314027904?*I, -0.500000000000000000000000000000000000000000000000000000000000? - 0.866025403784438646763723170752936183471402626905190314027904?*I] """ CIF = ComplexIntervalField(prec) CIFX = CIF['x'] ip = CIFX(p) ipd = CIFX(p.derivative()) irts = [] for rt in rts: irt = refine_root(ip, ipd, CIF(rt), CIF) if irt is None: return None irts.append(irt) return irts
def followstrand(f, x0, x1, y0a, prec=53): r""" Return a piecewise linear aproximation of the homotopy continuation of the root ``y0a`` from ``x0`` to ``x1``. INPUT: - ``f`` -- a polynomial in two variables - ``x0`` -- a complex value, where the homotopy starts - ``x1`` -- a complex value, where the homotopy ends - ``y0a`` -- 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 aproximation) - the piecewise linear path determined by the points has a tubular neighborhood where the actual homotopy continuation path lies, and no other root intersects it. EXAMPLES:: sage: from sage.schemes.curves.zariski_vankampen import followstrand # optional - sirocco 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) # optional - sirocco # abs tol 1e-15 [(0.0, -1.0, 0.0), (0.7500000000000001, -1.015090921153253, -0.24752813818386948), (1.0, -1.026166099551513, -0.32768940253604323)] """ CIF = ComplexIntervalField(prec) CC = ComplexField(prec) G = f.change_ring(QQbar).change_ring(CIF) (x, y) = G.variables() g = G.subs({x: (1-x)*CIF(x0) + x*CIF(x1)}) coefs = [] deg = g.total_degree() for d in range(deg + 1): for i in range(d + 1): c = CIF(g.coefficient({x: d-i, y: i})) cr = c.real() ci = c.imag() coefs += list(cr.endpoints()) coefs += list(ci.endpoints()) yr = CC(y0a).real() yi = CC(y0a).imag() from sage.libs.sirocco import contpath, contpath_mp try: if prec == 53: points = contpath(deg, coefs, yr, yi) else: points = contpath_mp(deg, coefs, yr, yi, prec) return points except Exception: return followstrand(f, x0, x1, y0a, 2*prec)
def exactly_checked_canonical_retriangulation(M, bits_prec, degree): """ Given a proto-canonical triangulation of a cusped (possibly non-orientable) manifold M, return its canonical retriangulation which is computed from exact shapes. The exact shapes are computed using snap (which uses the LLL-algorithm). The precision (in bits) and the maximal degree need to be specified (here 300 bits precision and polynomials of degree less than 4):: sage: from snappy import Manifold sage: M = Manifold("m412") sage: M.canonize() sage: K = exactly_checked_canonical_retriangulation(M, 300, 4) M's canonical cell decomposition has a cube, so non-tetrahedral:: sage: K.has_finite_vertices() True Has 12 tetrahedra after the retrianglation:: sage: K.num_tetrahedra() 12 Check that it fails on something which is not a proto-canonical triangulation:: sage: from snappy import Manifold sage: M = Manifold("m015") sage: exactly_checked_canonical_retriangulation(M, 500, 6) # doctest: +IGNORE_EXCEPTION_DETAIL Traceback (most recent call last): ... TiltProvenPositiveNumericalVerifyError: Numerical verification that tilt is negative has failed, tilt is actually positive. This is provably not the proto-canonical triangulation: 0.1645421638874662848910671879? <= 0 """ # Convert to decimal precision dec_prec = prec_bits_to_dec(bits_prec) # Try to find the exact shapes shapes = find_shapes_as_complex_sqrt_lin_combinations(M, dec_prec, degree) if not shapes: raise FindExactShapesError() # Build the cusp cross section c = RealCuspCrossSection.fromManifoldAndShapes(M, shapes) # Check that the exact solutions form a complete hyperbolic structure # We convert to intervals to check that the shapes are positive and # the angles add up to 2pi and not some other multiple of 2pi. c.check_polynomial_edge_equations_exactly() c.check_cusp_development_exactly() CIF = ComplexIntervalField(bits_prec) c.check_logarithmic_edge_equations_and_positivity(CIF) # Normalize cusp area. This is not needed when only 1 cusp if M.num_cusps() > 1: c.normalize_cusps() # Compute tilts c.compute_tilts() # Get the opacity of a face in the proto-canonical triangulation def get_opacity(tilt): # Get the tilt of the sign. The sign method is implemented # to use exact arithmetic to certify that the sign is 0 and # to use interval arithmetic (of increasing precision until a decision # can be made) to certify the sign otherwise. sign, interval = tilt.sign_with_interval() # Tilt is negative, return True if sign < 0: return True # Tilt is zero, return False if sign == 0: return False # Tilt is positive, raise exception if sign > 0: raise exceptions.TiltProvenPositiveNumericalVerifyError(interval) def index_of_face_corner(corner): face_index = t3m.simplex.comp(corner.Subsimplex).bit_length() - 1 return 4 * corner.Tetrahedron.Index + face_index # Opacities of all four faces of each tetrahedron, initialize with None. # The format is opacity of face 0, 1, 2, 3 of the first tetrahedron, # ... of second tetrahedron, ... opacities = (4 * len(c.mcomplex.Tetrahedra)) * [None] # For each face of the triangulation for face in c.mcomplex.Faces: opacity = get_opacity(face.Tilt) for corner in face.Corners: opacities[index_of_face_corner(corner)] = opacity if None in opacities: raise Exception("Mismatch with opacities") # If there are transparent faces, this triangulation is just the # proto-canonical triangulation. We need to call into the SnapPea # kernel to retriangulate (introduces finite vertices) if False in opacities: return M._canonical_retriangulation(opacities) # No transparent faces, this triangulation is the canonical triangulation. # Return it without introducing finite vertices. return M
def __init__(self, M, initial_shapes, bits_prec=None, dec_prec=None): """ Initializes the KrawczykShapesEngine given an orientable SnapPy Manifold M, approximated solutions initial_shapes to the gluing equations (e.g., as returned by M.tetrahedra_shapes('rect')) and the precision to be used for the desired computations in either bits bits_prec or decimal digits dec_prec. This requires Sage since it uses Sage's ComplexIntervalField for its computations. Note that this will choose an independent set of edge equations and one equation per cusp. It is known that a solution to such a subset of rectangular gluing equations is also a solution to the full set of rectangular gluing equations:: sage: from snappy import Manifold sage: M = Manifold("m019") sage: C = KrawczykShapesEngine(M, M.tetrahedra_shapes('rect'), bits_prec = 53) sage: C.expand_until_certified() True sage: C.certified_shapes # doctest: +ELLIPSIS (0.780552527850...? + 0.914473662967...?*I, 0.780552527850...? + 0.91447366296773?*I, 0.4600211755737...? + 0.6326241936052...?*I) Does not work with non-orientable manifolds:: sage: M = Manifold("m000") sage: KrawczykShapesEngine(M, M.tetrahedra_shapes('rect'), bits_prec = 53) Traceback (most recent call last): ... Exception: Manifold needs to be orientable Or some non-hyperbolic manifolds:: sage: Manifold("t02333(1,0)").tetrahedra_shapes(intervals = True) Traceback (most recent call last): ... RuntimeError: Could not certify shape intervals, either there are degenerate shapes or the precision must be increased. """ # Require sage if not _within_sage: raise SageNotAvailable( "Sorry, the verify module can only be used within Sage") # Convert to precision in bits if necessary if dec_prec: self.prec = prec_dec_to_bits(dec_prec) elif bits_prec: self.prec = bits_prec else: raise Exception("Need dec_prec or bits_prec") # Setup interval types of desired precision self.CIF = ComplexIntervalField(self.prec) self.RIF = RealIntervalField(self.prec) # Verify that manifold is orientable if not M.is_orientable(): raise Exception("Manifold needs to be orientable") # Initialize the shape intervals, they have zero length self.initial_shapes = vector( [self.CIF(shape) for shape in initial_shapes]) # Get an independent set of gluing equations from snap self.equations = snap.shapes.enough_gluing_equations(M) self._make_sparse_equations() self.identity = matrix.identity(self.CIF, len(self.initial_shapes)) CDF = ComplexDoubleField() # Could be sparse approx_deriv = self.log_gluing_LHS_derivatives( [CDF(shape) for shape in initial_shapes]) approx_inverse_double = approx_deriv.inverse() self.approx_inverse = approx_inverse_double.change_ring(self.CIF) # Compute the term z0 - c * f(z0) in the formula for # the Krawczyk interval K(z0, [z], f) value_at_initial_shapes = self.log_gluing_LHSs(self.initial_shapes) self.first_term = (self.initial_shapes - self.approx_inverse * value_at_initial_shapes) # Shapes have not been certified yet self.certified_shapes = None
def covariant_z0(F, z0_cov=False, prec=53, emb=None, error_limit=0.000001): r""" Return the covariant and Julia invariant from Cremona-Stoll [CS2003]_. In [CS2003]_ and [HS2018]_ the Julia invariant is denoted as `\Theta(F)` or `R(F, z(F))`. Note that you may get faster convergence if you first move `z_0(F)` to the fundamental domain before computing the true covariant INPUT: - ``F`` -- binary form of degree at least 3 with no multiple roots - ``z0_cov`` -- boolean, compute only the `z_0` invariant. Otherwise, solve the minimization problem - ``prec``-- positive integer. precision to use in CC - ``emb`` -- embedding into CC - ``error_limit`` -- sets the error tolerance (default:0.000001) OUTPUT: a complex number, a real number EXAMPLES:: sage: from sage.rings.polynomial.binary_form_reduce import covariant_z0 sage: R.<x,y> = QQ[] sage: F = 19*x^8 - 262*x^7*y + 1507*x^6*y^2 - 4784*x^5*y^3 + 9202*x^4*y^4\ ....: - 10962*x^3*y^5 + 7844*x^2*y^6 - 3040*x*y^7 + 475*y^8 sage: covariant_z0(F, prec=80, z0_cov=True) (1.3832330115323681438175 + 0.31233552177413614978744*I, 3358.4074848663492819259) sage: F = -x^8 + 6*x^7*y - 7*x^6*y^2 - 12*x^5*y^3 + 27*x^4*y^4\ ....: - 4*x^3*y^5 - 19*x^2*y^6 + 10*x*y^7 - 5*y^8 sage: covariant_z0(F, prec=80) (0.64189877107807122203366 + 1.1852516565091601348355*I, 3134.5148284344627168276) :: sage: R.<x,y> = QQ[] sage: covariant_z0(x^3 + 2*x^2*y - 3*x*y^2, z0_cov=True)[0] 0.230769230769231 + 0.799408065031789*I sage: -1/covariant_z0(-y^3 + 2*y^2*x + 3*y*x^2, z0_cov=True)[0] 0.230769230769231 + 0.799408065031789*I :: sage: R.<x,y> = QQ[] sage: covariant_z0(2*x^2*y - 3*x*y^2, z0_cov=True)[0] 0.750000000000000 + 1.29903810567666*I sage: -1/covariant_z0(-x^3 - x^2*y + 2*x*y^2, z0_cov=True)[0] + 1 0.750000000000000 + 1.29903810567666*I :: sage: R.<x,y> = QQ[] sage: covariant_z0(x^2*y - x*y^2, prec=100) # tol 1e-28 (0.50000000000000000000000000003 + 0.86602540378443864676372317076*I, 1.5396007178390020386910634147) TESTS:: sage: R.<x,y>=QQ[] sage: covariant_z0(x^2 + 24*x*y + y^2) Traceback (most recent call last): ... ValueError: must be at least degree 3 sage: covariant_z0((x+y)^3, z0_cov=True) Traceback (most recent call last): ... ValueError: cannot have multiple roots for z0 invariant sage: covariant_z0(x^3 + 3*x*y + y) Traceback (most recent call last): ... TypeError: must be a binary form sage: covariant_z0(-2*x^2*y^3 + 3*x*y^4 + 127*y^5) Traceback (most recent call last): ... ValueError: cannot have a root with multiplicity >= 5/2 sage: covariant_z0((x^2+2*y^2)^2) Traceback (most recent call last): ... ValueError: must have at least 3 distinct roots """ R = F.parent() d = ZZ(F.degree()) if R.ngens() != 2 or any(sum(t) != d for t in F.exponents()): raise TypeError('must be a binary form') if d < 3: raise ValueError('must be at least degree 3') f = F.subs({R.gen(1): 1}).univariate_polynomial() if f.degree() < d: # we have a root at infinity if f.constant_coefficient() != 0: # invert so we find all roots! mat = matrix(ZZ, 2, 2, [0, -1, 1, 0]) else: t = 0 while f(t) == 0: t += 1 mat = matrix(ZZ, 2, 2, [t, -1, 1, 0]) else: mat = matrix(ZZ, 2, 2, [1, 0, 0, 1]) f = F(list(mat * vector(R.gens()))).subs({ R.gen(1): 1 }).univariate_polynomial() # now we have a single variable polynomial with all the roots of F K = ComplexField(prec=prec) if f.base_ring() != K: if emb is None: f = f.change_ring(K) else: f = f.change_ring(emb) roots = f.roots() if max(ex for _, ex in roots) > 1 or f.degree() < d - 1: if z0_cov: raise ValueError('cannot have multiple roots for z0 invariant') else: # just need a starting point for Newton's method f = f.lc() * prod([p for p, ex in f.factor() ]) # removes multiple roots if f.degree() < 3: raise ValueError('must have at least 3 distinct roots') roots = f.roots() roots = [p for p, _ in roots] # finding quadratic Q_0, gives us our covariant, z_0 dF = f.derivative() n = ZZ(f.degree()) PR = PolynomialRing(K, 'x,y') x, y = PR.gens() # finds Stoll and Cremona's Q_0 q = sum([(1/(dF(r).abs()**(2/(n-2)))) * ((x-(r*y)) * (x-(r.conjugate()*y)))\ for r in roots]) # this is Q_0 , always positive def as long as F has distinct roots A = q.monomial_coefficient(x**2) B = q.monomial_coefficient(x * y) C = q.monomial_coefficient(y**2) # need positive root try: z = ((-B + ((B**2) - (4 * A * C)).sqrt()) / (2 * A)) except ValueError: raise ValueError("not enough precision") if z.imag() < 0: z = (-B - ((B**2) - (4 * A * C)).sqrt()) / (2 * A) if z0_cov: FM = f # for Julia's invariant else: # solve the minimization problem for 'true' covariant CF = ComplexIntervalField( prec=prec) # keeps trac of our precision error z = CF(z) FM = F(list(mat * vector(R.gens()))).subs({ R.gen(1): 1 }).univariate_polynomial() from sage.rings.polynomial.complex_roots import complex_roots L1 = complex_roots(FM, min_prec=prec) L = [] err = z.diameter() # making sure multiplicity isn't too large using convergence conditions in paper for p, e in L1: if e >= d / 2: raise ValueError( 'cannot have a root with multiplicity >= %s/2' % d) for _ in range(e): L.append(p) RCF = PolynomialRing(CF, 'u,t') a = RCF(0) c = RCF(0) u, t = RCF.gens() for l in L: a += u**2 / ((t - l) * (t - l.conjugate()) + u**2) c += (t - l.real()) / ((t - l) * (t - l.conjugate()) + u**2) # Newton's Method, to find solutions. Error bound is less than diameter of our z err = z.diameter() zz = z.diameter() g1 = a.numerator() - d / 2 * a.denominator() g2 = c.numerator() G = vector([g1, g2]) J = jacobian(G, [u, t]) v0 = vector([z.imag(), z.real()]) # z0 as starting point # finds our correct z while err <= zz: NJ = J.subs({u: v0[0], t: v0[1]}) NJinv = NJ.inverse() # inverse for CIF matrix seems to return fractions not CIF elements, fix them if NJinv.base_ring() != CF: NJinv = matrix(CF, 2, 2, [ CF(zw.numerator() / zw.denominator()) for zw in NJinv.list() ]) w = z v0 = v0 - NJinv * G.subs({u: v0[0], t: v0[1]}) z = v0[1].constant_coefficient( ) + v0[0].constant_coefficient() * CF.gen(0) err = z.diameter() # precision zz = (w - z).abs() # difference in w and z else: if err > error_limit or err.is_NaN(): raise ValueError( "accuracy of Newton's root not within tolerance(%s > %s), increase precision" % (err, error_limit)) if z.imag() <= z.diameter(): raise ArithmeticError( "Newton's method converged to z not in the upper half plane") z = z.center() # Julia's invariant if FM.base_ring() != ComplexField(prec=prec): FM = FM.change_ring(ComplexField(prec=prec)) tF = z.real() uF = z.imag() th = FM.lc().abs()**2 for r, ex in FM.roots(): for _ in range(ex): th = th * ((((r - tF).abs())**2 + uF**2) / uF) # undo shift and invert (if needed) # since F \cdot m ~ m^(-1)\cdot z # we apply m to z to undo m acting on F l = mat * vector([z, 1]) return l[0] / l[1], th
def followstrand(f, x0, x1, y0a, prec=53): """ Return a piecewise linear aproximation of the homotopy continuation of the root y0a from x0 to x1 INPUT: - ``f`` -- a polynomial in two variables - ``x0`` -- a complex value, where the homotopy starts - ``x1`` -- a complex value, where the homotopy ends - ``y0a`` -- an approximate solution of the polynomial $F(y)=f(x_0,y)$ - ``prec`` -- the precission to use OUTPUT: A list of values (t, ytr, yti) 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 aproximation) - the piecewise linear path determined by the points has a tubular neighborhood where the actual homotopy continuation path lies, and no other root intersects it. EXAMPLES:: 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) # optional [(0.0, -1.0, 0.0), (0.063125, -1.0001106593601545, -0.02104011456212618), (0.12230468750000001, -1.0004151100003031, -0.04075695242737829), (0.17778564453125, -1.0008762007617709, -0.059227299979315154), (0.28181243896484376, -1.0021948141328698, -0.09380038023464156), (0.3728358840942383, -1.0038270754728402, -0.123962953123039), (0.4524813985824585, -1.005613540368227, -0.15026634875747985), (0.5221712237596512, -1.0074443351354077, -0.17320066690539515), (0.5831498207896948, -1.009246118007067, -0.1931978258913501), (0.636506093190983, -1.0109719597443307, -0.21063630045261386), (0.6831928315421101, -1.0125937110987449, -0.2258465242053158), (0.7648946236565826, -1.0156754074572174, -0.2523480191006915), (0.8261709677424369, -1.0181837235093538, -0.2721208327884435), (0.8721282258068277, -1.0201720738092597, -0.2868892148703381), (0.9410641129034139, -1.0233210275568283, -0.3089391941950436), (1.0, -1.026166099551513, -0.3276894025360433)] """ CIF = ComplexIntervalField(prec) CC = ComplexField(prec) G = f.change_ring(QQbar).change_ring(CIF) (x, y) = G.variables() g = G.subs({x: (1 - x) * CIF(x0) + x * CIF(x1)}) coefs = [] deg = g.total_degree() for d in range(deg + 1): for i in range(d + 1): c = CIF(g.coefficient({x: d - i, y: i})) cr = c.real() ci = c.imag() coefs += list(cr.endpoints()) coefs += list(ci.endpoints()) yr = CC(y0a).real() yi = CC(y0a).imag() try: if prec == 53: points = contpath(deg, coefs, yr, yi) else: points = contpath_mp(deg, coefs, yr, yi, prec) return points except: return followstrand(f, x0, x1, y0a, 2 * prec)
def field_containing_real_and_imaginary_part_of_number_field(number_field): """ Given a Sage number field number_field with a complex embedding z, return (real_number_field, real_part, imag_part). The number field real_number_field is the smallest number field containing the real part and imaginary part of every element in number_field. real_part and imag_part are elements in real_number_field which comes with a real embedding such that under this embedding, we have z = real_part + imag_part * I. sage: from sage.rings.complex_field import ComplexField sage: CF = ComplexField() sage: x = var('x') sage: nf = NumberField(x**2 + 1, 'x', embedding = CF(1.0j)) sage: field_containing_real_and_imaginary_part_of_number_field(nf) (Number Field in x with defining polynomial x with x = 0, 0, 1) sage: nf = NumberField(x**2 + 7, 'x', embedding = CF(2.64575j)) sage: field_containing_real_and_imaginary_part_of_number_field(nf) (Number Field in x with defining polynomial x^2 - 7 with x = 2.645751311064591?, 0, x) sage: nf = NumberField(x**3 + x**2 + 23, 'x', embedding = CF(1.1096 + 2.4317j)) sage: field_containing_real_and_imaginary_part_of_number_field(nf) (Number Field in x with defining polynomial x^6 + 2*x^5 + 2*x^4 - 113/2*x^3 - 229/4*x^2 - 115/4*x - 575/8 with x = 3.541338405550421?, -20/14377*x^5 + 382/14377*x^4 + 526/14377*x^3 + 1533/14377*x^2 - 18262/14377*x - 10902/14377, 20/14377*x^5 - 382/14377*x^4 - 526/14377*x^3 - 1533/14377*x^2 + 32639/14377*x + 10902/14377) """ # Let p be the defining polynomial of the given number field. # Given the one complex equation p(z) = 0, translate it into two # real equations Re(p(x+y*I)) = 0, Im(p(x+y*I)) = 0. # equations are sage symbolic expressions in x and y. equations = [ _real_or_imaginary_part_for_polynomial_in_complex_variable( number_field.defining_polynomial(), start) for start in [0, 1] ] # In _solve_two_equations, we implemented a method that can solve # a system of two polynomial equations in two variables x and y # provided that y is in the number field generated by x. # If we are lucky, this is the case. # If we are unlucky, the number field containing x and y is generated # by x' = x + k * y where k is some small natural number. # The k mentioned above. We start with 0 and increase until we # succeed. k = 0 # The amount of extra precision beyond double precision we are working # with. We increase it if one of the above methods fails to find the right # factor of the above polynomials. extra_prec = 0 # Initialize the intervals for x and y with double prescision intervals CIF = ComplexIntervalField() z_val = CIF(number_field.gen_embedding()) x_val = z_val.real() y_val = z_val.imag() # Keep trying to find a k or increase precision until we succeed # Give up if k reaches 100 or we are at precision 16 times greater than # that of a double while k < 100 and extra_prec < 5: # Compute the interval for x' xprime_val = x_val + k * y_val # From the equations for x and y, get the equations for x' and y # where x' = x + k * y as abover equations_for_xprime = [ eqn.substitute(x=var('x') - k * var('y')) for eqn in equations ] try: # Try to find a solution to the two equations solution = _solve_two_equations(equations_for_xprime[0], equations_for_xprime[1], xprime_val, y_val) if solution: # We succeeded. We have a solution for the equations in # x' and y, thus, we need to do x = x' - k * y real_number_field, y_expression = solution x_expression = real_number_field.gen() - k * y_expression return real_number_field, x_expression, y_expression else: # No solution found. This means that y is not in the # number field generated by x'. Try a higher k in the # next iteration k += 1 except _IsolateFactorError: # We did not use enough precision. The given intervals for # x and y that are supposed to isolate a solution to the # system of two equations did not have enough precision to # succeed and give a unique answer. # Double the precision we will use from now on extra_prec += 1 # Recompute the intervals for x and y with the new precision. CIF = ComplexIntervalField(53 * 2**extra_prec) z_val = CIF(number_field.gen_embedding()) x_val = z_val.real() y_val = z_val.imag() # Give up return None
def refine_root(ip, ipd, irt, fld): """ We are given a polynomial and its derivative (with complex interval coefficients), an estimated root, and a complex interval field to use in computations. We use interval arithmetic to refine the root and prove that we have in fact isolated a unique root. If we succeed, we return the isolated root; if we fail, we return None. EXAMPLES:: sage: from sage.rings.polynomial.complex_roots import * sage: x = polygen(ZZ) sage: p = x^9 - 1 sage: ip = CIF['x'](p); ip x^9 - 1 sage: ipd = CIF['x'](p.derivative()); ipd 9*x^8 sage: irt = CIF(CC(cos(2*pi/9), sin(2*pi/9))); irt 0.76604444311897802? + 0.64278760968653926?*I sage: ip(irt) 0.?e-14 + 0.?e-14*I sage: ipd(irt) 6.89439998807080? - 5.78508848717885?*I sage: refine_root(ip, ipd, irt, CIF) 0.766044443118978? + 0.642787609686540?*I """ # There has got to be a better way to do this, but I don't know # what it is... # We start with a basic fact: if we do an interval Newton-Raphson # step, and the refined interval is contained in the original interval, # then the refined interval contains exactly one root. # Unfortunately, our initial estimated root almost certainly does not # contain the actual root (our initial interval is a point, which # is exactly equal to whatever floating-point estimate we got from # the external solver). So we need to do multiple Newton-Raphson # steps, and check this inclusion property on each step. # After a few steps of refinement, if the initial root estimate was # close to a root, we should have an essentially perfect interval # bound on the root (since Newton-Raphson has quadratic convergence), # unless either the real or imaginary component of the root is zero. # If the real or imaginary component is zero, then we could spend # a long time computing closer and closer approximations to that # component. (This doesn't happen for non-zero components, because # of the imprecision of floating-point numbers combined with the # outward interval rounding; but close to zero, MPFI provides # extremely precise numbers.) # If the root is actually a real root, but we start with an imaginary # component, we can bounce back and forth between having a positive # and negative imaginary component, without ever hitting zero. # To deal with this, on every other Newton-Raphson step, instead of # replacing the old interval with the new one, we take the union. # If the containment check continues to fail many times in a row, # we give up and return None; we also return None if we detect # that the slope in our current interval is not bounded away # from zero at any step. # After every refinement step, we check to see if the real or # imaginary component of our interval includes zero. If so, we # try setting it to exactly zero. This gives us a good chance of # detecting real roots. However, we do this replacement at most # once per component. refinement_steps = 10 smashed_real = False smashed_imag = False for i in range(refinement_steps): slope = ipd(irt) if slope.contains_zero(): return None center = fld(irt.center()) val = ip(center) nirt = center - val / slope # print irt, nirt, (nirt in irt), nirt.diameter(), irt.diameter(), center, val, slope if nirt in irt and (nirt.diameter() >= irt.diameter() >> 3 or i >= 8): # If the new diameter is much less than the original diameter, # then we have not yet converged. (Perhaps we were asked # for a particularly high-precision result.) So we don't # return yet. return nirt if i & 1: irt = nirt else: irt = irt.union(nirt) # If we don't find a root after a while, try (approximately) # tripling the size of the region. if i >= 6: rD = irt.real().absolute_diameter() iD = irt.imag().absolute_diameter() md = max(rD, iD) md_intv = RealIntervalField(rD.prec())(-md, md) md_cintv = ComplexIntervalField(rD.prec())(md_intv, md_intv) irt = irt + md_cintv if not smashed_real and irt.real().contains_zero(): irt = irt.parent()(0, irt.imag()) smashed_real = True if not smashed_imag and irt.imag().contains_zero(): irt = irt.parent()(irt.real(), 0) smashed_imag = True return None
def orient_circuit(circuit): r""" Reverses a circuit if it goes clockwise; otherwise leaves it unchanged. INPUT: - ``circuit`` -- a circuit in the graph of a Voronoi Diagram, given by a list of edges OUTPUT: The same circuit if it goes counterclockwise, and its reverse otherwise 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() [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: orient_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)] 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) [(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)] """ prec = 53 vectors = [v[1].vector() - v[0].vector() for v in circuit] while True: CIF = ComplexIntervalField(prec) totalangle = sum((CIF(*vectors[i]) / CIF(*vectors[i - 1])).argument() for i in range(len(vectors))) if totalangle < 0: return list(reversed([(c[1], c[0]) + c[2:] for c in circuit])) elif totalangle > 0: return circuit else: prec *= 2
def braid_in_segment(g, x0, x1): """ Return the braid formed by the `y` roots of ``f`` when `x` moves from ``x0`` to ``x1``. INPUT: - ``g`` -- a polynomial factorization in two variables - ``x0`` -- a complex number - ``x1`` -- a complex number OUTPUT: A braid. EXAMPLES:: sage: from sage.schemes.curves.zariski_vankampen import braid_in_segment # optional - sirocco sage: R.<x,y> = QQ[] sage: f = x^2 + y^3 sage: x0 = CC(1,0) sage: x1 = CC(1, 0.5) sage: braid_in_segment(f.factor(), x0, x1) # optional - sirocco s1 TESTS: Check that :trac:`26503` is fixed:: sage: wp = QQ['t']([1, 1, 1]).roots(QQbar)[0][0] sage: Kw.<wp> = NumberField(wp.minpoly(), embedding=wp) sage: R.<x, y> = Kw[] sage: z = -wp - 1 sage: f = y*(y + z)*x*(x - 1)*(x - y)*(x + z*y - 1)*(x + z*y + wp) sage: from sage.schemes.curves import zariski_vankampen as zvk sage: g = f.subs({x: x + 2*y}) sage: p1 = QQbar(sqrt(-1/3)) sage: p2 = QQbar(1/2+sqrt(-1/3)/2) sage: B = zvk.braid_in_segment(g.factor(),CC(p1),CC(p2)) # optional - sirocco sage: B # optional - sirocco s5*s3^-1 """ (x, y) = g.value().parent().gens() I = QQbar.gen() X0 = QQ(x0.real()) + I * QQ(x0.imag()) X1 = QQ(x1.real()) + I * QQ(x1.imag()) intervals = {} precision = {} y0s = [] for (f, naux) in g: if f.variables() == (y, ): F0 = QQbar[y](f.base_ring()[y](f)) else: F0 = QQbar[y](f(X0, y)) y0sf = F0.roots(multiplicities=False) y0s += list(y0sf) precision[f] = 53 while True: CIFp = ComplexIntervalField(precision[f]) intervals[f] = [r.interval(CIFp) for r in y0sf] if not any( a.overlaps(b) for a, b in itertools.combinations(intervals[f], 2)): break precision[f] *= 2 strands = [ followstrand(f[0], [p[0] for p in g if p[0] != f[0]], x0, x1, i.center(), precision[f[0]]) for f in g for i in intervals[f[0]] ] complexstrands = [[(QQ(a[0]), QQ(a[1]), QQ(a[2])) for a in b] for b in strands] centralbraid = braid_from_piecewise(complexstrands) initialstrands = [] finalstrands = [] initialintervals = roots_interval_cached(g.value(), X0) finalintervals = roots_interval_cached(g.value(), X1) for cs in complexstrands: ip = cs[0][1] + I * cs[0][2] fp = cs[-1][1] + I * cs[-1][2] matched = 0 for center, interval in initialintervals.items(): if ip in interval: initialstrands.append([(0, center.real(), center.imag()), (1, cs[0][1], cs[0][2])]) matched += 1 if matched == 0: raise ValueError("unable to match braid endpoint with root") if matched > 1: raise ValueError("braid endpoint mathes more than one root") matched = 0 for center, interval in finalintervals.items(): if fp in interval: finalstrands.append([(0, cs[-1][1], cs[-1][2]), (1, center.real(), center.imag())]) matched += 1 if matched == 0: raise ValueError("unable to match braid endpoint with root") if matched > 1: raise ValueError("braid endpoint mathes more than one root") initialbraid = braid_from_piecewise(initialstrands) finalbraid = braid_from_piecewise(finalstrands) return initialbraid * centralbraid * finalbraid
def roots_interval(f, x0): """ Find disjoint intervals that isolate the roots of a polynomial for a fixed value of the first variable. INPUT: - ``f`` -- a bivariate squarefree polynomial - ``x0`` -- a value where the first coordinate will be fixed The intervals are taken as big as possible to be able to detect when two approximate roots of `f(x_0, y)` correspond to the same exact root. The result is given as a dictionary, where the keys are approximations to the roots with rational real and imaginary parts, and the values are intervals containing them. EXAMPLES:: sage: from sage.schemes.curves.zariski_vankampen import roots_interval sage: R.<x,y> = QQ[] sage: f = y^3 - x^2 sage: ri = roots_interval(f, 1) sage: ri {-138907099/160396102*I - 1/2: -1.? - 1.?*I, 138907099/160396102*I - 1/2: -1.? + 1.?*I, 1: 1.? + 0.?*I} sage: [r.endpoints() for r in ri.values()] [(0.566987298107781 - 0.433012701892219*I, 1.43301270189222 + 0.433012701892219*I, 0.566987298107781 + 0.433012701892219*I, 1.43301270189222 - 0.433012701892219*I), (-0.933012701892219 - 1.29903810567666*I, -0.0669872981077806 - 0.433012701892219*I, -0.933012701892219 - 0.433012701892219*I, -0.0669872981077806 - 1.29903810567666*I), (-0.933012701892219 + 0.433012701892219*I, -0.0669872981077806 + 1.29903810567666*I, -0.933012701892219 + 1.29903810567666*I, -0.0669872981077806 + 0.433012701892219*I)] """ x, y = f.parent().gens() I = QQbar.gen() fx = QQbar[y](f.subs({x: QQ(x0.real()) + I * QQ(x0.imag())})) roots = fx.roots(QQbar, multiplicities=False) result = {} for i in range(len(roots)): r = roots[i] prec = 53 IF = ComplexIntervalField(prec) CF = ComplexField(prec) divisor = 4 diam = min((CF(r) - CF(r0)).abs() for r0 in roots[:i] + roots[i + 1:]) / divisor envelop = IF(diam) * IF((-1, 1), (-1, 1)) while not newton(fx, r, r + envelop) in r + envelop: prec += 53 IF = ComplexIntervalField(prec) CF = ComplexField(prec) divisor *= 2 diam = min([(CF(r) - CF(r0)).abs() for r0 in roots[:i] + roots[i + 1:]]) / divisor envelop = IF(diam) * IF((-1, 1), (-1, 1)) qapr = QQ(CF(r).real()) + QQbar.gen() * QQ(CF(r).imag()) if qapr not in r + envelop: raise ValueError("Could not approximate roots with exact values") result[qapr] = r + envelop return result
def followstrand(f, factors, x0, x1, y0a, prec=53): r""" Return a piecewise linear approximation of the homotopy continuation of the root ``y0a`` from ``x0`` to ``x1``. INPUT: - ``f`` -- an irreducible polynomial in two variables - ``factors`` -- a list of irreducible polynomials in two variables - ``x0`` -- a complex value, where the homotopy starts - ``x1`` -- a complex value, where the homotopy ends - ``y0a`` -- 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 in ``factors``, intersects it. EXAMPLES:: sage: from sage.schemes.curves.zariski_vankampen import followstrand # optional - sirocco 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) # optional - sirocco # 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) # optional - sirocco # 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)] """ if f.degree() == 1: CF = ComplexField(prec) g = f.change_ring(CF) (x, y) = g.parent().gens() y0 = CF[y](g.subs({x: x0})).roots()[0][0] y1 = CF[y](g.subs({x: x1})).roots()[0][0] res = [(0.0, y0.real(), y0.imag()), (1.0, y1.real(), y1.imag())] return res CIF = ComplexIntervalField(prec) CC = ComplexField(prec) G = f.change_ring(QQbar).change_ring(CIF) (x, y) = G.parent().gens() g = G.subs({x: (1 - x) * CIF(x0) + x * CIF(x1)}) coefs = [] deg = g.total_degree() for d in range(deg + 1): for i in range(d + 1): c = CIF(g.coefficient({x: d - i, y: i})) cr = c.real() ci = c.imag() coefs += list(cr.endpoints()) coefs += list(ci.endpoints()) yr = CC(y0a).real() yi = CC(y0a).imag() coefsfactors = [] degsfactors = [] for fc in factors: degfc = fc.degree() degsfactors.append(degfc) G = fc.change_ring(QQbar).change_ring(CIF) g = G.subs({x: (1 - x) * CIF(x0) + x * CIF(x1)}) for d in range(degfc + 1): for i in range(d + 1): c = CIF(g.coefficient({x: d - i, y: i})) cr = c.real() ci = c.imag() coefsfactors += list(cr.endpoints()) coefsfactors += list(ci.endpoints()) from sage.libs.sirocco import contpath, contpath_mp, contpath_comps, contpath_mp_comps try: if prec == 53: if factors: points = contpath_comps(deg, coefs, yr, yi, degsfactors, coefsfactors) else: points = contpath(deg, coefs, yr, yi) else: if factors: points = contpath_mp_comps(deg, coefs, yr, yi, prec, degsfactors, coefsfactors) else: points = contpath_mp(deg, coefs, yr, yi, prec) return points except Exception: return followstrand(f, factors, x0, x1, y0a, 2 * prec)