Local data for elliptic curves over number fields (including \QQ) at primes.

AUTHORS:

  • John Cremona: First version 2008-09-21 (refactoring code from ell_number_field.py and ell_rational_field.py)
class sage.schemes.elliptic_curves.ell_local_data.EllipticCurveLocalData(E, P, proof=None, algorithm='pari')

The class for the local reduction data of an elliptic curve.

Currently supported are elliptic curves defined over \QQ, and elliptic curves defined over a number field, at an arbitrary prime or prime ideal.

__init__(E, P, proof=None, algorithm='pari')

Initializes the reduction data for the elliptic curve E at the prime P.

INPUT:

  • E – an elliptic curve defined over a number field, or \QQ.
  • P – a prime ideal of the field, or a prime integer if the field is \QQ.
  • proof (bool)– if True, only use provably correct methods (default controlled by global proof module). Note that the proof module is number_field, not elliptic_curves, since the functions that actually need the flag are in number fields.
  • algorithm (string, default: “pari”) – Ignored unless the base field is \QQ. If “pari”, use the PARI C-library ellglobalred implementation of Tate’s algorithm over \QQ. If “generic”, use the general number field implementation.

Note

This function is not normally called directly by users, who may access the data via methods of the EllipticCurve classes.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve('14a1')
sage: EllipticCurveLocalData(E,2)
Local data at Principal ideal (2) of Integer Ring:
Reduction type: bad non-split multiplicative
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
Minimal discriminant valuation: 6
Conductor exponent: 1
Kodaira Symbol: I6
Tamagawa Number: 2
sage: EllipticCurveLocalData(E,2,algorithm="generic")
Local data at Principal ideal (2) of Integer Ring:
Reduction type: bad non-split multiplicative
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
Minimal discriminant valuation: 6
Conductor exponent: 1
Kodaira Symbol: I6
Tamagawa Number: 2
sage: EllipticCurveLocalData(E,2,algorithm="pari")
Local data at Principal ideal (2) of Integer Ring:
Reduction type: bad non-split multiplicative
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
Minimal discriminant valuation: 6
Conductor exponent: 1
Kodaira Symbol: I6
Tamagawa Number: 2
sage: EllipticCurveLocalData(E,2,algorithm="unknown")
...
ValueError: algorithm must be one of 'pari', 'generic'
sage: EllipticCurveLocalData(E,3)
Local data at Principal ideal (3) of Integer Ring:
Reduction type: good
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
Minimal discriminant valuation: 0
Conductor exponent: 0
Kodaira Symbol: I0
Tamagawa Number: 1
sage: EllipticCurveLocalData(E,7)            
Local data at Principal ideal (7) of Integer Ring:
Reduction type: bad split multiplicative
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
Minimal discriminant valuation: 3
Conductor exponent: 1
Kodaira Symbol: I3
Tamagawa Number: 3
__repr__()

Returns the string representation of this reduction data.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve('14a1')
sage: EllipticCurveLocalData(E,2).__repr__()
'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'
sage: EllipticCurveLocalData(E,3).__repr__()
'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'
sage: EllipticCurveLocalData(E,7).__repr__()
'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'
__weakref__
list of weak references to the object (if defined)
_tate(proof=None)

Tate’s algorithm for an elliptic curve over a number field.

Computes both local reduction data at a prime ideal and a local minimal model.

The model is not required to be integral on input. If P is principal, uses a generator as uniformizer, so it will not affect integrality or minimality at other primes. If P is not principal, the minimal model returned will preserve integrality at other primes, but not minimality.

Note

Called only by EllipticCurveLocalData.__init__().

OUTPUT:

(tuple) (Emin, p, val_disc, fp, KS, cp) where:

  • Emin (EllipticCurve) is a model (integral and) minimal at P
  • p (int) is the residue characteristic
  • val_disc (int) is the valuation of the local minimal discriminant
  • fp (int) is the valuation of the conductor
  • KS (string) is the Kodaira symbol
  • cp (int) is the Tamagawa number
bad_reduction_type()

