Path: blob/master/sage/schemes/elliptic_curves/ell_local_data.py
4107 views
# -*- coding: utf-8 -*-1r"""2Local data for elliptic curves over number fields34Let `E` be an elliptic curve over a number field `K` (including `\QQ`).5There are several local invariants at a finite place `v` that6can be computed via Tate's algorithm (see [Sil2] IV.9.4 or [Ta]).78These include the type of reduction (good, additive, multiplicative),9a minimal equation of `E` over `K_v`,10the Tamagawa number `c_v`, defined to be the index `[E(K_v):E^0(K_v)]`11of the points with good reduction among the local points, and the12exponent of the conductor `f_v`.1314The functions in this file will typically be called by using ``local_data``.1516EXAMPLES::1718sage: K.<i> = NumberField(x^2+1)19sage: E = EllipticCurve([(2+i)^2,(2+i)^7])20sage: pp = K.fractional_ideal(2+i)21sage: da = E.local_data(pp)22sage: da.has_bad_reduction()23True24sage: da.has_multiplicative_reduction()25False26sage: da.kodaira_symbol()27I0*28sage: da.tamagawa_number()29430sage: da.minimal_model()31Elliptic Curve defined by y^2 = x^3 + (4*i+3)*x + (-29*i-278) over Number Field in i with defining polynomial x^2 + 13233An example to show how the Neron model can change as one extends the field::3435sage: E = EllipticCurve([0,-1])36sage: E.local_data(2)37Local data at Principal ideal (2) of Integer Ring:38Reduction type: bad additive39Local minimal model: Elliptic Curve defined by y^2 = x^3 - 1 over Rational Field40Minimal discriminant valuation: 441Conductor exponent: 442Kodaira Symbol: II43Tamagawa Number: 14445sage: EK = E.base_extend(K)46sage: EK.local_data(1+i)47Local data at Fractional ideal (i + 1):48Reduction type: bad additive49Local minimal model: Elliptic Curve defined by y^2 = x^3 + (-1) over Number Field in i with defining polynomial x^2 + 150Minimal discriminant valuation: 851Conductor exponent: 252Kodaira Symbol: IV*53Tamagawa Number: 35455Or how the minimal equation changes::5657sage: E = EllipticCurve([0,8])58sage: E.is_minimal()59True60sage: EK = E.base_extend(K)61sage: da = EK.local_data(1+i)62sage: da.minimal_model()63Elliptic Curve defined by y^2 = x^3 + i over Number Field in i with defining polynomial x^2 + 16465REFERENCES:6667- [Sil2] Silverman, Joseph H., Advanced topics in the arithmetic of elliptic curves.68Graduate Texts in Mathematics, 151. Springer-Verlag, New York, 1994.6970- [Ta] Tate, John, Algorithm for determining the type of a singular fiber in an elliptic pencil.71Modular functions of one variable, IV, pp. 33--52. Lecture Notes in Math., Vol. 476,72Springer, Berlin, 1975.7374AUTHORS:7576- John Cremona: First version 2008-09-21 (refactoring code from77``ell_number_field.py`` and ``ell_rational_field.py``)7879- Chris Wuthrich: more documentation 2010-018081"""8283#*****************************************************************************84# Copyright (C) 2005 William Stein <[email protected]>85#86# Distributed under the terms of the GNU General Public License (GPL)87#88# This code is distributed in the hope that it will be useful,89# but WITHOUT ANY WARRANTY; without even the implied warranty of90# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU91# General Public License for more details.92#93# The full text of the GPL is available at:94#95# http://www.gnu.org/licenses/96#*****************************************************************************979899from sage.structure.sage_object import SageObject100from sage.misc.misc import verbose101102from sage.rings.all import PolynomialRing, QQ, ZZ, Integer, is_Ideal, is_NumberFieldElement, is_NumberFieldFractionalIdeal, is_NumberField103from constructor import EllipticCurve104from kodaira_symbol import KodairaSymbol105106class EllipticCurveLocalData(SageObject):107r"""108The class for the local reduction data of an elliptic curve.109110Currently supported are elliptic curves defined over `\QQ`, and111elliptic curves defined over a number field, at an arbitrary prime112or prime ideal.113114INPUT:115116- ``E`` -- an elliptic curve defined over a number field, or `\QQ`.117118- ``P`` -- a prime ideal of the field, or a prime integer if the field is `\QQ`.119120- ``proof`` (bool)-- if True, only use provably correct121methods (default controlled by global proof module). Note122that the proof module is number_field, not elliptic_curves,123since the functions that actually need the flag are in124number fields.125126- ``algorithm`` (string, default: "pari") -- Ignored unless the127base field is `\QQ`. If "pari", use the PARI C-library128``ellglobalred`` implementation of Tate's algorithm over129`\QQ`. If "generic", use the general number field130implementation.131132.. note::133134This function is not normally called directly by users, who135may access the data via methods of the EllipticCurve136classes.137138EXAMPLES::139140sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData141sage: E = EllipticCurve('14a1')142sage: EllipticCurveLocalData(E,2)143Local data at Principal ideal (2) of Integer Ring:144Reduction type: bad non-split multiplicative145Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field146Minimal discriminant valuation: 6147Conductor exponent: 1148Kodaira Symbol: I6149Tamagawa Number: 2150151"""152153def __init__(self, E, P, proof=None, algorithm="pari"):154r"""155Initializes the reduction data for the elliptic curve `E` at the prime `P`.156157INPUT:158159- ``E`` -- an elliptic curve defined over a number field, or `\QQ`.160161- ``P`` -- a prime ideal of the field, or a prime integer if the field is `\QQ`.162163- ``proof`` (bool)-- if True, only use provably correct164methods (default controlled by global proof module). Note165that the proof module is number_field, not elliptic_curves,166since the functions that actually need the flag are in167number fields.168169- ``algorithm`` (string, default: "pari") -- Ignored unless the170base field is `\QQ`. If "pari", use the PARI C-library171``ellglobalred`` implementation of Tate's algorithm over172`\QQ`. If "generic", use the general number field173implementation.174175.. note::176177This function is not normally called directly by users, who178may access the data via methods of the EllipticCurve179classes.180181EXAMPLES::182183sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData184sage: E = EllipticCurve('14a1')185sage: EllipticCurveLocalData(E,2)186Local data at Principal ideal (2) of Integer Ring:187Reduction type: bad non-split multiplicative188Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field189Minimal discriminant valuation: 6190Conductor exponent: 1191Kodaira Symbol: I6192Tamagawa Number: 2193194::195196sage: EllipticCurveLocalData(E,2,algorithm="generic")197Local data at Principal ideal (2) of Integer Ring:198Reduction type: bad non-split multiplicative199Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field200Minimal discriminant valuation: 6201Conductor exponent: 1202Kodaira Symbol: I6203Tamagawa Number: 2204205::206207sage: EllipticCurveLocalData(E,2,algorithm="pari")208Local data at Principal ideal (2) of Integer Ring:209Reduction type: bad non-split multiplicative210Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field211Minimal discriminant valuation: 6212Conductor exponent: 1213Kodaira Symbol: I6214Tamagawa Number: 2215216::217218sage: EllipticCurveLocalData(E,2,algorithm="unknown")219Traceback (most recent call last):220...221ValueError: algorithm must be one of 'pari', 'generic'222223::224225sage: EllipticCurveLocalData(E,3)226Local data at Principal ideal (3) of Integer Ring:227Reduction type: good228Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field229Minimal discriminant valuation: 0230Conductor exponent: 0231Kodaira Symbol: I0232Tamagawa Number: 1233234::235236sage: EllipticCurveLocalData(E,7)237Local data at Principal ideal (7) of Integer Ring:238Reduction type: bad split multiplicative239Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field240Minimal discriminant valuation: 3241Conductor exponent: 1242Kodaira Symbol: I3243Tamagawa Number: 3244"""245self._curve = E246K = E.base_field()247p = check_prime(K,P) # error handling done in that function248if algorithm != "pari" and algorithm != "generic":249raise ValueError, "algorithm must be one of 'pari', 'generic'"250251self._reduction_type = None252if K is QQ:253self._prime = ZZ.ideal(p)254else:255self._prime = p256257if algorithm=="pari" and K is QQ:258Eint = E.integral_model()259data = Eint.pari_curve().elllocalred(p)260self._fp = data[0].python()261self._KS = KodairaSymbol(data[1].python())262self._cp = data[3].python()263# We use a global minimal model since we can:264self._Emin_reduced = Eint.minimal_model()265self._val_disc = self._Emin_reduced.discriminant().valuation(p)266if self._fp>0:267self._reduction_type = Eint.ap(p) # = 0,-1 or +1268else:269self._Emin, ch, self._val_disc, self._fp, self._KS, self._cp, self._split = self._tate(proof)270if self._fp>0:271if self._Emin.c4().valuation(p)>0:272self._reduction_type = 0273elif self._split:274self._reduction_type = +1275else:276self._reduction_type = -1277278def __repr__(self):279r"""280Returns the string representation of this reduction data.281282EXAMPLES::283284sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData285sage: E = EllipticCurve('14a1')286sage: EllipticCurveLocalData(E,2).__repr__()287'Local data at Principal ideal (2) of Integer Ring:\nReduction type: bad non-split multiplicative\nLocal minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field\nMinimal discriminant valuation: 6\nConductor exponent: 1\nKodaira Symbol: I6\nTamagawa Number: 2'288sage: EllipticCurveLocalData(E,3).__repr__()289'Local data at Principal ideal (3) of Integer Ring:\nReduction type: good\nLocal minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field\nMinimal discriminant valuation: 0\nConductor exponent: 0\nKodaira Symbol: I0\nTamagawa Number: 1'290sage: EllipticCurveLocalData(E,7).__repr__()291'Local data at Principal ideal (7) of Integer Ring:\nReduction type: bad split multiplicative\nLocal minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field\nMinimal discriminant valuation: 3\nConductor exponent: 1\nKodaira Symbol: I3\nTamagawa Number: 3'292"""293red_type = "good"294if not self._reduction_type is None:295red_type = ["bad non-split multiplicative","bad additive","bad split multiplicative"][1+self._reduction_type]296return "Local data at %s:\nReduction type: %s\nLocal minimal model: %s\nMinimal discriminant valuation: %s\nConductor exponent: %s\nKodaira Symbol: %s\nTamagawa Number: %s"%(self._prime,red_type,self.minimal_model(),self._val_disc,self._fp,self._KS,self._cp)297298def minimal_model(self, reduce=True):299"""300Return the (local) minimal model from this local reduction data.301302INPUT:303304- ``reduce`` -- (default: True) if set to True the EC returned305by Tate's algorithm will be306"reduced" as specified in _reduce_model() for curves over307number fields.308309EXAMPLES::310311sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData312sage: E = EllipticCurve([0,0,0,0,64]); E313Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field314sage: data = EllipticCurveLocalData(E,2)315sage: data.minimal_model()316Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field317sage: data.minimal_model() == E.local_minimal_model(2)318True319320To demonstrate the behaviour of the parameter ``reduce``::321322sage: K.<a> = NumberField(x^3+x+1)323sage: E = EllipticCurve(K, [0, 0, a, 0, 1])324sage: E.local_data(K.ideal(a-1)).minimal_model()325Elliptic Curve defined by y^2 + a*y = x^3 + 1 over Number Field in a with defining polynomial x^3 + x + 1326sage: E.local_data(K.ideal(a-1)).minimal_model(reduce=False)327Elliptic Curve defined by y^2 + (a+2)*y = x^3 + 3*x^2 + 3*x + (-a+1) over Number Field in a with defining polynomial x^3 + x + 1328329sage: E = EllipticCurve([2, 1, 0, -2, -1])330sage: E.local_data(ZZ.ideal(2), algorithm="generic").minimal_model(reduce=False)331Elliptic Curve defined by y^2 + 2*x*y + 2*y = x^3 + x^2 - 4*x - 2 over Rational Field332sage: E.local_data(ZZ.ideal(2), algorithm="pari").minimal_model(reduce=False)333Traceback (most recent call last):334...335ValueError: the argument reduce must not be False if algorithm=pari is used336sage: E.local_data(ZZ.ideal(2), algorithm="generic").minimal_model()337Elliptic Curve defined by y^2 = x^3 - x^2 - 3*x + 2 over Rational Field338sage: E.local_data(ZZ.ideal(2), algorithm="pari").minimal_model()339Elliptic Curve defined by y^2 = x^3 - x^2 - 3*x + 2 over Rational Field340"""341if reduce:342try:343return self._Emin_reduced344except AttributeError:345pass346self._Emin_reduced = self._Emin._reduce_model()347return self._Emin_reduced348else:349try:350return self._Emin351except AttributeError:352raise ValueError, "the argument reduce must not be False if algorithm=pari is used"353354def prime(self):355"""356Return the prime ideal associated with this local reduction data.357358EXAMPLES::359360sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData361sage: E = EllipticCurve([0,0,0,0,64]); E362Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field363sage: data = EllipticCurveLocalData(E,2)364sage: data.prime()365Principal ideal (2) of Integer Ring366"""367return self._prime368369def conductor_valuation(self):370"""371Return the valuation of the conductor from this local reduction data.372373EXAMPLES::374375sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData376sage: E = EllipticCurve([0,0,0,0,64]); E377Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field378sage: data = EllipticCurveLocalData(E,2)379sage: data.conductor_valuation()3802381"""382return self._fp383384def kodaira_symbol(self):385r"""386Return the Kodaira symbol from this local reduction data.387388EXAMPLES::389390sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData391sage: E = EllipticCurve([0,0,0,0,64]); E392Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field393sage: data = EllipticCurveLocalData(E,2)394sage: data.kodaira_symbol()395IV396"""397return self._KS398399def tamagawa_number(self):400r"""401Return the Tamagawa number from this local reduction data.402403This is the index `[E(K_v):E^0(K_v)]`.404405EXAMPLES::406407sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData408sage: E = EllipticCurve([0,0,0,0,64]); E409Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field410sage: data = EllipticCurveLocalData(E,2)411sage: data.tamagawa_number()4123413"""414return self._cp415416def tamagawa_exponent(self):417r"""418Return the Tamagawa index from this local reduction data.419420This is the exponent of `E(K_v)/E^0(K_v)`; in most cases it is421the same as the Tamagawa index.422423EXAMPLES::424425sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData426sage: E = EllipticCurve('816a1')427sage: data = EllipticCurveLocalData(E,2)428sage: data.kodaira_symbol()429I2*430sage: data.tamagawa_number()4314432sage: data.tamagawa_exponent()4332434435sage: E = EllipticCurve('200c4')436sage: data = EllipticCurveLocalData(E,5)437sage: data.kodaira_symbol()438I4*439sage: data.tamagawa_number()4404441sage: data.tamagawa_exponent()4422443"""444cp = self._cp445if cp!=4:446return cp447ks = self._KS448if ks._roman==1 and ks._n%2==0 and ks._starred:449return 2450return 4451452def bad_reduction_type(self):453r"""454Return the type of bad reduction of this reduction data.455456OUTPUT:457458(int or ``None``):459460- +1 for split multiplicative reduction461- -1 for non-split multiplicative reduction462- 0 for additive reduction463- ``None`` for good reduction464465EXAMPLES::466467sage: E=EllipticCurve('14a1')468sage: [(p,E.local_data(p).bad_reduction_type()) for p in prime_range(15)]469[(2, -1), (3, None), (5, None), (7, 1), (11, None), (13, None)]470471sage: K.<a>=NumberField(x^3-2)472sage: P17a, P17b = [P for P,e in K.factor(17)]473sage: E = EllipticCurve([0,0,0,0,2*a+1])474sage: [(p,E.local_data(p).bad_reduction_type()) for p in [P17a,P17b]]475[(Fractional ideal (4*a^2 - 2*a + 1), None), (Fractional ideal (2*a + 1), 0)]476"""477return self._reduction_type478479def has_good_reduction(self):480r"""481Return True if there is good reduction.482483EXAMPLES::484485sage: E = EllipticCurve('14a1')486sage: [(p,E.local_data(p).has_good_reduction()) for p in prime_range(15)]487[(2, False), (3, True), (5, True), (7, False), (11, True), (13, True)]488489sage: K.<a> = NumberField(x^3-2)490sage: P17a, P17b = [P for P,e in K.factor(17)]491sage: E = EllipticCurve([0,0,0,0,2*a+1])492sage: [(p,E.local_data(p).has_good_reduction()) for p in [P17a,P17b]]493[(Fractional ideal (4*a^2 - 2*a + 1), True),494(Fractional ideal (2*a + 1), False)]495"""496return self._reduction_type is None497498def has_bad_reduction(self):499r"""500Return True if there is bad reduction.501502EXAMPLES::503504sage: E = EllipticCurve('14a1')505sage: [(p,E.local_data(p).has_bad_reduction()) for p in prime_range(15)]506[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]507508::509510sage: K.<a> = NumberField(x^3-2)511sage: P17a, P17b = [P for P,e in K.factor(17)]512sage: E = EllipticCurve([0,0,0,0,2*a+1])513sage: [(p,E.local_data(p).has_bad_reduction()) for p in [P17a,P17b]]514[(Fractional ideal (4*a^2 - 2*a + 1), False),515(Fractional ideal (2*a + 1), True)]516"""517return not self._reduction_type is None518519def has_multiplicative_reduction(self):520r"""521Return True if there is multiplicative reduction.522523.. note::524525See also ``has_split_multiplicative_reduction()`` and526``has_nonsplit_multiplicative_reduction()``.527528EXAMPLES::529530sage: E = EllipticCurve('14a1')531sage: [(p,E.local_data(p).has_multiplicative_reduction()) for p in prime_range(15)]532[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]533534::535536sage: K.<a> = NumberField(x^3-2)537sage: P17a, P17b = [P for P,e in K.factor(17)]538sage: E = EllipticCurve([0,0,0,0,2*a+1])539sage: [(p,E.local_data(p).has_multiplicative_reduction()) for p in [P17a,P17b]]540[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]541"""542return self._reduction_type in (-1,+1)543544def has_split_multiplicative_reduction(self):545r"""546Return True if there is split multiplicative reduction.547548EXAMPLES::549550sage: E = EllipticCurve('14a1')551sage: [(p,E.local_data(p).has_split_multiplicative_reduction()) for p in prime_range(15)]552[(2, False), (3, False), (5, False), (7, True), (11, False), (13, False)]553554::555556sage: K.<a> = NumberField(x^3-2)557sage: P17a, P17b = [P for P,e in K.factor(17)]558sage: E = EllipticCurve([0,0,0,0,2*a+1])559sage: [(p,E.local_data(p).has_split_multiplicative_reduction()) for p in [P17a,P17b]]560[(Fractional ideal (4*a^2 - 2*a + 1), False),561(Fractional ideal (2*a + 1), False)]562"""563return self._reduction_type == +1564565def has_nonsplit_multiplicative_reduction(self):566r"""567Return True if there is non-split multiplicative reduction.568569EXAMPLES::570571sage: E = EllipticCurve('14a1')572sage: [(p,E.local_data(p).has_nonsplit_multiplicative_reduction()) for p in prime_range(15)]573[(2, True), (3, False), (5, False), (7, False), (11, False), (13, False)]574575::576577sage: K.<a> = NumberField(x^3-2)578sage: P17a, P17b = [P for P,e in K.factor(17)]579sage: E = EllipticCurve([0,0,0,0,2*a+1])580sage: [(p,E.local_data(p).has_nonsplit_multiplicative_reduction()) for p in [P17a,P17b]]581[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]582"""583return self._reduction_type == -1584585def has_additive_reduction(self):586r"""587Return True if there is additive reduction.588589EXAMPLES::590591sage: E = EllipticCurve('27a1')592sage: [(p,E.local_data(p).has_additive_reduction()) for p in prime_range(15)]593[(2, False), (3, True), (5, False), (7, False), (11, False), (13, False)]594595::596597sage: K.<a> = NumberField(x^3-2)598sage: P17a, P17b = [P for P,e in K.factor(17)]599sage: E = EllipticCurve([0,0,0,0,2*a+1])600sage: [(p,E.local_data(p).has_additive_reduction()) for p in [P17a,P17b]]601[(Fractional ideal (4*a^2 - 2*a + 1), False),602(Fractional ideal (2*a + 1), True)]603"""604return self._reduction_type == 0605606def _tate(self, proof = None):607r"""608Tate's algorithm for an elliptic curve over a number field.609610Computes both local reduction data at a prime ideal and a611local minimal model.612613The model is not required to be integral on input. If `P` is614principal, uses a generator as uniformizer, so it will not615affect integrality or minimality at other primes. If `P` is not616principal, the minimal model returned will preserve617integrality at other primes, but not minimality.618619.. note::620621Called only by ``EllipticCurveLocalData.__init__()``.622623OUTPUT:624625(tuple) ``(Emin, p, val_disc, fp, KS, cp)`` where:626627- ``Emin`` (EllipticCurve) is a model (integral and) minimal at P628- ``p`` (int) is the residue characteristic629- ``val_disc`` (int) is the valuation of the local minimal discriminant630- ``fp`` (int) is the valuation of the conductor631- ``KS`` (string) is the Kodaira symbol632- ``cp`` (int) is the Tamagawa number633634635EXAMPLES (this raised a type error in sage prior to 4.4.4, see ticket #7930) ::636637sage: E = EllipticCurve('99d1')638639sage: R.<X> = QQ[]640sage: K.<t> = NumberField(X^3 + X^2 - 2*X - 1)641sage: L.<s> = NumberField(X^3 + X^2 - 36*X - 4)642643sage: EK = E.base_extend(K)644sage: toK = EK.torsion_order()645sage: da = EK.local_data() # indirect doctest646647sage: EL = E.base_extend(L)648sage: da = EL.local_data() # indirect doctest649650EXAMPLES:651652The following example shows that the bug at #9324 is fixed::653654sage: K.<a> = NumberField(x^2-x+6)655sage: E = EllipticCurve([0,0,0,-53160*a-43995,-5067640*a+19402006])656sage: E.conductor() # indirect doctest657Fractional ideal (18, 6*a)658659The following example shows that the bug at #9417 is fixed::660661sage: K.<a> = NumberField(x^2+18*x+1)662sage: E = EllipticCurve(K, [0, -36, 0, 320, 0])663sage: E.tamagawa_number(K.ideal(2))6644665666"""667E = self._curve668P = self._prime669K = E.base_ring()670OK = K.maximal_order()671t = verbose("Running Tate's algorithm with P = %s"%P, level=1)672F = OK.residue_field(P)673p = F.characteristic()674675# In case P is not principal we mostly use a uniformiser which676# is globally integral (with positive valuation at some other677# primes); for this to work, it is essential that we can678# reduce (mod P) elements of K which are not integral (but are679# P-integral). However, if the model is non-minimal and we680# end up dividing a_i by pi^i then at that point we use a681# uniformiser pi which has non-positive valuation at all other682# primes, so that we can divide by it without losing683# integrality at other primes.684685principal_flag = P.is_principal()686if principal_flag:687pi = P.gens_reduced()[0]688verbose("P is principal, generator pi = %s"%pi, t, 1)689else:690pi = K.uniformizer(P, 'positive')691verbose("P is not principal, uniformizer pi = %s"%pi, t, 1)692pi2 = pi*pi; pi3 = pi*pi2; pi4 = pi*pi3693pi_neg = None694prime = pi if K is QQ else P695696pval = lambda x: x.valuation(prime)697pdiv = lambda x: x.is_zero() or pval(x) > 0698# Since ResidueField is cached in a way that699# does not care much about embeddings of number700# fields, it can happen that F.p.ring() is different701# from K. This is a problem: If F.p.ring() has no702# embedding but K has, then there is no coercion703# from F.p.ring().maximal_order() to K. But it is704# no problem to do an explicit conversion in that705# case (Simon King, trac ticket #8800).706707from sage.categories.pushout import pushout, CoercionException708try:709if hasattr(F.p.ring(), 'maximal_order'): # it is not ZZ710_tmp_ = pushout(F.p.ring().maximal_order(),K)711pinv = lambda x: F.lift(~F(x))712proot = lambda x,e: F.lift(F(x).nth_root(e, extend = False, all = True)[0])713preduce = lambda x: F.lift(F(x))714except CoercionException: # the pushout does not exist, we need conversion715pinv = lambda x: K(F.lift(~F(x)))716proot = lambda x,e: K(F.lift(F(x).nth_root(e, extend = False, all = True)[0]))717preduce = lambda x: K(F.lift(F(x)))718719def _pquadroots(a, b, c):720r"""721Local function returning True iff `ax^2 + bx + c` has roots modulo `P`722"""723(a, b, c) = (F(a), F(b), F(c))724if a == 0:725return (b != 0) or (c == 0)726elif p == 2:727return len(PolynomialRing(F, "x")([c,b,a]).roots()) > 0728else:729return (b**2 - 4*a*c).is_square()730def _pcubicroots(b, c, d):731r"""732Local function returning the number of roots of `x^3 +733b*x^2 + c*x + d` modulo `P`, counting multiplicities734"""735736return sum([rr[1] for rr in PolynomialRing(F, 'x')([F(d), F(c), F(b), F(1)]).roots()],0)737738if p == 2:739halfmodp = OK(Integer(0))740else:741halfmodp = pinv(Integer(2))742743A = E.a_invariants()744A = [0, A[0], A[1], A[2], A[3], 0, A[4]]745indices = [1,2,3,4,6]746if min([pval(a) for a in A if a != 0]) < 0:747verbose("Non-integral model at P: valuations are %s; making integral"%([pval(a) for a in A if a != 0]), t, 1)748e = 0749for i in range(7):750if A[i] != 0:751e = max(e, (-pval(A[i])/i).ceil())752pie = pi**e753for i in range(7):754if A[i] != 0:755A[i] *= pie**i756verbose("P-integral model is %s, with valuations %s"%([A[i] for i in indices], [pval(A[i]) for i in indices]), t, 1)757758split = None # only relevant for multiplicative reduction759760(a1, a2, a3, a4, a6) = (A[1], A[2], A[3], A[4], A[6])761while True:762C = EllipticCurve([a1, a2, a3, a4, a6]);763(b2, b4, b6, b8) = C.b_invariants()764(c4, c6) = C.c_invariants()765delta = C.discriminant()766val_disc = pval(delta)767768if val_disc == 0:769## Good reduction already770cp = 1771fp = 0772KS = KodairaSymbol("I0")773break #return774775# Otherwise, we change coordinates so that p | a3, a4, a6776if p == 2:777if pdiv(b2):778r = proot(a4, 2)779t = proot(((r + a2)*r + a4)*r + a6, 2)780else:781temp = pinv(a1)782r = temp * a3783t = temp * (a4 + r*r)784elif p == 3:785if pdiv(b2):786r = proot(-b6, 3)787else:788r = -pinv(b2) * b4789t = a1 * r + a3790else:791if pdiv(c4):792r = -pinv(12) * b2793else:794r = -pinv(12*c4) * (c6 + b2 * c4)795t = -halfmodp * (a1 * r + a3)796r = preduce(r)797t = preduce(t)798verbose("Before first transform C = %s"%C)799verbose("[a1,a2,a3,a4,a6] = %s"%([a1, a2, a3, a4, a6]))800C = C.rst_transform(r, 0, t)801(a1, a2, a3, a4, a6) = C.a_invariants()802(b2, b4, b6, b8) = C.b_invariants()803if min([pval(a) for a in (a1, a2, a3, a4, a6) if a != 0]) < 0:804raise RuntimeError, "Non-integral model after first transform!"805verbose("After first transform %s\n, [a1,a2,a3,a4,a6] = %s\n, valuations = %s"%([r, 0, t], [a1, a2, a3, a4, a6], [pval(a1), pval(a2), pval(a3), pval(a4), pval(a6)]), t, 2)806if pval(a3) == 0:807raise RuntimeError, "p does not divide a3 after first transform!"808if pval(a4) == 0:809raise RuntimeError, "p does not divide a4 after first transform!"810if pval(a6) == 0:811raise RuntimeError, "p does not divide a6 after first transform!"812813# Now we test for Types In, II, III, IV814# NB the c invariants never change.815816if not pdiv(c4):817# Multiplicative reduction: Type In (n = val_disc)818split = False819if _pquadroots(1, a1, -a2):820cp = val_disc821split = True822elif Integer(2).divides(val_disc):823cp = 2824else:825cp = 1826KS = KodairaSymbol("I%s"%val_disc)827fp = 1828break #return829830# Additive reduction831832if pval(a6) < 2:833## Type II834KS = KodairaSymbol("II")835fp = val_disc836cp = 1837break #return838if pval(b8) < 3:839## Type III840KS = KodairaSymbol("III")841fp = val_disc - 1842cp = 2843break #return844if pval(b6) < 3:845## Type IV846cp = 1847a3t = preduce(a3/pi)848a6t = preduce(a6/pi2)849if _pquadroots(1, a3t, -a6t): cp = 3850KS = KodairaSymbol("IV")851fp = val_disc - 2852break #return853854# If our curve is none of these types, we change coords so that855# p | a1, a2; p^2 | a3, a4; p^3 | a6856if p == 2:857s = proot(a2, 2) # so s^2=a2 (mod pi)858t = pi*proot(a6/pi2, 2) # so t^2=a6 (mod pi^3)859elif p == 3:860s = a1 # so a1'=2s+a1=3a1=0 (mod pi)861t = a3 # so a3'=2t+a3=3a3=0 (mod pi^2)862else:863s = -a1*halfmodp # so a1'=2s+a1=0 (mod pi)864t = -a3*halfmodp # so a3'=2t+a3=0 (mod pi^2)865C = C.rst_transform(0, s, t)866(a1, a2, a3, a4, a6) = C.a_invariants()867(b2, b4, b6, b8) = C.b_invariants()868verbose("After second transform %s\n[a1, a2, a3, a4, a6] = %s\nValuations: %s"%([0, s, t], [a1,a2,a3,a4,a6],[pval(a1),pval(a2),pval(a3),pval(a4),pval(a6)]), t, 2)869if pval(a1) == 0:870raise RuntimeError, "p does not divide a1 after second transform!"871if pval(a2) == 0:872raise RuntimeError, "p does not divide a2 after second transform!"873if pval(a3) < 2:874raise RuntimeError, "p^2 does not divide a3 after second transform!"875if pval(a4) < 2:876raise RuntimeError, "p^2 does not divide a4 after second transform!"877if pval(a6) < 3:878raise RuntimeError, "p^3 does not divide a6 after second transform!"879if min(pval(a1), pval(a2), pval(a3), pval(a4), pval(a6)) < 0:880raise RuntimeError, "Non-integral model after second transform!"881882# Analyze roots of the cubic T^3 + bT^2 + cT + d = 0 mod P, where883# b = a2/p, c = a4/p^2, d = a6/p^3884b = preduce(a2/pi)885c = preduce(a4/pi2)886d = preduce(a6/pi3)887bb = b*b888cc = c*c889bc = b*c890w = 27*d*d - bb*cc + 4*b*bb*d - 18*bc*d + 4*c*cc891x = 3*c - bb892if pdiv(w):893if pdiv(x):894sw = 3895else:896sw = 2897else:898sw = 1899verbose("Analyzing roots of cubic T^3 + %s*T^2 + %s*T + %s, case %s"%(b, c, d, sw), t, 1)900if sw == 1:901## Three distinct roots - Type I*0902verbose("Distinct roots", t, 1)903KS = KodairaSymbol("I0*")904cp = 1 + _pcubicroots(b, c, d)905fp = val_disc - 4906break #return907elif sw == 2:908## One double root - Type I*m for some m909verbose("One double root", t, 1)910## Change coords so that the double root is T = 0 mod p911if p == 2:912r = proot(c, 2)913elif p == 3:914r = c * pinv(b)915else:916r = (bc - 9*d)*pinv(2*x)917r = pi * preduce(r)918C = C.rst_transform(r, 0, 0)919(a1, a2, a3, a4, a6) = C.a_invariants()920(b2, b4, b6, b8) = C.b_invariants()921# The rest of this branch is just to compute cp, fp, KS.922# We use pi to keep transforms integral.923ix = 3; iy = 3; mx = pi2; my = mx924while True:925a2t = preduce(a2 / pi)926a3t = preduce(a3 / my)927a4t = preduce(a4 / (pi*mx))928a6t = preduce(a6 / (mx*my))929if pdiv(a3t*a3t + 4*a6t):930if p == 2:931t = my*proot(a6t, 2)932else:933t = my*preduce(-a3t*halfmodp)934C = C.rst_transform(0, 0, t)935(a1, a2, a3, a4, a6) = C.a_invariants()936(b2, b4, b6, b8) = C.b_invariants()937my *= pi938iy += 1939a2t = preduce(a2 / pi)940a3t = preduce(a3/my)941a4t = preduce(a4/(pi*mx))942a6t = preduce(a6/(mx*my))943if pdiv(a4t*a4t - 4*a6t*a2t):944if p == 2:945r = mx*proot(a6t*pinv(a2t), 2)946else:947r = mx*preduce(-a4t*pinv(2*a2t))948C = C.rst_transform(r, 0, 0)949(a1, a2, a3, a4, a6) = C.a_invariants()950(b2, b4, b6, b8) = C.b_invariants()951mx *= pi952ix += 1 # and stay in loop953else:954if _pquadroots(a2t, a4t, a6t):955cp = 4956else:957cp = 2958break # exit loop959else:960if _pquadroots(1, a3t, -a6t):961cp = 4962else:963cp = 2964break965KS = KodairaSymbol("I%s*"%(ix+iy-5))966fp = val_disc - ix - iy + 1967break #return968else: # sw == 3969## The cubic has a triple root970verbose("Triple root", t, 1)971## First we change coordinates so that T = 0 mod p972if p == 2:973r = b974elif p == 3:975r = proot(-d, 3)976else:977r = -b * pinv(3)978r = pi*preduce(r)979C = C.rst_transform(r, 0, 0)980(a1, a2, a3, a4, a6) = C.a_invariants()981(b2, b4, b6, b8) = C.b_invariants()982verbose("After third transform %s\n[a1,a2,a3,a4,a6] = %s\nValuations: %s"%([r,0,0],[a1,a2,a3,a4,a6],[pval(ai) for ai in [a1,a2,a3,a4,a6]]), t, 2)983if min(pval(ai) for ai in [a1,a2,a3,a4,a6]) < 0:984raise RuntimeError, "Non-integral model after third transform!"985if pval(a2) < 2 or pval(a4) < 3 or pval(a6) < 4:986raise RuntimeError, "Cubic after transform does not have a triple root at 0"987a3t = preduce(a3/pi2)988a6t = preduce(a6/pi4)989# We test for Type IV*990if not pdiv(a3t*a3t + 4*a6t):991cp = 3 if _pquadroots(1, a3t, -a6t) else 1992KS = KodairaSymbol("IV*")993fp = val_disc - 6994break #return995# Now change coordinates so that p^3|a3, p^5|a6996if p==2:997t = -pi2*proot(a6t, 2)998else:999t = pi2*preduce(-a3t*halfmodp)1000C = C.rst_transform(0, 0, t)1001(a1, a2, a3, a4, a6) = C.a_invariants()1002(b2, b4, b6, b8) = C.b_invariants()1003# We test for types III* and II*1004if pval(a4) < 4:1005## Type III*1006KS = KodairaSymbol("III*")1007fp = val_disc - 71008cp = 21009break #return1010if pval(a6) < 6:1011## Type II*1012KS = KodairaSymbol("II*")1013fp = val_disc - 81014cp = 11015break #return1016if pi_neg is None:1017if principal_flag:1018pi_neg = pi1019else:1020pi_neg = K.uniformizer(P, 'negative')1021pi_neg2 = pi_neg*pi_neg1022pi_neg3 = pi_neg*pi_neg21023pi_neg4 = pi_neg*pi_neg31024pi_neg6 = pi_neg4*pi_neg21025a1 /= pi_neg1026a2 /= pi_neg21027a3 /= pi_neg31028a4 /= pi_neg41029a6 /= pi_neg61030verbose("Non-minimal equation, dividing out...\nNew model is %s"%([a1, a2, a3, a4, a6]), t, 1)1031return (C, p, val_disc, fp, KS, cp, split)103210331034def check_prime(K,P):1035r"""1036Function to check that `P` determines a prime of `K`, and return that ideal.10371038INPUT:10391040- ``K`` -- a number field (including `\QQ`).10411042- ``P`` -- an element of ``K`` or a (fractional) ideal of ``K``.10431044OUTPUT:10451046- If ``K`` is `\QQ`: the prime integer equal to or which generates `P`.10471048- If ``K`` is not `\QQ`: the prime ideal equal to or generated by `P`.10491050.. note::10511052If `P` is not a prime and does not generate a prime, a TypeError is raised.10531054EXAMPLES::10551056sage: from sage.schemes.elliptic_curves.ell_local_data import check_prime1057sage: check_prime(QQ,3)105831059sage: check_prime(QQ,ZZ.ideal(31))1060311061sage: K.<a>=NumberField(x^2-5)1062sage: check_prime(K,a)1063Fractional ideal (a)1064sage: check_prime(K,a+1)1065Fractional ideal (a + 1)1066sage: [check_prime(K,P) for P in K.primes_above(31)]1067[Fractional ideal (5/2*a + 1/2), Fractional ideal (5/2*a - 1/2)]1068"""1069if K is QQ:1070if isinstance(P, (int,long,Integer)):1071P = Integer(P)1072if P.is_prime():1073return P1074else:1075raise TypeError, "%s is not prime"%P1076else:1077if is_Ideal(P) and P.base_ring() is ZZ and P.is_prime():1078return P.gen()1079raise TypeError, "%s is not a prime ideal of %s"%(P,ZZ)10801081if not is_NumberField(K):1082raise TypeError, "%s is not a number field"%K10831084if is_NumberFieldFractionalIdeal(P):1085if P.is_prime():1086return P1087else:1088raise TypeError, "%s is not a prime ideal of %s"%(P,K)10891090if is_NumberFieldElement(P):1091if P in K:1092P = K.ideal(P)1093else:1094raise TypeError, "%s is not an element of %s"%(P,K)1095if P.is_prime():1096return P1097else:1098raise TypeError, "%s is not a prime ideal of %s"%(P,K)10991100raise TypeError, "%s is not a valid prime of %s"%(P,K)110111021103