Canonical heights for elliptic curves over number fields¶
Also, rigorous lower bounds for the canonical height of non-torsion points, implementing the algorithms in [CS2006] (over \(\QQ\)) and [Tho2010], which also refer to [CPS2006].
AUTHORS:
- Robert Bradshaw (2010): initial version
- John Cremona (2014): added many docstrings and doctests
-
class
sage.schemes.elliptic_curves.height.
EllipticCurveCanonicalHeight
(E)¶ Class for computing canonical heights of points on elliptic curves defined over number fields, including rigorous lower bounds for the canonical height of non-torsion points.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import EllipticCurveCanonicalHeight sage: E = EllipticCurve([0,0,0,0,1]) sage: EllipticCurveCanonicalHeight(E) EllipticCurveCanonicalHeight object associated to Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field
Normally this object would be created like this:
sage: E.height_function() EllipticCurveCanonicalHeight object associated to Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field
-
B
(n, mu)¶ Return the value \(B_n(\mu)\).
INPUT:
n
(int) - a positive integermu
(real) - a positive real number
OUTPUT:
The real value \(B_n(\mu)\) as defined in [Tho2010], section 5.
EXAMPLES:
Example 10.2 from [Tho2010]:
sage: K.<i>=QuadraticField(-1) sage: E = EllipticCurve([0,1-i,i,-i,0]) sage: H = E.height_function()
In [Tho2010] the value is given as 0.772:
sage: RealField(12)( H.B(5, 0.01) ) 0.777
-
DE
(n)¶ Return the value \(D_E(n)\).
INPUT:
n
(int) - a positive integer
OUTPUT:
The value \(D_E(n)\) as defined in [Tho2010], section 4.
EXAMPLES:
sage: K.<i>=QuadraticField(-1) sage: E = EllipticCurve([0,0,0,1+5*i,3+i]) sage: H = E.height_function() sage: [H.DE(n) for n in srange(1,6)] [0, 2*log(5) + 2*log(2), 0, 2*log(13) + 2*log(5) + 4*log(2), 0]
-
ME
()¶ Return the norm of the ideal \(M_E\).
OUTPUT:
The norm of the ideal \(M_E\) as defined in [Tho2010], section 3.1. This is \(1\) if \(E\) is a global minimal model, and in general measures the non-minimality of \(E\).
EXAMPLES:
sage: K.<i>=QuadraticField(-1) sage: E = EllipticCurve([0,0,0,1+5*i,3+i]) sage: H = E.height_function() sage: H.ME() 1 sage: E = EllipticCurve([0,0,0,0,1]) sage: E.height_function().ME() 1 sage: E = EllipticCurve([0,0,0,0,64]) sage: E.height_function().ME() 4096 sage: E.discriminant()/E.minimal_model().discriminant() 4096
-
S
(xi1, xi2, v)¶ Return the union of intervals \(S^{(v)}(\xi_1,\xi_2)\).
INPUT:
xi1, xi2
(real) - real numbers with \(\xi_1\le\xi_2\).v
(embedding) - a real embedding of the field.
OUTPUT:
The union of intervals \(S^{(v)}(\xi_1,\xi_2)\) defined in [Tho2010] section 6.1.
EXAMPLES:
An example over \(\QQ\):
sage: E = EllipticCurve('389a') sage: v = QQ.places()[0] sage: H = E.height_function() sage: H.S(2,3,v) ([0.224512677391895, 0.274544821597130] U [0.725455178402870, 0.775487322608105])
An example over a number field:
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,0,a]) sage: v = K.real_places()[0] sage: H = E.height_function() sage: H.S(9,10,v) ([0.0781194447253472, 0.0823423732016403] U [0.917657626798360, 0.921880555274653])
-
Sn
(xi1, xi2, n, v)¶ Return the union of intervals \(S_n^{(v)}(\xi_1,\xi_2)\).
INPUT:
xi1, xi2
(real) - real numbers with \(\xi_1\le\xi_2\).n
(integer) - a positive integer.v
(embedding) - a real embedding of the field.
OUTPUT:
The union of intervals \(S_n^{(v)}(\xi_1,\xi_2)\) defined in [Tho2010] (Lemma 6.1).
EXAMPLES:
An example over \(\QQ\):
sage: E = EllipticCurve('389a') sage: v = QQ.places()[0] sage: H = E.height_function() sage: H.S(2,3,v) , H.Sn(2,3,1,v) (([0.224512677391895, 0.274544821597130] U [0.725455178402870, 0.775487322608105]), ([0.224512677391895, 0.274544821597130] U [0.725455178402870, 0.775487322608105])) sage: H.Sn(2,3,6,v) ([0.0374187795653158, 0.0457574702661884] U [0.120909196400478, 0.129247887101351] U [0.204085446231982, 0.212424136932855] U [0.287575863067145, 0.295914553768017] U [0.370752112898649, 0.379090803599522] U [0.454242529733812, 0.462581220434684] U [0.537418779565316, 0.545757470266188] U [0.620909196400478, 0.629247887101351] U [0.704085446231982, 0.712424136932855] U [0.787575863067145, 0.795914553768017] U [0.870752112898649, 0.879090803599522] U [0.954242529733812, 0.962581220434684])
An example over a number field:
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,0,a]) sage: v = K.real_places()[0] sage: H = E.height_function() sage: H.S(2,3,v) , H.Sn(2,3,1,v) (([0.142172065860075, 0.172845716928584] U [0.827154283071416, 0.857827934139925]), ([0.142172065860075, 0.172845716928584] U [0.827154283071416, 0.857827934139925])) sage: H.Sn(2,3,6,v) ([0.0236953443100124, 0.0288076194880974] U [0.137859047178569, 0.142971322356654] U [0.190362010976679, 0.195474286154764] U [0.304525713845236, 0.309637989023321] U [0.357028677643346, 0.362140952821431] U [0.471192380511903, 0.476304655689988] U [0.523695344310012, 0.528807619488097] U [0.637859047178569, 0.642971322356654] U [0.690362010976679, 0.695474286154764] U [0.804525713845236, 0.809637989023321] U [0.857028677643346, 0.862140952821431] U [0.971192380511903, 0.976304655689988])
-
alpha
(v, tol=0.01)¶ Return the constant \(\alpha_v\) associated to the embedding
v
.INPUT:
v
– an embedding of the base field into \(\RR\) or \(\CC\)
OUTPUT:
The constant \(\alpha_v\). In the notation of [CPS2006] and [Tho2010] (section 3.2), \(\alpha_v^3=\epsilon_v\). The result is cached since it only depends on the curve.
EXAMPLES:
Example 1 from [CPS2006]:
sage: K.<i>=QuadraticField(-1) sage: E = EllipticCurve([0,0,0,1+5*i,3+i]) sage: H = E.height_function() sage: alpha = H.alpha(K.places()[0]) sage: alpha 1.12272013439355
Compare with \(\log(\epsilon_v)=0.344562...\) in [CPS2006]:
sage: 3*alpha.log() 0.347263296676126
-
base_field
()¶ Return the base field.
EXAMPLES:
sage: E = EllipticCurve([0,0,0,0,1]) sage: H = E.height_function() sage: H.base_field() Rational Field
-
complex_intersection_is_empty
(Bk, v, verbose=False, use_half=True)¶ Returns True iff an intersection of \(T_n^{(v)}\) sets is empty.
INPUT:
Bk
(list) - a list of reals.v
(embedding) - a complex embedding of the number field.verbose
(boolean, default False) - verbosity flag.use_half
(boolean, default False) - if True, use only half the fundamental region.
OUTPUT:
True or False, according as the intersection of the unions of intervals \(T_n^{(v)}(-b,b)\) for \(b\) in the list
Bk
(see [Tho2010], section 7) is empty or not. WhenBk
is the list of \(b=\sqrt{B_n(\mu)}\) for \(n=1,2,3,\dots\) for some \(\mu>0\) this means that all non-torsion points on \(E\) with everywhere good reduction have canonical height strictly greater than \(\mu\), by [Tho2010], Proposition 7.8.EXAMPLES:
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,0,a]) sage: v = K.complex_embeddings()[0] sage: H = E.height_function()
The following two lines prove that the heights of non-torsion points on \(E\) with everywhere good reduction have canonical height strictly greater than 0.02, but fail to prove the same for 0.03. For the first proof, using only \(n=1,2,3\) is not sufficient:
sage: H.complex_intersection_is_empty([H.B(n,0.02) for n in [1,2,3]],v) # long time (~6s) False sage: H.complex_intersection_is_empty([H.B(n,0.02) for n in [1,2,3,4]],v) True sage: H.complex_intersection_is_empty([H.B(n,0.03) for n in [1,2,3,4]],v) # long time (4s) False
Using \(n\le6\) enables us to prove the lower bound 0.03. Note that it takes longer when the result is
False
than when it isTrue
:sage: H.complex_intersection_is_empty([H.B(n,0.03) for n in [1..6]],v) True
-
curve
()¶ Return the elliptic curve.
EXAMPLES:
sage: E = EllipticCurve([0,0,0,0,1]) sage: H = E.height_function() sage: H.curve() Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field
-
e_p
(p)¶ Return the exponent of the group over the residue field at
p
.INPUT:
p
- a prime ideal of \(K\) (or a prime number if \(K=\QQ\)).
OUTPUT:
A positive integer \(e_p\), the exponent of the group of nonsingular points on the reduction of the elliptic curve modulo \(p\). The result is cached.
EXAMPLES:
sage: K.<i>=QuadraticField(-1) sage: E = EllipticCurve([0,0,0,1+5*i,3+i]) sage: H = E.height_function() sage: H.e_p(K.prime_above(2)) 2 sage: H.e_p(K.prime_above(3)) 10 sage: H.e_p(K.prime_above(5)) 9 sage: E.conductor().norm().factor() 2^10 * 20921 sage: p1, p2 = K.primes_above(20921) sage: E.local_data(p1) Local data at Fractional ideal (-40*i + 139): Reduction type: bad split multiplicative ... sage: H.e_p(p1) 20920 sage: E.local_data(p2) Local data at Fractional ideal (40*i + 139): Reduction type: good ... sage: H.e_p(p2) 20815
-
fk_intervals
(v=None, N=20, domain=Complex Interval Field with 53 bits of precision)¶ Return a function approximating the Weierstrass function, with error.
INPUT:
v
(embedding) - an embedding of the number field. If None (default) use the real embedding if the field is \(\QQ\) and raise an error for other fields.N
(int) - The number of terms to use in the \(q\)-expansion of \(\wp\).domain
(complex field) - the model of \(\CC\) to use, for exampleCDF
ofCIF
(default).
OUTPUT:
A pair of functions fk, err which can be evaluated at complex numbers \(z\) (in the correct
domain
) to give an approximation to \(\wp(z)\) and an upper bound on the error, respectively. The Weierstrass function returned is with respect to the normalised lattice \([1,\tau]\) associated to the given embedding.EXAMPLES:
sage: E = EllipticCurve('37a') sage: L = E.period_lattice() sage: w1, w2 = L.normalised_basis() sage: z = CDF(0.3, 0.4)
Compare the value give by the standard elliptic exponential (scaled since
fk
is with respect to the normalised lattice):sage: L.elliptic_exponential(z*w2, to_curve=False)[0] * w2 ** 2 -1.82543539306049 - 2.49336319992847*I
to the value given by this function, and see the error:
sage: fk, err = E.height_function().fk_intervals(N=10) sage: fk(CIF(z)) -1.82543539306049? - 2.49336319992847?*I sage: err(CIF(z)) 2.71750621458744e-31
The same, but in the domain
CDF
instad ofCIF
:sage: fk, err = E.height_function().fk_intervals(N=10, domain=CDF) sage: fk(z) -1.8254353930604... - 2.493363199928...*I
-
min
(tol, n_max, verbose=False)¶ Returns a lower bound for all points of infinite order.
INPUT:
tol
- tolerance in output (see below).n_max
- how many multiples to use in iteration.verbose
(boolean, default False) - verbosity flag.
OUTPUT:
A positive real \(\mu\) for which it has been established rigorously that every point of infinite order on the elliptic curve (defined over its ground field) has canonical height greater than \(\mu\), and such that it is not possible (at least without increasing
n_max
) to prove the same for \(\mu\cdot\text{tol}\).EXAMPLES:
Example 1 from [CS2006] (where the same lower bound of 0.1126 was given):
sage: E = EllipticCurve([1, 0, 1, 421152067, 105484554028056]) # 60490d1 sage: E.height_function().min(.0001, 5) 0.0011263287309893311
Example 10.1 from [Tho2010] (where a lower bound of 0.18 was given):
sage: K.<i> = QuadraticField(-1) sage: E = EllipticCurve([0,0,0,91-26*i,-144-323*i]) sage: H = E.height_function() sage: H.min(0.1,4) # long time (8.1s) 0.1621049443313762
Example 10.2 from [Tho2010]:
sage: K.<i> = QuadraticField(-1) sage: E = EllipticCurve([0,1-i,i,-i,0]) sage: H = E.height_function() sage: H.min(0.01,5) # long time (4s) 0.020153685521979152
In this example the point \(P=(0,0)\) has height 0.023 so our lower bound is quite good:
sage: P = E((0,0)) sage: P.height() 0.0230242154471211
Example 10.3 from [Tho2010] (where the same bound of 0.0625 is given):
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,-3*a-a^2,a^2]) sage: H = E.height_function() sage: H.min(0.1,5) # long time (7s) 0.0625
More examples over \(\QQ\):
sage: E = EllipticCurve('37a') sage: h = E.height_function() sage: h.min(.01, 5) 0.03987318057488725 sage: E.gen(0).height() 0.0511114082399688
After base change the lower bound can decrease:
sage: K.<a> = QuadraticField(-5) sage: E.change_ring(K).height_function().min(0.5, 10) # long time (8s) 0.04419417382415922 sage: E = EllipticCurve('389a') sage: h = E.height_function() sage: h.min(0.1, 5) 0.05731275270029196 sage: [P.height() for P in E.gens()] [0.686667083305587, 0.327000773651605]
-
min_gr
(tol, n_max, verbose=False)¶ Returns a lower bound for points of infinite order with good reduction.
INPUT:
tol
- tolerance in output (see below).n_max
- how many multiples to use in iteration.verbose
(boolean, default False) - verbosity flag.
OUTPUT:
A positive real \(\mu\) for which it has been established rigorously that every point of infinite order on the elliptic curve (defined over its ground field), which has good reduction at all primes, has canonical height greater than \(\mu\), and such that it is not possible (at least without increasing
n_max
) to prove the same for \(\mu\cdot\text{tol}\).EXAMPLES:
Example 1 from [CS2006] (where a lower bound of 1.9865 was given):
sage: E = EllipticCurve([1, 0, 1, 421152067, 105484554028056]) # 60490d1 sage: E.height_function().min_gr(.0001, 5) 1.98684388146518
Example 10.1 from [Tho2010] (where a lower bound of 0.18 was given):
sage: K.<i> = QuadraticField(-1) sage: E = EllipticCurve([0,0,0,91-26*i,-144-323*i]) sage: H = E.height_function() sage: H.min_gr(0.1,4) # long time (8.1s) 0.1621049443313762
Example 10.2 from [Tho2010]:
sage: K.<i> = QuadraticField(-1) sage: E = EllipticCurve([0,1-i,i,-i,0]) sage: H = E.height_function() sage: H.min_gr(0.01, 5) # long time 0.020153685521979152
In this example the point \(P=(0,0)\) has height 0.023 so our lower bound is quite good:
sage: P = E((0,0)) sage: P.has_good_reduction() True sage: P.height() 0.0230242154471211
Example 10.3 from [Tho2010] (where the same bound of 0.25 is given):
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,-3*a-a^2,a^2]) sage: H = E.height_function() sage: H.min_gr(0.1,5) # long time (7.2s) 0.25
-
psi
(xi, v)¶ Return the normalised elliptic log of a point with this x-coordinate.
INPUT:
xi
(real) - the real x-coordinate of a point on the curve in the connected component with respect to a real embedding.v
(embedding) - a real embedding of the number field.
OUTPUT:
A real number in the interval [0.5,1] giving the elliptic logarithm of a point on \(E\) with \(x\)-coordinate
xi
, on the connected component with respect to the embedding \(v\), scaled by the real period.EXAMPLES:
An example over \(\QQ\):
sage: E = EllipticCurve('389a') sage: v = QQ.places()[0] sage: L = E.period_lattice(v) sage: P = E.lift_x(10/9) sage: L(P) 1.53151606047462 sage: L(P) / L.real_period() 0.615014189772115 sage: H = E.height_function() sage: H.psi(10/9,v) 0.615014189772115
An example over a number field:
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,0,a]) sage: P = E.lift_x(1/3*a^2 + a + 5/3) sage: v = K.real_places()[0] sage: L = E.period_lattice(v) sage: L(P) 3.51086196882538 sage: L(P) / L.real_period() 0.867385122699931 sage: xP = v(P.xy()[0]) sage: H = E.height_function() sage: H.psi(xP,v) 0.867385122699931 sage: H.psi(1.23,v) 0.785854718241495
-
real_intersection_is_empty
(Bk, v)¶ Returns True iff an intersection of \(S_n^{(v)}\) sets is empty.
INPUT:
Bk
(list) - a list of reals.v
(embedding) - a real embedding of the number field.
OUTPUT:
True or False, according as the intersection of the unions of intervals \(S_n^{(v)}(-b,b)\) for \(b\) in the list
Bk
is empty or not. WhenBk
is the list of \(b=B_n(\mu)\) for \(n=1,2,3,\dots\) for some \(\mu>0\) this means that all non-torsion points on \(E\) with everywhere good reduction have canonical height strictly greater than \(\mu\), by [Tho2010], Proposition 6.2.EXAMPLES:
An example over \(\QQ\):
sage: E = EllipticCurve('389a') sage: v = QQ.places()[0] sage: H = E.height_function()
The following two lines prove that the heights of non-torsion points on \(E\) with everywhere good reduction have canonical height strictly greater than 0.2, but fail to prove the same for 0.3:
sage: H.real_intersection_is_empty([H.B(n,0.2) for n in srange(1,10)],v) True sage: H.real_intersection_is_empty([H.B(n,0.3) for n in srange(1,10)],v) False
An example over a number field:
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,0,a]) sage: v = K.real_places()[0] sage: H = E.height_function()
The following two lines prove that the heights of non-torsion points on \(E\) with everywhere good reduction have canonical height strictly greater than 0.07, but fail to prove the same for 0.08:
sage: H.real_intersection_is_empty([H.B(n,0.07) for n in srange(1,5)],v) # long time (3.3s) True sage: H.real_intersection_is_empty([H.B(n,0.08) for n in srange(1,5)],v) False
-
tau
(v)¶ Return the normalised upper half-plane parameter \(\tau\) for the period lattice with respect to the embedding \(v\).
INPUT:
v
(embedding) - a real or complex embedding of the number field.
OUTPUT:
(Complex) \(\tau = \omega_1/\omega_2\) in the fundamental region of the upper half-plane.
EXAMPLES:
sage: E = EllipticCurve('37a') sage: H = E.height_function() sage: H.tau(QQ.places()[0]) 1.22112736076463*I
-
test_mu
(mu, N, verbose=True)¶ Return
True
if we can prove that \(\mu\) is a lower bound.INPUT:
mu
(real) - a positive real numberN
(integer) - upper bound on the multiples to be used.verbose
(boolean, default True) - verbosity flag.
OUTPUT:
True
orFalse
, according to whether we succeed in proving that \(\mu\) is a lower bound for the canonical heights of points of infinite order with everywhere good reduction.Note
A
True
result is rigorous;False
only means that the attempt failed: trying again with larger \(N\) may yieldTrue
.EXAMPLES:
sage: K.<a> = NumberField(x^3-2) sage: E = EllipticCurve([0,0,0,0,a]) sage: H = E.height_function()
This curve does have a point of good reduction whose canonical point is approximately 1.68:
sage: P = E.gens(lim3=5)[0]; P (1/3*a^2 + a + 5/3 : -2*a^2 - 4/3*a - 5/3 : 1) sage: P.height() 1.68038085233673 sage: P.has_good_reduction() True
Using \(N=5\) we can prove that 0.1 is a lower bound (in fact we only need \(N=2\)), but not that 0.2 is:
sage: H.test_mu(0.1, 5) B_1(0.100000000000000) = 1.51580969677387 B_2(0.100000000000000) = 0.932072561526720 True sage: H.test_mu(0.2, 5) B_1(0.200000000000000) = 2.04612906979932 B_2(0.200000000000000) = 3.09458988474327 B_3(0.200000000000000) = 27.6251108409484 B_4(0.200000000000000) = 1036.24722370223 B_5(0.200000000000000) = 3.67090854562318e6 False
Since 0.1 is a lower bound we can deduce that the point \(P\) is either primitive or divisible by either 2 or 3. In fact it is primitive:
sage: (P.height()/0.1).sqrt() 4.09924487233530 sage: P.division_points(2) [] sage: P.division_points(3) []
-
wp_c
(v)¶ Return a bound for the Weierstrass \(\wp\)-function.
INPUT:
v
(embedding) - a real or complex embedding of the number field.
OUTPUT:
(Real) \(c>0\) such that
\[|\wp(z) - z^-2| \le \frac{c^2|z|^2}{1-c|z|^2}\]whenever \(c|z|^2<1\). Given the recurrence relations for the Laurent series expansion of \(\wp\), it is easy to see that there is such a constant \(c\). [Reference?]
EXAMPLES:
sage: E = EllipticCurve('37a') sage: H = E.height_function() sage: H.wp_c(QQ.places()[0]) 2.68744508779950 sage: K.<i>=QuadraticField(-1) sage: E = EllipticCurve([0,0,0,1+5*i,3+i]) sage: H = E.height_function() sage: H.wp_c(K.places()[0]) 2.66213425640096
-
wp_intervals
(v=None, N=20, abs_only=False)¶ Return a function approximating the Weierstrass function.
INPUT:
v
(embedding) - an embedding of the number field. If None (default) use the real embedding if the field is \(\QQ\) and raise an error for other fields.N
(int, default 20) - The number of terms to use in the \(q\)-expansion of \(\wp\).abs_only
(boolean, default False) - flag to determine whether (if True) the error adjustment should use the absolute value or (if False) the real and imaginary parts.
OUTPUT:
A function wp which can be evaluated at complex numbers \(z\) to give an approximation to \(\wp(z)\). The Weierstrass function returned is with respect to the normalised lattice \([1,\tau]\) associated to the given embedding. For \(z\) which are not near a lattice point the function
fk
is used, otherwise a better approximation is used.EXAMPLES:
sage: E = EllipticCurve('37a') sage: wp = E.height_function().wp_intervals() sage: z = CDF(0.3, 0.4) sage: wp(CIF(z)) -1.82543539306049? - 2.4933631999285?*I sage: L = E.period_lattice() sage: w1, w2 = L.normalised_basis() sage: L.elliptic_exponential(z*w2, to_curve=False)[0] * w2^2 -1.82543539306049 - 2.49336319992847*I sage: z = CDF(0.3, 0.1) sage: wp(CIF(z)) 8.5918243572165? - 5.4751982004351?*I sage: L.elliptic_exponential(z*w2, to_curve=False)[0] * w2^2 8.59182435721650 - 5.47519820043503*I
-
wp_on_grid
(v, N, half=False)¶ Return an array of the values of \(\wp\) on an \(N\times N\) grid.
INPUT:
v
(embedding) - an embedding of the number field.N
(int) - The number of terms to use in the \(q\)-expansion of \(\wp\).half
(boolean, default False) - if True, use an array of size \(N\times N/2\) instead of \(N\times N\).
OUTPUT:
An array of size either \(N\times N/2\) or \(N\times N\) whose \((i,j)\) entry is the value of the Weierstrass \(\wp\)-function at \((i+.5)/N + (j+.5)*\tau/N\), a grid of points in the fundamental region for the lattice \([1,\tau]\).
EXAMPLES:
sage: E = EllipticCurve('37a') sage: H = E.height_function() sage: v = QQ.places()[0]
The array of values on the grid shows symmetry, since \(\wp\) is even:
sage: H.wp_on_grid(v,4) array([[25.43920182, 5.28760943, 5.28760943, 25.43920182], [ 6.05099485, 1.83757786, 1.83757786, 6.05099485], [ 6.05099485, 1.83757786, 1.83757786, 6.05099485], [25.43920182, 5.28760943, 5.28760943, 25.43920182]])
The array of values on the half-grid:
sage: H.wp_on_grid(v,4,True) array([[25.43920182, 5.28760943], [ 6.05099485, 1.83757786], [ 6.05099485, 1.83757786], [25.43920182, 5.28760943]])
-
-
class
sage.schemes.elliptic_curves.height.
UnionOfIntervals
(endpoints)¶ A class representing a finite union of closed intervals in \(\RR\) which can be scaled, shifted, intersected, etc.
The intervals are represented as an ordered list of their endpoints, which may include \(-\infty\) and \(+\infty\).
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import UnionOfIntervals sage: R = UnionOfIntervals([1,2,3,infinity]); R ([1, 2] U [3, +Infinity]) sage: R + 5 ([6, 7] U [8, +Infinity]) sage: ~R ([-Infinity, 1] U [2, 3]) sage: ~R | (10*R + 100) ([-Infinity, 1] U [2, 3] U [110, 120] U [130, +Infinity])
Todo
Unify
UnionOfIntervals
with the classRealSet
introduced by trac ticket #13125; see trac ticket #16063.-
finite_endpoints
()¶ Returns the finite endpoints of this union of intervals.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import UnionOfIntervals sage: UnionOfIntervals([0,1]).finite_endpoints() [0, 1] sage: UnionOfIntervals([-infinity, 0, 1, infinity]).finite_endpoints() [0, 1]
-
classmethod
intersection
(L)¶ Return the intersection of a list of UnionOfIntervals.
INPUT:
L
(list) – a list of UnionOfIntervals instances
OUTPUT:
A new UnionOfIntervals instance representing the intersection of the UnionOfIntervals in the list.
Note
This is a class method.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import UnionOfIntervals sage: A = UnionOfIntervals([1,3,5,7]); A ([1, 3] U [5, 7]) sage: B = A+1; B ([2, 4] U [6, 8]) sage: A.intersection([A,B]) ([2, 3] U [6, 7])
-
intervals
()¶ Returns the intervals in self, as a list of 2-tuples.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import UnionOfIntervals sage: UnionOfIntervals(list(range(10))).intervals() [(0, 1), (2, 3), (4, 5), (6, 7), (8, 9)] sage: UnionOfIntervals([-infinity, pi, 17, infinity]).intervals() [(-Infinity, pi), (17, +Infinity)]
-
is_empty
()¶ Returns whether self is empty.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import UnionOfIntervals sage: UnionOfIntervals([3,4]).is_empty() False sage: all = UnionOfIntervals([-infinity, infinity]) sage: all.is_empty() False sage: (~all).is_empty() True sage: A = UnionOfIntervals([0,1]) & UnionOfIntervals([2,3]) sage: A.is_empty() True
-
static
join
(L, condition)¶ Utility function to form the union or intersection of a list of UnionOfIntervals.
INPUT:
L
(list) – a list of UnionOfIntervals instancescondition
(function) – eitherany
orall
, or some other boolean function of a list of boolean values.
OUTPUT:
A new UnionOfIntervals instance representing the subset of ‘RR’ equal to those reals in any/all/condition of the UnionOfIntervals in the list.
Note
This is a static method for the class.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import UnionOfIntervals sage: A = UnionOfIntervals([1,3,5,7]); A ([1, 3] U [5, 7]) sage: B = A+1; B ([2, 4] U [6, 8]) sage: A.join([A,B],any) # union ([1, 4] U [5, 8]) sage: A.join([A,B],all) # intersection ([2, 3] U [6, 7]) sage: A.join([A,B],sum) # symmetric difference ([1, 2] U [3, 4] U [5, 6] U [7, 8])
-
classmethod
union
(L)¶ Return the union of a list of UnionOfIntervals.
INPUT:
L
(list) – a list of UnionOfIntervals instances
OUTPUT:
A new UnionOfIntervals instance representing the union of the UnionOfIntervals in the list.
Note
This is a class method.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import UnionOfIntervals sage: A = UnionOfIntervals([1,3,5,7]); A ([1, 3] U [5, 7]) sage: B = A+1; B ([2, 4] U [6, 8]) sage: A.union([A,B]) ([1, 4] U [5, 8])
-
-
sage.schemes.elliptic_curves.height.
eps
(err, is_real)¶ Return a Real or Complex interval centered on 0 with radius err.
INPUT:
err
(real) – a positive real number, the radius of the intervalis_real
(boolean) – if True, returns a real interval in RIF, else a complex interval in CIF
OUTPUT:
An element of RIF or CIF (as specified), centered on 0, with given radius.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import eps sage: eps(0.01, True) 0.0? sage: eps(0.01, False) 0.0? + 0.0?*I
-
sage.schemes.elliptic_curves.height.
inf_max_abs
(f, g, D)¶ Returns \(\inf_D(\max(|f|, |g|))\).
INPUT:
f
,g
(polynomials) – real univariate polynomialsD
(UnionOfIntervals) – a subset of \(\RR\)
OUTPUT:
A real number approximating the value of \(\inf_D(\max(|f|, |g|))\).
ALGORITHM:
The extreme values must occur at an endpoint of a subinterval of \(D\) or at a point where one of \(f\), \(f'\), \(g\), \(g'\), \(f\pm g\) is zero.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import inf_max_abs, UnionOfIntervals sage: x = polygen(RR) sage: f = (x-10)^4+1 sage: g = 2*x^3+100 sage: inf_max_abs(f,g,UnionOfIntervals([1,2,3,4,5,6])) 425.638201706391 sage: r0 = (f-g).roots()[0][0] sage: r0 5.46053402234697 sage: max(abs(f(r0)),abs(g(r0))) 425.638201706391
-
sage.schemes.elliptic_curves.height.
min_on_disk
(f, tol, max_iter=10000)¶ Returns the minimum of a real-valued complex function on a square.
INPUT:
f
– a function from CIF to RIFtol
(real) – a positive real numbermax_iter
(integer, default 10000) – a positive integer bounding the number of iterations to be used
OUTPUT:
A 2-tuple \((s,t)\), where \(t=f(s)\) and \(s\) is a CIF element contained in the disk \(|z|\le1\), at which \(f\) takes its minumum value.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import min_on_disk sage: f = lambda x: (x^2+100).abs() sage: s, t = min_on_disk(f, 0.0001) sage: s, f(s), t (0.01? + 1.00?*I, 99.01?, 99.0000000000000)
-
sage.schemes.elliptic_curves.height.
nonneg_region
(f)¶ Returns the UnionOfIntervals representing the region where
f
is non-negative.INPUT:
f
(polynomial) – a univariate polynomial over \(\RR\).
OUTPUT:
A UnionOfIntervals representing the set \(\{x \in\RR mid f(x) \ge 0\}\).
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import nonneg_region sage: x = polygen(RR) sage: nonneg_region(x^2-1) ([-Infinity, -1.00000000000000] U [1.00000000000000, +Infinity]) sage: nonneg_region(1-x^2) ([-1.00000000000000, 1.00000000000000]) sage: nonneg_region(1-x^3) ([-Infinity, 1.00000000000000]) sage: nonneg_region(x^3-1) ([1.00000000000000, +Infinity]) sage: nonneg_region((x-1)*(x-2)) ([-Infinity, 1.00000000000000] U [2.00000000000000, +Infinity]) sage: nonneg_region(-(x-1)*(x-2)) ([1.00000000000000, 2.00000000000000]) sage: nonneg_region((x-1)*(x-2)*(x-3)) ([1.00000000000000, 2.00000000000000] U [3.00000000000000, +Infinity]) sage: nonneg_region(-(x-1)*(x-2)*(x-3)) ([-Infinity, 1.00000000000000] U [2.00000000000000, 3.00000000000000]) sage: nonneg_region(x^4+1) ([-Infinity, +Infinity]) sage: nonneg_region(-x^4-1) ()
-
sage.schemes.elliptic_curves.height.
rat_term_CIF
(z, try_strict=True)¶ Compute the value of \(u/(1-u)^2\) in
CIF
, where \(u=\exp(2\pi i z)\).INPUT:
z
(complex) – a CIF elementtry_strict
(bool) – flag
EXAMPLES:
sage: from sage.schemes.elliptic_curves.height import rat_term_CIF sage: z = CIF(0.5,0.2) sage: rat_term_CIF(z) -0.172467461182437? + 0.?e-16*I sage: rat_term_CIF(z, False) -0.172467461182437? + 0.?e-16*I