Return the type of bad reduction of this reduction data.

OUTPUT:

(int or None):

  • +1 for split multiplicative reduction
  • -1 for non-split multiplicative reduction
  • 0 for additive reduction
  • None for good reduction

EXAMPLES:

sage: E=EllipticCurve('14a1')
sage: [(p,E.local_data(p).bad_reduction_type()) for p in prime_range(15)]
[(2, -1), (3, None), (5, None), (7, 1), (11, None), (13, None)]

sage: K.<a>=NumberField(x^3-2)
sage: P17a, P17b = [P for P,e in K.factor(17)]
sage: E = EllipticCurve([0,0,0,0,2*a+1])
sage: [(p,E.local_data(p).bad_reduction_type()) for p in [P17a,P17b]]           
[(Fractional ideal (4*a^2 - 2*a + 1), None), (Fractional ideal (2*a + 1), 0)]
conductor_valuation()

Return the valuation of the conductor from this local reduction data.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve([0,0,0,0,64]); E
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
sage: data = EllipticCurveLocalData(E,2)
sage: data.conductor_valuation()
2
has_additive_reduction()

Return True if there is additive reduction.

EXAMPLES:

sage: E = EllipticCurve('27a1')
sage: [(p,E.local_data(p).has_additive_reduction()) for p in prime_range(15)]
[(2, False), (3, True), (5, False), (7, False), (11, False), (13, False)]
sage: K.<a> = NumberField(x^3-2)
sage: P17a, P17b = [P for P,e in K.factor(17)]
sage: E = EllipticCurve([0,0,0,0,2*a+1])
sage: [(p,E.local_data(p).has_additive_reduction()) for p in [P17a,P17b]]           
[(Fractional ideal (4*a^2 - 2*a + 1), False),
(Fractional ideal (2*a + 1), True)]
has_bad_reduction()

Return True if there is bad reduction.

EXAMPLES:

sage: E = EllipticCurve('14a1')
sage: [(p,E.local_data(p).has_bad_reduction()) for p in prime_range(15)]
[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]
sage: K.<a> = NumberField(x^3-2)
sage: P17a, P17b = [P for P,e in K.factor(17)]
sage: E = EllipticCurve([0,0,0,0,2*a+1])
sage: [(p,E.local_data(p).has_bad_reduction()) for p in [P17a,P17b]]           
[(Fractional ideal (4*a^2 - 2*a + 1), False),
(Fractional ideal (2*a + 1), True)]
has_good_reduction()

Return True if there is good reduction.

EXAMPLES:

sage: E = EllipticCurve('14a1')
sage: [(p,E.local_data(p).has_good_reduction()) for p in prime_range(15)]
[(2, False), (3, True), (5, True), (7, False), (11, True), (13, True)]

sage: K.<a> = NumberField(x^3-2)
sage: P17a, P17b = [P for P,e in K.factor(17)]
sage: E = EllipticCurve([0,0,0,0,2*a+1])
sage: [(p,E.local_data(p).has_good_reduction()) for p in [P17a,P17b]]           
[(Fractional ideal (4*a^2 - 2*a + 1), True),
(Fractional ideal (2*a + 1), False)]
has_multiplicative_reduction()

Return True if there is multiplicative reduction.

Note

See also has_split_multiplicative_reduction() and has_nonsplit_multiplicative_reduction().

EXAMPLES:

sage: E = EllipticCurve('14a1')
sage: [(p,E.local_data(p).has_multiplicative_reduction()) for p in prime_range(15)]
[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]
sage: K.<a> = NumberField(x^3-2)
sage: P17a, P17b = [P for P,e in K.factor(17)]
sage: E = EllipticCurve([0,0,0,0,2*a+1])
sage: [(p,E.local_data(p).has_multiplicative_reduction()) for p in [P17a,P17b]]           
[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]
has_nonsplit_multiplicative_reduction()

Return True if there is non-split multiplicative reduction.

EXAMPLES:

