Path: blob/master/sage/schemes/elliptic_curves/period_lattice.py
4159 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, is_RealField, is_ComplexField, QQbar, AA98from sage.rings.real_mpfr import RealNumber as RealNumber99from sage.rings.complex_number import ComplexNumber as ComplexNumber100from sage.rings.number_field.number_field import refine_embedding101from sage.rings.infinity import Infinity102from sage.schemes.elliptic_curves.constructor import EllipticCurve103from sage.misc.cachefunc import cached_method104105class PeriodLattice(FreeModule_generic_pid):106"""107The class for the period lattice of an algebraic variety.108"""109pass110111class PeriodLattice_ell(PeriodLattice):112r"""113The class for the period lattice of an elliptic curve.114115Currently supported are elliptic curves defined over `\QQ`, and116elliptic curves defined over a number field with a real or complex117embedding, where the lattice constructed depends on that118embedding.119"""120121def __init__(self, E, embedding=None):122r"""123Initialises the period lattice by storing the elliptic curve and the embedding.124125INPUT:126127- ``E`` -- an elliptic curve128129- ``embedding`` (defult: ``None``) -- an embedding of the base130field `K` of ``E`` into a real or complex field. If131``None``:132133- use the built-in coercion to `\RR` for `K=\QQ`;134135- use the first embedding into `\RR` given by136``K.embeddings(RealField())``, if there are any;137138- use the first embedding into `\CC` given by139``K.embeddings(ComplexField())``, if `K` is totally complex.140141.. note::142143No periods are computed on creation of the lattice; see the144functions ``basis()``, ``normalised_basis()`` and145``real_period()`` for precision setting.146147EXAMPLES:148149This function is not normally called directly, but will be150called by the period_lattice() function of classes151ell_number_field and ell_rational_field::152153sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell154sage: E = EllipticCurve('37a')155sage: PeriodLattice_ell(E)156Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field157158::159160sage: K.<a> = NumberField(x^3-2)161sage: emb = K.embeddings(RealField())[0]162sage: E = EllipticCurve([0,1,0,a,a])163sage: L = PeriodLattice_ell(E,emb); L164Period 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:165From: Number Field in a with defining polynomial x^3 - 2166To: Algebraic Real Field167Defn: a |--> 1.259921049894873?168169sage: emb = K.embeddings(ComplexField())[0]170sage: L = PeriodLattice_ell(E,emb); L171Period 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:172From: Number Field in a with defining polynomial x^3 - 2173To: Algebraic Field174Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I175176TESTS::177178sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell179sage: K.<a> = NumberField(x^3-2)180sage: emb = K.embeddings(RealField())[0]181sage: E = EllipticCurve([0,1,0,a,a])182sage: L = PeriodLattice_ell(E,emb)183sage: L == loads(dumps(L))184True185"""186# First we cache the elliptic curve with this period lattice:187188self.E = E189190# Next we cache the embedding into QQbar or AA which extends191# the given embedding:192193K = E.base_field()194if embedding is None:195embs = K.embeddings(AA)196real = len(embs)>0197if not real:198embs = K.embeddings(QQbar)199embedding = embs[0]200else:201embedding = refine_embedding(embedding,Infinity)202real = embedding(K.gen()).imag().is_zero()203204self.embedding = embedding205206# Next we compute and cache (in self.real_flag) the type of207# the lattice: +1 for real rectangular, -1 for real208# non-rectangular, 0 for non-real:209210self.real_flag = 0211if real:212self.real_flag = +1213if embedding(E.discriminant())<0:214self.real_flag = -1215216# The following algebraic data associated to E and the217# embedding is cached:218#219# Ebar: the curve E base-changed to QQbar (or AA)220# f2: the 2-division polynomial of Ebar221# ei: the roots e1, e2, e3 of f2, as elements of QQbar (or AA)222#223# The ei are used both for period computation and elliptic224# logarithms.225226self.Ebar = self.E.change_ring(self.embedding)227self.f2 = self.Ebar.two_division_polynomial()228if self.real_flag == 1: # positive discriminant229self._ei = self.f2.roots(AA,multiplicities=False)230self._ei.sort() # e1 < e2 < e3231e1, e2, e3 = self._ei232elif self.real_flag == -1: # negative discriminant233self._ei = self.f2.roots(QQbar,multiplicities=False)234self._ei = list(sorted(self._ei,key=lambda z: z.imag()))235e1, e3, e2 = self._ei # so e3 is real236e3 = AA(e3)237self._ei = [e1, e2, e3]238else:239self._ei = self.f2.roots(QQbar,multiplicities=False)240e1, e2, e3 = self._ei241242# The quantities sqrt(e_i-e_j) are cached (as elements of243# QQbar) to be used in period computations:244245self._abc = (e3-e1).sqrt(), (e3-e2).sqrt(), (e2-e1).sqrt()246247PeriodLattice.__init__(self, base_ring=ZZ, rank=2, degree=1, sparse=False)248249def __cmp__(self, other):250r"""251Comparison function for period lattices252253TESTS::254255sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell256sage: K.<a> = NumberField(x^3-2)257sage: E = EllipticCurve([0,1,0,a,a])258sage: embs = K.embeddings(ComplexField())259sage: L1,L2,L3 = [PeriodLattice_ell(E,e) for e in embs]260sage: L1 < L2 < L3261True262263"""264if not isinstance(other, PeriodLattice_ell): return -1265t = cmp(self.E, other.E)266if t: return t267a = self.E.base_field().gen()268return cmp(self.embedding(a), other.embedding(a))269270def __repr__(self):271"""272Returns the string representation of this period lattice.273274EXAMPLES::275276sage: E = EllipticCurve('37a')277sage: E.period_lattice()278Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field279280::281282sage: K.<a> = NumberField(x^3-2)283sage: emb = K.embeddings(RealField())[0]284sage: E = EllipticCurve([0,1,0,a,a])285sage: L = E.period_lattice(emb); L286Period 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:287From: Number Field in a with defining polynomial x^3 - 2288To: Algebraic Real Field289Defn: a |--> 1.259921049894873?290"""291if self.E.base_field() is QQ:292return "Period lattice associated to %s"%(self.E)293else:294return "Period lattice associated to %s with respect to the embedding %s"%(self.E, self.embedding)295296def __call__(self, P, prec=None):297r"""298Return the elliptic logarithm of a point `P`.299300INPUT:301302- ``P`` (point) -- a point on the elliptic curve associated303with this period lattice.304305- ``prec`` (default: ``None``) -- precision in bits (default306precision if ``None``).307308OUTPUT:309310(complex number) The elliptic logarithm of the point `P` with311respect to this period lattice. If `E` is the elliptic curve312and `\sigma:K\to\CC` the embedding, then the returned value `z`313is such that `z\pmod{L}` maps to `\sigma(P)` under the314standard Weierstrass isomorphism from `\CC/L` to `\sigma(E)`.315316EXAMPLES::317318sage: E = EllipticCurve('389a')319sage: L = E.period_lattice()320sage: E.discriminant() > 0321True322sage: L.real_flag3231324sage: P = E([-1,1])325sage: P.is_on_identity_component ()326False327sage: L(P, prec=96)3280.4793482501902193161295330101 + 0.985868850775824102211203849...*I329sage: Q=E([3,5])330sage: Q.is_on_identity_component()331True332sage: L(Q, prec=96)3331.931128271542559442488585220334335Note that this is actually the inverse of the Weierstrass isomorphism::336337sage: L.elliptic_exponential(L(Q))338(3.00000000000000 : 5.00000000000000 : 1.00000000000000)339340An example with negative discriminant, and a torsion point::341342sage: E = EllipticCurve('11a1')343sage: L = E.period_lattice()344sage: E.discriminant() < 0345True346sage: L.real_flag347-1348sage: P = E([16,-61])349sage: L(P)3500.253841860855911351sage: L.real_period() / L(P)3525.00000000000000353"""354return self.elliptic_logarithm(P,prec)355356@cached_method357def basis(self, prec=None, algorithm='sage'):358r"""359Return a basis for this period lattice as a 2-tuple.360361INPUT:362363- ``prec`` (default: ``None``) -- precision in bits (default364precision if ``None``).365366- ``algorithm`` (string, default 'sage') -- choice of367implementation (for real embeddings only) between 'sage'368(native Sage implementation) or 'pari' (use the PARI369library: only available for real embeddings).370371OUTPUT:372373(tuple of Complex) `(\omega_1,\omega_2)` where the lattice is374`\ZZ\omega_1 + \ZZ\omega_2`. If the lattice is real then375`\omega_1` is real and positive, `\Im(\omega_2)>0` and376`\Re(\omega_1/\omega_2)` is either `0` (for rectangular377lattices) or `\frac{1}{2}` (for non-rectangular lattices).378Otherwise, `\omega_1/\omega_2` is in the fundamental region of379the upper half-plane. If the latter normalisation is required380for real lattices, use the function ``normalised_basis()``381instead.382383EXAMPLES::384385sage: E = EllipticCurve('37a')386sage: E.period_lattice().basis()387(2.99345864623196, 2.45138938198679*I)388389This shows that the issue reported at trac \#3954 is fixed::390391sage: E = EllipticCurve('37a')392sage: b1 = E.period_lattice().basis(prec=30)393sage: b2 = E.period_lattice().basis(prec=30)394sage: b1 == b2395True396397This shows that the issue reported at trac \#4064 is fixed::398399sage: E = EllipticCurve('37a')400sage: E.period_lattice().basis(prec=30)[0].parent()401Real Field with 30 bits of precision402sage: E.period_lattice().basis(prec=100)[0].parent()403Real Field with 100 bits of precision404405::406407sage: K.<a> = NumberField(x^3-2)408sage: emb = K.embeddings(RealField())[0]409sage: E = EllipticCurve([0,1,0,a,a])410sage: L = E.period_lattice(emb)411sage: L.basis(64)412(3.81452977217854509, 1.90726488608927255 + 1.34047785962440202*I)413414sage: emb = K.embeddings(ComplexField())[0]415sage: L = E.period_lattice(emb)416sage: w1,w2 = L.basis(); w1,w2417(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)418sage: L.is_real()419False420sage: tau = w1/w2; tau4210.387694505032876 + 1.30821088214407*I422"""423# We divide into two cases: (1) Q, or a number field with a424# real embedding; (2) a number field with a complex embedding.425# In each case the periods are computed by a different426# internal function.427428if self.is_real():429return self._compute_periods_real(prec=prec, algorithm=algorithm)430else:431return self._compute_periods_complex(prec=prec)432433@cached_method434def normalised_basis(self, prec=None, algorithm='sage'):435r"""436Return a normalised basis for this period lattice as a 2-tuple.437438INPUT:439440- ``prec`` (default: ``None``) -- precision in bits (default441precision if ``None``).442443- ``algorithm`` (string, default 'sage') -- choice of444implementation (for real embeddings only) between 'sage'445(native Sage implementation) or 'pari' (use the PARI446library: only available for real embeddings).447448OUTPUT:449450(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has451the form `\ZZ\omega_1 + \ZZ\omega_2`. The basis is normalised452so that `\omega_1/\omega_2` is in the fundamental region of453the upper half-plane. For an alternative normalisation for454real lattices (with the first period real), use the function455basis() instead.456457EXAMPLES::458459sage: E = EllipticCurve('37a')460sage: E.period_lattice().normalised_basis()461(2.99345864623196, -2.45138938198679*I)462463::464465sage: K.<a> = NumberField(x^3-2)466sage: emb = K.embeddings(RealField())[0]467sage: E = EllipticCurve([0,1,0,a,a])468sage: L = E.period_lattice(emb)469sage: L.normalised_basis(64)470(1.90726488608927255 - 1.34047785962440202*I, -1.90726488608927255 - 1.34047785962440202*I)471472sage: emb = K.embeddings(ComplexField())[0]473sage: L = E.period_lattice(emb)474sage: w1,w2 = L.normalised_basis(); w1,w2475(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)476sage: L.is_real()477False478sage: tau = w1/w2; tau4790.387694505032876 + 1.30821088214407*I480"""481w1, w2 = periods = self.basis(prec=prec, algorithm=algorithm)482periods, mat = normalise_periods(w1,w2)483return periods484485@cached_method486def _compute_periods_real(self, prec=None, algorithm='sage'):487r"""488Internal function to compute the periods (real embedding case).489490INPUT:491492493- `prec` (int or ``None`` (default)) -- floating point494precision (in bits); if None, use the default precision.495496- `algorithm` (string, default 'sage') -- choice of implementation between497- `pari`: use the PARI library498499- `sage`: use a native Sage implementation (with the same underlying algorithm).500501502OUTPUT:503504(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has505the form `\ZZ\omega_1 + \ZZ\omega_2`, `\omega_1` is real and506`\omega_1/\omega_2` has real part either `0` or `frac{1}{2}`.507508EXAMPLES::509510sage: K.<a> = NumberField(x^3-2)511sage: E = EllipticCurve([0,1,0,a,a])512sage: embs = K.embeddings(CC)513sage: Ls = [E.period_lattice(e) for e in embs]514sage: [L.is_real() for L in Ls]515[False, False, True]516sage: Ls[2]._compute_periods_real(100)517(3.8145297721785450936365098936,5181.9072648860892725468182549468 + 1.3404778596244020196600112394*I)519sage: Ls[2]._compute_periods_real(100, algorithm='pari')520(3.8145297721785450936365098936,5211.9072648860892725468182549468 - 1.3404778596244020196600112394*I)522"""523if prec is None:524prec = RealField().precision()525R = RealField(prec)526C = ComplexField(prec)527528if algorithm=='pari':529if self.E.base_field() is QQ:530periods = self.E.pari_curve(prec).omega().python()531return (R(periods[0]), C(periods[1]))532533from sage.libs.pari.all import pari534E_pari = pari([R(self.embedding(ai).real()) for ai in self.E.a_invariants()]).ellinit(precision=prec)535periods = E_pari.omega().python()536return (R(periods[0]), C(periods[1]))537538if algorithm!='sage':539raise ValueError, "invalid value of 'algorithm' parameter"540541pi = R.pi()542# Up to now everything has been exact in AA or QQbar, but now543# we must go transcendental. Only now is the desired544# precision used!545if self.real_flag == 1: # positive discriminant546a, b, c = (R(x) for x in self._abc)547w1 = R(pi/a.agm(b)) # least real period548w2 = C(0,pi/a.agm(c)) # least pure imaginary period549else:550a = C(self._abc[0])551x, y, r = a.real().abs(), a.imag().abs(), a.abs()552w1 = R(pi/r.agm(x)) # least real period553w2 = R(pi/r.agm(y)) # least pure imaginary period /i554w2 = C(w1,w2)/2555556return (w1,w2)557558@cached_method559def _compute_periods_complex(self, prec=None, normalise=True):560r"""561Internal function to compute the periods (complex embedding case).562563INPUT:564565- `prec` (int or ``None`` (default)) -- floating point precision (in bits); if None,566use the default precision.567568- `normalise` (bool, default True) -- whether to normalise the569basis after computation.570571OUTPUT:572573(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has574the form `\ZZ\omega_1 + \ZZ\omega_2`. If `normalise` is575`True`, the basis is normalised so that `(\omega_1/\omega_2)`576is in the fundamental region of the upper half plane.577578EXAMPLES::579580sage: K.<a> = NumberField(x^3-2)581sage: E = EllipticCurve([0,1,0,a,a])582sage: embs = K.embeddings(CC)583sage: Ls = [E.period_lattice(e) for e in embs]584sage: [L.is_real() for L in Ls]585[False, False, True]586sage: L = Ls[0]587sage: w1,w2 = L._compute_periods_complex(100); w1,w2588(-1.3758860416607626645495991458 - 2.5856094662444337042877901304*I, -2.1033990784735587243397865076 + 0.42837877646062187766760569686*I)589sage: tau = w1/w2; tau5900.38769450503287609349437509561 + 1.3082108821440725664008561928*I591sage: tau.real()5920.38769450503287609349437509561593sage: tau.abs()5941.3644496111593345713923386773595596Without normalisation::597598sage: w1,w2 = L._compute_periods_complex(normalise=False); w1,w2599(2.10339907847356 - 0.428378776460622*I, 0.727513036812796 - 3.01398824270506*I)600sage: tau = w1/w2; tau6010.293483964608883 + 0.627038168678760*I602sage: tau.real()6030.293483964608883604sage: tau.abs() # > 16050.692321964451917606"""607if prec is None:608prec = RealField().precision()609C = ComplexField(prec)610611# Up to now everything has been exact in AA, but now we612# must go transcendental. Only now is the desired613# precision used!614pi = C.pi()615a, b, c = (C(x) for x in self._abc)616if (a+b).abs() < (a-b).abs(): b=-b617if (a+c).abs() < (a-c).abs(): c=-c618w1 = pi/a.agm(b)619w2 = pi*C.gen()/a.agm(c)620if (w1/w2).imag()<0: w2=-w2621if normalise:622w1w2, mat = normalise_periods(w1,w2)623return w1w2624return (w1,w2)625626def is_real(self):627r"""628Return True if this period lattice is real.629630EXAMPLES::631632sage: f = EllipticCurve('11a')633sage: f.period_lattice().is_real()634True635636::637638sage: K.<i> = QuadraticField(-1)639sage: E = EllipticCurve(K,[0,0,0,i,2*i])640sage: emb = K.embeddings(ComplexField())[0]641sage: L = E.period_lattice(emb)642sage: L.is_real()643False644645::646647sage: K.<a> = NumberField(x^3-2)648sage: E = EllipticCurve([0,1,0,a,a])649sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)]650[False, False, True]651652653ALGORITHM:654655The lattice is real if it is associated to a real embedding;656such lattices are stable under conjugation.657"""658return self.real_flag!=0659660def is_rectangular(self):661r"""662Return True if this period lattice is rectangular.663664.. note::665666Only defined for real lattices; a RuntimeError is raised for667non-real lattices.668669EXAMPLES::670671sage: f = EllipticCurve('11a')672sage: f.period_lattice().basis()673(1.26920930427955, 0.634604652139777 + 1.45881661693850*I)674sage: f.period_lattice().is_rectangular()675False676677::678679sage: f = EllipticCurve('37b')680sage: f.period_lattice().basis()681(1.08852159290423, 1.76761067023379*I)682sage: f.period_lattice().is_rectangular()683True684685ALGORITHM:686687The period lattice is rectangular precisely if the688discriminant of the Weierstrass equation is positive, or689equivalently if the number of real components is 2.690"""691if self.is_real():692return self.real_flag == +1693raise RuntimeError, "Not defined for non-real lattices."694695def real_period(self, prec = None, algorithm='sage'):696"""697Returns the real period of this period lattice.698699INPUT:700701- ``prec`` (int or ``None`` (default)) -- real precision in702bits (default real precision if ``None``)703704- ``algorithm`` (string, default 'sage') -- choice of705implementation (for real embeddings only) between 'sage'706(native Sage implementation) or 'pari' (use the PARI707library: only available for real embeddings).708709.. note::710711Only defined for real lattices; a RuntimeError is raised for712non-real lattices.713714EXAMPLES::715716sage: E = EllipticCurve('37a')717sage: E.period_lattice().real_period()7182.99345864623196719720::721722sage: K.<a> = NumberField(x^3-2)723sage: emb = K.embeddings(RealField())[0]724sage: E = EllipticCurve([0,1,0,a,a])725sage: L = E.period_lattice(emb)726sage: L.real_period(64)7273.81452977217854509728"""729if self.is_real():730return self.basis(prec,algorithm)[0]731raise RuntimeError, "Not defined for non-real lattices."732733def omega(self, prec = None):734r"""735Returns the real or complex volume of this period lattice.736737INPUT:738739- ``prec`` (int or ``None``(default)) -- real precision in740bits (default real precision if ``None``)741742OUTPUT:743744(real) For real lattices, this is the real period times the745number of connected components. For non-real lattices it is746the complex area.747748.. note::749750If the curve is defined over `\QQ` and is given by a751*minimal* Weierstrass equation, then this is the correct752period in the BSD conjecture, i.e., it is the least real753period * 2 when the period lattice is rectangular. More754generally the product of this quantity over all embeddings755appears in the generalised BSD formula.756757758EXAMPLES::759760sage: E = EllipticCurve('37a')761sage: E.period_lattice().omega()7625.98691729246392763764This is not a minimal model::765766sage: E = EllipticCurve([0,-432*6^2])767sage: E.period_lattice().omega()7680.486109385710056769770If you were to plug the above omega into the BSD conjecture, you771would get nonsense. The following works though::772773sage: F = E.minimal_model()774sage: F.period_lattice().omega()7750.972218771420113776777::778779sage: K.<a> = NumberField(x^3-2)780sage: emb = K.embeddings(RealField())[0]781sage: E = EllipticCurve([0,1,0,a,a])782sage: L = E.period_lattice(emb)783sage: L.omega(64)7843.81452977217854509785786A complex example (taken from J.E.Cremona and E.Whitley,787*Periods of cusp forms and elliptic curves over imaginary788quadratic fields*, Mathematics of Computation 62 No. 205789(1994), 407-429)::790791sage: K.<i> = QuadraticField(-1)792sage: E = EllipticCurve([0,1-i,i,-i,0])793sage: L = E.period_lattice(K.embeddings(CC)[0])794sage: L.omega()7958.80694160502647796"""797if self.is_real():798n_components = (self.real_flag+3)//2799return self.real_period(prec) * n_components800else:801return self.complex_area()802803@cached_method804def basis_matrix(self, prec=None, normalised=False):805r"""806Return the basis matrix of this period lattice.807808INPUT:809810- ``prec`` (int or ``None``(default)) -- real precision in811bits (default real precision if ``None``).812813- ``normalised`` (bool, default None) -- if True and the814embedding is real, use the normalised basis (see815``normalised_basis()``) instead of the default.816817OUTPUT:818819A 2x2 real matrix whose rows are the lattice basis vectors,820after identifying `\CC` with `\RR^2`.821822EXAMPLES::823824sage: E = EllipticCurve('37a')825sage: E.period_lattice().basis_matrix()826[ 2.99345864623196 0.000000000000000]827[0.000000000000000 2.45138938198679]828829::830831sage: K.<a> = NumberField(x^3-2)832sage: emb = K.embeddings(RealField())[0]833sage: E = EllipticCurve([0,1,0,a,a])834sage: L = E.period_lattice(emb)835sage: L.basis_matrix(64)836[ 3.81452977217854509 0.000000000000000000]837[ 1.90726488608927255 1.34047785962440202]838839See \#4388::840841sage: L = EllipticCurve('11a1').period_lattice()842sage: L.basis_matrix()843[ 1.26920930427955 0.000000000000000]844[0.634604652139777 1.45881661693850]845sage: L.basis_matrix(normalised=True)846[0.634604652139777 -1.45881661693850]847[-1.26920930427955 0.000000000000000]848849::850851sage: L = EllipticCurve('389a1').period_lattice()852sage: L.basis_matrix()853[ 2.49021256085505 0.000000000000000]854[0.000000000000000 1.97173770155165]855sage: L.basis_matrix(normalised=True)856[ 2.49021256085505 0.000000000000000]857[0.000000000000000 -1.97173770155165]858"""859from sage.matrix.all import Matrix860861if normalised:862return Matrix([list(w) for w in self.normalised_basis(prec)])863864w1,w2 = self.basis(prec)865if self.is_real():866return Matrix([[w1,0],list(w2)])867else:868return Matrix([list(w) for w in w1,w2])869870def complex_area(self, prec=None):871"""872Return the area of a fundamental domain for the period lattice873of the elliptic curve.874875INPUT:876877- ``prec`` (int or ``None``(default)) -- real precision in878bits (default real precision if ``None``).879880EXAMPLES::881882sage: E = EllipticCurve('37a')883sage: E.period_lattice().complex_area()8847.33813274078958885886::887888sage: K.<a> = NumberField(x^3-2)889sage: embs = K.embeddings(ComplexField())890sage: E = EllipticCurve([0,1,0,a,a])891sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)]892[False, False, True]893sage: [E.period_lattice(emb).complex_area() for emb in embs]894[6.02796894766694, 6.02796894766694, 5.11329270448345]895"""896w1,w2 = self.basis(prec)897return (w1*w2.conjugate()).imag().abs()898899def sigma(self, z, prec = None, flag=0):900r"""901Returns the value of the Weierstrass sigma function for this elliptic curve period lattice.902903INPUT:904905- ``z`` -- a complex number906907- ``prec`` (default: ``None``) -- real precision in bits908(default real precision if None).909910- ``flag`` --9119120: (default) ???;9139141: computes an arbitrary determination of log(sigma(z))9159162, 3: same using the product expansion instead of theta series. ???917918.. note::919920The reason for the ???'s above, is that the PARI921documentation for ellsigma is very vague. Also this is922only implemented for curves defined over `\QQ`.923924TODO:925926This function does not use any of the PeriodLattice functions927and so should be moved to ell_rational_field.928929EXAMPLES::930931sage: EllipticCurve('389a1').period_lattice().sigma(CC(2,1))9322.60912163570108 - 0.200865080824587*I933"""934if prec is None:935prec = RealField().precision()936try:937return self.E.pari_curve(prec).ellsigma(z, flag)938except AttributeError:939raise NotImplementedError, "sigma function not yet implemented for period lattices of curves not defined over Q."940941def curve(self):942r"""943Return the elliptic curve associated with this period lattice.944945EXAMPLES::946947sage: E = EllipticCurve('37a')948sage: L = E.period_lattice()949sage: L.curve() is E950True951952::953954sage: K.<a> = NumberField(x^3-2)955sage: E = EllipticCurve([0,1,0,a,a])956sage: L = E.period_lattice(K.embeddings(RealField())[0])957sage: L.curve() is E958True959960sage: L = E.period_lattice(K.embeddings(ComplexField())[0])961sage: L.curve() is E962True963"""964return self.E965966def ei(self):967r"""968Return the x-coordinates of the 2-division points of the elliptic curve associated with this period lattice, as elements of QQbar.969970EXAMPLES::971972sage: E = EllipticCurve('37a')973sage: L = E.period_lattice()974sage: L.ei()975[-1.107159871688768?, 0.2695944364054446?, 0.8375654352833230?]976977::978979sage: K.<a> = NumberField(x^3-2)980sage: E = EllipticCurve([0,1,0,a,a])981sage: L = E.period_lattice(K.embeddings(RealField())[0])982sage: L.ei()983[0.?e-19 - 1.122462048309373?*I, 0.?e-19 + 1.122462048309373?*I, -1]984985sage: L = E.period_lattice(K.embeddings(ComplexField())[0])986sage: L.ei()987[-1.000000000000000? + 0.?e-1...*I,988-0.9720806486198328? - 0.561231024154687?*I,9890.9720806486198328? + 0.561231024154687?*I]990"""991return self._ei992993def coordinates(self, z, rounding=None):994r"""995Returns the coordinates of a complex number w.r.t. the lattice basis996997INPUT:998999- ``z`` (complex) -- A complex number.10001001- ``rounding`` (default ``None``) -- whether and how to round the1002output (see below).10031004OUTPUT:10051006When ``rounding`` is ``None`` (the default), returns a tuple1007of reals `x`, `y` such that `z=xw_1+yw_2` where `w_1`, `w_2`1008are a basis for the lattice (normalised in the case of complex1009embeddings).10101011When ``rounding`` is 'round', returns a tuple of integers `n_1`,1012`n_2` which are the closest integers to the `x`, `y` defined1013above. If `z` is in the lattice these are the coordinates of1014`z` with respect to the lattice basis.10151016When ``rounding`` is 'floor', returns a tuple of integers1017`n_1`, `n_2` which are the integer parts to the `x`, `y`1018defined above. These are used in :meth:``.reduce``10191020EXAMPLES::10211022sage: E = EllipticCurve('389a')1023sage: L = E.period_lattice()1024sage: w1, w2 = L.basis(prec=100)1025sage: P = E([-1,1])1026sage: zP = P.elliptic_logarithm(precision=100); zP10270.47934825019021931612953301006 + 0.98586885077582410221120384908*I1028sage: L.coordinates(zP)1029(0.19249290511394227352563996419, 0.50000000000000000000000000000)1030sage: sum([x*w for x,w in zip(L.coordinates(zP), L.basis(prec=100))])10310.47934825019021931612953301006 + 0.98586885077582410221120384908*I10321033sage: L.coordinates(12*w1+23*w2)1034(12.000000000000000000000000000, 23.000000000000000000000000000)1035sage: L.coordinates(12*w1+23*w2, rounding='floor')1036(11, 22)1037sage: L.coordinates(12*w1+23*w2, rounding='round')1038(12, 23)1039"""1040C = z.parent()1041if is_RealField(C):1042C = ComplexField(C.precision())1043z = C(z)1044else:1045if is_ComplexField(C):1046pass1047else:1048try:1049C = ComplexField()1050z = C(z)1051except TypeError:1052raise TypeError, "%s is not a complex number"%z1053prec = C.precision()1054from sage.matrix.all import Matrix1055from sage.modules.all import vector1056if self.real_flag:1057w1,w2 = self.basis(prec)1058M = Matrix([[w1,0], list(w2)])**(-1)1059else:1060w1,w2 = self.normalised_basis(prec)1061M = Matrix([list(w1), list(w2)])**(-1)1062u,v = vector(z)*M1063# Now z = u*w1+v*w21064if rounding=='round':1065return u.round(), v.round()1066if rounding=='floor':1067return u.floor(), v.floor()1068return u,v10691070def reduce(self, z):1071r"""1072Reduce a complex number modulo the lattice10731074INPUT:10751076- ``z`` (complex) -- A complex number.10771078OUTPUT:10791080(complex) the reduction of `z` modulo the lattice, lying in1081the fundamental period parallelogram with respect to the1082lattice basis. For curves defined over the reals (i.e. real1083embeddings) the output will be real when possible.10841085EXAMPLES::10861087sage: E = EllipticCurve('389a')1088sage: L = E.period_lattice()1089sage: w1, w2 = L.basis(prec=100)1090sage: P = E([-1,1])1091sage: zP = P.elliptic_logarithm(precision=100); zP10920.47934825019021931612953301006 + 0.98586885077582410221120384908*I1093sage: z = zP+10*w1-20*w2; z109425.381473858740770069343110929 - 38.448885180257139986236950114*I1095sage: L.reduce(z)10960.47934825019021931612953301006 + 0.98586885077582410221120384908*I1097sage: L.elliptic_logarithm(2*P)10980.9586965003804391099sage: L.reduce(L.elliptic_logarithm(2*P))11000.9586965003804391101sage: L.reduce(L.elliptic_logarithm(2*P)+10*w1-20*w2)11020.9586965003804441103"""1104C = z.parent()1105z_is_real = False1106if is_RealField(C):1107z_is_real = True1108C = ComplexField(C.precision())1109z = C(z)1110else:1111if is_ComplexField(C):1112z_is_real = z.is_real()1113else:1114try:1115C = ComplexField()1116z = C(z)1117z_is_real = z.is_real()1118except TypeError:1119raise TypeError, "%s is not a complex number"%z1120prec = C.precision()1121if self.real_flag:1122w1,w2 = self.basis(prec) # w1 real1123else:1124w1,w2 = self.normalised_basis(prec)1125# print "z = ",z1126# print "w1,w2 = ",w1,w21127u,v = self.coordinates(z, rounding='floor')1128# print "u,v = ",u,v1129z = z-u*w1-v*w211301131# Final adjustments for the real case.11321133# NB We assume here that when the embedding is real then the1134# point is also real!11351136if self.real_flag == 0: return z1137if self.real_flag == -1:1138k = (z.imag()/w2.imag()).round()1139z = z-k*w21140return C(z.real(),0)11411142if ((2*z.imag()/w2.imag()).round())%2:1143return C(z.real(),w2.imag()/2)1144else:1145return C(z.real(),0)11461147def elliptic_logarithm(self, P, prec=None, reduce=True):1148r"""1149Return the elliptic logarithm of a point.11501151INPUT:11521153- ``P`` (point) -- A point on the elliptic curve associated1154with this period lattice.11551156- ``prec`` (default: ``None``) -- real precision in bits1157(default real precision if None).11581159- ``reduce`` (default: ``True``) -- if ``True``, the result1160is reduced with respect to the period lattice basis.11611162OUTPUT:11631164(complex number) The elliptic logarithm of the point `P` with1165respect to this period lattice. If `E` is the elliptic curve1166and `\sigma:K\to\CC` the embedding, the the returned value `z`1167is such that `z\pmod{L}` maps to `\sigma(P)` under the1168standard Weierstrass isomorphism from `\CC/L` to `\sigma(E)`.1169If ``reduce`` is ``True``, the output is reduced so that it is1170in the fundamental period parallelogram with respect to the1171normalised lattice basis.11721173ALGORITHM: Uses the complex AGM. See [Cremona2010]_ for details.11741175.. [Cremona2010] J. E. Cremona and T. Thongjunthug, The1176Complex AGM, periods of elliptic curves over $\CC$ and1177complex elliptic logarithms. Preprint 2010.11781179EXAMPLES::11801181sage: E = EllipticCurve('389a')1182sage: L = E.period_lattice()1183sage: E.discriminant() > 01184True1185sage: L.real_flag118611187sage: P = E([-1,1])1188sage: P.is_on_identity_component ()1189False1190sage: L.elliptic_logarithm(P, prec=96)11910.4793482501902193161295330101 + 0.9858688507758241022112038491*I1192sage: Q=E([3,5])1193sage: Q.is_on_identity_component()1194True1195sage: L.elliptic_logarithm(Q, prec=96)11961.93112827154255944248858522011971198Note that this is actually the inverse of the Weierstrass isomorphism::11991200sage: L.elliptic_exponential(_)1201(3.00000000000000000000000000... : 5.00000000000000000000000000... : 1.000000000000000000000000000)12021203An example with negative discriminant, and a torsion point::12041205sage: E = EllipticCurve('11a1')1206sage: L = E.period_lattice()1207sage: E.discriminant() < 01208True1209sage: L.real_flag1210-11211sage: P = E([16,-61])1212sage: L.elliptic_logarithm(P)12130.2538418608559111214sage: L.real_period() / L.elliptic_logarithm(P)12155.0000000000000012161217An example where precision is problematic::12181219sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g11220sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1])1221sage: L = E.period_lattice()1222sage: L.ei()1223[5334.003952567705? - 1.964393150436?e-6*I, 5334.003952567705? + 1.964393150436?e-6*I, -10668.25790513541?]1224sage: L.elliptic_logarithm(P,prec=100)12250.2765620401410706146407620309712261227Some complex examples, taken from the paper by Cremona and Thongjunthug::12281229sage: K.<i> = QuadraticField(-1)1230sage: a4 = 9*i-101231sage: a6 = 21-i1232sage: E = EllipticCurve([0,0,0,a4,a6])1233sage: e1 = 3-2*i; e2 = 1+i; e3 = -4+i1234sage: emb = K.embeddings(CC)[1]1235sage: L = E.period_lattice(emb)1236sage: P = E(2-i,4+2*i)12371238By default, the output is reduced with respect to the1239normalised lattice basis, so that its coordinates with respect1240to that basis lie in the interval [0,1)::12411242sage: z = L.elliptic_logarithm(P,prec=100); z12430.70448375537782208460499649302 - 0.79246725643650979858266018068*I1244sage: L.coordinates(z)1245(0.46247636364807931766105406092, 0.79497588726808704200760395829)12461247Using ``reduce=False`` this step can be omitted. In this case1248the coordinates are usually in the interval [-0.5,0.5), but1249this is not guaranteed. This option is mainly for testing1250purposes::12511252sage: z = L.elliptic_logarithm(P,prec=100, reduce=False); z12530.57002153834710752778063503023 + 0.46476340520469798857457031393*I1254sage: L.coordinates(z)1255(0.46247636364807931766105406092, -0.20502411273191295799239604171)12561257The elliptic logs of the 2-torsion points are half-periods::12581259sage: L.elliptic_logarithm(E(e1,0),prec=100)12600.64607575874356525952487867052 + 0.22379609053909448304176885364*I1261sage: L.elliptic_logarithm(E(e2,0),prec=100)12620.71330686725892253793705940192 - 0.40481924028150941053684639367*I1263sage: L.elliptic_logarithm(E(e3,0),prec=100)12640.067231108515357278412180731396 - 0.62861533082060389357861524731*I12651266We check this by doubling and seeing that the resulting1267coordinates are integers::12681269sage: L.coordinates(2*L.elliptic_logarithm(E(e1,0),prec=100))1270(1.0000000000000000000000000000, 0.00000000000000000000000000000)1271sage: L.coordinates(2*L.elliptic_logarithm(E(e2,0),prec=100))1272(1.0000000000000000000000000000, 1.0000000000000000000000000000)1273sage: L.coordinates(2*L.elliptic_logarithm(E(e3,0),prec=100))1274(0.00000000000000000000000000000, 1.0000000000000000000000000000)12751276::12771278sage: a4 = -78*i + 1041279sage: a6 = -216*i - 3121280sage: E = EllipticCurve([0,0,0,a4,a6])1281sage: emb = K.embeddings(CC)[1]1282sage: L = E.period_lattice(emb)1283sage: P = E(3+2*i,14-7*i)1284sage: L.elliptic_logarithm(P)12850.297147783912228 - 0.546125549639461*I1286sage: L.coordinates(L.elliptic_logarithm(P))1287(0.628653378040238, 0.371417754610223)1288sage: e1 = 1+3*i; e2 = -4-12*i; e3=-e1-e21289sage: L.coordinates(L.elliptic_logarithm(E(e1,0)))1290(0.500000000000000, 0.500000000000000)1291sage: L.coordinates(L.elliptic_logarithm(E(e2,0)))1292(1.00000000000000, 0.500000000000000)1293sage: L.coordinates(L.elliptic_logarithm(E(e3,0)))1294(0.500000000000000, 0.000000000000000)12951296TESTS (see #10026 and #11767)::12971298sage: K.<w> = QuadraticField(2)1299sage: E = EllipticCurve([ 0, -1, 1, -3*w -4, 3*w + 4 ])1300sage: T = E.simon_two_descent()1301sage: P,Q = T[2]1302sage: embs = K.embeddings(CC)1303sage: Lambda = E.period_lattice(embs[0])1304sage: Lambda.elliptic_logarithm(P,100)13054.71001311261996727669736009981306sage: R.<x> = QQ[]1307sage: K.<a> = NumberField(x^2 + x + 5)1308sage: E = EllipticCurve(K, [0,0,1,-3,-5])1309sage: P = E([0,a])1310sage: Lambda = P.curve().period_lattice(K.embeddings(ComplexField(600))[0])1311sage: Lambda.elliptic_logarithm(P, prec=600)1312-0.842248166487739393375018008381693990800588864069506187033873183845246233548058477561706400464057832396643843146464236956684557207157300006542470428493573195030603817094900751609464 - 0.571366031453267388121279381354098224265947866751130917440598461117775339240176310729173301979590106474259885638797913383502735083088736326391919063211421189027226502851390118943491*I1313sage: K.<a> = QuadraticField(-5)1314sage: E = EllipticCurve([1,1,a,a,0])1315sage: P = E(0,0)1316sage: L = P.curve().period_lattice(K.embeddings(ComplexField())[0])1317sage: L.elliptic_logarithm(P, prec=500)13181.17058357737548897849026170185581196033579563441850967539191867385734983296504066660506637438866628981886518901958717288150400849746892393771983141354 - 1.13513899565966043682474529757126359416758251309237866586896869548539516543734207347695898664875799307727928332953834601460994992792519799260968053875*I1319sage: L.elliptic_logarithm(P, prec=1000)13201.17058357737548897849026170185581196033579563441850967539191867385734983296504066660506637438866628981886518901958717288150400849746892393771983141354014895386251320571643977497740116710952913769943240797618468987304985625823413440999754037939123032233879499904283600304184828809773650066658885672885 - 1.13513899565966043682474529757126359416758251309237866586896869548539516543734207347695898664875799307727928332953834601460994992792519799260968053875387282656993476491590607092182964878750169490985439873220720963653658829712494879003124071110818175013453207439440032582917366703476398880865439217473*I1321"""1322if not P.curve() is self.E:1323raise ValueError, "Point is on the wrong curve"1324if prec is None:1325prec = RealField().precision()1326if P.is_zero():1327return ComplexField(prec)(0)1328# Note: using log2(prec) + 3 guard bits is usually enough.1329# To avoid computing a logarithm, we use 40 guard bits which1330# should be largely enough in practice.1331prec2 = prec + 401332R = RealField(prec2)1333C = ComplexField(prec2)1334pi = R.pi()1335e1,e2,e3 = self._ei1336a1,a2,a3 = [self.embedding(a) for a in self.E.ainvs()[:3]]1337xP, yP = [self.embedding(coord) for coord in P.xy()]1338wP = 2*yP+a1*xP+a313391340# We treat the case of 2-torsion points separately. (Note1341# that Cohen's algorithm does not handle these properly.)13421343if wP.is_zero(): # 2-torsion treated separately1344w1,w2 = self._compute_periods_complex(prec,normalise=False)1345if xP==e1:1346z = w2/21347else:1348if xP==e3:1349z = w1/21350else:1351z = (w1+w2)/21352if reduce:1353z = self.reduce(z)1354return z13551356# NB The first block of code works fine for real embeddings as1357# well as complex embeddings. The special code for real1358# embeddings uses only real arithmetic in the iteration, and is1359# based on Cremona and Thongjunthug.13601361# An older version, based on Cohen's Algorithm 7.4.8 also uses1362# only real arithmetic, and gives different normalisations,1363# but also causes problems (see #10026). It is left in but1364# commented out below.13651366if self.real_flag==0: # complex case13671368a = C((e1-e3).sqrt())1369b = C((e1-e2).sqrt())1370if (a+b).abs() < (a-b).abs(): b=-b1371r = C(((xP-e3)/(xP-e2)).sqrt())1372if r.real()<0: r=-r1373t = -C(wP)/(2*r*(xP-e2))1374# eps controls the end of the loop. Since we aim at a target1375# precision of prec bits, eps = 2^(-prec) is enough.1376eps = R(1) >> prec1377while True:1378s = b*r+a1379a, b = (a+b)/2, (a*b).sqrt()1380if (a+b).abs() < (a-b).abs(): b=-b1381r = (a*(r+1)/s).sqrt()1382if (r.abs()-1).abs() < eps: break1383if r.real()<0: r=-r1384t *= r1385z = ((a/t).arctan())/a1386z = ComplexField(prec)(z)1387if reduce:1388z = self.reduce(z)1389return z13901391if self.real_flag==-1: # real, connected case1392z = C(self._abc[0]) # sqrt(e3-e1)1393a, y, b = z.real(), z.imag(), z.abs()1394uv = (xP-e1).sqrt()1395u, v = uv.real().abs(), uv.imag().abs()1396r = (u*a/(u*a+v*y)).sqrt()1397t = -r*R(wP)/(2*(u**2+v**2))1398on_egg = False1399else: # real, disconnected case1400a = R(e3-e1).sqrt()1401b = R(e3-e2).sqrt()1402if (a+b).abs() < (a-b).abs(): b=-b1403on_egg = (xP<e3)1404if on_egg:1405r = a/R(e3-xP).sqrt()1406t = r*R(wP)/(2*R(xP-e1))1407else:1408r = R((xP-e1)/(xP-e2)).sqrt()1409t = -R(wP)/(2*r*R(xP-e2))14101411# eps controls the end of the loop. Since we aim at a target1412# precision of prec bits, eps = 2^(-prec) is enough.1413eps = R(1) >> prec1414while True:1415s = b*r+a1416a, b = (a+b)/2, (a*b).sqrt()1417r = (a*(r+1)/s).sqrt()1418if (r-1).abs() < eps: break1419t *= r1420z = ((a/t).arctan())/a1421if on_egg:1422w1,w2 = self._compute_periods_real(prec)1423z += w2/21424z = ComplexField(prec)(z)1425if reduce:1426z = self.reduce(z)1427return z142814291430# if self.real_flag==-1:1431# z = self._abc[1] # sqrt(e3-e2)1432# beta = z.norm()1433# alpha = 2*(e3-e2).real()1434# a = beta.sqrt()*21435# b = (alpha+2*beta).sqrt()1436# c = (xP-e3+beta)/(xP-e3).sqrt()1437# else:1438# on_egg = (xP<e3)1439# if on_egg:1440# y3 = -(a1*e1+a3)/21441# lam = (yP-y3)/(xP-e1)1442# x3 = lam*(lam+a1)-a2-xP-e11443# yP = lam*(xP-x3)-yP-a1*x3-a31444# xP = x31445# wP = 2*yP+a1*xP+a31446# a = self._abc[0] # sqrt(e3-e1)1447# b = self._abc[1] # sqrt(e3-e2)1448# c = (xP-e1).sqrt()1449# # So far the values of a,b,c are algebraic (in AA)1450# a = R(a)1451# b = R(b)1452# c = R(c)1453# a,b,c = extended_agm_iteration(a,b,c)1454# if self.real_flag==-1:1455# z = (a/c).arcsin()1456# if wP*((xP-e3)*(xP-e3)-beta*beta) >= 0:1457# z = pi - z1458# if wP > 0:1459# z += pi1460# z /= a1461# elif self.real_flag==+1:1462# z = (a/c).arcsin()/a1463# w1 = w2 = None1464# if wP > 0:1465# if w1 is None:1466# w1, w2 = self.basis(prec)1467# z = w1 - z1468# if on_egg:1469# if w2 is None:1470# w1, w2 = self.basis(prec)1471# z += w2/21472# z = ComplexField(prec)(z)1473# if reduce:1474# z = self.reduce(z)1475# return z14761477def elliptic_exponential(self, z, to_curve=True):1478r"""1479Return the elliptic exponential of a complex number.14801481INPUT:14821483- ``z`` (complex) -- A complex number (viewed modulo this period lattice).14841485- ``to_curve`` (bool, default True): see below.14861487OUTPUT:14881489- If ``to_curve`` is False, a 2-tuple of real or complex1490numbers representing the point `(x,y) = (\wp(z),\wp'(z))`1491where `\wp` denotes the Weierstrass `\wp`-function with1492respect to this lattice.14931494- If ``to_curve`` is True, the point `(X,Y) =1495(x-b_2/12,y-(a_1(x-b_2/12)-a_3)/2)` as a point in `E(\RR)`1496or `E(\CC)`, with `(x,y) = (\wp(z),\wp'(z))` as above, where1497`E` is the elliptic curve over `\RR` or `\CC` whose period1498lattice this is.14991500- If the lattice is real and `z` is also real then the output1501is a pair of real numbers if ``to_curve`` is True, or a1502point in `E(\RR)` if ``to_curve`` is False.15031504.. note::15051506The precision is taken from that of the input ``z``.15071508EXAMPLES::15091510sage: E = EllipticCurve([1,1,1,-8,6])1511sage: P = E(1,-2)1512sage: L = E.period_lattice()1513sage: z = L(P); z15141.170447572400901515sage: L.elliptic_exponential(z)1516(0.999999999999999 : -2.00000000000000 : 1.00000000000000)1517sage: _.curve()1518Elliptic 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 precision1519sage: L.elliptic_exponential(z,to_curve=False)1520(1.41666666666667, -1.00000000000000)1521sage: z = L(P,prec=201); z15221.170447572400895922989921884823714935044725616774510079941891523sage: L.elliptic_exponential(z)1524(1.00000000000000000000000000000000000000000000000000000000000 : -2.00000000000000000000000000000000000000000000000000000000000 : 1.00000000000000000000000000000000000000000000000000000000000)15251526Examples over number fields::15271528sage: x = polygen(QQ)1529sage: K.<a> = NumberField(x^3-2)1530sage: embs = K.embeddings(CC)1531sage: E = EllipticCurve('37a')1532sage: EK = E.change_ring(K)1533sage: Li = [EK.period_lattice(e) for e in embs]1534sage: P = EK(-1,-1)1535sage: Q = EK(a-1,1-a^2)1536sage: zi = [L.elliptic_logarithm(P) for L in Li]1537sage: [c.real() for c in Li[0].elliptic_exponential(zi[0])]1538[-1.00000000000000, -1.00000000000000, 1.00000000000000]1539sage: [c.real() for c in Li[0].elliptic_exponential(zi[1])]1540[-1.00000000000000, -1.00000000000000, 1.00000000000000]1541sage: [c.real() for c in Li[0].elliptic_exponential(zi[2])]1542[-1.00000000000000, -1.00000000000000, 1.00000000000000]15431544sage: zi = [L.elliptic_logarithm(Q) for L in Li]1545sage: Li[0].elliptic_exponential(zi[0])1546(-1.62996052494744 - 1.09112363597172*I : 1.79370052598410 - 1.37472963699860*I : 1.00000000000000)1547sage: [embs[0](c) for c in Q]1548[-1.62996052494744 - 1.09112363597172*I, 1.79370052598410 - 1.37472963699860*I, 1.00000000000000]1549sage: Li[1].elliptic_exponential(zi[1])1550(-1.62996052494744 + 1.09112363597172*I : 1.79370052598410 + 1.37472963699860*I : 1.00000000000000)1551sage: [embs[1](c) for c in Q]1552[-1.62996052494744 + 1.09112363597172*I, 1.79370052598410 + 1.37472963699860*I, 1.00000000000000]1553sage: [c.real() for c in Li[2].elliptic_exponential(zi[2])]1554[0.259921049894873, -0.587401051968199, 1.00000000000000]1555sage: [embs[2](c) for c in Q]1556[0.259921049894873, -0.587401051968200, 1.00000000000000]15571558Test to show that #8820 is fixed::15591560sage: E = EllipticCurve('37a')1561sage: K.<a> = QuadraticField(-5)1562sage: L = E.change_ring(K).period_lattice(K.places()[0])1563sage: L.elliptic_exponential(CDF(.1,.1))1564(0.0000142854026029... - 49.9960001066650*I : 249.520141250950 + 250.019855549131*I : 1.00000000000000)1565sage: L.elliptic_exponential(CDF(.1,.1), to_curve=False)1566(0.0000142854026029... - 49.9960001066650*I, 250.020141250950 + 250.019855549131*I)15671568`z=0` is treated as a special case::15691570sage: E = EllipticCurve([1,1,1,-8,6])1571sage: L = E.period_lattice()1572sage: L.elliptic_exponential(0)1573(0.000000000000000 : 1.00000000000000 : 0.000000000000000)1574sage: L.elliptic_exponential(0, to_curve=False)1575(+infinity, +infinity)15761577::15781579sage: E = EllipticCurve('37a')1580sage: K.<a> = QuadraticField(-5)1581sage: L = E.change_ring(K).period_lattice(K.places()[0])1582sage: P = L.elliptic_exponential(0); P1583(0.000000000000000 : 1.00000000000000 : 0.000000000000000)1584sage: P.parent()1585Abelian 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 precision15861587Very small `z` are handled properly (see #8820)::15881589sage: K.<a> = QuadraticField(-1)1590sage: E = EllipticCurve([0,0,0,a,0])1591sage: L = E.period_lattice(K.complex_embeddings()[0])1592sage: L.elliptic_exponential(1e-100)1593(0.000000000000000 : 1.00000000000000 : 0.000000000000000)15941595The elliptic exponential of `z` is returned as (0 : 1 : 0) if1596the coordinates of z with respect to the period lattice are1597approximately integral::15981599sage: (100/log(2.0,10))/0.81600415.2410118609201601sage: L.elliptic_exponential((RealField(415)(1e-100))).is_zero()1602True1603sage: L.elliptic_exponential((RealField(420)(1e-100))).is_zero()1604False1605"""1606C = z.parent()1607z_is_real = False1608if is_RealField(C):1609z_is_real = True1610C = ComplexField(C.precision())1611z = C(z)1612else:1613if is_ComplexField(C):1614z_is_real = z.is_real()1615else:1616try:1617C = ComplexField()1618z = C(z)1619z_is_real = z.is_real()1620except TypeError:1621raise TypeError, "%s is not a complex number"%z1622prec = C.precision()16231624# test for the point at infinity:16251626eps = (C(2)**(-0.8*prec)).real() ## to test integrality w.r.t. lattice within 20%1627if all([(t.round()-t).abs() < eps for t in self.coordinates(z)]):1628K = z.parent()1629if to_curve:1630return self.curve().change_ring(K)(0)1631else:1632return (K('+infinity'),K('+infinity'))16331634# general number field code (including QQ):16351636# We do not use PARI's ellztopoint function since it is only1637# defined for curves over the reals (note that PARI only1638# computes the period lattice basis in that case). But Sage1639# can compute the period lattice basis over CC, and then1640# PARI's ellwp function works fine.16411642# NB converting the PARI values to Sage values might land up1643# in real/complex fields of spuriously higher precision than1644# the input, since PARI's precision is in word-size chunks.1645# So we force the results back into the real/complex fields of1646# the same precision as the input.16471648from sage.libs.all import pari16491650x,y = pari(self.basis(prec=prec)).ellwp(pari(z),flag=1)1651x,y = [C(t) for t in (x,y)]16521653if self.real_flag and z_is_real:1654x = x.real()1655y = y.real()16561657if to_curve:1658a1,a2,a3,a4,a6 = [self.embedding(a) for a in self.E.ainvs()]1659b2 = self.embedding(self.E.b2())1660x = x - (b2/12)1661y = y - (a1*x+a3)/21662K = x.parent()1663EK = EllipticCurve(K,[a1,a2,a3,a4,a6])1664return EK.point((x,y,K(1)), check=False)1665else:1666return (x,y)16671668def reduce_tau(tau):1669r"""1670Transform a point in the upper half plane to the fundamental region.16711672INPUT:16731674- ``tau`` (complex) -- a complex number with positive imaginary part16751676OUTPUT:16771678(tuple) `(\tau',[a,b,c,d])` where `a,b,c,d` are integers such that16791680- `ad-bc=1`;1681- `\tau`=(a\tau+b)/(c\tau+d)`;1682- `|\tau'|\ge1`;1683- `|\Re(\tau')|\le\frac{1}{2}`.16841685EXAMPLES::16861687sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau1688sage: reduce_tau(CC(1.23,3.45))1689(0.230000000000000 + 3.45000000000000*I, [1, -1, 0, 1])1690sage: reduce_tau(CC(1.23,0.0345))1691(-0.463960069171512 + 1.35591888067914*I, [-5, 6, 4, -5])1692sage: reduce_tau(CC(1.23,0.0000345))1693(0.130000000001761 + 2.89855072463768*I, [13, -16, 100, -123])1694"""1695assert tau.imag()>01696tau_orig = tau1697a, b = ZZ(1), ZZ(0)1698c, d = b, a1699k = tau.real().round()1700tau -= k1701a -= k*c1702b -= k*d1703while tau.abs()<0.999:1704tau = -1/tau1705a, b, c, d = c, d, -a, -b1706k = tau.real().round()1707tau -= k1708a -= k*c1709b -= k*d1710assert a*d-b*c==11711assert tau.abs()>=0.999 and tau.real().abs() <= 0.51712return tau, [a,b,c,d]17131714def normalise_periods(w1,w2):1715r"""1716Normalise the period basis `(w_1,w_2)` so that `w_1/w_2` is in the fundamental region.17171718INPUT:17191720- ``w1,w2`` (complex) -- two complex numbers with non-real ratio17211722OUTPUT:17231724(tuple) `((\omega_1',\omega_2'),[a,b,c,d])` where `a,b,c,d` are1725integers such that17261727- `ad-bc=\pm1`;1728- `(\omega_1',\omega_2') = (a\omega_1+b\omega_2,c\omega_1+d\omega_2)`;1729- `\tau=\omega_1'/\omega_2'` is in the upper half plane;1730- `|\tau|\ge1` and `|\Re(\tau)|\le\frac{1}{2}`.17311732EXAMPLES::17331734sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau, normalise_periods1735sage: w1 = CC(1.234, 3.456)1736sage: w2 = CC(1.234, 3.456000001)1737sage: w1/w2 # in lower half plane!17380.999999999743367 - 9.16334785827644e-11*I1739sage: w1w2, abcd = normalise_periods(w1,w2)1740sage: a,b,c,d = abcd1741sage: w1w2 == (a*w1+b*w2, c*w1+d*w2)1742True1743sage: w1w2[0]/w1w2[1]17441.23400010389203e9*I1745sage: a*d-b*c # note change of orientation1746-117471748"""1749tau = w1/w21750s = +11751if tau.imag()<0:1752w2 = -w21753tau = -tau1754s = -11755tau, abcd = reduce_tau(tau)1756a, b, c, d = abcd1757if s<0:1758abcd = (a,-b,c,-d)1759return (a*w1+b*w2,c*w1+d*w2), abcd176017611762def extended_agm_iteration(a,b,c):1763r"""1764Internal function for the extended AGM used in elliptic logarithm computation.1765INPUT:17661767- ``a``, ``b``, ``c`` (real or complex) -- three real or complex numbers.17681769OUTPUT:17701771(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)`.17721773EXAMPLES::17741775sage: from sage.schemes.elliptic_curves.period_lattice import extended_agm_iteration1776sage: extended_agm_iteration(RR(1),RR(2),RR(3))1777(1.45679103104691, 1.45679103104691, 3.21245294970054)1778sage: extended_agm_iteration(CC(1,2),CC(2,3),CC(3,4))1779(1.46242448156430 + 2.47791311676267*I,17801.46242448156430 + 2.47791311676267*I,17813.22202144343535 + 4.28383734262540*I)17821783TESTS::17841785sage: extended_agm_iteration(1,2,3)1786Traceback (most recent call last):1787...1788ValueError: values must be real or complex numbers17891790"""1791if not isinstance(a, (RealNumber,ComplexNumber)):1792raise ValueError, "values must be real or complex numbers"1793eps = a.parent().one().real()>>(a.parent().precision()-10)1794while True:1795a1 = (a + b)/21796b1 = (a*b).sqrt()1797delta = (b**2 - a**2)/c**21798f = (1 + (1 + delta).sqrt())/21799if (f.abs()-1).abs() < eps:1800return a,b,c1801c*=f1802a,b = a1,b1180318041805