Path: blob/master/src/sage/schemes/elliptic_curves/gal_reps.py
8820 views
# -*- coding: utf-8 -*-1r"""2Galois representations attached to elliptic curves34Given an elliptic curve `E` over `\QQ`5and a rational prime number `p`, the `p^n`-torsion6`E[p^n]` points of `E` is a representation of the7absolute Galois group `G_{\QQ}` of `\QQ`. As `n` varies8we obtain the Tate module `T_p E` which is a9a representation of `G_\QQ` on a free `\ZZ_p`-module10of rank `2`. As `p` varies the representations11are compatible.1213Currently sage can decide whether the Galois module14`E[p]` is reducible, i.e., if `E` admits an isogeny15of degree `p`, and whether the image of16the representation on `E[p]` is surjective onto17`\text{Aut}(E[p]) = GL_2(\mathbb{F}_p)`.1819The following are the most useful functions for the class ``GaloisRepresentation``.2021For the reducibility:2223- ``is_reducible(p)``24- ``is_irreducible(p)``25- ``reducible_primes()``2627For the image:2829- ``is_surjective(p)``30- ``non_surjective()``31- ``image_type(p)``3233For the classification of the representation3435- ``is_semistable(p)``36- ``is_unramified(p, ell)``37- ``is_crystalline(p)``3839EXAMPLES::4041sage: E = EllipticCurve('196a1')42sage: rho = E.galois_representation()43sage: rho.is_irreducible(7)44True45sage: rho.is_reducible(3)46True47sage: rho.is_irreducible(2)48True49sage: rho.is_surjective(2)50False51sage: rho.is_surjective(3)52False53sage: rho.is_surjective(5)54True55sage: rho.reducible_primes()56[3]57sage: rho.non_surjective()58[2, 3]59sage: rho.image_type(2)60'The image is cyclic of order 3.'61sage: rho.image_type(3)62'The image is contained in a Borel subgroup as there is a 3-isogeny.'63sage: rho.image_type(5)64'The image is all of GL_2(F_5).'6566For semi-stable curve it is known that the representation is67surjective if and only if it is irreducible::6869sage: E = EllipticCurve('11a1')70sage: rho = E.galois_representation()71sage: rho.non_surjective()72[5]73sage: rho.reducible_primes()74[5]7576For cm curves it is not true that there are only finitely many primes for which the77Galois representation mod p is surjective onto `GL_2(\mathbb{F}_p)`::7879sage: E = EllipticCurve('27a1')80sage: rho = E.galois_representation()81sage: rho.non_surjective()82[0]83sage: rho.reducible_primes()84[3]85sage: E.has_cm()86True87sage: rho.image_type(11)88'The image is contained in the normalizer of a non-split Cartan group. (cm)'8990REFERENCES:9192.. [Se1] Jean-Pierre Serre,93Propriétés galoisiennes des points d'ordre fini94des courbes elliptiques.95Invent. Math. 15 (1972), no. 4, 259--331.96.. [Se2] Jean-Pierre Serre,97Sur les représentations modulaires de degré982 de `\text{Gal}(\bar\QQ/\QQ)`.99Duke Math. J. 54 (1987), no. 1, 179--230.100.. [Co] Alina Carmen Cojocaru,101On the surjectivity of the Galois representations102associated to non-CM elliptic curves.103With an appendix by Ernst Kani.104Canad. Math. Bull. 48 (2005), no. 1, 16--31.105106AUTHORS:107108- chris wuthrich (02/10) - moved from ell_rational_field.py.109110"""111112######################################################################113# Copyright (C) 2010 William Stein <[email protected]>114#115# Distributed under the terms of the GNU General Public License (GPL)116#117# This code is distributed in the hope that it will be useful,118# but WITHOUT ANY WARRANTY; without even the implied warranty of119# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU120# General Public License for more details.121#122# The full text of the GPL is available at:123#124# http://www.gnu.org/licenses/125######################################################################126127from sage.structure.sage_object import SageObject128import sage.rings.arith as arith129import sage.misc.misc as misc130import sage.rings.all as rings131from sage.rings.all import RealField, GF132from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing133134from math import sqrt135from sage.libs.pari.all import pari136137def _ex_set(p):138"""139Gives the set of the only values of trace^2/det140that appear in a exceptional subgroup in PGL_2(F_p).141142EXAMPLES::143144sage: from sage.schemes.elliptic_curves.gal_reps import _ex_set145sage: for p in prime_range(3,30): print p, _ex_set(p)1463 [0, 1, 2, 1]1475 [0, 1, 2, 4]1487 [0, 1, 2, 4]14911 [0, 1, 2, 4, 9, 5]15013 [0, 1, 2, 4]15117 [0, 1, 2, 4]15219 [0, 1, 2, 4, 16, 6]15323 [0, 1, 2, 4]15429 [0, 1, 2, 4, 25, 7]155"""156k = GF(p)157res = [ k(0), k(1), k(2), k(4) ]158R = k['X']159f = R([1,-3,1]) #(X**2 - 3*X+1)160ro = f.roots()161for a in ro:162if a[0] not in res:163res.append(a[0])164return res165166167class GaloisRepresentation(SageObject):168r"""169The compatible family of Galois representation170attached to an elliptic curve over the rational numbers.171172Given an elliptic curve `E` over `\QQ` and a rational173prime number `p`, the `p^n`-torsion `E[p^n]` points of174`E` is a representation of the absolute Galois group.175As `n` varies we obtain the Tate module `T_p E` which is176a representation of the absolute Galois group on a free177`\ZZ_p`-module of rank `2`. As `p` varies the178representations are compatible.179180EXAMPLES::181182sage: rho = EllipticCurve('11a1').galois_representation()183sage: rho184Compatible family of Galois representations associated to the Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field185186"""187188def __init__(self, E):189r"""190see ``GaloisRepresentation`` for documentation191192EXAMPLES::193194sage: rho = EllipticCurve('11a1').galois_representation()195sage: loads(rho.dumps()) == rho196True197198"""199self.__image_type = {}200self._E = E201202def __repr__(self):203r"""204string representation of the class205206EXAMPLES::207208sage: rho = EllipticCurve([0,1]).galois_representation()209sage: rho210Compatible family of Galois representations associated to the Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field211212"""213return "Compatible family of Galois representations associated to the " + repr(self._E)214215def __eq__(self,other):216r"""217Compares two Galois representations.218We define tho compatible families of representations219attached to elliptic curves to be isomorphic if the curves are equal.220221EXAMPLES::222223sage: rho = EllipticCurve('11a1').galois_representation()224sage: rho2 = EllipticCurve('11a2').galois_representation()225sage: rho == rho226True227sage: rho == rho2228False229sage: rho == 34230False231232"""233# if rho_E = rho_E' then the L-functions agree,234# so E and E' are isogenous235# except for p=2236# the mod p representations will be different237# for p dividing the degree of the isogeny238# anyway, there should not be a _compatible_239# isomorphism between rho and rho' unless E240# is isomorphic to E'241# Note that rho can not depend on the Weierstrass model242if not type(self) == type(other):243return False244return self._E.is_isomorphic(other._E)245246def elliptic_curve(self):247r"""248The elliptic curve associated to this representation.249250EXAMPLES::251252sage: E = EllipticCurve('11a1')253sage: rho = E.galois_representation()254sage: rho.elliptic_curve() == E255True256257"""258from copy import copy259return copy(self._E)260261#####################################################################262# reducibility263#####################################################################264265def is_reducible(self, p):266r"""267Return True if the mod-p representation is268reducible. This is equivalent to the existence of an269isogeny defined over `\QQ` of degree `p` from the270elliptic curve.271272INPUT:273274- ``p`` - a prime number275276OUTPUT:277278- a boolean279280The answer is cached.281282EXAMPLES::283284sage: rho = EllipticCurve('121a').galois_representation()285sage: rho.is_reducible(7)286False287sage: rho.is_reducible(11)288True289sage: EllipticCurve('11a').galois_representation().is_reducible(5)290True291sage: rho = EllipticCurve('11a2').galois_representation()292sage: rho.is_reducible(5)293True294sage: EllipticCurve('11a2').torsion_order()2951296297"""298try:299return self.__is_reducible[p]300except AttributeError:301self.__is_reducible = {}302except KeyError:303pass304305if not arith.is_prime(p):306raise ValueError('p (=%s) must be prime' % p)307# we do is_surjective first, since this is308# much easier than computing isogeny_class309t = self._is_surjective(p, A=-1)310if t:311self.__is_reducible[p] = False312return False # definitely not reducible313isogeny_matrix = self._E.isogeny_class().matrix(fill=True)314v = isogeny_matrix.row(0) # first row315for a in v:316if a != 0 and a % p == 0:317self.__is_reducible[p] = True318return True319self.__is_reducible[p] = False320return False321322def is_irreducible(self, p):323r"""324Return True if the mod p representation is irreducible.325326INPUT:327328- ``p`` - a prime number329330OUTPUT:331332- a boolean333334EXAMPLES::335336sage: rho = EllipticCurve('37b').galois_representation()337sage: rho.is_irreducible(2)338True339sage: rho.is_irreducible(3)340False341sage: rho.is_reducible(2)342False343sage: rho.is_reducible(3)344True345"""346return not self.is_reducible(p)347348def reducible_primes(self):349r"""350Returns a list of the primes `p` such that the mod-`p`351representation is reducible. For all other primes the352representation is irreducible.353354EXAMPLES::355356sage: rho = EllipticCurve('225a').galois_representation()357sage: rho.reducible_primes()358[3]359"""360try:361return self.__reducible_primes362except AttributeError:363pass364isocls = self._E.isogeny_class()365X = set(isocls.matrix().list())366R = [p for p in X if arith.is_prime(p)]367self.__reducible_primes = R368return R369370#####################################################################371# image372#####################################################################373374def is_surjective(self, p, A=1000):375r"""376Return True if the mod-p representation is377surjective onto `Aut(E[p]) = GL_2(\mathbb{F}_p)`.378379False if it is not, or None if we were unable to380determine whether it is or not.381382INPUT:383384- ``p`` - int (a prime number)385386- ``A`` - int (a bound on the number of a_p to use)387388OUTPUT:389390- boolean. True if the mod-p representation is surjective391and False if not.392393The answer is cached.394395EXAMPLES::396397sage: rho = EllipticCurve('37b').galois_representation()398sage: rho.is_surjective(2)399True400sage: rho.is_surjective(3)401False402403::404405sage: rho = EllipticCurve('121a1').galois_representation()406sage: rho.non_surjective()407[11]408sage: rho.is_surjective(5)409True410sage: rho.is_surjective(11)411False412413sage: rho = EllipticCurve('121d1').galois_representation()414sage: rho.is_surjective(5)415False416sage: rho.is_surjective(11)417True418419Here is a case, in which the algorithm does not return an answer::420421sage: rho = EllipticCurve([0,0,1,2580,549326]).galois_representation()422sage: rho.is_surjective(7)423424In these cases, one can use image_type to get more information about the image::425426sage: rho.image_type(7)427'The image is contained in the normalizer of a split Cartan group.'428429430REMARKS:4314321. If `p \geq 5` then the mod-p representation is433surjective if and only if the p-adic representation is434surjective. When `p = 2, 3` there are435counterexamples. See papers of Dokchitsers and Elkies436for more details.4374382. For the primes `p=2` and 3, this will always answer either439True or False. For larger primes it might give None.440441"""442if not arith.is_prime(p):443raise TypeError("p (=%s) must be prime." % p)444A = int(A)445key = (p, A)446try:447return self.__is_surjective[key]448except KeyError:449pass450except AttributeError:451self.__is_surjective = {}452453ans = self._is_surjective(p, A)454self.__is_surjective[key] = ans455return ans456457458def _is_surjective(self, p, A):459r"""460helper function for ``is_surjective``461462The value of `A` is as before, except that463`A=-1` is a special code to stop before464testing reducibility (to avoid an infinite loop).465466EXAMPLES::467468sage: rho = EllipticCurve('37b').galois_representation()469sage: rho._is_surjective(7,100)470True471472TEST for :trac:`8451`::473474sage: E = EllipticCurve('648a1')475sage: rho = E.galois_representation()476sage: rho._is_surjective(5,1000)477478"""479T = self._E.torsion_subgroup().order()480if T % p == 0 and p != 2 :481# we could probably determine the group structure directly482self.__image_type[p] = "The image is meta-cyclic inside a Borel subgroup as there is a %s-torsion point on the curve."%p483return False484485R = rings.PolynomialRing(self._E.base_ring(), 'x')486x = R.gen()487488if p == 2:489# E is isomorphic to [0,b2,0,8*b4,16*b6]490b2,b4,b6,b8=self._E.b_invariants()491f = x**3 + b2*x**2 + 8*b4*x + 16*b6492if not f.is_irreducible():493if len(f.roots()) > 2:494self.__image_type[p] = "The image is trivial as all 2-torsion points are rational."495else:496self.__image_type[p] = "The image is cyclic of order 2 as there is exactly one rational 2-torsion point."497return False #, '2-torsion'498if arith.is_square(f.discriminant()):499self.__image_type[p] = "The image is cyclic of order 3."500return False #, "A3"501self.__image_type[p] = "The image is all of GL_2(F_2), i.e. a symmetric group of order 6."502return True #, None503504if p == 3:505# Algorithm: Let f be the 3-division polynomial, which is506# a polynomial of degree 4. Then I claim that this507# polynomial has Galois group S_4 if and only if the508# representation rhobar_{E,3} is surjective. If the group509# is S_4, then S_4 is a quotient of the image of510# rhobar_{E,3}. Since S_4 has order 24 and GL_2(F_3)511# has order 48, the only possibility we have to consider512# is that the image of rhobar is isomorphic to S_4.513# But this is not the case because S_4 is not a subgroup514# of GL_2(F_3). If it were, it would be normal, since515# it would have index 2. But there is a *unique* normal516# subgroup of GL_2(F_3) of index 2, namely SL_2(F_3),517# and SL_2(F_3) is not isomorphic to S_4 (S_4 has a normal518# subgroup of index 2 and SL_2(F_3) does not.)519# (What's a simple way to see that SL_2(F_3) is the520# unique index-2 normal subgroup? I didn't see an obvious521# reason, so just used the NormalSubgroups command in MAGMA522# and it output exactly one of index 2.)523524#sage: G = SymmetricGroup(4)525#sage: [H.group_id() for H in G.conjugacy_classes_subgroups()]526#[[1, 1], [2, 1], [2, 1], [3, 1], [4, 2], [4, 2], [4, 1], [6, 1], [8, 3], [12, 3], [24, 12]]527#sage: G = GL(2,GF(3)).as_matrix_group().as_permutation_group()528#sage: [H.group_id() for H in G.conjugacy_classes_subgroups()]529#[[1, 1], [2, 1], [2, 1], [3, 1], [4, 2], [4, 1], [6, 2], [6, 1], [6, 1], [8, 4], [8, 1], [8, 3], [12, 4], [16, 8], [24, 3], [48, 29]]530531# Here's Noam Elkies proof for the other direction:532533#> Let E be an elliptic curve over Q. Is the mod-3534#> representation E[3] surjective if and only if the535#> (degree 4) division polynomial has Galois group S_4? I536#> can see why the group being S_4 implies the537#> representation is surjective, but the converse is not538#> clear to me.539# I would have thought that this is the easier part: to540# say that E[3] is surjective is to say the 3-torsion541# field Q(E[3]) has Galois group GL_2(Z/3) over Q. Let542# E[3]+ be the subfield fixed by the element -1 of543# GL_2(Z/3). Then E[3] has Galois group PGL_2(Z/3), which544# is identified with S_4 by its action on the four545# 3-element subgroups of E[3]. Each such subgroup is in546# turn determined by the x-coordinate shared by its two547# nonzero points. So, if E[3] is surjective then any548# permutation of those x-coordinates is realized by some549# element of Gal(E[3]+/Q). Thus the Galois group of the550# division polynomial (whose roots are those551# x-coordinates) maps surjectively to S_4, which means it552# equals S_4.553554f = self._E.division_polynomial(3)555if not f.is_irreducible():556return False #, "reducible_3-divpoly"557n = pari(f).polgalois()[0]558if n == 24:559self.__image_type[p] = "The image is all of GL_2(F_3)."560return True #, None561else:562return False #, "3-divpoly_galgroup_order_%s"%n563564if self._E.has_cm():565return False #, "CM"566567# Now we try to prove that the rep IS surjective.568569Np = self._E.conductor() * p570signs = []571# this follows Proposition 19 in Serre.572# there was a bug in the original implementation,573# counter-examples were 324b1 and 324d1, 648a1 and 648c1 for p=5574exclude_exceptional_image = False575ex_setp = _ex_set(p)576ell = 4577k = GF(p)578579if A == -1:580Am = 1000581else:582Am = A583while ell < Am:584ell = arith.next_prime(ell)585if Np % ell != 0:586a_ell = self._E.ap(ell)587if a_ell % p != 0:588if not exclude_exceptional_image:589u = k(a_ell)**2 * k(ell)**(-1)590if u not in ex_setp:591exclude_exceptional_image = True592s = arith.kronecker(a_ell**2 - 4*ell, p)593if s != 0 and s not in signs:594signs.append(s)595if len(signs) == 2 and exclude_exceptional_image:596self.__image_type[p] = "The image is all of GL_2(F_%s)."%p597return True #,None598599if A == -1: # we came in from is reducible. Now go out with False600return False601602if self.is_reducible(p):603return False #, Borel604605# if we reach this, then we do not know if it is surjective. Most likely606# not but we can't be certain. See trac 11271.607misc.verbose("We can not conclude if the representation is surjective or not. Increasing the parameter A may help.")608return None609610def non_surjective(self, A=1000):611r"""612Returns a list of primes p such that the mod-p representation613*might* not be surjective. If `p` is not in the returned list,614then the mod-p representation is provably surjective.615616By a theorem of Serre, there are only finitely617many primes in this list, except when the curve has618complex multiplication.619620If the curve has CM, we simply return the621sequence [0] and do no further computation.622623INPUT:624625- ``A`` - an integer (default 1000). By increasing this parameter626the resulting set might get smaller.627628OUTPUT:629630- ``list`` - if the curve has CM, returns [0].631Otherwise, returns a list of primes where mod-p representation is632very likely not surjective. At any prime not in this list, the633representation is definitely surjective.634635EXAMPLES::636637sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A638sage: E.galois_representation().non_surjective() # CM curve639[0]640641::642643sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11)644sage: E.galois_representation().non_surjective()645[5]646647sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A648sage: E.galois_representation().non_surjective()649[]650651sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C652sage: E.galois_representation().non_surjective()653[13]654655::656657sage: E = EllipticCurve([1,-1,1,-9965,385220]) # 9999a1658sage: rho = E.galois_representation()659sage: rho.non_surjective()660[2]661662sage: E = EllipticCurve('324b1')663sage: rho = E.galois_representation()664sage: rho.non_surjective()665[3, 5]666667ALGORITHM:668We first find an upper bound `B` on the possible primes. If `E`669is semi-stable, we can take `B=11` by a result of Mazur. There is670a bound by Serre in the case that the `j`-invariant is not integral671in terms of the smallest prime of good reduction. Finally672there is an unconditional bound by Cojocaru, but which depends673on the conductor of `E`.674For the prime below that bound we call ``is_surjective``.675676"""677if self._E.has_cm():678misc.verbose("cm curve")679return [0]680N = self._E.conductor()681if self._E.is_semistable():682# Mazur's bound683C = 11684misc.verbose("semistable -- so bound is 11")685elif not self._E.j_invariant().is_integral():686# prop 24 in Serre687vs = self._E.j_invariant().denominator().prime_factors()688C1 = arith.gcd([-arith.valuation(self._E.j_invariant(),v) for v in vs])689p0 = 2690while self._E.has_bad_reduction(p0):691p0 = arith.next_prime(p0+1)692C2 = (sqrt(p0)+1)**8693C = max(C1,C2)694misc.verbose("j is not integral -- Serre's bound is %s"%C)695C3 = 1 + 4*sqrt(6)*int(N)/3 * sqrt(misc.mul([1+1.0/int(p) for p,_ in arith.factor(N)]))696C = min(C,C3)697misc.verbose("conductor = %s, and bound is %s"%(N,C))698else:699# Cojocaru's bound (depends on the conductor)700C = 1 + 4*sqrt(6)*int(N)/3 * sqrt(misc.mul([1+1.0/int(p) for p,_ in arith.factor(N)]))701misc.verbose("conductor = %s, and bound is %s"%(N,C))702B = []703p = 2704while p <= C:705t = self.is_surjective(p, A=A)706misc.verbose("(%s,%s)"%(p,t))707# both False and None will be appended here.708if not t:709B.append(p)710p = arith.next_prime(p)711return B712713def image_type(self, p):714r"""715Returns a string describing the image of the716mod-p representation.717The result is provably correct, but only718indicates what sort of an image we have. If719one wishes to determine the exact group one720needs to work a bit harder. The probabilistic721method of image_classes or Sutherland's galrep722package can give a very good guess what the723image should be.724725INPUT:726727- ``p`` a prime number728729OUTPUT:730731- a string.732733EXAMPLES ::734735sage: E = EllipticCurve('14a1')736sage: rho = E.galois_representation()737sage: rho.image_type(5)738'The image is all of GL_2(F_5).'739740sage: E = EllipticCurve('11a1')741sage: rho = E.galois_representation()742sage: rho.image_type(5)743'The image is meta-cyclic inside a Borel subgroup as there is a 5-torsion point on the curve.'744745sage: EllipticCurve('27a1').galois_representation().image_type(5)746'The image is contained in the normalizer of a non-split Cartan group. (cm)'747sage: EllipticCurve('30a1').galois_representation().image_type(5)748'The image is all of GL_2(F_5).'749sage: EllipticCurve("324b1").galois_representation().image_type(5)750'The image in PGL_2(F_5) is the exceptional group S_4.'751752sage: E = EllipticCurve([0,0,0,-56,4848])753sage: rho = E.galois_representation()754755sage: rho.image_type(5)756'The image is contained in the normalizer of a split Cartan group.'757758sage: EllipticCurve('49a1').galois_representation().image_type(7)759'The image is contained in a Borel subgroup as there is a 7-isogeny.'760761sage: EllipticCurve('121c1').galois_representation().image_type(11)762'The image is contained in a Borel subgroup as there is a 11-isogeny.'763sage: EllipticCurve('121d1').galois_representation().image_type(11)764'The image is all of GL_2(F_11).'765sage: EllipticCurve('441f1').galois_representation().image_type(13)766'The image is contained in a Borel subgroup as there is a 13-isogeny.'767768sage: EllipticCurve([1,-1,1,-5,2]).galois_representation().image_type(5)769'The image is contained in the normalizer of a non-split Cartan group.'770sage: EllipticCurve([0,0,1,-25650,1570826]).galois_representation().image_type(5)771'The image is contained in the normalizer of a split Cartan group.'772sage: EllipticCurve([1,-1,1,-2680,-50053]).galois_representation().image_type(7) # the dots (...) in the output fix #11937 (installed 'Kash' may give additional output); long time (2s on sage.math, 2014)773'The image is a... group of order 18.'774sage: EllipticCurve([1,-1,0,-107,-379]).galois_representation().image_type(7) # the dots (...) in the output fix #11937 (installed 'Kash' may give additional output); long time (1s on sage.math, 2014)775'The image is a... group of order 36.'776sage: EllipticCurve([0,0,1,2580,549326]).galois_representation().image_type(7)777'The image is contained in the normalizer of a split Cartan group.'778779Test :trac:`14577`::780781sage: EllipticCurve([0, 1, 0, -4788, 109188]).galois_representation().image_type(13)782'The image in PGL_2(F_13) is the exceptional group S_4.'783784Test :trac:`14752`::785786sage: EllipticCurve([0, 0, 0, -1129345880,-86028258620304]).galois_representation().image_type(11)787'The image is contained in the normalizer of a non-split Cartan group.'788789For `p=2`::790791sage: E = EllipticCurve('11a1')792sage: rho = E.galois_representation()793sage: rho.image_type(2)794'The image is all of GL_2(F_2), i.e. a symmetric group of order 6.'795796sage: rho = EllipticCurve('14a1').galois_representation()797sage: rho.image_type(2)798'The image is cyclic of order 2 as there is exactly one rational 2-torsion point.'799800sage: rho = EllipticCurve('15a1').galois_representation()801sage: rho.image_type(2)802'The image is trivial as all 2-torsion points are rational.'803804sage: rho = EllipticCurve('196a1').galois_representation()805sage: rho.image_type(2)806'The image is cyclic of order 3.'807808`p=3`::809810sage: rho = EllipticCurve('33a1').galois_representation()811sage: rho.image_type(3)812'The image is all of GL_2(F_3).'813814sage: rho = EllipticCurve('30a1').galois_representation()815sage: rho.image_type(3)816'The image is meta-cyclic inside a Borel subgroup as there is a 3-torsion point on the curve.'817818sage: rho = EllipticCurve('50b1').galois_representation()819sage: rho.image_type(3)820'The image is contained in a Borel subgroup as there is a 3-isogeny.'821822sage: rho = EllipticCurve('3840h1').galois_representation()823sage: rho.image_type(3)824'The image is contained in a dihedral group of order 8.'825826sage: rho = EllipticCurve('32a1').galois_representation()827sage: rho.image_type(3)828'The image is a semi-dihedral group of order 16, gap.SmallGroup([16,8]).'829830831ALGORITHM: Mainly based on Serre's paper.832"""833if not arith.is_prime(p):834raise TypeError("p (=%s) must be prime." % p)835try:836return self.__image_type[p]837except KeyError:838pass839except AttributeError:840self.__image_type = {}841842# we eliminate step by step possibilities.843# checks if the image is all of GL_2, while doing so, it may detect certain other classes already844# see _is_surjective. It will have set __image_type845846ans = self._is_surjective(p, A=1000)847try:848return self.__image_type[p]849except KeyError:850pass851852# check if the rep is reducible853854if self.is_reducible(p):855self.__image_type[p] = "The image is contained in a Borel subgroup as there is a %s-isogeny."%p856return self.__image_type[p]857858# if we are then the image of rho is not surjective and not contained in a Borel subgroup859# there are three cases left:860# it could be in a normalizer of a split Cartan,861# normalizer of a non-split Cartan,862# or the image in PGL_2 is one of the three exceptional groups A_4 S_4 A_5863864non_split_str = "The image is contained in the normalizer of a non-split Cartan group."865split_str = "The image is contained in the normalizer of a split Cartan group."866s4_str = "The image in PGL_2(F_%s) is the exceptional group S_4."%p867a4_str = "The image in PGL_2(F_%s) is the exceptional group A_4."%p868a5_str = "The image in PGL_2(F_%s) is the exceptional group A_5."%p869870# we first treat p=3 and 5 seperately. p=2 has already been done.871872if p == 3:873# this implies that the image of rhobar in PGL_2 = S_4874# determines completely the image of rho875f = self._E.division_polynomial(3)876if not f.is_irreducible():877# must be a product of two polynomials of degree 2878self.__image_type[p] = "The image is contained in a dihedral group of order 8."879return self.__image_type[p]880n = pari(f).polgalois()[0]881# the following is due to a simple classification of all subgroups of GL_2(F_3)882if n == 2:883self.__image_type[p] = "The image is a cyclic group of order 4."884elif n == 4:885for ell in prime_range(5,1000):886if ell % 3 == 2 and self._E.ap(ell) % 3 != 0:887# there is an element of order 8 in the image888self.__image_type[p] = "The image is a cyclic group of order 8."889return self.__image_type[p]890self.__image_type[p] = "The image is a group of order 8, most likely a quaternion group."891elif n == 6:892self.__image_type[p] = "The image is a dihedral group of order 12."893elif n == 8:894self.__image_type[p] = "The image is a semi-dihedral group of order 16, gap.SmallGroup([16,8])."895elif n == 12:896self.__image_type[p] = "The image is SL_2(F_3)."897else:898raise RuntimeError("Bug in image_type for p = 3.")899return self.__image_type[p]900901# we also eliminate cm curves902903if self._E.has_cm():904if self._E.is_good(p) and self._E.is_ordinary(p):905self.__image_type[p] = split_str + " (cm)"906return self.__image_type[p]907if self._E.is_supersingular(p):908self.__image_type[p] = non_split_str + " (cm)"909return self.__image_type[p]910else:911# if the reduction is bad (additive nec.) then912# the image is in a Borel subgroup and we should have found this before913raise NotImplementedError("image_type is not implemented for cm-curves at bad prime.")914915# now to p=5 where a lot of non-standard thing happen916# we run through primes if we hit an element of order 6 in PGL_2, we know that it is the normaliser of a NON-split Cartan917# if we find both an element of order 3 and one of order 4, we know that we have a S_4 in PGL_2918919if p == 5:920# we filter here a few cases and leave the rest to the computation of the galois group later921ell = 1922k = GF(p)923Np = self._E.conductor()*p924has_an_el_order_4 = False925has_an_el_order_3 = False926while ell < 10000:927ell = arith.next_prime(ell)928if Np % ell != 0:929a_ell = self._E.ap(ell)930u = k(a_ell)**2 * k(ell)**(-1)931if u == 3:932misc.verbose("found an element of order 6",2)933# found an element of order 6:934self.__image_type[p] = non_split_str935return self.__image_type[p]936937if u == 2 and not has_an_el_order_4:938# found an element of order 4939misc.verbose("found an element of order 4",2)940has_an_el_order_4 = True941if has_an_el_order_3:942self.__image_type[p] = s4_str943return self.__image_type[p]944945if u == 1 and not has_an_el_order_3:946# found an element of order 3947misc.verbose("found an element of order 3",2)948has_an_el_order_3 = True949if has_an_el_order_4:950self.__image_type[p] = s4_str951return self.__image_type[p]952953misc.verbose("p=5 and we could not determine the image, yet", 2)954# we have not yet determined the image, there are only the following possible subgroups of PGL_2955# (unless we were unlucky and none of the elements of order 6 showed up above, for instance)956# A_4 of order 12 with elements of order 2 and 3957# D_4 of order 8 with elements of order 2 and 4 (normalizer of the SPLIT)958# Z/2 x Z/2 of order 4 with elements of order 2 (in the normalizer of the SPLIT)959# S_3 of order 6 with elements of order 2 and 3 (inside the normalizer of the NON-split)960961# we compute the splitting field of the 5-division polynomial. Its degree is equal to the above order or the double of it.962# That allows us to determine almost all cases.963964f = self._E.division_polynomial(5)965from sage.rings.number_field.splitting_field import SplittingFieldAbort966try:967K = f.splitting_field('x', degree_multiple=240, abort_degree=24)968except SplittingFieldAbort:969pass970else:971if K.degree() in (4,8,16):972self.__image_type[p] = split_str973return self.__image_type[p]974if K.degree() == 24:975self.__image_type[p] = a4_str976return self.__image_type[p]977if K.degree() == 6:978self.__image_type[p] = non_split_str979return self.__image_type[p]980if K.degree() == 12:981# PGL - image could be a S_3 in the normalizer of the split or A4982self.__image_type[p] = "The image is of order 6 or 12. Probably contained in the normalizer of the split Cartan g."983return self.__image_type[p]984985## now E has no cm, is not semi-stable,986## p > 5,987## rho_p it is (most probably) not surjective , it is not contained in a Borel subgroup.988# trying to detect that the image is not exceptional in PGL_2989# this uses Serre 2.6.iii and Prop 19990# the existence of certain classes in the image rules out certain cases.991# we run over the small prime until we are left with only one case992# for p = 5 this could never distinguish any from an exceptional S_4 or A_4,993# that is why the case 5 is treated a part before994995else:996ex_setp = _ex_set(p)997ell = 1998k = GF(p)999Np = self._E.conductor()*p1000could_be_exc = 11001could_be_split = 11002could_be_non_split = 11003# loops over primes as long as we still have two options left1004while ell < 10000 and (could_be_exc + could_be_split + could_be_non_split > 1):1005ell = arith.next_prime(ell)1006if Np % ell != 0:1007a_ell = self._E.ap(ell)1008u = k(a_ell)**2 * k(ell)**(-1)1009if (u not in ex_setp) and could_be_exc == 1:1010# it can not be in the exceptional1011misc.verbose("the image cannot be exceptional, found u=%s"%u,2)1012could_be_exc = 01013if a_ell != 0 and arith.kronecker(a_ell**2 - 4*ell,p) == 1 and could_be_non_split == 1:1014# it can not be in the noramlizer of the non-split Cartan1015misc.verbose("the image cannot be non-split, found u=%s"%u,2)1016could_be_non_split = 01017if a_ell != 0 and arith.kronecker(a_ell**2 - 4*ell,p) == -1 and could_be_split == 1:1018# it can not be in the noramlizer of the split Cartan1019misc.verbose("the image cannot be split, found u=%s"%u,2)1020could_be_split = 010211022assert could_be_exc + could_be_split + could_be_non_split > 0, "bug in image_type."10231024if could_be_exc + could_be_split + could_be_non_split == 1:1025# it is only one of the three cases:1026if could_be_split == 1 :1027self.__image_type[p] = split_str1028return self.__image_type[p]1029if could_be_non_split == 1 :1030self.__image_type[p] = non_split_str1031return self.__image_type[p]1032if could_be_exc == 1:1033# here we can distinguish further1034could_be_a4 = 11035could_be_s4 = 11036could_be_a5 = 11037if p % 5 != 1 and p % 5 != 4 :1038could_be_a5 = 01039# elements of order 5 # bug corrected see trac 145771040R = k['X']1041f = R([1,-3,1]) #(X**2 - 3*X+1)1042el5 = f.roots()1043# loops over primes as long as we still have two options left1044while ell < 10000 and (could_be_s4 + could_be_a4 + could_be_a5 > 1):1045ell = arith.next_prime(ell)1046if Np % ell != 0:1047a_ell = self._E.ap(ell)1048u = k(a_ell)**2 * k(ell)**(-1)1049if u == 2:1050# it can not be A4 not A5 as they have no elements of order 41051could_be_a4 = 01052could_be_a5 = 01053if u in el5 :1054# it can not be A4 or S4 as they have no elements of order 51055could_be_a4 = 01056could_be_s4 = 010571058assert (could_be_s4 + could_be_a4 + could_be_a5 > 0), "bug in image_type."10591060if could_be_s4 + could_be_a4 + could_be_a5 == 1:1061if could_be_s4 == 1:1062self.__image_type[p] = s4_str1063return self.__image_type[p]1064if could_be_a4 == 1:1065self.__image_type[p] = a4_str1066return self.__image_type[p]1067if could_be_a5 == 1:1068self.__image_type[p] = a5_str1069return self.__image_type[p]10701071else:1072self.__image_type[p] = "The image in PGL_2(F_%s) is an exceptional group A_4, S_4 or A_5, but we could not determine which one."%p1073return self.__image_type[p]10741075# If all fails, we probably have a fairly small group and we can try to detect it using the galois_group1076if p <= 13:1077K = self._E.division_field(p, 'z')1078d = K.absolute_degree()10791080misc.verbose("field of degree %s. try to compute galois group"%(d),2)1081try:1082G = K.galois_group()1083except Exception:1084self.__image_type[p] = "The image is a group of order %s."%d1085return self.__image_type[p]10861087else:1088if G.is_abelian():1089ab = ""1090else:1091ab = "non-"1092self.__image_type[p] = "The image is a " + ab + "abelian group of order %s."%G.order()1093return self.__image_type[p]10941095## everything failed :10961097self.__image_type[p] = "The image could not be determined. Sorry."1098return self.__image_type[p]109911001101def image_classes(self,p,bound=10000):1102r"""1103This function returns, given the representation `\rho`1104a list of `p` values that add up to 1, representing the1105frequency of the conjugacy classes of the projective image1106of `\rho` in `PGL_2(\mathbb{F}_p)`.11071108Let `M` be a matrix in `GL_2(\mathbb{F}_p)`, then define1109`u(M) = \text{tr}(M)^2/\det(M)`, which only depends on the1110conjugacy class of `M` in `PGL_2(\mathbb{F}_p)`. Hence this defines1111a map `u: PGL_2(\mathbb{F}_p) \to \mathbb{F}_p`, which is almost1112a bijection between conjugacy classes of the source1113and `\mathbb{F}_p` (the elements of order `p` and the identity1114map to `4` and both classes of elements of order 2 map to 0).11151116This function returns the frequency with which the values of1117`u` appeared among the images of the Frobenius elements1118`a_{\ell}`at `\ell` for good primes `\ell\neq p` below a1119given ``bound``.11201121INPUT:11221123- a prime ``p``11241125- a natural number ``bound`` (optional, default=10000)11261127OUTPUT:11281129- a list of `p` real numbers in the interval `[0,1]` adding up to 111301131EXAMPLES::11321133sage: E = EllipticCurve('14a1')1134sage: rho = E.galois_representation()1135sage: rho.image_classes(5)1136[0.2095, 0.1516, 0.2445, 0.1728, 0.2217]11371138sage: E = EllipticCurve('11a1')1139sage: rho = E.galois_representation()1140sage: rho.image_classes(5)1141[0.2467, 0.0000, 0.5049, 0.0000, 0.2484]11421143::11441145sage: EllipticCurve('27a1').galois_representation().image_classes(5)1146[0.5839, 0.1645, 0.0000, 0.1702, 0.08143]1147sage: EllipticCurve('30a1').galois_representation().image_classes(5)1148[0.1956, 0.1801, 0.2543, 0.1728, 0.1972]1149sage: EllipticCurve('32a1').galois_representation().image_classes(5)1150[0.6319, 0.0000, 0.2492, 0.0000, 0.1189]1151sage: EllipticCurve('900a1').galois_representation().image_classes(5)1152[0.5852, 0.1679, 0.0000, 0.1687, 0.07824]1153sage: EllipticCurve('441a1').galois_representation().image_classes(5)1154[0.5860, 0.1646, 0.0000, 0.1679, 0.08150]1155sage: EllipticCurve('648a1').galois_representation().image_classes(5)1156[0.3945, 0.3293, 0.2388, 0.0000, 0.03749]11571158::11591160sage: EllipticCurve('784h1').galois_representation().image_classes(7)1161[0.5049, 0.0000, 0.0000, 0.0000, 0.4951, 0.0000, 0.0000]1162sage: EllipticCurve('49a1').galois_representation().image_classes(7)1163[0.5045, 0.0000, 0.0000, 0.0000, 0.4955, 0.0000, 0.0000]11641165sage: EllipticCurve('121c1').galois_representation().image_classes(11)1166[0.1001, 0.0000, 0.0000, 0.0000, 0.1017, 0.1953, 0.1993, 0.0000, 0.0000, 0.2010, 0.2026]1167sage: EllipticCurve('121d1').galois_representation().image_classes(11)1168[0.08869, 0.07974, 0.08706, 0.08137, 0.1001, 0.09439, 0.09764, 0.08218, 0.08625, 0.1017, 0.1009]11691170sage: EllipticCurve('441f1').galois_representation().image_classes(13)1171[0.08232, 0.1663, 0.1663, 0.1663, 0.08232, 0.0000, 0.1549, 0.0000, 0.0000, 0.0000, 0.0000, 0.1817, 0.0000]11721173REMARKS:11741175Conjugacy classes of subgroups of `PGL_2(\mathbb{F}_5)`11761177For the case `p=5`, the order of an element determines almost the value of `u`:11781179+-------+---+---+---+---+--------+1180|`u` | 0 | 1 | 2 | 3 | 4 |1181+-------+---+---+---+---+--------+1182|orders | 2 | 3 | 4 | 6 | 1 or 5 |1183+-------+---+---+---+---+--------+11841185Here we give here the full table of all conjugacy classes of subgroups with the values1186that ``image_classes`` should give (as ``bound`` tends to `\infty`). Comparing with the output1187of the above examples, it is now easy to guess what the image is.11881189+---------+-----+------------------------------------------+1190|subgroup |order| frequencies of values of `u` |1191+=========+=====+==========================================+1192| trivial | 1 | [0.0000, 0.0000, 0.0000, 0.0000, 1.000] |1193+---------+-----+------------------------------------------+1194| cyclic | 2 | [0.5000, 0.0000, 0.0000, 0.0000, 0.5000] |1195+---------+-----+------------------------------------------+1196| cyclic | 2 | [0.5000, 0.0000, 0.0000, 0.0000, 0.5000]|1197+---------+-----+------------------------------------------+1198| cyclic | 3 | [0.0000, 0.6667, 0.0000, 0.0000, 0.3333]|1199+---------+-----+------------------------------------------+1200| Klein | 4 | [0.7500, 0.0000, 0.0000, 0.0000, 0.2500]|1201+---------+-----+------------------------------------------+1202| cyclic | 4 | [0.2500, 0.0000, 0.5000, 0.0000, 0.2500]|1203+---------+-----+------------------------------------------+1204| Klein | 4 | [0.7500, 0.0000, 0.0000, 0.0000, 0.2500]|1205+---------+-----+------------------------------------------+1206| cyclic | 5 | [0.0000, 0.0000, 0.0000, 0.0000, 1.000]|1207+---------+-----+------------------------------------------+1208| cyclic | 6 | [0.1667, 0.3333, 0.0000, 0.3333, 0.1667]|1209+---------+-----+------------------------------------------+1210| `S_3` | 6 | [0.5000, 0.3333, 0.0000, 0.0000, 0.1667]|1211+---------+-----+------------------------------------------+1212| `S_3` | 6 | [0.5000, 0.3333, 0.0000, 0.0000, 0.1667] |1213+---------+-----+------------------------------------------+1214| `D_4` | 8 | [0.6250, 0.0000, 0.2500, 0.0000, 0.1250]|1215+---------+-----+------------------------------------------+1216| `D_5` | 10 | [0.5000, 0.0000, 0.0000, 0.0000, 0.5000]|1217+---------+-----+------------------------------------------+1218| `A_4` | 12 | [0.2500, 0.6667, 0.0000, 0.0000, 0.08333]|1219+---------+-----+------------------------------------------+1220| `D_6` | 12 | [0.5833, 0.1667, 0.0000, 0.1667, 0.08333]|1221+---------+-----+------------------------------------------+1222| Borel | 20 | [0.2500, 0.0000, 0.5000, 0.0000, 0.2500]|1223+---------+-----+------------------------------------------+1224| `S_4` | 24 | [0.3750, 0.3333, 0.2500, 0.0000, 0.04167]|1225+---------+-----+------------------------------------------+1226| `PSL_2` | 60 | [0.2500, 0.3333, 0.0000, 0.0000, 0.4167]|1227+---------+-----+------------------------------------------+1228| `PGL_2` | 120| [0.2083, 0.1667, 0.2500, 0.1667, 0.2083] |1229+---------+-----+------------------------------------------+1230"""1231res = [0 for i in range(p)]1232co = 01233ell = 11234while ell <= bound:1235ell = arith.next_prime(ell)1236if ell != p and self._E.is_good(ell):1237d = (self._E.ap(ell)**2 * ell.inverse_mod(p)) % p1238res[d] += 11239co += 11240# print res1241Rt = RealField(16)1242res = [Rt(x)/Rt(co) for x in res]1243return res12441245#####################################################################1246# classification of ell and p-adic reps1247#####################################################################12481249# ell-adic reps12501251def is_unramified(self,p,ell):1252r"""1253Returns true if the Galois representation to `GL_2(\ZZ_p)` is unramified at `\ell`, i.e.1254if the inertia group at a place above `\ell` in `\text{Gal}(\bar\QQ/\QQ)` has trivial image in1255`GL_2(\ZZ_p)`.12561257For a Galois representation attached to an elliptic curve `E`, this returns True if `\ell\neq p`1258and `E` has good reduction at `\ell`.12591260INPUT:12611262- ``p`` a prime1263- ``ell`` another prime12641265OUTPUT:12661267- Boolean12681269EXAMPLES::12701271sage: rho = EllipticCurve('20a3').galois_representation()1272sage: rho.is_unramified(5,7)1273True1274sage: rho.is_unramified(5,5)1275False1276sage: rho.is_unramified(7,5)1277False12781279This says that the 5-adic representation is unramified at 7, but the 7-adic representation is ramified at 5.1280"""1281if not arith.is_prime(p):1282raise ValueError('p (=%s) must be prime' % p)1283if not arith.is_prime(ell):1284raise ValueError('ell (=%s) must be prime' % ell)1285return (ell != p) and self._E.has_good_reduction(ell)12861287def is_unipotent(self,p,ell):1288r"""1289Returns true if the Galois representation to `GL_2(\ZZ_p)` is unipotent at `\ell\neq p`, i.e.1290if the inertia group at a place above `\ell` in `\text{Gal}(\bar\QQ/\QQ)` maps into a Borel subgroup.12911292For a Galois representation attached to an elliptic curve `E`, this returns True if1293`E` has semi-stable reduction at `\ell`.12941295INPUT:12961297- ``p`` a prime1298- ``ell`` a different prime12991300OUTPUT:13011302- Boolean13031304EXAMPLES::13051306sage: rho = EllipticCurve('120a1').galois_representation()1307sage: rho.is_unipotent(2,5)1308True1309sage: rho.is_unipotent(5,2)1310False1311sage: rho.is_unipotent(5,7)1312True1313sage: rho.is_unipotent(5,3)1314True1315sage: rho.is_unipotent(5,5)1316Traceback (most recent call last):1317...1318ValueError: unipotent is not defined for l = p, use semistable instead.1319"""1320if not arith.is_prime(p):1321raise ValueError('p (=%s) must be prime' % p)1322if not arith.is_prime(ell):1323raise ValueError('ell (=%s) must be prime' % ell)1324if ell == p:1325raise ValueError("unipotent is not defined for l = p, use semistable instead.")1326return not self._E.has_additive_reduction(ell)13271328def is_quasi_unipotent(self,p,ell):1329r"""1330Returns true if the Galois representation to `GL_2(\ZZ_p)` is quasi-unipotent at `\ell\neq p`, i.e. if there is a fintie extension `K/\QQ` such that the inertia group at a place above `\ell` in `\text{Gal}(\bar\QQ/K)` maps into a Borel subgroup.13311332For a Galois representation attached to an elliptic curve `E`, this returns always True.13331334INPUT:13351336- ``p`` a prime1337- ``ell`` a different prime13381339OUTPUT:13401341- Boolean13421343EXAMPLES::13441345sage: rho = EllipticCurve('11a3').galois_representation()1346sage: rho.is_quasi_unipotent(11,13)1347True1348"""1349if not arith.is_prime(p):1350raise ValueError('p (=%s) must be prime' % p)1351if not arith.is_prime(ell):1352raise ValueError('ell (=%s) must be prime' % ell)1353if ell == p:1354raise ValueError("quasi unipotent is not defined for l = p, use semistable instead.")1355return True13561357# p-adic reps13581359def is_ordinary(self,p):1360r"""1361Returns true if the `p`-adic Galois representation to `GL_2(\ZZ_p)` is ordinary, i.e.1362if the image of the decomposition group in `\text{Gal}(\bar\QQ/\QQ)` above he prime `p` maps into1363a Borel subgroup.13641365For an elliptic curve `E`, this is to ask whether `E` is ordinary at `p`, i.e. good ordinary or multiplicative.13661367INPUT:13681369- ``p`` a prime13701371OUTPUT:13721373- a Boolean13741375EXAMPLES::13761377sage: rho = EllipticCurve('11a3').galois_representation()1378sage: rho.is_ordinary(11)1379True1380sage: rho.is_ordinary(5)1381True1382sage: rho.is_ordinary(19)1383False1384"""1385if not arith.is_prime(p):1386raise ValueError('p (=%s) must be prime' % p)1387if self._E.has_additive_reduction(p):1388raise NotImplementedError('is_ordinary is only implemented for semi-stable representations')1389return self._E.has_multiplicative_reduction(p) or (self._E.has_good_reduction(p) and self._E.ap(p) % p != 0)13901391def is_crystalline(self,p):1392r"""1393Returns true is the `p`-adic Galois representation to `GL_2(\ZZ_p)` is crystalline.13941395For an elliptic curve `E`, this is to ask whether `E` has good reduction at `p`.13961397INPUT:13981399- ``p`` a prime14001401OUTPUT:14021403- a Boolean14041405EXAMPLES::14061407sage: rho = EllipticCurve('64a1').galois_representation()1408sage: rho.is_crystalline(5)1409True1410sage: rho.is_crystalline(2)1411False1412"""1413if not arith.is_prime(p):1414raise ValueError('p (=%s) must be prime' % p)1415return self._E.has_good_reduction(p)14161417def is_potentially_crystalline(self,p):1418r"""1419Returns true is the `p`-adic Galois representation to `GL_2(\ZZ_p)` is potentially crystalline, i.e.1420if there is a finite extension `K/\QQ_p` such that the `p`-adic representation becomes crystalline.14211422For an elliptic curve `E`, this is to ask whether `E` has potentially good reduction at `p`.14231424INPUT:14251426- ``p`` a prime14271428OUTPUT:14291430- a Boolean14311432EXAMPLES::14331434sage: rho = EllipticCurve('37b1').galois_representation()1435sage: rho.is_potentially_crystalline(37)1436False1437sage: rho.is_potentially_crystalline(7)1438True1439"""1440if not arith.is_prime(p):1441raise ValueError('p (=%s) must be prime' % p)1442return self._E.j_invariant().valuation(p) >= 0144314441445def is_semistable(self,p):1446r"""1447Returns true if the `p`-adic Galois representation to `GL_2(\ZZ_p)` is semistable.14481449For an elliptic curve `E`, this is to ask whether `E` has semistable reduction at `p`.14501451INPUT:14521453- ``p`` a prime14541455OUTPUT:14561457- a Boolean14581459EXAMPLES::14601461sage: rho = EllipticCurve('20a3').galois_representation()1462sage: rho.is_semistable(2)1463False1464sage: rho.is_semistable(3)1465True1466sage: rho.is_semistable(5)1467True1468"""1469if not arith.is_prime(p):1470raise ValueError('p (=%s) must be prime' % p)1471return not self._E.has_additive_reduction(p)14721473def is_potentially_semistable(self,p):1474r"""1475Returns true if the `p`-adic Galois representation to `GL_2(\ZZ_p)` is potentially semistable.14761477For an elliptic curve `E`, this returns True always14781479INPUT:14801481- ``p`` a prime14821483OUTPUT:14841485- a Boolean14861487EXAMPLES::14881489sage: rho = EllipticCurve('27a2').galois_representation()1490sage: rho.is_potentially_semistable(3)1491True1492"""1493if not arith.is_prime(p):1494raise ValueError('p (=%s) must be prime' % p)1495return True149614971498