sage: E = EllipticCurve('14a1')
sage: [(p,E.local_data(p).has_nonsplit_multiplicative_reduction()) for p in prime_range(15)]
[(2, True), (3, False), (5, False), (7, False), (11, False), (13, False)]
sage: K.<a> = NumberField(x^3-2)
sage: P17a, P17b = [P for P,e in K.factor(17)]
sage: E = EllipticCurve([0,0,0,0,2*a+1])
sage: [(p,E.local_data(p).has_nonsplit_multiplicative_reduction()) for p in [P17a,P17b]]           
[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]
has_split_multiplicative_reduction()

Return True if there is split multiplicative reduction.

EXAMPLES:

sage: E = EllipticCurve('14a1')
sage: [(p,E.local_data(p).has_split_multiplicative_reduction()) for p in prime_range(15)]
[(2, False), (3, False), (5, False), (7, True), (11, False), (13, False)]
sage: K.<a> = NumberField(x^3-2)
sage: P17a, P17b = [P for P,e in K.factor(17)]
sage: E = EllipticCurve([0,0,0,0,2*a+1])
sage: [(p,E.local_data(p).has_split_multiplicative_reduction()) for p in [P17a,P17b]]           
[(Fractional ideal (4*a^2 - 2*a + 1), False),
(Fractional ideal (2*a + 1), False)]
kodaira_symbol()

Return the Kodaira symbol from this local reduction data.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve([0,0,0,0,64]); E
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
sage: data = EllipticCurveLocalData(E,2)
sage: data.kodaira_symbol()
IV
minimal_model()

Return the (local) minimal model from this local reduction data.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve([0,0,0,0,64]); E
Elliptic Curve defined by y^2  = x^3 + 64 over Rational Field
sage: data = EllipticCurveLocalData(E,2)
sage: data.minimal_model()
Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field
sage: data.minimal_model() == E.local_minimal_model(2)
True
prime()

Return the prime ideal associated with this local reduction data.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve([0,0,0,0,64]); E
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
sage: data = EllipticCurveLocalData(E,2)
sage: data.prime()
Principal ideal (2) of Integer Ring
tamagawa_exponent()

Return the Tamagawa index from this local reduction data.

This is the exponent of E(K_v)/E^0(K_v); in most cases it is the same as the Tamagawa index.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve('816a1')
sage: data = EllipticCurveLocalData(E,2)
sage: data.kodaira_symbol()
I2*
sage: data.tamagawa_number()
4
sage: data.tamagawa_exponent()
2

sage: E = EllipticCurve('200c4')
sage: data = EllipticCurveLocalData(E,5)
sage: data.kodaira_symbol()
I4*
sage: data.tamagawa_number()
4
sage: data.tamagawa_exponent()
2
tamagawa_number()

Return the Tamagawa number from this local reduction data.

This is the index [E(K_v):E^0(K_v)].

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
sage: E = EllipticCurve([0,0,0,0,64]); E
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
sage: data = EllipticCurveLocalData(E,2)
sage: data.tamagawa_number()
3
sage.schemes.elliptic_curves.ell_local_data.check_prime(K, P)

Function to check that P determines a prime of K, and return that ideal.

INPUT:

  • K – a number field (including \QQ).
  • P – an element of K or a (fractional) ideal of K.

OUTPUT:

  • If K is \QQ: the prime integer equal to or which generates P.
  • If K is not \QQ`: the prime ideal equal to or generated by P.

Note

If P is not a prime and does not generate a prime, a TypeError is raised.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.ell_local_data import check_prime
sage: check_prime(QQ,3)
3
sage: check_prime(QQ,ZZ.ideal(31))
31
sage: K.<a>=NumberField(x^2-5)
sage: check_prime(K,a)
Fractional ideal (a)
sage: check_prime(K,a+1)
Fractional ideal (a + 1)
sage: [check_prime(K,P) for P in K.primes_above(31)]
[Fractional ideal (-5/2*a - 1/2), Fractional ideal (-5/2*a + 1/2)]

Previous topic

Torsion subgroups of elliptic curves over number fields (including \QQ).

Next topic

Kodaira symbols.

This Page