ainvs, [extra=None]) |
an,
analytic_rank,
anlist,
ap,
complex_area,
conductor,
CPS_height_bound,
cremona_label,
database_curve,
eval_modular_form,
formal_group,
formal_inverse,
formal_log,
formal_mult,
formal_n_isogeny,
formal_w,
formal_x,
formal_y,
gens,
has_cm,
heegner_discriminants,
heegner_discriminants_list,
heegner_index,
heegner_index_bound,
heegner_point_height,
is_integral,
is_minimal,
is_semistable,
is_surjective,
isogeny_class,
kodaira_type,
L1_vanishes,
L_ratio,
Lambda,
Lseries,
Lseries_at1,
Lseries_deriv_at1,
Lseries_extended,
minimal_model,
modular_degree,
modular_parametrization,
mwrank_curve,
newform,
ngens,
non_surjective,
Np,
omega,
padic_sigma,
pari_curve,
pari_mincurve,
period_lattice,
point_search,
quadratic_twist,
rank,
real_components,
regulator,
root_number,
satisfies_heegner_hypothesis,
saturation,
selmer_rank_bound,
sha_an,
shabound,
shabound_kato,
shabound_kolyvagin,
sigma,
silverman_height_bound,
tamagawa_number,
tamagawa_product,
torsion_order,
torsion_subgroup,
two_descent,
two_selmer_shabound,
two_torsion_rank,
watkins_data
These methods are defined as follows:
n) |
The n-th Fourier coefficient of the modular form corresponding to this elliptic curve, where n is a positive integer.
) |
n) |
The Fourier coefficients up to and including a_n of the modular form attached to this elliptic curve. The ith element of the return list is a[i].
sage: E = EllipticCurve([0, -1, 1, -10, -20]) sage: E.anlist(3) [0, 1, -2, -1] sage: E = EllipticCurve([0,1]) sage: E.anlist(20) [0, 1, 0, 0, 0, 0, 0, -4, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 8, 0]
p) |
The p-th Fourier coefficient of the modular form corresponding to this elliptic curve, where p is prime.
) |
Return the area of a fundamental domain for the period lattice of the elliptic curve.
[algorithm=pari]) |
Returns the conductor of the elliptic curve.
INPUT: algorithm -- str, (default: "pari") "pari" -- use the PARI C-library ellglobalred implementation of Tate's algorithm "mwrank" -- use Cremona's mwrank implementation of Tate's algorithm; can be faster if the curve has integer coefficients (TODO: limited to small conductor until mwrank gets integer factorization) "all" -- use both implementations, verify that the results are the same (or raise an error), and output the common value.
sage: E = EllipticCurve([1, -1, 1, -29372, -1932937]) sage: E.conductor(algorithm="pari") 3006 sage: E.conductor(algorithm="mwrank") 3006 sage: E.conductor(algorithm="all") 3006
NOTE: The conductor computed using each algorithm is cached separately. Thus calling E.conductor("pari"), then E.conductor("mwrank") and getting the same result checks that both systems compute the same answer.
) |
Return the Cremona-Prickett-Siksek height bound. This is a floating point number B such that if P is a point on the curve, then the naive logarithmetic height of P is off from the canonical height by at most B.
[space=False]) |
Return the Cremona label associated to (the minimal model) of this curve, if it is known. If not, raise a RuntimeError exception.
) |
Return the curve in the elliptic curve database isomorphic to this curve, if possible. Otherwise raise a RuntimeError exception.
The following examples assumes that the elliptic curve database is installed.
sage: E = EllipticCurve([0,1,2,3,4]) sage: E.database_curve() Elliptic Curve defined by y^2 = x^3 + x^2 + 3*x + 5 over Rational Field
NOTES: The model of the curve in the database can be different than the Weierstrass model for this curve, since database models are always minimal.
points, prec) |
[prec=10]) |
[prec=20]) |
[prec=20]) |
n, [prec=10]) |
n, [prec=20]) |
[prec=20]) |
Return the formal power series
to precision
of Proposition IV.1.1 of [Silverman
AEC1]. This is the formal expansion of
about the
formal parameter
at
.
We compute w using the recursive procedure (4.1) on page 18 of Bluher's `A leisurely introduction to formal groups and elliptic curves', which I downloaded from http://www.math.uiuc.edu/Algebraic-Number-Theory/0076/
[prec=20]) |
Return the formal power series x(t) in terms of the local parameter t = -x/y at infinity.
[prec=20]) |
Return the formal power series y(t) in terms of the local parameter t = -x/y at infinity.
[verbose=10], [rank1_search=False]) |
Compute and return generators for the Mordell-Weil group E(Q) modulo torsion.
INPUT: verbose -- (default: None), if specified changes the verbosity of mwrank computations. rank1_search -- (default: 16), if the curve has analytic rank 1, try to find a generator by a direct search up to this logarithmic height. If this fails the usual mwrank procedure is called. OUTPUT: generators -- List of generators for the Mordell-Weil group.
IMPLEMENTATION: Uses Cremona's mwrank C library.
) |
bound) |
n) |
List of the first n Heegner discriminants for self.
D, [min_p=True], [prec=5], [verbose=3]) |
Assume self has rank 0.
Return (an interval that contains) the square of the odd part of the index of the Heegner point in the group of K-rational points *modulo torsion* on the twist of the elliptic curve by D.
If 0 is in the interval of the height of the Heegner point computed to the given prec, then this function returns 0.
INPUT: D (int) -- Heegner discriminant min_p (int) -- (default: 3) only rule out primes >= min_p dividing the index. verbose (bool) -- (default: True) prec (int) -- (default: 5), use prec*sqrt(N) + 20 terms of L-series in computations, where N is the conductor. OUTPUT: Interval.
[D=21], [prec=True], [verbose=5], [max_height=0]) |
Assume self has rank 0.
Return a list v of primes such that if an odd prime p divides the index of the the Heegner point in the group of rational points *modulo torsion*, then p is in v.
If 0 is in the interval of the height of the Heegner point computed to the given prec, then this function returns v = 0. This does not mean that the Heegner point is torsion, just that it is very likely torsion.
If we obtain no information from a search up to max_height, e.g., if the Siksek et al. bound is bigger than max_height, then we return v = -1.
INPUT: D (int) -- (deault: 0) Heegner discriminant; if 0, use the first discriminant < -4 that satisfies the Heegner hypothesis verbose (bool) -- (default: True) prec (int) -- (default: 5), use prec*sqrt(N) + 20 terms of L-series in computations, where N is the conductor. max_height (float) -- should be <= 21; bound on logarithmic naive height used in point searches. Make smaller to make this function faster, at the expense of possibly obtaining a worse answer. A good range is between 13 and 21. OUTPUT: v -- list or int (bad primes or 0 or -1) D -- the discriminant that was used (this is useful if D was automatically selected).
D, [prec=2]) |
Use the Gross-Zagier formula to compute the Neron-Tate canonical height over K of the Heegner point corresponding to D, as an Interval (since it's computed to some precision using L-functions).
INPUT: D (int) -- fundamental discriminant prec (int) -- (default: 2), use prec*sqrt(N) + 20 terms of L-series in computations, where N is the conductor. OUTPUT: Interval that contains the height of the Heegner point.
) |
) |
) |
p, [A=1000]) |
Return True if the mod-p representation attached to E is surjective, False if it is not, or None if we were unable to determine whether it is or not.
INPUT: p -- int (a prime number) A -- int (a bound on the number of a_p to use) OUTPUT: a 2-tuple: -- surjective or (probably) not -- information about what it is if not surjective
REMARKS:
1. If p >= 5 then the mod-p representation is surjective if and only if the p-adic representation is surjective. When p = 2, 3 there are counterexamples. See a very recent paper of Elkies for more details when p=3.
2. When p <= 3 this function always gives the correct result irregardless of A, since it explicitly determines the p-division polynomial.
) |
Return all curves over the rational field in the isogeny class of this elliptic curve.
NOTE: Currently implemented only using the database.
p) |
Local Kodaira type of the elliptic curve at p.
1 means good reduction (type
), 2, 3 and 4 mean types II,
III and IV, respectively,
with
means
type
; finally the opposite values -1, -2,
etc. refer to the starred types
,
, etc.
) |
Returns whether or not L(E,1) = 0. The result is provably correct if the Manin constant of the associated optimal quotient is <= 2. This hypothesis on the Manin constant is true for all curves of conductor <= 40000 (by Cremona) and all semistable curves (i.e., squarefree conductor).
sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11) sage: E.L1_vanishes() False sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11) sage: E.L1_vanishes() False sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1) sage: E.L1_vanishes() True sage: E = EllipticCurve([0, 1, 1, -2, 0]) # 389A (rank 2) sage: E.L1_vanishes() True sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A (CM curve)) sage: E.L1_vanishes() True sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C (13-isogeny) sage: E.L1_vanishes() False
WARNING: It's conceivable that machine floats are not large enough precision for the computation; if this could be the case a RuntimeError is raised. The curve's real period would have to be very small for this to occur.
ALGORITHM: Compute the root number. If it is -1 then L(E,s) vanishes to odd order at 1, hence vanishes. If it is +1, use a result about modular symbols and Mazur's "Rational Isogenies" paper to determine a provably correct bound (assuming Manin constant is <= 2) so that we can determine whether L(E,1) = 0.
AUTHOR: William Stein, 2005-04-20.
) |
Returns the ratio L(E,1)/Omega as an exact rational number. The result is *provably* correct if the Manin constant of the associated optimal quotient is <= 2. This hypothesis on the Manin constant is true for all semistable curves (i.e., squarefree conductor), by a theorem of Mazur from his "Rational Isogenies of Prime Degree" paper.
sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11) sage: E.L_ratio() 1/5 sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11) sage: E.L_ratio() 1/25 sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1) sage: E.L_ratio() 0 sage: E = EllipticCurve([0, 1, 1, -2, 0]) # 389A (rank 2) sage: E.L_ratio() 0 sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A (CM curve)) sage: E.L_ratio() 0 sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C (13-isogeny) sage: E.L_ratio() 1
WARNING: It's conceivable that machine floats are not large enough precision for the computation; if this could be the case a RuntimeError is raised. The curve's real period would have to be very small for this to occur.
ALGORITHM: Compute the root number. If it is -1 then L(E,s) vanishes to odd order at 1, hence vanishes. If it is +1, use a result about modular symbols and Mazur's "Rational Isogenies" paper to determine a provably correct bound (assuming Manin constant is <= 2) so that we can determine whether L(E,1) = 0.
AUTHOR: William Stein, 2005-04-20.
s, prec) |
Returns the value of the Lambda-series of the elliptic curve E at s can be any complex number. Use Lambda(s) if s is real, since it will be much faster.
Warning: This function uses sage.functions.transcendental.gamma_inc, which currently calls mathematica.
s) |
Returns the value of the L-series of the elliptic curve E at s, where s must be a real number.
NOTES:
(1) If the conductor of the curve is large, say >
, then this function
will take a long time, since it uses an O(sqrt(N)) algorithm.
(2) There is an algorithm to compute L(E,s) for any complex s, but it involves computing the incomplete Gamma function, which is not implemented in SAGE.
sage: E = EllipticCurve([1,2,3,4,5]) sage: E.Lseries(1) 0.0
(The following may have low order bits that are different on a 64-bit machine.)
sage: E.Lseries(1.1) 0.2854910076781481284428923431
[k=0]) |
Compute L(E,1) using k terms of the series for L(E,1) form page 406 of Henri Cohen's book "A Course in Computational Algebraic Number Theory". If the k argument is not specified, then it defaults to sqrt(N), where N is the conductor.
The real precision of the computation is the precision of Python floats.
INPUT: k -- (optional) an integer, defaults to sqrt(N). OUTPUT: float -- L(E,1) float -- a bound on the error in the approximation; this is a proveably correct upper bound on the sum of the tail end of the series used to compute L(E,1).
ALGORITHM: 1. Compute the root number eps. If it is -1, return 0.
2. Compute the Fourier coefficients a_n, for n up to and including k.
3. Compute the sum
where N is the conductor of E.
4. Compute a bound on the tail end of the series, which is
For a proof see [Grigov-Jorza-Patrascu-Patrikis-Stein].
[k=0]) |
Compute L'(E,1) using k terms of the series for L'(E,1) form page 406 of Henri Cohen's book"A Course in Computational Algebraic Number Theory".
The real precision of the computation is the precision of Python floats.
INPUT: OUTPUT:
ALGORITHM: 1. Compute the root number eps. If it is 1, return 0.
2. Compute the Fourier coefficients a_n, for n up to and including k.
3. Compute the sum
where N is the conductor of E, and E_1 is the exponential integral function.
4. Compute a bound on the tail end of the series, which is
For a proof see [Grigov-Jorza-Patrascu-Patrikis-Stein]. This is exactly the same as the bound for the approximation to
s, prec) |
Returns the value of the L-series of the elliptic curve E at s can be any complex number using prec terms of the power series expansion.
Warning: This function uses sage.functions.transcendental.gamma_inc, which currently calls mathematica.
) |
Return the unique minimal Weierstrass equation for this
elliptic curve. This is the model with minimal discriminant
and
.
) |
Return the modular degree of this elliptic curve, computing using Mark Watkins's C implementation of his analytic symmetric square L-function algorithm.
The correctness of this function is subject to the following three hypothesis:
NOTE: Mark has a much more robust implementation available in MAGMA.
) |
Computes and returns ...
[verbose=False]) |
) |
Returns the newform attached to this elliptic curve.
) |
[A=1000]) |
Returns a list of primes p such that the mod-p representation rho_E,p *might* not be surjective (this list usually contains 2, because of shortcomings of the algorithm). If p is not in the returned list, then rho_E,p is provably surjective (see A. Cojocaru's paper). If the curve has CM then infinitely many representations are not surjective, so we simply return the sequence [(0,"cm")] and do no further computation.
INPUT: A -- an integer OUTPUT: list -- if curve has CM, returns [(0,"cm")]. Otherwise, returns a list of primes where mod-p representation very likely not surjective. At any prime not in this list, the representation is definitely surjective.
sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A sage: E.non_surjective() # CM curve [(0, 'cm')]
sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11) sage: E.non_surjective() [(5, '5-torsion')]
sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A sage: E.non_surjective() []
sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C sage: E.non_surjective() [(13, [1])]
ALGORITHM: When p<=3 use division polynomials. For 5 <= p <= B, where B is Cojocaru's bound, use the results in Section 2 of Serre's inventiones paper"Sur Les Representations Modulaires Deg Degre 2 de Galqbar Over Q."
p) |
The number of points on E modulo p, where p is a prime of good reduction.
sage: E = EllipticCurve([0, -1, 1, -10, -20]) sage: E.Np(2) 5 sage: E.Np(3) 5 sage: E.conductor() 11 sage: E.Np(11) Traceback (most recent call last): ... ArithmeticError: p (=11) must be a prime of good reduction
) |
Returns the real period. This is the correct period in the BSD conjecture, i.e., it is the least real period * 2 when the period lattice is rectangular.
p, [prec=20]) |
Returns the p-adic sigma function of the elliptic curve as a power series in t to precision prec.
) |
Return the PARI curve corresponding to this elliptic curve.
sage: E = EllipticCurve([0, 0,1,-1,0]) sage: e = E.pari_curve() sage: type(e) <type 'gen.gen'> sage: e.type() 't_VEC' sage: e.ellan(10) [1, -2, -3, 2, -2, 6, -1, 0, 6, 4]
sage: E = EllipticCurve(RationalField(), ['1/3', '2/3']) sage: e = E.pari_curve() sage: e.type() 't_VEC' sage: e[:5] [0, 0, 0, 1/3, 2/3]
) |
Return the PARI curve corresponding to a minimal model for this elliptic curve.
sage: E = EllipticCurve(RationalField(), ['1/3', '2/3']) sage: e = E.pari_mincurve() sage: e[:5] [0, 0, 0, 27, 486] sage: E.conductor() 47232 sage: e.ellglobalred() [47232, [1, 0, 0, 0], 2]
) |
Return a basis for the period lattice of the elliptic curve over Q as a 2-tuple.
The basis has the form [omega_1, omega_2], where Im(omega_1/omega_2) > 0 and omega_1 is real.
INPUT: -- an elliptic curve OUTPUT: omega_1 -- complex number omega_2 -- complex number
height_limit, [verbose=True]) |
Search for points on a curve up to an input bound on the naive logarithmic height.
INPUT: height_limit (float) -- bound on naive height (at most 21, or mwrank overflows) verbose (bool) -- (default: True) OUTPUT: points (list) -- list of points found
IMPLEMENTATION: Uses Cremona's mwrank package.
D) |
Return the global minimal model of the quadratic twist of this curve by D.
[use_database=False], [verbose=False]) |
Return the rank of this elliptic curve, assuming no conjectures.
If we fail to provably compute the rank, raises a RuntimeError exception.
INPUT: use_database (bool) -- (default: False), if True, try to look up the regulator in the Cremona database. proof (bool) -- (default: True), if True then rank is proven correct; otherwise, computation of rank is not certain. verbose -- (default: None), if specified changes the verbosity of mwrank computations. OUTPUT: rank (int) -- the rank of the elliptic curve.
IMPLEMENTATION: Uses L-functions and mwrank.
) |
Returns 1 if there is 1 real component and 2 if there are 2.
[use_database=None], [verbose=False]) |
Returns the regulator of this curve, which must be defined over Q.
INPUT: use_database -- bool (default: False), if True, try to look up the regulator in the Cremona database. verbose -- (default: None), if specified changes the verbosity of mwrank computations.
sage: E = EllipticCurve([0, 0, 1, -1, 0]) sage: E.regulator() 0.05111140823996884
) |
Returns the root number of this elliptic curve.
This is 1 if the order of vanishing of the L-function L(E,s) at 1 is even, and -1 if it is odd.
D) |
Returns True precisely when D is a fundamental discriminant that satisfies the Heegner hypothesis for this elliptic curve.
points, [verbose=False], [max_prime=0], [odd_primes_only=False]) |
Given a list of rational points on E, compute the saturation in E(Q) of the subgroup they generate.
INPUT: points (list) -- list of points on E verbose (bool) -- (default: False), if True, give verbose output max_prime (int) -- (default: 0), saturation is performed for all primes up to max_prime. If max_prime==0, perform saturation at *all* primes, i.e., compute the true saturation. odd_primes_only (bool) -- only do saturation at odd primes OUTPUT: saturation (list) -- points that form a basis for the saturation index (int) -- the index of the group generated by points in their saturation regulator (float) -- regulator of saturated points.
IMPLEMENTATION: Uses Cremona's mwrank package. With max_prime=0, we call mwrank with successively larger prime bounds until the full saturation is provably found. The results of saturation at the previous primes is stored in each case, so this should be reasonably fast.
) |
Bound on the rank of the curve, computed using the 2-selmer group. This is the rank of the curve minus the rank of the 2-torsion, minus a number determined by whatever mwrank was able to determine related to Sha[2]. Thus in many cases, this is the actual rank of the curve.
EXAMPLE: The following is the curve 960D1, which has rank 0, but Sha of order 4.
sage: E = EllipticCurve([0, -1, 0, -900, -10098]) sage: E.selmer_rank_bound() 0 It gives 0 instead of 2, because it knows Sha is nontrivial. In contrast, for the curve 571A, also with rank 0 and Sha of order 4, we get a worse bound: sage: E = EllipticCurve([0, -1, 1, -929, -10595]) sage: E.selmer_rank_bound() 2 sage: E.rank() 0
[use_database=False]) |
Returns the Birch and Swinnerton-Dyer conjectural order of Sha, unless the analytic rank is > 1, in which case this function returns 0.
This result is proved correct if the order of vanishing is 0 and the Manin constant is <= 2.
If the optional parameter use_database is True (default: False), this function returns the analytic order of Sha as listed in Cremona's tables, if this curve appears in Cremona's tables.
sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11) sage: E.sha_an() 1 sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11) sage: E.sha_an() 1
The smallest conductor curve with nontrivial Sha:
sage: E = EllipticCurve([1,1,1,-352,-2689]) # 66B3 sage: E.sha_an() 4
The four optimal quotients with nontrivial Sha and conductor <= 1000:
sage: E = EllipticCurve([0, -1, 1, -929, -10595]) # 571A sage: E.sha_an() 4 sage: E = EllipticCurve([1, 1, 0, -1154, -15345]) # 681B sage: E.sha_an() 9 sage: E = EllipticCurve([0, -1, 0, -900, -10098]) # 960D sage: E.sha_an() 4 sage: E = EllipticCurve([0, 1, 0, -20, -42]) # 960N sage: E.sha_an() 4
The smallest conductor curve of rank > 1:
sage: E = EllipticCurve([0, 1, 1, -2, 0]) # 389A (rank 2) sage: E.sha_an() 0
The following are examples that require computation of the Mordell-Weil group and regulator:
sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1) sage: E.sha_an() 1
sage: E = EllipticCurve("1610F3") sage: E.sha_an() 4
sage: e = EllipticCurve([1, 0, 0, -19491080, -33122512122]) # 15834T2 sage: e.sha_an() # takes about 22 seconds 25
) |
Compute a provably correct bound on the order of the Shafarevich-Tate group of this curve. The bound is a either False (no bound) or a list B of primes such that any divisor of Sha is in this list.
) |
Returns a list p of primes such tha theorems of Kato's and
others (e.g., as explained in a paper/thesis of Grigor Grigorov)
imply that if p divides
then
is in the list.
If L(E,1) = 0, then Kato's theorem gives no information, so this function returns False.
THEOREM (Kato): Suppose p >= 5 is a prime so the p-adic
representation rho_E,p is surjective. Then
divides
.
sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11) sage: E.shabound_kato() [2, 3, 5] sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11) sage: E.shabound_kato() [2, 3, 5] sage: E = EllipticCurve([1,1,1,-352,-2689]) # 66B3 sage: E.shabound_kato() [2, 3]
For the following curve one really has 25 |
(by Grigorov-Stein paper):
sage: E = EllipticCurve([1, -1, 0, -332311, -73733731]) # 1058D1 sage: E.shabound_kato() [2, 3, 5] sage: E.non_surjective() []
For this one, Sha is divisible by 7.
sage: E = EllipticCurve([0, 0, 0, -4062871, -3152083138]) # 3364C1 sage: E.shabound_kato() [2, 3, 7]
No information about curves of rank > 0:
sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1) sage: E.shabound_kato() False
[D=False], [regulator=None], [ignore_nonsurj_hypothesis=0]) |
Given a fundamental discriminant D (=-3,-4) that satisfies the
Heegner hypothesis, return a list of primes so that
Kolyvagin's theorem (as in Gross's paper) implies that any
prime divisor of
is in this list.
INPUT: D -- (optional) a fundamental discriminant < -4 that satisfies the Heegner hypothesis for E; if not given, use the first such D regulator -- (optional) regulator of E(K); if not given, will be computed (which could take a long time) ignore_nonsurj_hypothesis (optional: default False) -- If True, then gives the bound coming from Heegner point index, but without any hypothesis on surjectivity of the mod-p representation. OUTPUT: bound and index
More precisely:
0 - if E/K has complex multiplication or analytic rank >= 2 or B - list of primes such that if p divides Sha(E/K), then p is in B.
and
I - the odd part of the index of the Heegner point in the full group of K-rational points on E. (If E has CM, returns 0.)
REMARKS: (1) We do not have to assume that the Manin constant is 1 (or a power of 2). If the Manin constant were divisible by a prime, that prime would get included in the list of bad primes.
(2) We assume the Gross-Zagier theorem is True under the hypothesis that gcd(N,D) = 1, instead of the stronger hypothesis gcd(2*N,D)=1 that is in the original Gross-Zagier paper. That Gross-Zagier is true when gcd(N,D)=1 is"well-known" to the experts, but doesn't seem to written up well in the literature.
(3) Correctness of the computation is guaranteed using
interval arithmetic, under the assumption that the
regulator, square root, and period lattice are
computed to precision at least
, i.e., they are
correct up to addition or a real number with absolute
value less than
.
z, [flag=0]) |
Returns the value of the Weierstrass sigma function of the lattice associated to this elliptic curve E.
INPUT: z -- a complex number flag -- 0 - default ??? 1 - computes an arbitrary determination of log(sigma(z)) 2, 3 - same using the product expansion instead of theta series. ??? OUTPUT: a complex number
NOTE: The reason for the ???'s above, is that the PARI documentation for ellsigma is very vague.
) |
Return the Silverman height bound. This is a floating point number B such that if P is a point on the curve, then the naive logarithmetic height of P is off from the canonical height by at most B.
Note that the CPS_height_bound is typically much better than the Silverman bound.
p) |
The Tamagawa number of the elliptic curve at
.
) |
Returns the product of the Tamagawa numbers.
) |
Return the order of the torsion subgroup.
[flag=0]) |
Returns the torsion subgroup of this elliptic curve.
The flag is passed onto PARI and has the same algorithm meaning as there. So: If flag = 0, use Doud's algorithm; if flag = 1, use Lutz-Nagell.
[verbose=1], [sel=-1], [first_limit=5], [second_limit=20], [n_aux=False], [second_descent=True]) |
Compute 2-descent data for this curve. The actual data is part of self.mwrank_curve().
INPUT: verbose -- (default: True) print what mwrank is doing sel -- (default: False) ?? first_limit -- (default: 20) ?? second_limit-- (default: 5) ?? n_aux -- (default: -1) ?? second_descent -- (default: 1) ??
) |
Returns a bound on the dimension of Sha(E)[2], computed using a 2-descent.
) |
) |
Return a dict of the data computed by Mark Watkins's ec program applied to this elliptic curve.
See About this document... for information on suggesting changes.