def zeta_symmetric(s): r""" Completed function `\xi(s)` that satisfies `\xi(s) = \xi(1-s)` and has zeros at the same points as the Riemann zeta function. INPUT: - ``s`` - real or complex number If s is a real number the computation is done using the MPFR library. When the input is not real, the computation is done using the PARI C library. More precisely, .. MATH:: xi(s) = \gamma(s/2 + 1) * (s-1) * \pi^{-s/2} * \zeta(s). EXAMPLES:: sage: zeta_symmetric(0.7) 0.497580414651127 sage: zeta_symmetric(1-0.7) 0.497580414651127 sage: RR = RealField(200) sage: zeta_symmetric(RR(0.7)) 0.49758041465112690357779107525638385212657443284080589766062 sage: C.<i> = ComplexField() sage: zeta_symmetric(0.5 + i*14.0) 0.000201294444235258 + 1.49077798716757e-19*I sage: zeta_symmetric(0.5 + i*14.1) 0.0000489893483255687 + 4.40457132572236e-20*I sage: zeta_symmetric(0.5 + i*14.2) -0.0000868931282620101 + 7.11507675693612e-20*I REFERENCE: - I copied the definition of xi from http://web.viu.ca/pughg/RiemannZeta/RiemannZetaLong.html """ if not (is_ComplexNumber(s) or is_RealNumber(s)): s = ComplexField()(s) R = s.parent() if s == 1: # deal with poles, hopefully return R(0.5) return (s / 2 + 1).gamma() * (s - 1) * (R.pi()**(-s / 2)) * s.zeta()
def get_bound_poly(F, prec=53, norm_type='norm', emb=None): """ The hyperbolic distance from `j` which must contain the smallest poly. This defines the maximum possible distance from `j` to the `z_0` covariant in the hyperbolic 3-space for which the associated `F` could have smaller coefficients. INPUT: - ``F`` -- binary form of degree at least 3 with no multiple roots - ``prec``-- positive integer. precision to use in CC - ``norm_type`` -- string, either norm or height - ``emb`` -- embedding into CC OUTPUT: a positive real number EXAMPLES:: sage: from sage.rings.polynomial.binary_form_reduce import get_bound_poly sage: R.<x,y> = QQ[] sage: F = -2*x^3 + 2*x^2*y + 3*x*y^2 + 127*y^3 sage: get_bound_poly(F) # tol 1e-12 28.0049336543295 sage: get_bound_poly(F, norm_type='height') # tol 1e-11 111.890642019092 """ if F.base_ring() != ComplexField(prec=prec): if emb is None: compF = F.change_ring(ComplexField(prec=prec)) else: compF = F.change_ring(emb) else: compF = F n = F.degree() assert (n > 2), "degree 2 polynomial" z0F, thetaF = covariant_z0(compF, prec=prec, emb=emb) if norm_type == 'norm': # euclidean norm squared normF = (sum([abs(i)**2 for i in compF.coefficients()])) target = (2**(n - 1)) * normF / thetaF elif norm_type == 'height': hF = exp(max([c.global_height(prec=prec) for c in F.coefficients()])) # height target = (2**(n - 1)) * (n + 1) * (hF**2) / thetaF else: raise ValueError('type must be norm or height') return cosh(epsinv(F, target, prec=prec))
def check_integrals_of_motion(self, affine_parameter, solution_key=None): r""" Check the constancy of the four integrals of motion INPUT: - ``affine_parameter`` -- value of the affine parameter `\lambda` - ``solution_key`` -- (default: ``None``) string denoting the numerical solution to use for evaluating the various integrals of motion; if ``None``, the latest solution computed by :meth:`integrate` is used. OUTPUT: - a `SageMath table <https://doc.sagemath.org/html/en/reference/misc/sage/misc/table.html>`_ with the absolute and relative differences with respect to the initial values. """ CF = ComplexField(16) lambda_min = self.domain().lower_bound() res = [[ "quantity", "value", "initial value", "diff.", "relative diff." ]] mu2 = self.evaluate_mu2(affine_parameter, solution_key=solution_key) mu20 = self.evaluate_mu2(lambda_min, solution_key=solution_key) diff = mu2 - mu20 rel_diff = CF(diff / mu20) if mu20 != 0 else "-" res.append([r"$\mu^2$", mu2, mu20, CF(diff), rel_diff]) E = self.evaluate_E(affine_parameter, solution_key=solution_key) E0 = self.evaluate_E(lambda_min, solution_key=solution_key) diff = E - E0 rel_diff = CF(diff / E0) if E0 != 0 else "-" res.append([r"$E$", E, E0, CF(diff), rel_diff]) L = self.evaluate_L(affine_parameter, solution_key=solution_key) L0 = self.evaluate_L(lambda_min, solution_key=solution_key) diff = L - L0 rel_diff = CF(diff / L0) if L0 != 0 else "-" res.append([r"$L$", L, L0, CF(diff), rel_diff]) Q = self.evaluate_Q(affine_parameter, solution_key=solution_key) Q0 = self.evaluate_Q(lambda_min, solution_key=solution_key) diff = Q - Q0 rel_diff = CF(diff / Q0) if Q0 != 0 else "-" res.append([r"$Q$", Q, Q0, CF(diff), rel_diff]) return table(res, align="center")
def q_log(q, u): """ Determines, if possible, an integer n such that q^n = u. Requires that both q and u belong to either QQ or some rational function field over QQ. q must not be zero or a root of unity. A ValueError is thrown if no n exists. """ if q in QQ and u in QQ: qq, uu = q, u else: q, u = canonical_coercion(q, u) ev = dict((y, hash(y)) for y in u.parent().gens_dict_recursive()) qq, uu = q(**ev), u(**ev) n = ComplexField(53)(uu.n().log() / qq.n().log()).real_part().round() if q**n == u: return n else: raise ValueError
def get_bound_dynamical(F, f, m=1, dynatomic=True, prec=53, emb=None): """ The hyperbolic distance from `j` which must contain the smallest map. This defines the maximum possible distance from `j` to the `z_0` covariant of the associated binary form `F` in the hyperbolic 3-space for which the map `f`` could have smaller coefficients. INPUT: - ``F`` -- binary form of degree at least 3 with no multiple roots associated to ``f`` - ``f`` -- a dynamical system on `P^1` - ``m`` - positive integer. the period used to create ``F`` - ``dynatomic`` -- boolean. whether ``F`` is the periodic points or the formal periodic points of period ``m`` for ``f`` - ``prec``-- positive integer. precision to use in CC - ``emb`` -- embedding into CC OUTPUT: a positive real number EXAMPLES:: sage: from sage.dynamics.arithmetic_dynamics.endPN_minimal_model import get_bound_dynamical sage: P.<x,y> = ProjectiveSpace(QQ,1) sage: f = DynamicalSystem([50*x^2 + 795*x*y + 2120*y^2, 265*x^2 + 106*y^2]) sage: get_bound_dynamical(f.dynatomic_polynomial(1), f) 35.5546923182219 """ def coshdelta(z): #The cosh of the hyperbolic distance from z = t+uj to j return (z.norm() + 1) / (2 * z.imag()) if F.base_ring() != ComplexField(prec=prec): if emb is None: compF = F.change_ring(ComplexField(prec=prec)) else: compF = F.change_ring(emb) else: compF = F n = F.degree() z0F, thetaF = covariant_z0(compF, prec=prec, emb=emb) d = f.degree() hF = e**f.global_height(prec=prec) #get precomputed constants C,k if m == 1: C = 4 * d + 2 k = 2 else: Ck_values = { (False, 2, 2): (322, 6), (False, 2, 3): (385034, 14), (False, 2, 4): (4088003923454, 30), (False, 3, 2): (18044, 8), (False, 4, 2): (1761410, 10), (False, 5, 2): (269283820, 12), (True, 2, 2): (43, 4), (True, 2, 3): (106459, 12), (True, 2, 4): (39216735905, 24), (True, 3, 2): (1604, 6), (True, 4, 2): (114675, 8), (True, 5, 2): (14158456, 10) } try: C, k = Ck_values[(dynatomic, d, m)] except KeyError: raise ValueError("constants not computed for this (m,d) pair") if n == 2 and d == 2: #bound with epsilonF = 1 bound = 2 * ((2 * C * (hF**k)) / (thetaF)) else: bound = cosh(epsinv(F, (2**(n - 1)) * C * (hF**k) / thetaF, prec=prec)) return bound
def elliptic_j(z, prec=53): r""" Returns the elliptic modular `j`-function evaluated at `z`. INPUT: - ``z`` (complex) -- a complex number with positive imaginary part. - ``prec`` (default: 53) -- precision in bits for the complex field. OUTPUT: (complex) The value of `j(z)`. ALGORITHM: Calls the ``pari`` function ``ellj()``. AUTHOR: John Cremona EXAMPLES:: sage: elliptic_j(CC(i)) 1728.00000000000 sage: elliptic_j(sqrt(-2.0)) 8000.00000000000 sage: z = ComplexField(100)(1,sqrt(11))/2 sage: elliptic_j(z) -32768.000... sage: elliptic_j(z).real().round() -32768 :: sage: tau = (1 + sqrt(-163))/2 sage: (-elliptic_j(tau.n(100)).real().round())^(1/3) 640320 This example shows the need for higher precision than the default one of the `ComplexField`, see :trac:`28355`:: sage: -elliptic_j(tau) # rel tol 1e-2 2.62537412640767e17 - 732.558854258998*I sage: -elliptic_j(tau,75) # rel tol 1e-2 2.625374126407680000000e17 - 0.0001309913593909879441262*I sage: -elliptic_j(tau,100) # rel tol 1e-2 2.6253741264076799999999999999e17 - 1.3012822400356887122945119790e-12*I sage: (-elliptic_j(tau, 100).real().round())^(1/3) 640320 """ CC = z.parent() if not isinstance(CC, sage.rings.abc.ComplexField): from sage.rings.complex_mpfr import ComplexField CC = ComplexField(prec) try: z = CC(z) except ValueError: raise ValueError("elliptic_j only defined for complex arguments.") from sage.libs.all import pari return CC(pari(z).ellj())
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 = [] # 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.zero() c = RCF.zero() u, t = RCF.gens() for l in L: denom = ((t - l) * (t - l.conjugate()) + u**2) a += u**2 / denom c += (t - l.real()) / denom # 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().lower() # difference in w and z else: # despite there is no break, this happens 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().upper() <= 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 epsinv(F, target, prec=53, target_tol=0.001, z=None, emb=None): """ Compute a bound on the hyperbolic distance. The true minimum will be within the computed bound. It is computed as the inverse of epsilon_F from [HS2018]_. INPUT: - ``F`` -- binary form of degree at least 3 with no multiple roots - ``target`` -- positive real number. The value we want to attain, i.e., the value we are taking the inverse of - ``prec``-- positive integer. precision to use in CC - ``target_tol`` -- positive real number. The tolerance with which we attain the target value. - ``z`` -- complex number. ``z_0`` covariant for F. - ``emb`` -- embedding into CC OUTPUT: a real number delta satisfying target + target_tol > eps_F(delta) > target. EXAMPLES:: sage: from sage.rings.polynomial.binary_form_reduce import epsinv sage: R.<x,y> = QQ[] sage: epsinv(-2*x^3 + 2*x^2*y + 3*x*y^2 + 127*y^3, 31.5022020249597) # tol 1e-12 4.02520895942207 """ def RQ(delta): # this is the quotient R(F_0,z)/R(F_0,z(F)) for a generic z # at distance delta from j. See Lemma 4.2 in [HS2018]. cd = cosh(delta).n(prec=prec) sd = sinh(delta).n(prec=prec) return prod([cd + (cost * phi[0] + sint * phi[1]) * sd for phi in phis]) def epsF(delta): pol = RQ(delta) # get R quotient in terms of z S = PolynomialRing(C, 'v') g = S([(i - d) * pol[i - d] for i in range(2 * d + 1)]) # take derivative drts = [e for e in g.roots(ring=C, multiplicities=False) if (e.norm() - 1).abs() < 0.1] # find min return min([pol(r / r.abs()).real() for r in drts]) C = ComplexField(prec=prec) R = F.parent() d = F.degree() if z is None: z, th = covariant_z0(F, prec=prec, emb=emb) else: # need to do our own input checking 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() # now we have a single variable polynomial if (max(ex for p, ex in f.roots(ring=C)) >= QQ(d)/2 or f.degree() < QQ(d)/2): raise ValueError('cannot have root with multiplicity >= deg(F)/2') R = RealField(prec=prec) PR = PolynomialRing(R, 't') t = PR.gen(0) # compute phi_1, ..., phi_k # first find F_0 and its roots # this change of variables on f moves z(f) to j, i.e. produces F_0 rts = f(z.imag()*t + z.real()).roots(ring=C) phis = [] # stereographic projection of roots for r, e in rts: phis.extend([[2*r.real()/(r.norm()+1), (r.norm()-1)/(r.norm()+1)]]) if d != f.degree(): # include roots at infinity phis.extend([(d - f.degree()) * [0, 1]]) # for writing RQ in terms of generic z to minimize LC = LaurentSeriesRing(C, 'u', default_prec=2 * d + 2) u = LC.gen(0) cost = (u + u**(-1)) / 2 sint = (u - u**(-1)) / (2 * C.gen(0)) # first find an interval containing the desired value # then use regula falsi on log eps_F # d -> delta value in interval [0,1] # v in value in interval [1,epsF(1)] dl = R(0.0) vl = R(1.0) du = R(1.0) vu = epsF(du) while vu < target: # compute the next value of epsF for delta = 2*delta dl = du vl = vu du *= 2 vu = epsF(du) # now dl < delta <= du logt = target.log() l2 = (vu.log() - logt).n(prec=prec) l1 = (vl.log() - logt).n(prec=prec) dn = (dl*l2 - du*l1)/(l2 - l1) vn = epsF(dn) dl = du vl = vu du = dn vu = vn while (du - dl).abs() >= target_tol or max(vl, vu) < target: l2 = (vu.log() - logt).n(prec=prec) l1 = (vl.log() - logt).n(prec=prec) dn = (dl * l2 - du * l1) / (l2 - l1) vn = epsF(dn) dl = du vl = vu du = dn vu = vn return max(dl, du)
def complex_roots(p, skip_squarefree=False, retval='interval', min_prec=0): """ Compute the complex roots of a given polynomial with exact coefficients (integer, rational, Gaussian rational, and algebraic coefficients are supported). Returns a list of pairs of a root and its multiplicity. Roots are returned as a ComplexIntervalFieldElement; each interval includes exactly one root, and the intervals are disjoint. By default, the algorithm will do a squarefree decomposition to get squarefree polynomials. The skip_squarefree parameter lets you skip this step. (If this step is skipped, and the polynomial has a repeated root, then the algorithm will loop forever!) You can specify retval='interval' (the default) to get roots as complex intervals. The other options are retval='algebraic' to get elements of QQbar, or retval='algebraic_real' to get only the real roots, and to get them as elements of AA. EXAMPLES:: sage: from sage.rings.polynomial.complex_roots import complex_roots sage: x = polygen(ZZ) sage: complex_roots(x^5 - x - 1) [(1.167303978261419?, 1), (-0.764884433600585? - 0.352471546031727?*I, 1), (-0.764884433600585? + 0.352471546031727?*I, 1), (0.181232444469876? - 1.083954101317711?*I, 1), (0.181232444469876? + 1.083954101317711?*I, 1)] sage: v=complex_roots(x^2 + 27*x + 181) Unfortunately due to numerical noise there can be a small imaginary part to each root depending on CPU, compiler, etc, and that affects the printing order. So we verify the real part of each root and check that the imaginary part is small in both cases:: sage: v # random [(-14.61803398874990?..., 1), (-12.3819660112501...? + 0.?e-27*I, 1)] sage: sorted((v[0][0].real(),v[1][0].real())) [-14.61803398874989?, -12.3819660112501...?] sage: v[0][0].imag().upper() < 1e25 True sage: v[1][0].imag().upper() < 1e25 True sage: K.<im> = QuadraticField(-1) sage: eps = 1/2^100 sage: x = polygen(K) sage: p = (x-1)*(x-1-eps)*(x-1+eps)*(x-1-eps*im)*(x-1+eps*im) This polynomial actually has all-real coefficients, and is very, very close to (x-1)^5:: sage: [RR(QQ(a)) for a in list(p - (x-1)^5)] [3.87259191484932e-121, -3.87259191484932e-121] sage: rts = complex_roots(p) sage: [ComplexIntervalField(10)(rt[0] - 1) for rt in rts] [-7.8887?e-31, 0, 7.8887?e-31, -7.8887?e-31*I, 7.8887?e-31*I] We can get roots either as intervals, or as elements of QQbar or AA. :: sage: p = (x^2 + x - 1) sage: p = p * p(x*im) sage: p -x^4 + (im - 1)*x^3 + im*x^2 + (-im - 1)*x + 1 Two of the roots have a zero real component; two have a zero imaginary component. These zero components will be found slightly inaccurately, and the exact values returned are very sensitive to the (non-portable) results of NumPy. So we post-process the roots for printing, to get predictable doctest results. :: sage: def tiny(x): ....: return x.contains_zero() and x.absolute_diameter() < 1e-14 sage: def smash(x): ....: x = CIF(x[0]) # discard multiplicity ....: if tiny(x.imag()): return x.real() ....: if tiny(x.real()): return CIF(0, x.imag()) sage: rts = complex_roots(p); type(rts[0][0]), sorted(map(smash, rts)) (<class 'sage.rings.complex_interval.ComplexIntervalFieldElement'>, [-1.618033988749895?, -0.618033988749895?*I, 1.618033988749895?*I, 0.618033988749895?]) sage: rts = complex_roots(p, retval='algebraic'); type(rts[0][0]), sorted(map(smash, rts)) (<class 'sage.rings.qqbar.AlgebraicNumber'>, [-1.618033988749895?, -0.618033988749895?*I, 1.618033988749895?*I, 0.618033988749895?]) sage: rts = complex_roots(p, retval='algebraic_real'); type(rts[0][0]), rts (<class 'sage.rings.qqbar.AlgebraicReal'>, [(-1.618033988749895?, 1), (0.618033988749895?, 1)]) TESTS: Verify that :trac:`12026` is fixed:: sage: f = matrix(QQ, 8, lambda i, j: 1/(i + j + 1)).charpoly() sage: from sage.rings.polynomial.complex_roots import complex_roots sage: len(complex_roots(f)) 8 """ if skip_squarefree: factors = [(p, 1)] else: factors = p.squarefree_decomposition() prec = 53 while True: CC = ComplexField(prec) CCX = CC['x'] all_rts = [] ok = True for (factor, exp) in factors: cfac = CCX(factor) rts = cfac.roots(multiplicities=False) # Make sure the number of roots we found is the degree. If # we don't find that many roots, it's because the # precision isn't big enough and though the (possibly # exact) polynomial "factor" is squarefree, it is not # squarefree as an element of CCX. if len(rts) < factor.degree(): ok = False break irts = interval_roots(factor, rts, max(prec, min_prec)) if irts is None: ok = False break if retval != 'interval': factor = QQbar.common_polynomial(factor) for irt in irts: all_rts.append((irt, factor, exp)) if ok and intervals_disjoint([rt for (rt, fac, mult) in all_rts]): all_rts = sort_complex_numbers_for_display(all_rts) if retval == 'interval': return [(rt, mult) for (rt, fac, mult) in all_rts] elif retval == 'algebraic': return [(QQbar.polynomial_root(fac, rt), mult) for (rt, fac, mult) in all_rts] elif retval == 'algebraic_real': rts = [] for (rt, fac, mult) in all_rts: qqbar_rt = QQbar.polynomial_root(fac, rt) if qqbar_rt.imag().is_zero(): rts.append((AA(qqbar_rt), mult)) return rts prec = prec * 2
def extract(cls, obj): """ Takes an object extracted by the json parser and decodes the special-formating dictionaries used to store special types. """ if isinstance(obj, dict) and "data" in obj: if len(obj) == 2 and "__ComplexList__" in obj: return [complex(*v) for v in obj["data"]] elif len(obj) == 2 and "__QQList__" in obj: assert SAGE_MODE return [QQ(tuple(v)) for v in obj["data"]] elif len(obj) == 3 and "__NFList__" in obj and "base" in obj: assert SAGE_MODE base = cls.extract(obj["base"]) return [cls._extract(base, c) for c in obj["data"]] elif len(obj) == 2 and "__IntDict__" in obj: if SAGE_MODE: return {Integer(k): cls.extract(v) for k, v in obj["data"]} else: return {int(k): cls.extract(v) for k, v in obj["data"]} elif len(obj) == 3 and "__Vector__" in obj and "base" in obj: assert SAGE_MODE base = cls.extract(obj["base"]) return vector([cls._extract(base, v) for v in obj["data"]]) elif len(obj) == 2 and "__Rational__" in obj: assert SAGE_MODE return Rational(*obj["data"]) elif len(obj) == 3 and "__RealLiteral__" in obj and "prec" in obj: assert SAGE_MODE return LmfdbRealLiteral(RealField(obj["prec"]), obj["data"]) elif len(obj) == 2 and "__complex__" in obj: return complex(*obj["data"]) elif len(obj) == 3 and "__Complex__" in obj and "prec" in obj: assert SAGE_MODE return ComplexNumber(ComplexField(obj["prec"]), *obj["data"]) elif len(obj) == 3 and "__NFElt__" in obj and "parent" in obj: assert SAGE_MODE return cls._extract(cls.extract(obj["parent"]), obj["data"]) elif (len(obj) == 3 and ("__NFRelative__" in obj or "__NFAbsolute__" in obj) and "vname" in obj): assert SAGE_MODE poly = cls.extract(obj["data"]) return NumberField(poly, name=obj["vname"]) elif len(obj) == 2 and "__NFCyclotomic__" in obj: assert SAGE_MODE return CyclotomicField(obj["data"]) elif len(obj) == 2 and "__IntegerRing__" in obj: assert SAGE_MODE return ZZ elif len(obj) == 2 and "__RationalField__" in obj: assert SAGE_MODE return QQ elif len( obj) == 3 and "__RationalPoly__" in obj and "vname" in obj: assert SAGE_MODE return QQ[obj["vname"]]([QQ(tuple(v)) for v in obj["data"]]) elif (len(obj) == 4 and "__Poly__" in obj and "vname" in obj and "base" in obj): assert SAGE_MODE base = cls.extract(obj["base"]) return base[obj["vname"]]( [cls._extract(base, c) for c in obj["data"]]) elif (len(obj) == 5 and "__PowerSeries__" in obj and "vname" in obj and "base" in obj and "prec" in obj): assert SAGE_MODE base = cls.extract(obj["base"]) prec = infinity if obj["prec"] == "inf" else int(obj["prec"]) return base[[obj["vname"] ]]([cls._extract(base, c) for c in obj["data"]], prec=prec) elif len(obj) == 2 and "__date__" in obj: return datetime.datetime.strptime(obj["data"], "%Y-%m-%d").date() elif len(obj) == 2 and "__time__" in obj: return datetime.datetime.strptime(obj["data"], "%H:%M:%S.%f").time() elif len(obj) == 2 and "__datetime__" in obj: return datetime.datetime.strptime(obj["data"], "%Y-%m-%d %H:%M:%S.%f") return obj
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)
def braid_in_segment(f, x0, x1): """ Return the braid formed by the `y` roots of ``f`` when `x` moves from ``x0`` to ``x1``. INPUT: - ``f`` -- a polynomial 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, 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,CC(p1),CC(p2)) # optional - sirocco sage: B.left_normal_form() # optional - sirocco (1, s5) """ CC = ComplexField(64) (x, y) = f.variables() I = QQbar.gen() X0 = QQ(x0.real()) + I * QQ(x0.imag()) X1 = QQ(x1.real()) + I * QQ(x1.imag()) F0 = QQbar[y](f(X0, y)) y0s = F0.roots(multiplicities=False) strands = [followstrand(f, x0, x1, CC(a)) for a in y0s] complexstrands = [[(a[0], CC(a[1], a[2])) for a in b] for b in strands] centralbraid = braid_from_piecewise(complexstrands) initialstrands = [] y0aps = [c[0][1] for c in complexstrands] used = [] for y0ap in y0aps: distances = [((y0ap - y0).norm(), y0) for y0 in y0s] y0 = sorted(distances)[0][1] if y0 in used: raise ValueError("different roots are too close") used.append(y0) initialstrands.append([(0, y0), (1, y0ap)]) initialbraid = braid_from_piecewise(initialstrands) F1 = QQbar[y](f(X1, y)) y1s = F1.roots(multiplicities=False) finalstrands = [] y1aps = [c[-1][1] for c in complexstrands] used = [] for y1ap in y1aps: distances = [((y1ap - y1).norm(), y1) for y1 in y1s] y1 = sorted(distances)[0][1] if y1 in used: raise ValueError("different roots are too close") used.append(y1) finalstrands.append([(0, y1ap), (1, y1)]) finallbraid = braid_from_piecewise(finalstrands) return initialbraid * centralbraid * finallbraid
def followstrand(f, 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`` -- 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 approximation) - 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)