Path: blob/master/src/sage/schemes/elliptic_curves/period_lattice.py
8820 views
# -*- coding: utf-8 -*-1r"""2Period lattices of elliptic curves and related functions34Let `E` be an elliptic curve defined over a number field `K`5(including `\QQ`). We attach a period lattice (a discrete rank 26subgroup of `\CC`) to each embedding of `K` into `\CC`.78In the case of real embeddings, the lattice is stable under complex9conjugation and is called a real lattice. These have two types:10rectangular, (the real curve has two connected components and positive11discriminant) or non-rectangular (one connected component, negative12discriminant).1314The periods are computed to arbitrary precision using the AGM (Gauss's15Arithmetic-Geometric Mean).1617EXAMPLES::1819sage: K.<a> = NumberField(x^3-2)20sage: E = EllipticCurve([0,1,0,a,a])2122First we try a real embedding::2324sage: emb = K.embeddings(RealField())[0]25sage: L = E.period_lattice(emb); L26Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:27From: Number Field in a with defining polynomial x^3 - 228To: Algebraic Real Field29Defn: a |--> 1.259921049894873?3031The first basis period is real::3233sage: L.basis()34(3.81452977217855, 1.90726488608927 + 1.34047785962440*I)35sage: L.is_real()36True3738For a basis `\omega_1,\omega_2` normalised so that `\omega_1/\omega_2`39is in the fundamental region of the upper half-plane, use the function40``normalised_basis()`` instead::4142sage: L.normalised_basis()43(1.90726488608927 - 1.34047785962440*I, -1.90726488608927 - 1.34047785962440*I)4445Next a complex embedding::4647sage: emb = K.embeddings(ComplexField())[0]48sage: L = E.period_lattice(emb); L49Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:50From: Number Field in a with defining polynomial x^3 - 251To: Algebraic Field52Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I5354In this case, the basis `\omega_1`, `\omega_2` is always normalised so55that `\tau = \omega_1/\omega_2` is in the fundamental region in the56upper half plane::5758sage: w1,w2 = L.basis(); w1,w259(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)60sage: L.is_real()61False62sage: tau = w1/w2; tau630.387694505032876 + 1.30821088214407*I64sage: L.normalised_basis()65(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)6667We test that bug #8415 (caused by a PARI bug fixed in v2.3.5) is OK::6869sage: E = EllipticCurve('37a')70sage: K.<a> = QuadraticField(-7)71sage: EK = E.change_ring(K)72sage: EK.period_lattice(K.complex_embeddings()[0])73Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 + (-1)*x over Number Field in a with defining polynomial x^2 + 7 with respect to the embedding Ring morphism:74From: Number Field in a with defining polynomial x^2 + 775To: Algebraic Field76Defn: a |--> -2.645751311064591?*I77787980AUTHORS:8182- ?: initial version.8384- John Cremona:8586- Adapted to handle real embeddings of number fields, September 2008.8788- Added basis_matrix function, November 20088990- Added support for complex embeddings, May 2009.9192- Added complex elliptic logs, March 2010; enhanced, October 2010.9394"""9596from sage.modules.free_module import FreeModule_generic_pid97from sage.rings.all import ZZ, QQ, RealField, ComplexField, QQbar, AA98from sage.rings.real_mpfr import is_RealField99from sage.rings.complex_field import is_ComplexField100from sage.rings.real_mpfr import RealNumber as RealNumber101from sage.rings.complex_number import ComplexNumber as ComplexNumber102from sage.rings.number_field.number_field import refine_embedding103from sage.rings.infinity import Infinity104from sage.schemes.elliptic_curves.constructor import EllipticCurve105from sage.misc.cachefunc import cached_method106107class PeriodLattice(FreeModule_generic_pid):108"""109The class for the period lattice of an algebraic variety.110"""111pass112113class PeriodLattice_ell(PeriodLattice):114r"""115The class for the period lattice of an elliptic curve.116117Currently supported are elliptic curves defined over `\QQ`, and118elliptic curves defined over a number field with a real or complex119embedding, where the lattice constructed depends on that120embedding.121"""122123def __init__(self, E, embedding=None):124r"""125Initialises the period lattice by storing the elliptic curve and the embedding.126127INPUT:128129- ``E`` -- an elliptic curve130131- ``embedding`` (defult: ``None``) -- an embedding of the base132field `K` of ``E`` into a real or complex field. If133``None``:134135- use the built-in coercion to `\RR` for `K=\QQ`;136137- use the first embedding into `\RR` given by138``K.embeddings(RealField())``, if there are any;139140- use the first embedding into `\CC` given by141``K.embeddings(ComplexField())``, if `K` is totally complex.142143.. note::144145No periods are computed on creation of the lattice; see the146functions ``basis()``, ``normalised_basis()`` and147``real_period()`` for precision setting.148149EXAMPLES:150151This function is not normally called directly, but will be152called by the period_lattice() function of classes153ell_number_field and ell_rational_field::154155sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell156sage: E = EllipticCurve('37a')157sage: PeriodLattice_ell(E)158Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field159160::161162sage: K.<a> = NumberField(x^3-2)163sage: emb = K.embeddings(RealField())[0]164sage: E = EllipticCurve([0,1,0,a,a])165sage: L = PeriodLattice_ell(E,emb); L166Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:167From: Number Field in a with defining polynomial x^3 - 2168To: Algebraic Real Field169Defn: a |--> 1.259921049894873?170171sage: emb = K.embeddings(ComplexField())[0]172sage: L = PeriodLattice_ell(E,emb); L173Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:174From: Number Field in a with defining polynomial x^3 - 2175To: Algebraic Field176Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I177178TESTS::179180sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell181sage: K.<a> = NumberField(x^3-2)182sage: emb = K.embeddings(RealField())[0]183sage: E = EllipticCurve([0,1,0,a,a])184sage: L = PeriodLattice_ell(E,emb)185sage: L == loads(dumps(L))186True187"""188# First we cache the elliptic curve with this period lattice:189190self.E = E191192# Next we cache the embedding into QQbar or AA which extends193# the given embedding:194195K = E.base_field()196if embedding is None:197embs = K.embeddings(AA)198real = len(embs)>0199if not real:200embs = K.embeddings(QQbar)201embedding = embs[0]202else:203embedding = refine_embedding(embedding,Infinity)204real = embedding(K.gen()).imag().is_zero()205206self.embedding = embedding207208# Next we compute and cache (in self.real_flag) the type of209# the lattice: +1 for real rectangular, -1 for real210# non-rectangular, 0 for non-real:211212self.real_flag = 0213if real:214self.real_flag = +1215if embedding(E.discriminant())<0:216self.real_flag = -1217218# The following algebraic data associated to E and the219# embedding is cached:220#221# Ebar: the curve E base-changed to QQbar (or AA)222# f2: the 2-division polynomial of Ebar223# ei: the roots e1, e2, e3 of f2, as elements of QQbar (or AA)224#225# The ei are used both for period computation and elliptic226# logarithms.227228self.Ebar = self.E.change_ring(self.embedding)229self.f2 = self.Ebar.two_division_polynomial()230if self.real_flag == 1: # positive discriminant231self._ei = self.f2.roots(AA,multiplicities=False)232self._ei.sort() # e1 < e2 < e3233e1, e2, e3 = self._ei234elif self.real_flag == -1: # negative discriminant235self._ei = self.f2.roots(QQbar,multiplicities=False)236self._ei = list(sorted(self._ei,key=lambda z: z.imag()))237e1, e3, e2 = self._ei # so e3 is real238e3 = AA(e3)239self._ei = [e1, e2, e3]240else:241self._ei = self.f2.roots(QQbar,multiplicities=False)242e1, e2, e3 = self._ei243244# The quantities sqrt(e_i-e_j) are cached (as elements of245# QQbar) to be used in period computations:246247self._abc = (e3-e1).sqrt(), (e3-e2).sqrt(), (e2-e1).sqrt()248249PeriodLattice.__init__(self, base_ring=ZZ, rank=2, degree=1, sparse=False)250251def __cmp__(self, other):252r"""253Comparison function for period lattices254255TESTS::256257sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell258sage: K.<a> = NumberField(x^3-2)259sage: E = EllipticCurve([0,1,0,a,a])260sage: embs = K.embeddings(ComplexField())261sage: L1,L2,L3 = [PeriodLattice_ell(E,e) for e in embs]262sage: L1 < L2 < L3263True264265"""266if not isinstance(other, PeriodLattice_ell): return -1267t = cmp(self.E, other.E)268if t: return t269a = self.E.base_field().gen()270return cmp(self.embedding(a), other.embedding(a))271272def __repr__(self):273"""274Returns the string representation of this period lattice.275276EXAMPLES::277278sage: E = EllipticCurve('37a')279sage: E.period_lattice()280Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field281282::283284sage: K.<a> = NumberField(x^3-2)285sage: emb = K.embeddings(RealField())[0]286sage: E = EllipticCurve([0,1,0,a,a])287sage: L = E.period_lattice(emb); L288Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:289From: Number Field in a with defining polynomial x^3 - 2290To: Algebraic Real Field291Defn: a |--> 1.259921049894873?292"""293if self.E.base_field() is QQ:294return "Period lattice associated to %s"%(self.E)295else:296return "Period lattice associated to %s with respect to the embedding %s"%(self.E, self.embedding)297298def __call__(self, P, prec=None):299r"""300Return the elliptic logarithm of a point `P`.301302INPUT:303304- ``P`` (point) -- a point on the elliptic curve associated305with this period lattice.306307- ``prec`` (default: ``None``) -- precision in bits (default308precision if ``None``).309310OUTPUT:311312(complex number) The elliptic logarithm of the point `P` with313respect to this period lattice. If `E` is the elliptic curve314and `\sigma:K\to\CC` the embedding, then the returned value `z`315is such that `z\pmod{L}` maps to `\sigma(P)` under the316standard Weierstrass isomorphism from `\CC/L` to `\sigma(E)`.317318EXAMPLES::319320sage: E = EllipticCurve('389a')321sage: L = E.period_lattice()322sage: E.discriminant() > 0323True324sage: L.real_flag3251326sage: P = E([-1,1])327sage: P.is_on_identity_component ()328False329sage: L(P, prec=96)3300.4793482501902193161295330101 + 0.985868850775824102211203849...*I331sage: Q=E([3,5])332sage: Q.is_on_identity_component()333True334sage: L(Q, prec=96)3351.931128271542559442488585220336337Note that this is actually the inverse of the Weierstrass isomorphism::338339sage: L.elliptic_exponential(L(Q))340(3.00000000000000 : 5.00000000000000 : 1.00000000000000)341342An example with negative discriminant, and a torsion point::343344sage: E = EllipticCurve('11a1')345sage: L = E.period_lattice()346sage: E.discriminant() < 0347True348sage: L.real_flag349-1350sage: P = E([16,-61])351sage: L(P)3520.253841860855911353sage: L.real_period() / L(P)3545.00000000000000355"""356return self.elliptic_logarithm(P,prec)357358@cached_method359def basis(self, prec=None, algorithm='sage'):360r"""361Return a basis for this period lattice as a 2-tuple.362363INPUT:364365- ``prec`` (default: ``None``) -- precision in bits (default366precision if ``None``).367368- ``algorithm`` (string, default 'sage') -- choice of369implementation (for real embeddings only) between 'sage'370(native Sage implementation) or 'pari' (use the PARI371library: only available for real embeddings).372373OUTPUT:374375(tuple of Complex) `(\omega_1,\omega_2)` where the lattice is376`\ZZ\omega_1 + \ZZ\omega_2`. If the lattice is real then377`\omega_1` is real and positive, `\Im(\omega_2)>0` and378`\Re(\omega_1/\omega_2)` is either `0` (for rectangular379lattices) or `\frac{1}{2}` (for non-rectangular lattices).380Otherwise, `\omega_1/\omega_2` is in the fundamental region of381the upper half-plane. If the latter normalisation is required382for real lattices, use the function ``normalised_basis()``383instead.384385EXAMPLES::386387sage: E = EllipticCurve('37a')388sage: E.period_lattice().basis()389(2.99345864623196, 2.45138938198679*I)390391This shows that the issue reported at trac \#3954 is fixed::392393sage: E = EllipticCurve('37a')394sage: b1 = E.period_lattice().basis(prec=30)395sage: b2 = E.period_lattice().basis(prec=30)396sage: b1 == b2397True398399This shows that the issue reported at trac \#4064 is fixed::400401sage: E = EllipticCurve('37a')402sage: E.period_lattice().basis(prec=30)[0].parent()403Real Field with 30 bits of precision404sage: E.period_lattice().basis(prec=100)[0].parent()405Real Field with 100 bits of precision406407::408409sage: K.<a> = NumberField(x^3-2)410sage: emb = K.embeddings(RealField())[0]411sage: E = EllipticCurve([0,1,0,a,a])412sage: L = E.period_lattice(emb)413sage: L.basis(64)414(3.81452977217854509, 1.90726488608927255 + 1.34047785962440202*I)415416sage: emb = K.embeddings(ComplexField())[0]417sage: L = E.period_lattice(emb)418sage: w1,w2 = L.basis(); w1,w2419(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)420sage: L.is_real()421False422sage: tau = w1/w2; tau4230.387694505032876 + 1.30821088214407*I424"""425# We divide into two cases: (1) Q, or a number field with a426# real embedding; (2) a number field with a complex embedding.427# In each case the periods are computed by a different428# internal function.429430if self.is_real():431return self._compute_periods_real(prec=prec, algorithm=algorithm)432else:433return self._compute_periods_complex(prec=prec)434435@cached_method436def normalised_basis(self, prec=None, algorithm='sage'):437r"""438Return a normalised basis for this period lattice as a 2-tuple.439440INPUT:441442- ``prec`` (default: ``None``) -- precision in bits (default443precision if ``None``).444445- ``algorithm`` (string, default 'sage') -- choice of446implementation (for real embeddings only) between 'sage'447(native Sage implementation) or 'pari' (use the PARI448library: only available for real embeddings).449450OUTPUT:451452(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has453the form `\ZZ\omega_1 + \ZZ\omega_2`. The basis is normalised454so that `\omega_1/\omega_2` is in the fundamental region of455the upper half-plane. For an alternative normalisation for456real lattices (with the first period real), use the function457basis() instead.458459EXAMPLES::460461sage: E = EllipticCurve('37a')462sage: E.period_lattice().normalised_basis()463(2.99345864623196, -2.45138938198679*I)464465::466467sage: K.<a> = NumberField(x^3-2)468sage: emb = K.embeddings(RealField())[0]469sage: E = EllipticCurve([0,1,0,a,a])470sage: L = E.period_lattice(emb)471sage: L.normalised_basis(64)472(1.90726488608927255 - 1.34047785962440202*I, -1.90726488608927255 - 1.34047785962440202*I)473474sage: emb = K.embeddings(ComplexField())[0]475sage: L = E.period_lattice(emb)476sage: w1,w2 = L.normalised_basis(); w1,w2477(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)478sage: L.is_real()479False480sage: tau = w1/w2; tau4810.387694505032876 + 1.30821088214407*I482"""483w1, w2 = periods = self.basis(prec=prec, algorithm=algorithm)484periods, mat = normalise_periods(w1,w2)485return periods486487@cached_method488def _compute_periods_real(self, prec=None, algorithm='sage'):489r"""490Internal function to compute the periods (real embedding case).491492INPUT:493494495- `prec` (int or ``None`` (default)) -- floating point496precision (in bits); if None, use the default precision.497498- `algorithm` (string, default 'sage') -- choice of implementation between499- `pari`: use the PARI library500501- `sage`: use a native Sage implementation (with the same underlying algorithm).502503504OUTPUT:505506(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has507the form `\ZZ\omega_1 + \ZZ\omega_2`, `\omega_1` is real and508`\omega_1/\omega_2` has real part either `0` or `frac{1}{2}`.509510EXAMPLES::511512sage: K.<a> = NumberField(x^3-2)513sage: E = EllipticCurve([0,1,0,a,a])514sage: embs = K.embeddings(CC)515sage: Ls = [E.period_lattice(e) for e in embs]516sage: [L.is_real() for L in Ls]517[False, False, True]518sage: Ls[2]._compute_periods_real(100)519(3.8145297721785450936365098936,5201.9072648860892725468182549468 + 1.3404778596244020196600112394*I)521sage: Ls[2]._compute_periods_real(100, algorithm='pari')522(3.8145297721785450936365098936,5231.9072648860892725468182549468 - 1.3404778596244020196600112394*I)524"""525if prec is None:526prec = RealField().precision()527R = RealField(prec)528C = ComplexField(prec)529530if algorithm=='pari':531if self.E.base_field() is QQ:532periods = self.E.pari_curve(prec).omega().python()533return (R(periods[0]), C(periods[1]))534535from sage.libs.pari.all import pari536E_pari = pari([R(self.embedding(ai).real()) for ai in self.E.a_invariants()]).ellinit(precision=prec)537periods = E_pari.omega().python()538return (R(periods[0]), C(periods[1]))539540if algorithm!='sage':541raise ValueError, "invalid value of 'algorithm' parameter"542543pi = R.pi()544# Up to now everything has been exact in AA or QQbar, but now545# we must go transcendental. Only now is the desired546# precision used!547if self.real_flag == 1: # positive discriminant548a, b, c = (R(x) for x in self._abc)549w1 = R(pi/a.agm(b)) # least real period550w2 = C(0,pi/a.agm(c)) # least pure imaginary period551else:552a = C(self._abc[0])553x, y, r = a.real().abs(), a.imag().abs(), a.abs()554w1 = R(pi/r.agm(x)) # least real period555w2 = R(pi/r.agm(y)) # least pure imaginary period /i556w2 = C(w1,w2)/2557558return (w1,w2)559560@cached_method561def _compute_periods_complex(self, prec=None, normalise=True):562r"""563Internal function to compute the periods (complex embedding case).564565INPUT:566567- `prec` (int or ``None`` (default)) -- floating point precision (in bits); if None,568use the default precision.569570- `normalise` (bool, default True) -- whether to normalise the571basis after computation.572573OUTPUT:574575(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has576the form `\ZZ\omega_1 + \ZZ\omega_2`. If `normalise` is577`True`, the basis is normalised so that `(\omega_1/\omega_2)`578is in the fundamental region of the upper half plane.579580EXAMPLES::581582sage: K.<a> = NumberField(x^3-2)583sage: E = EllipticCurve([0,1,0,a,a])584sage: embs = K.embeddings(CC)585sage: Ls = [E.period_lattice(e) for e in embs]586sage: [L.is_real() for L in Ls]587[False, False, True]588sage: L = Ls[0]589sage: w1,w2 = L._compute_periods_complex(100); w1,w2590(-1.3758860416607626645495991458 - 2.5856094662444337042877901304*I, -2.1033990784735587243397865076 + 0.42837877646062187766760569686*I)591sage: tau = w1/w2; tau5920.38769450503287609349437509561 + 1.3082108821440725664008561928*I593sage: tau.real()5940.38769450503287609349437509561595sage: tau.abs()5961.3644496111593345713923386773597598Without normalisation::599600sage: w1,w2 = L._compute_periods_complex(normalise=False); w1,w2601(2.10339907847356 - 0.428378776460622*I, 0.727513036812796 - 3.01398824270506*I)602sage: tau = w1/w2; tau6030.293483964608883 + 0.627038168678760*I604sage: tau.real()6050.293483964608883606sage: tau.abs() # > 16070.692321964451917608"""609if prec is None:610prec = RealField().precision()611C = ComplexField(prec)612613# Up to now everything has been exact in AA, but now we614# must go transcendental. Only now is the desired615# precision used!616pi = C.pi()617a, b, c = (C(x) for x in self._abc)618if (a+b).abs() < (a-b).abs(): b=-b619if (a+c).abs() < (a-c).abs(): c=-c620w1 = pi/a.agm(b)621w2 = pi*C.gen()/a.agm(c)622if (w1/w2).imag()<0: w2=-w2623if normalise:624w1w2, mat = normalise_periods(w1,w2)625return w1w2626return (w1,w2)627628def is_real(self):629r"""630Return True if this period lattice is real.631632EXAMPLES::633634sage: f = EllipticCurve('11a')635sage: f.period_lattice().is_real()636True637638::639640sage: K.<i> = QuadraticField(-1)641sage: E = EllipticCurve(K,[0,0,0,i,2*i])642sage: emb = K.embeddings(ComplexField())[0]643sage: L = E.period_lattice(emb)644sage: L.is_real()645False646647::648649sage: K.<a> = NumberField(x^3-2)650sage: E = EllipticCurve([0,1,0,a,a])651sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)]652[False, False, True]653654655ALGORITHM:656657The lattice is real if it is associated to a real embedding;658such lattices are stable under conjugation.659"""660return self.real_flag!=0661662def is_rectangular(self):663r"""664Return True if this period lattice is rectangular.665666.. note::667668Only defined for real lattices; a RuntimeError is raised for669non-real lattices.670671EXAMPLES::672673sage: f = EllipticCurve('11a')674sage: f.period_lattice().basis()675(1.26920930427955, 0.634604652139777 + 1.45881661693850*I)676sage: f.period_lattice().is_rectangular()677False678679::680681sage: f = EllipticCurve('37b')682sage: f.period_lattice().basis()683(1.08852159290423, 1.76761067023379*I)684sage: f.period_lattice().is_rectangular()685True686687ALGORITHM:688689The period lattice is rectangular precisely if the690discriminant of the Weierstrass equation is positive, or691equivalently if the number of real components is 2.692"""693if self.is_real():694return self.real_flag == +1695raise RuntimeError, "Not defined for non-real lattices."696697def real_period(self, prec = None, algorithm='sage'):698"""699Returns the real period of this period lattice.700701INPUT:702703- ``prec`` (int or ``None`` (default)) -- real precision in704bits (default real precision if ``None``)705706- ``algorithm`` (string, default 'sage') -- choice of707implementation (for real embeddings only) between 'sage'708(native Sage implementation) or 'pari' (use the PARI709library: only available for real embeddings).710711.. note::712713Only defined for real lattices; a RuntimeError is raised for714non-real lattices.715716EXAMPLES::717718sage: E = EllipticCurve('37a')719sage: E.period_lattice().real_period()7202.99345864623196721722::723724sage: K.<a> = NumberField(x^3-2)725sage: emb = K.embeddings(RealField())[0]726sage: E = EllipticCurve([0,1,0,a,a])727sage: L = E.period_lattice(emb)728sage: L.real_period(64)7293.81452977217854509730"""731if self.is_real():732return self.basis(prec,algorithm)[0]733raise RuntimeError, "Not defined for non-real lattices."734735def omega(self, prec = None):736r"""737Returns the real or complex volume of this period lattice.738739INPUT:740741- ``prec`` (int or ``None``(default)) -- real precision in742bits (default real precision if ``None``)743744OUTPUT:745746(real) For real lattices, this is the real period times the747number of connected components. For non-real lattices it is748the complex area.749750.. note::751752If the curve is defined over `\QQ` and is given by a753*minimal* Weierstrass equation, then this is the correct754period in the BSD conjecture, i.e., it is the least real755period * 2 when the period lattice is rectangular. More756generally the product of this quantity over all embeddings757appears in the generalised BSD formula.758759760EXAMPLES::761762sage: E = EllipticCurve('37a')763sage: E.period_lattice().omega()7645.98691729246392765766This is not a minimal model::767768sage: E = EllipticCurve([0,-432*6^2])769sage: E.period_lattice().omega()7700.486109385710056771772If you were to plug the above omega into the BSD conjecture, you773would get nonsense. The following works though::774775sage: F = E.minimal_model()776sage: F.period_lattice().omega()7770.972218771420113778779::780781sage: K.<a> = NumberField(x^3-2)782sage: emb = K.embeddings(RealField())[0]783sage: E = EllipticCurve([0,1,0,a,a])784sage: L = E.period_lattice(emb)785sage: L.omega(64)7863.81452977217854509787788A complex example (taken from J.E.Cremona and E.Whitley,789*Periods of cusp forms and elliptic curves over imaginary790quadratic fields*, Mathematics of Computation 62 No. 205791(1994), 407-429)::792793sage: K.<i> = QuadraticField(-1)794sage: E = EllipticCurve([0,1-i,i,-i,0])795sage: L = E.period_lattice(K.embeddings(CC)[0])796sage: L.omega()7978.80694160502647798"""799if self.is_real():800n_components = (self.real_flag+3)//2801return self.real_period(prec) * n_components802else:803return self.complex_area()804805@cached_method806def basis_matrix(self, prec=None, normalised=False):807r"""808Return the basis matrix of this period lattice.809810INPUT:811812- ``prec`` (int or ``None``(default)) -- real precision in813bits (default real precision if ``None``).814815- ``normalised`` (bool, default None) -- if True and the816embedding is real, use the normalised basis (see817``normalised_basis()``) instead of the default.818819OUTPUT:820821A 2x2 real matrix whose rows are the lattice basis vectors,822after identifying `\CC` with `\RR^2`.823824EXAMPLES::825826sage: E = EllipticCurve('37a')827sage: E.period_lattice().basis_matrix()828[ 2.99345864623196 0.000000000000000]829[0.000000000000000 2.45138938198679]830831::832833sage: K.<a> = NumberField(x^3-2)834sage: emb = K.embeddings(RealField())[0]835sage: E = EllipticCurve([0,1,0,a,a])836sage: L = E.period_lattice(emb)837sage: L.basis_matrix(64)838[ 3.81452977217854509 0.000000000000000000]839[ 1.90726488608927255 1.34047785962440202]840841See \#4388::842843sage: L = EllipticCurve('11a1').period_lattice()844sage: L.basis_matrix()845[ 1.26920930427955 0.000000000000000]846[0.634604652139777 1.45881661693850]847sage: L.basis_matrix(normalised=True)848[0.634604652139777 -1.45881661693850]849[-1.26920930427955 0.000000000000000]850851::852853sage: L = EllipticCurve('389a1').period_lattice()854sage: L.basis_matrix()855[ 2.49021256085505 0.000000000000000]856[0.000000000000000 1.97173770155165]857sage: L.basis_matrix(normalised=True)858[ 2.49021256085505 0.000000000000000]859[0.000000000000000 -1.97173770155165]860"""861from sage.matrix.all import Matrix862863if normalised:864return Matrix([list(w) for w in self.normalised_basis(prec)])865866w1,w2 = self.basis(prec)867if self.is_real():868return Matrix([[w1,0],list(w2)])869else:870return Matrix([list(w) for w in w1,w2])871872def complex_area(self, prec=None):873"""874Return the area of a fundamental domain for the period lattice875of the elliptic curve.876877INPUT:878879- ``prec`` (int or ``None``(default)) -- real precision in880bits (default real precision if ``None``).881882EXAMPLES::883884sage: E = EllipticCurve('37a')885sage: E.period_lattice().complex_area()8867.33813274078958887888::889890sage: K.<a> = NumberField(x^3-2)891sage: embs = K.embeddings(ComplexField())892sage: E = EllipticCurve([0,1,0,a,a])893sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)]894[False, False, True]895sage: [E.period_lattice(emb).complex_area() for emb in embs]896[6.02796894766694, 6.02796894766694, 5.11329270448345]897"""898w1,w2 = self.basis(prec)899return (w1*w2.conjugate()).imag().abs()900901def sigma(self, z, prec = None, flag=0):902r"""903Returns the value of the Weierstrass sigma function for this elliptic curve period lattice.904905INPUT:906907- ``z`` -- a complex number908909- ``prec`` (default: ``None``) -- real precision in bits910(default real precision if None).911912- ``flag`` --9139140: (default) ???;9159161: computes an arbitrary determination of log(sigma(z))9179182, 3: same using the product expansion instead of theta series. ???919920.. note::921922The reason for the ???'s above, is that the PARI923documentation for ellsigma is very vague. Also this is924only implemented for curves defined over `\QQ`.925926TODO:927928This function does not use any of the PeriodLattice functions929and so should be moved to ell_rational_field.930931EXAMPLES::932933sage: EllipticCurve('389a1').period_lattice().sigma(CC(2,1))9342.60912163570108 - 0.200865080824587*I935"""936if prec is None:937prec = RealField().precision()938try:939return self.E.pari_curve(prec).ellsigma(z, flag)940except AttributeError:941raise NotImplementedError, "sigma function not yet implemented for period lattices of curves not defined over Q."942943def curve(self):944r"""945Return the elliptic curve associated with this period lattice.946947EXAMPLES::948949sage: E = EllipticCurve('37a')950sage: L = E.period_lattice()951sage: L.curve() is E952True953954::955956sage: K.<a> = NumberField(x^3-2)957sage: E = EllipticCurve([0,1,0,a,a])958sage: L = E.period_lattice(K.embeddings(RealField())[0])959sage: L.curve() is E960True961962sage: L = E.period_lattice(K.embeddings(ComplexField())[0])963sage: L.curve() is E964True965"""966return self.E967968def ei(self):969r"""970Return the x-coordinates of the 2-division points of the elliptic curve associated with this period lattice, as elements of QQbar.971972EXAMPLES::973974sage: E = EllipticCurve('37a')975sage: L = E.period_lattice()976sage: L.ei()977[-1.107159871688768?, 0.2695944364054446?, 0.8375654352833230?]978979::980981sage: K.<a> = NumberField(x^3-2)982sage: E = EllipticCurve([0,1,0,a,a])983sage: L = E.period_lattice(K.embeddings(RealField())[0])984sage: L.ei()985[0.?e-19 - 1.122462048309373?*I, 0.?e-19 + 1.122462048309373?*I, -1]986987sage: L = E.period_lattice(K.embeddings(ComplexField())[0])988sage: L.ei()989[-1.000000000000000? + 0.?e-1...*I,990-0.9720806486198328? - 0.561231024154687?*I,9910.9720806486198328? + 0.561231024154687?*I]992"""993return self._ei994995def coordinates(self, z, rounding=None):996r"""997Returns the coordinates of a complex number w.r.t. the lattice basis998999INPUT:10001001- ``z`` (complex) -- A complex number.10021003- ``rounding`` (default ``None``) -- whether and how to round the1004output (see below).10051006OUTPUT:10071008When ``rounding`` is ``None`` (the default), returns a tuple1009of reals `x`, `y` such that `z=xw_1+yw_2` where `w_1`, `w_2`1010are a basis for the lattice (normalised in the case of complex1011embeddings).10121013When ``rounding`` is 'round', returns a tuple of integers `n_1`,1014`n_2` which are the closest integers to the `x`, `y` defined1015above. If `z` is in the lattice these are the coordinates of1016`z` with respect to the lattice basis.10171018When ``rounding`` is 'floor', returns a tuple of integers1019`n_1`, `n_2` which are the integer parts to the `x`, `y`1020defined above. These are used in :meth:``.reduce``10211022EXAMPLES::10231024sage: E = EllipticCurve('389a')1025sage: L = E.period_lattice()1026sage: w1, w2 = L.basis(prec=100)1027sage: P = E([-1,1])1028sage: zP = P.elliptic_logarithm(precision=100); zP10290.47934825019021931612953301006 + 0.98586885077582410221120384908*I1030sage: L.coordinates(zP)1031(0.19249290511394227352563996419, 0.50000000000000000000000000000)1032sage: sum([x*w for x,w in zip(L.coordinates(zP), L.basis(prec=100))])10330.47934825019021931612953301006 + 0.98586885077582410221120384908*I10341035sage: L.coordinates(12*w1+23*w2)1036(12.000000000000000000000000000, 23.000000000000000000000000000)1037sage: L.coordinates(12*w1+23*w2, rounding='floor')1038(11, 22)1039sage: L.coordinates(12*w1+23*w2, rounding='round')1040(12, 23)1041"""1042C = z.parent()1043if is_RealField(C):1044C = ComplexField(C.precision())1045z = C(z)1046else:1047if is_ComplexField(C):1048pass1049else:1050try:1051C = ComplexField()1052z = C(z)1053except TypeError:1054raise TypeError, "%s is not a complex number"%z1055prec = C.precision()1056from sage.matrix.all import Matrix1057from sage.modules.all import vector1058if self.real_flag:1059w1,w2 = self.basis(prec)1060M = Matrix([[w1,0], list(w2)])**(-1)1061else:1062w1,w2 = self.normalised_basis(prec)1063M = Matrix([list(w1), list(w2)])**(-1)1064u,v = vector(z)*M1065# Now z = u*w1+v*w21066if rounding=='round':1067return u.round(), v.round()1068if rounding=='floor':1069return u.floor(), v.floor()1070return u,v10711072def reduce(self, z):1073r"""1074Reduce a complex number modulo the lattice10751076INPUT:10771078- ``z`` (complex) -- A complex number.10791080OUTPUT:10811082(complex) the reduction of `z` modulo the lattice, lying in1083the fundamental period parallelogram with respect to the1084lattice basis. For curves defined over the reals (i.e. real1085embeddings) the output will be real when possible.10861087EXAMPLES::10881089sage: E = EllipticCurve('389a')1090sage: L = E.period_lattice()1091sage: w1, w2 = L.basis(prec=100)1092sage: P = E([-1,1])1093sage: zP = P.elliptic_logarithm(precision=100); zP10940.47934825019021931612953301006 + 0.98586885077582410221120384908*I1095sage: z = zP+10*w1-20*w2; z109625.381473858740770069343110929 - 38.448885180257139986236950114*I1097sage: L.reduce(z)10980.47934825019021931612953301006 + 0.98586885077582410221120384908*I1099sage: L.elliptic_logarithm(2*P)11000.9586965003804391101sage: L.reduce(L.elliptic_logarithm(2*P))11020.9586965003804391103sage: L.reduce(L.elliptic_logarithm(2*P)+10*w1-20*w2)11040.9586965003804441105"""1106C = z.parent()1107z_is_real = False1108if is_RealField(C):1109z_is_real = True1110C = ComplexField(C.precision())1111z = C(z)1112else:1113if is_ComplexField(C):1114z_is_real = z.is_real()1115else:1116try:1117C = ComplexField()1118z = C(z)1119z_is_real = z.is_real()1120except TypeError:1121raise TypeError, "%s is not a complex number"%z1122prec = C.precision()1123if self.real_flag:1124w1,w2 = self.basis(prec) # w1 real1125else:1126w1,w2 = self.normalised_basis(prec)1127# print "z = ",z1128# print "w1,w2 = ",w1,w21129u,v = self.coordinates(z, rounding='floor')1130# print "u,v = ",u,v1131z = z-u*w1-v*w211321133# Final adjustments for the real case.11341135# NB We assume here that when the embedding is real then the1136# point is also real!11371138if self.real_flag == 0: return z1139if self.real_flag == -1:1140k = (z.imag()/w2.imag()).round()1141z = z-k*w21142return C(z.real(),0)11431144if ((2*z.imag()/w2.imag()).round())%2:1145return C(z.real(),w2.imag()/2)1146else:1147return C(z.real(),0)11481149def elliptic_logarithm(self, P, prec=None, reduce=True):1150r"""1151Return the elliptic logarithm of a point.11521153INPUT:11541155- ``P`` (point) -- A point on the elliptic curve associated1156with this period lattice.11571158- ``prec`` (default: ``None``) -- real precision in bits1159(default real precision if None).11601161- ``reduce`` (default: ``True``) -- if ``True``, the result1162is reduced with respect to the period lattice basis.11631164OUTPUT:11651166(complex number) The elliptic logarithm of the point `P` with1167respect to this period lattice. If `E` is the elliptic curve1168and `\sigma:K\to\CC` the embedding, the the returned value `z`1169is such that `z\pmod{L}` maps to `\sigma(P)` under the1170standard Weierstrass isomorphism from `\CC/L` to `\sigma(E)`.1171If ``reduce`` is ``True``, the output is reduced so that it is1172in the fundamental period parallelogram with respect to the1173normalised lattice basis.11741175ALGORITHM: Uses the complex AGM. See [Cremona2010]_ for details.11761177.. [Cremona2010] J. E. Cremona and T. Thongjunthug, The1178Complex AGM, periods of elliptic curves over $\CC$ and1179complex elliptic logarithms. Preprint 2010.11801181EXAMPLES::11821183sage: E = EllipticCurve('389a')1184sage: L = E.period_lattice()1185sage: E.discriminant() > 01186True1187sage: L.real_flag118811189sage: P = E([-1,1])1190sage: P.is_on_identity_component ()1191False1192sage: L.elliptic_logarithm(P, prec=96)11930.4793482501902193161295330101 + 0.9858688507758241022112038491*I1194sage: Q=E([3,5])1195sage: Q.is_on_identity_component()1196True1197sage: L.elliptic_logarithm(Q, prec=96)11981.93112827154255944248858522011991200Note that this is actually the inverse of the Weierstrass isomorphism::12011202sage: L.elliptic_exponential(_)1203(3.00000000000000000000000000... : 5.00000000000000000000000000... : 1.000000000000000000000000000)12041205An example with negative discriminant, and a torsion point::12061207sage: E = EllipticCurve('11a1')1208sage: L = E.period_lattice()1209sage: E.discriminant() < 01210True1211sage: L.real_flag1212-11213sage: P = E([16,-61])1214sage: L.elliptic_logarithm(P)12150.2538418608559111216sage: L.real_period() / L.elliptic_logarithm(P)12175.0000000000000012181219An example where precision is problematic::12201221sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g11222sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1])1223sage: L = E.period_lattice()1224sage: L.ei()1225[5334.003952567705? - 1.964393150436?e-6*I, 5334.003952567705? + 1.964393150436?e-6*I, -10668.25790513541?]1226sage: L.elliptic_logarithm(P,prec=100)12270.2765620401410706146407620309712281229Some complex examples, taken from the paper by Cremona and Thongjunthug::12301231sage: K.<i> = QuadraticField(-1)1232sage: a4 = 9*i-101233sage: a6 = 21-i1234sage: E = EllipticCurve([0,0,0,a4,a6])1235sage: e1 = 3-2*i; e2 = 1+i; e3 = -4+i1236sage: emb = K.embeddings(CC)[1]1237sage: L = E.period_lattice(emb)1238sage: P = E(2-i,4+2*i)12391240By default, the output is reduced with respect to the1241normalised lattice basis, so that its coordinates with respect1242to that basis lie in the interval [0,1)::12431244sage: z = L.elliptic_logarithm(P,prec=100); z12450.70448375537782208460499649302 - 0.79246725643650979858266018068*I1246sage: L.coordinates(z)1247(0.46247636364807931766105406092, 0.79497588726808704200760395829)12481249Using ``reduce=False`` this step can be omitted. In this case1250the coordinates are usually in the interval [-0.5,0.5), but1251this is not guaranteed. This option is mainly for testing1252purposes::12531254sage: z = L.elliptic_logarithm(P,prec=100, reduce=False); z12550.57002153834710752778063503023 + 0.46476340520469798857457031393*I1256sage: L.coordinates(z)1257(0.46247636364807931766105406092, -0.20502411273191295799239604171)12581259The elliptic logs of the 2-torsion points are half-periods::12601261sage: L.elliptic_logarithm(E(e1,0),prec=100)12620.64607575874356525952487867052 + 0.22379609053909448304176885364*I1263sage: L.elliptic_logarithm(E(e2,0),prec=100)12640.71330686725892253793705940192 - 0.40481924028150941053684639367*I1265sage: L.elliptic_logarithm(E(e3,0),prec=100)12660.067231108515357278412180731396 - 0.62861533082060389357861524731*I12671268We check this by doubling and seeing that the resulting1269coordinates are integers::12701271sage: L.coordinates(2*L.elliptic_logarithm(E(e1,0),prec=100))1272(1.0000000000000000000000000000, 0.00000000000000000000000000000)1273sage: L.coordinates(2*L.elliptic_logarithm(E(e2,0),prec=100))1274(1.0000000000000000000000000000, 1.0000000000000000000000000000)1275sage: L.coordinates(2*L.elliptic_logarithm(E(e3,0),prec=100))1276(0.00000000000000000000000000000, 1.0000000000000000000000000000)12771278::12791280sage: a4 = -78*i + 1041281sage: a6 = -216*i - 3121282sage: E = EllipticCurve([0,0,0,a4,a6])1283sage: emb = K.embeddings(CC)[1]1284sage: L = E.period_lattice(emb)1285sage: P = E(3+2*i,14-7*i)1286sage: L.elliptic_logarithm(P)12870.297147783912228 - 0.546125549639461*I1288sage: L.coordinates(L.elliptic_logarithm(P))1289(0.628653378040238, 0.371417754610223)1290sage: e1 = 1+3*i; e2 = -4-12*i; e3=-e1-e21291sage: L.coordinates(L.elliptic_logarithm(E(e1,0)))1292(0.500000000000000, 0.500000000000000)1293sage: L.coordinates(L.elliptic_logarithm(E(e2,0)))1294(1.00000000000000, 0.500000000000000)1295sage: L.coordinates(L.elliptic_logarithm(E(e3,0)))1296(0.500000000000000, 0.000000000000000)12971298TESTS (see #10026 and #11767)::12991300sage: K.<w> = QuadraticField(2)1301sage: E = EllipticCurve([ 0, -1, 1, -3*w -4, 3*w + 4 ])1302sage: T = E.simon_two_descent()1303sage: P,Q = T[2]1304sage: embs = K.embeddings(CC)1305sage: Lambda = E.period_lattice(embs[0])1306sage: Lambda.elliptic_logarithm(P,100)13074.71001311261996727669736009981308sage: R.<x> = QQ[]1309sage: K.<a> = NumberField(x^2 + x + 5)1310sage: E = EllipticCurve(K, [0,0,1,-3,-5])1311sage: P = E([0,a])1312sage: Lambda = P.curve().period_lattice(K.embeddings(ComplexField(600))[0])1313sage: Lambda.elliptic_logarithm(P, prec=600)1314-0.842248166487739393375018008381693990800588864069506187033873183845246233548058477561706400464057832396643843146464236956684557207157300006542470428493573195030603817094900751609464 - 0.571366031453267388121279381354098224265947866751130917440598461117775339240176310729173301979590106474259885638797913383502735083088736326391919063211421189027226502851390118943491*I1315sage: K.<a> = QuadraticField(-5)1316sage: E = EllipticCurve([1,1,a,a,0])1317sage: P = E(0,0)1318sage: L = P.curve().period_lattice(K.embeddings(ComplexField())[0])1319sage: L.elliptic_logarithm(P, prec=500)13201.17058357737548897849026170185581196033579563441850967539191867385734983296504066660506637438866628981886518901958717288150400849746892393771983141354 - 1.13513899565966043682474529757126359416758251309237866586896869548539516543734207347695898664875799307727928332953834601460994992792519799260968053875*I1321sage: L.elliptic_logarithm(P, prec=1000)13221.17058357737548897849026170185581196033579563441850967539191867385734983296504066660506637438866628981886518901958717288150400849746892393771983141354014895386251320571643977497740116710952913769943240797618468987304985625823413440999754037939123032233879499904283600304184828809773650066658885672885 - 1.13513899565966043682474529757126359416758251309237866586896869548539516543734207347695898664875799307727928332953834601460994992792519799260968053875387282656993476491590607092182964878750169490985439873220720963653658829712494879003124071110818175013453207439440032582917366703476398880865439217473*I1323"""1324if not P.curve() is self.E:1325raise ValueError, "Point is on the wrong curve"1326if prec is None:1327prec = RealField().precision()1328if P.is_zero():1329return ComplexField(prec)(0)1330# Note: using log2(prec) + 3 guard bits is usually enough.1331# To avoid computing a logarithm, we use 40 guard bits which1332# should be largely enough in practice.1333prec2 = prec + 401334R = RealField(prec2)1335C = ComplexField(prec2)1336pi = R.pi()1337e1,e2,e3 = self._ei1338a1,a2,a3 = [self.embedding(a) for a in self.E.ainvs()[:3]]1339xP, yP = [self.embedding(coord) for coord in P.xy()]1340wP = 2*yP+a1*xP+a313411342# We treat the case of 2-torsion points separately. (Note1343# that Cohen's algorithm does not handle these properly.)13441345if wP.is_zero(): # 2-torsion treated separately1346w1,w2 = self._compute_periods_complex(prec,normalise=False)1347if xP==e1:1348z = w2/21349else:1350if xP==e3:1351z = w1/21352else:1353z = (w1+w2)/21354if reduce:1355z = self.reduce(z)1356return z13571358# NB The first block of code works fine for real embeddings as1359# well as complex embeddings. The special code for real1360# embeddings uses only real arithmetic in the iteration, and is1361# based on Cremona and Thongjunthug.13621363# An older version, based on Cohen's Algorithm 7.4.8 also uses1364# only real arithmetic, and gives different normalisations,1365# but also causes problems (see #10026). It is left in but1366# commented out below.13671368if self.real_flag==0: # complex case13691370a = C((e1-e3).sqrt())1371b = C((e1-e2).sqrt())1372if (a+b).abs() < (a-b).abs(): b=-b1373r = C(((xP-e3)/(xP-e2)).sqrt())1374if r.real()<0: r=-r1375t = -C(wP)/(2*r*(xP-e2))1376# eps controls the end of the loop. Since we aim at a target1377# precision of prec bits, eps = 2^(-prec) is enough.1378eps = R(1) >> prec1379while True:1380s = b*r+a1381a, b = (a+b)/2, (a*b).sqrt()1382if (a+b).abs() < (a-b).abs(): b=-b1383r = (a*(r+1)/s).sqrt()1384if (r.abs()-1).abs() < eps: break1385if r.real()<0: r=-r1386t *= r1387z = ((a/t).arctan())/a1388z = ComplexField(prec)(z)1389if reduce:1390z = self.reduce(z)1391return z13921393if self.real_flag==-1: # real, connected case1394z = C(self._abc[0]) # sqrt(e3-e1)1395a, y, b = z.real(), z.imag(), z.abs()1396uv = (xP-e1).sqrt()1397u, v = uv.real().abs(), uv.imag().abs()1398r = (u*a/(u*a+v*y)).sqrt()1399t = -r*R(wP)/(2*(u**2+v**2))1400on_egg = False1401else: # real, disconnected case1402a = R(e3-e1).sqrt()1403b = R(e3-e2).sqrt()1404if (a+b).abs() < (a-b).abs(): b=-b1405on_egg = (xP<e3)1406if on_egg:1407r = a/R(e3-xP).sqrt()1408t = r*R(wP)/(2*R(xP-e1))1409else:1410r = R((xP-e1)/(xP-e2)).sqrt()1411t = -R(wP)/(2*r*R(xP-e2))14121413# eps controls the end of the loop. Since we aim at a target1414# precision of prec bits, eps = 2^(-prec) is enough.1415eps = R(1) >> prec1416while True:1417s = b*r+a1418a, b = (a+b)/2, (a*b).sqrt()1419r = (a*(r+1)/s).sqrt()1420if (r-1).abs() < eps: break1421t *= r1422z = ((a/t).arctan())/a1423if on_egg:1424w1,w2 = self._compute_periods_real(prec)1425z += w2/21426z = ComplexField(prec)(z)1427if reduce:1428z = self.reduce(z)1429return z143014311432# if self.real_flag==-1:1433# z = self._abc[1] # sqrt(e3-e2)1434# beta = z.norm()1435# alpha = 2*(e3-e2).real()1436# a = beta.sqrt()*21437# b = (alpha+2*beta).sqrt()1438# c = (xP-e3+beta)/(xP-e3).sqrt()1439# else:1440# on_egg = (xP<e3)1441# if on_egg:1442# y3 = -(a1*e1+a3)/21443# lam = (yP-y3)/(xP-e1)1444# x3 = lam*(lam+a1)-a2-xP-e11445# yP = lam*(xP-x3)-yP-a1*x3-a31446# xP = x31447# wP = 2*yP+a1*xP+a31448# a = self._abc[0] # sqrt(e3-e1)1449# b = self._abc[1] # sqrt(e3-e2)1450# c = (xP-e1).sqrt()1451# # So far the values of a,b,c are algebraic (in AA)1452# a = R(a)1453# b = R(b)1454# c = R(c)1455# a,b,c = extended_agm_iteration(a,b,c)1456# if self.real_flag==-1:1457# z = (a/c).arcsin()1458# if wP*((xP-e3)*(xP-e3)-beta*beta) >= 0:1459# z = pi - z1460# if wP > 0:1461# z += pi1462# z /= a1463# elif self.real_flag==+1:1464# z = (a/c).arcsin()/a1465# w1 = w2 = None1466# if wP > 0:1467# if w1 is None:1468# w1, w2 = self.basis(prec)1469# z = w1 - z1470# if on_egg:1471# if w2 is None:1472# w1, w2 = self.basis(prec)1473# z += w2/21474# z = ComplexField(prec)(z)1475# if reduce:1476# z = self.reduce(z)1477# return z14781479def elliptic_exponential(self, z, to_curve=True):1480r"""1481Return the elliptic exponential of a complex number.14821483INPUT:14841485- ``z`` (complex) -- A complex number (viewed modulo this period lattice).14861487- ``to_curve`` (bool, default True): see below.14881489OUTPUT:14901491- If ``to_curve`` is False, a 2-tuple of real or complex1492numbers representing the point `(x,y) = (\wp(z),\wp'(z))`1493where `\wp` denotes the Weierstrass `\wp`-function with1494respect to this lattice.14951496- If ``to_curve`` is True, the point `(X,Y) =1497(x-b_2/12,y-(a_1(x-b_2/12)-a_3)/2)` as a point in `E(\RR)`1498or `E(\CC)`, with `(x,y) = (\wp(z),\wp'(z))` as above, where1499`E` is the elliptic curve over `\RR` or `\CC` whose period1500lattice this is.15011502- If the lattice is real and `z` is also real then the output1503is a pair of real numbers if ``to_curve`` is True, or a1504point in `E(\RR)` if ``to_curve`` is False.15051506.. note::15071508The precision is taken from that of the input ``z``.15091510EXAMPLES::15111512sage: E = EllipticCurve([1,1,1,-8,6])1513sage: P = E(1,-2)1514sage: L = E.period_lattice()1515sage: z = L(P); z15161.170447572400901517sage: L.elliptic_exponential(z)1518(0.999999999999999 : -2.00000000000000 : 1.00000000000000)1519sage: _.curve()1520Elliptic Curve defined by y^2 + 1.00000000000000*x*y + 1.00000000000000*y = x^3 + 1.00000000000000*x^2 - 8.00000000000000*x + 6.00000000000000 over Real Field with 53 bits of precision1521sage: L.elliptic_exponential(z,to_curve=False)1522(1.41666666666667, -1.00000000000000)1523sage: z = L(P,prec=201); z15241.170447572400895922989921884823714935044725616774510079941891525sage: L.elliptic_exponential(z)1526(1.00000000000000000000000000000000000000000000000000000000000 : -2.00000000000000000000000000000000000000000000000000000000000 : 1.00000000000000000000000000000000000000000000000000000000000)15271528Examples over number fields::15291530sage: x = polygen(QQ)1531sage: K.<a> = NumberField(x^3-2)1532sage: embs = K.embeddings(CC)1533sage: E = EllipticCurve('37a')1534sage: EK = E.change_ring(K)1535sage: Li = [EK.period_lattice(e) for e in embs]1536sage: P = EK(-1,-1)1537sage: Q = EK(a-1,1-a^2)1538sage: zi = [L.elliptic_logarithm(P) for L in Li]1539sage: [c.real() for c in Li[0].elliptic_exponential(zi[0])]1540[-1.00000000000000, -1.00000000000000, 1.00000000000000]1541sage: [c.real() for c in Li[0].elliptic_exponential(zi[1])]1542[-1.00000000000000, -1.00000000000000, 1.00000000000000]1543sage: [c.real() for c in Li[0].elliptic_exponential(zi[2])]1544[-1.00000000000000, -1.00000000000000, 1.00000000000000]15451546sage: zi = [L.elliptic_logarithm(Q) for L in Li]1547sage: Li[0].elliptic_exponential(zi[0])1548(-1.62996052494744 - 1.09112363597172*I : 1.79370052598410 - 1.37472963699860*I : 1.00000000000000)1549sage: [embs[0](c) for c in Q]1550[-1.62996052494744 - 1.09112363597172*I, 1.79370052598410 - 1.37472963699860*I, 1.00000000000000]1551sage: Li[1].elliptic_exponential(zi[1])1552(-1.62996052494744 + 1.09112363597172*I : 1.79370052598410 + 1.37472963699860*I : 1.00000000000000)1553sage: [embs[1](c) for c in Q]1554[-1.62996052494744 + 1.09112363597172*I, 1.79370052598410 + 1.37472963699860*I, 1.00000000000000]1555sage: [c.real() for c in Li[2].elliptic_exponential(zi[2])]1556[0.259921049894873, -0.587401051968199, 1.00000000000000]1557sage: [embs[2](c) for c in Q]1558[0.259921049894873, -0.587401051968200, 1.00000000000000]15591560Test to show that #8820 is fixed::15611562sage: E = EllipticCurve('37a')1563sage: K.<a> = QuadraticField(-5)1564sage: L = E.change_ring(K).period_lattice(K.places()[0])1565sage: L.elliptic_exponential(CDF(.1,.1))1566(0.0000142854026029... - 49.9960001066650*I : 249.520141250950 + 250.019855549131*I : 1.00000000000000)1567sage: L.elliptic_exponential(CDF(.1,.1), to_curve=False)1568(0.0000142854026029... - 49.9960001066650*I, 250.020141250950 + 250.019855549131*I)15691570`z=0` is treated as a special case::15711572sage: E = EllipticCurve([1,1,1,-8,6])1573sage: L = E.period_lattice()1574sage: L.elliptic_exponential(0)1575(0.000000000000000 : 1.00000000000000 : 0.000000000000000)1576sage: L.elliptic_exponential(0, to_curve=False)1577(+infinity, +infinity)15781579::15801581sage: E = EllipticCurve('37a')1582sage: K.<a> = QuadraticField(-5)1583sage: L = E.change_ring(K).period_lattice(K.places()[0])1584sage: P = L.elliptic_exponential(0); P1585(0.000000000000000 : 1.00000000000000 : 0.000000000000000)1586sage: P.parent()1587Abelian group of points on Elliptic Curve defined by y^2 + 1.00000000000000*y = x^3 + (-1.00000000000000)*x over Complex Field with 53 bits of precision15881589Very small `z` are handled properly (see #8820)::15901591sage: K.<a> = QuadraticField(-1)1592sage: E = EllipticCurve([0,0,0,a,0])1593sage: L = E.period_lattice(K.complex_embeddings()[0])1594sage: L.elliptic_exponential(1e-100)1595(0.000000000000000 : 1.00000000000000 : 0.000000000000000)15961597The elliptic exponential of `z` is returned as (0 : 1 : 0) if1598the coordinates of z with respect to the period lattice are1599approximately integral::16001601sage: (100/log(2.0,10))/0.81602415.2410118609201603sage: L.elliptic_exponential((RealField(415)(1e-100))).is_zero()1604True1605sage: L.elliptic_exponential((RealField(420)(1e-100))).is_zero()1606False1607"""1608C = z.parent()1609z_is_real = False1610if is_RealField(C):1611z_is_real = True1612C = ComplexField(C.precision())1613z = C(z)1614else:1615if is_ComplexField(C):1616z_is_real = z.is_real()1617else:1618try:1619C = ComplexField()1620z = C(z)1621z_is_real = z.is_real()1622except TypeError:1623raise TypeError, "%s is not a complex number"%z1624prec = C.precision()16251626# test for the point at infinity:16271628eps = (C(2)**(-0.8*prec)).real() ## to test integrality w.r.t. lattice within 20%1629if all([(t.round()-t).abs() < eps for t in self.coordinates(z)]):1630K = z.parent()1631if to_curve:1632return self.curve().change_ring(K)(0)1633else:1634return (K('+infinity'),K('+infinity'))16351636# general number field code (including QQ):16371638# We do not use PARI's ellztopoint function since it is only1639# defined for curves over the reals (note that PARI only1640# computes the period lattice basis in that case). But Sage1641# can compute the period lattice basis over CC, and then1642# PARI's ellwp function works fine.16431644# NB converting the PARI values to Sage values might land up1645# in real/complex fields of spuriously higher precision than1646# the input, since PARI's precision is in word-size chunks.1647# So we force the results back into the real/complex fields of1648# the same precision as the input.16491650from sage.libs.all import pari16511652x,y = pari(self.basis(prec=prec)).ellwp(pari(z),flag=1)1653x,y = [C(t) for t in (x,y)]16541655if self.real_flag and z_is_real:1656x = x.real()1657y = y.real()16581659if to_curve:1660a1,a2,a3,a4,a6 = [self.embedding(a) for a in self.E.ainvs()]1661b2 = self.embedding(self.E.b2())1662x = x - (b2/12)1663y = y - (a1*x+a3)/21664K = x.parent()1665EK = EllipticCurve(K,[a1,a2,a3,a4,a6])1666return EK.point((x,y,K(1)), check=False)1667else:1668return (x,y)16691670def reduce_tau(tau):1671r"""1672Transform a point in the upper half plane to the fundamental region.16731674INPUT:16751676- ``tau`` (complex) -- a complex number with positive imaginary part16771678OUTPUT:16791680(tuple) `(\tau',[a,b,c,d])` where `a,b,c,d` are integers such that16811682- `ad-bc=1`;1683- `\tau`=(a\tau+b)/(c\tau+d)`;1684- `|\tau'|\ge1`;1685- `|\Re(\tau')|\le\frac{1}{2}`.16861687EXAMPLES::16881689sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau1690sage: reduce_tau(CC(1.23,3.45))1691(0.230000000000000 + 3.45000000000000*I, [1, -1, 0, 1])1692sage: reduce_tau(CC(1.23,0.0345))1693(-0.463960069171512 + 1.35591888067914*I, [-5, 6, 4, -5])1694sage: reduce_tau(CC(1.23,0.0000345))1695(0.130000000001761 + 2.89855072463768*I, [13, -16, 100, -123])1696"""1697assert tau.imag()>01698tau_orig = tau1699a, b = ZZ(1), ZZ(0)1700c, d = b, a1701k = tau.real().round()1702tau -= k1703a -= k*c1704b -= k*d1705while tau.abs()<0.999:1706tau = -1/tau1707a, b, c, d = c, d, -a, -b1708k = tau.real().round()1709tau -= k1710a -= k*c1711b -= k*d1712assert a*d-b*c==11713assert tau.abs()>=0.999 and tau.real().abs() <= 0.51714return tau, [a,b,c,d]17151716def normalise_periods(w1,w2):1717r"""1718Normalise the period basis `(w_1,w_2)` so that `w_1/w_2` is in the fundamental region.17191720INPUT:17211722- ``w1,w2`` (complex) -- two complex numbers with non-real ratio17231724OUTPUT:17251726(tuple) `((\omega_1',\omega_2'),[a,b,c,d])` where `a,b,c,d` are1727integers such that17281729- `ad-bc=\pm1`;1730- `(\omega_1',\omega_2') = (a\omega_1+b\omega_2,c\omega_1+d\omega_2)`;1731- `\tau=\omega_1'/\omega_2'` is in the upper half plane;1732- `|\tau|\ge1` and `|\Re(\tau)|\le\frac{1}{2}`.17331734EXAMPLES::17351736sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau, normalise_periods1737sage: w1 = CC(1.234, 3.456)1738sage: w2 = CC(1.234, 3.456000001)1739sage: w1/w2 # in lower half plane!17400.999999999743367 - 9.16334785827644e-11*I1741sage: w1w2, abcd = normalise_periods(w1,w2)1742sage: a,b,c,d = abcd1743sage: w1w2 == (a*w1+b*w2, c*w1+d*w2)1744True1745sage: w1w2[0]/w1w2[1]17461.23400010389203e9*I1747sage: a*d-b*c # note change of orientation1748-117491750"""1751tau = w1/w21752s = +11753if tau.imag()<0:1754w2 = -w21755tau = -tau1756s = -11757tau, abcd = reduce_tau(tau)1758a, b, c, d = abcd1759if s<0:1760abcd = (a,-b,c,-d)1761return (a*w1+b*w2,c*w1+d*w2), abcd176217631764def extended_agm_iteration(a,b,c):1765r"""1766Internal function for the extended AGM used in elliptic logarithm computation.1767INPUT:17681769- ``a``, ``b``, ``c`` (real or complex) -- three real or complex numbers.17701771OUTPUT:17721773(3-tuple) `(a_0,b_0,c_0)`, the limit of the iteration `(a,b,c) \mapsto ((a+b)/2,\sqrt{ab},(c+\sqrt(c^2+b^2-a^2))/2)`.17741775EXAMPLES::17761777sage: from sage.schemes.elliptic_curves.period_lattice import extended_agm_iteration1778sage: extended_agm_iteration(RR(1),RR(2),RR(3))1779(1.45679103104691, 1.45679103104691, 3.21245294970054)1780sage: extended_agm_iteration(CC(1,2),CC(2,3),CC(3,4))1781(1.46242448156430 + 2.47791311676267*I,17821.46242448156430 + 2.47791311676267*I,17833.22202144343535 + 4.28383734262540*I)17841785TESTS::17861787sage: extended_agm_iteration(1,2,3)1788Traceback (most recent call last):1789...1790ValueError: values must be real or complex numbers17911792"""1793if not isinstance(a, (RealNumber,ComplexNumber)):1794raise ValueError, "values must be real or complex numbers"1795eps = a.parent().one().real()>>(a.parent().precision()-10)1796while True:1797a1 = (a + b)/21798b1 = (a*b).sqrt()1799delta = (b**2 - a**2)/c**21800f = (1 + (1 + delta).sqrt())/21801if (f.abs()-1).abs() < eps:1802return a,b,c1803c*=f1804a,b = a1,b1180518061807