Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

# -*- coding: utf-8 -*- 

""" 

Construct elliptic curves as Jacobians 

 

An elliptic curve is a genus one curve with a designated point. The 

Jacobian of a genus-one curve can be defined as the set of line 

bundles on the curve, and is isomorphic to the original genus-one 

curve. It is also an elliptic curve with the trivial line bundle as 

designated point. The utility of this construction is that we can 

construct elliptic curves without having to specify which point we 

take as the origin. 

 

EXAMPLES:: 

 

sage: R.<u,v,w> = QQ[] 

sage: Jacobian(u^3+v^3+w^3) 

Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

sage: Jacobian(u^4+v^4+w^2) 

Elliptic Curve defined by y^2 = x^3 - 4*x over Rational Field 

 

sage: C = Curve(u^3+v^3+w^3) 

sage: Jacobian(C) 

Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

 

sage: P2.<u,v,w> = ProjectiveSpace(2, QQ) 

sage: C = P2.subscheme(u^3+v^3+w^3) 

sage: Jacobian(C) 

Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

 

One can also define Jacobians of varieties that are not genus-one 

curves. These are not implemented in this module, but we call the 

relevant functionality:: 

 

sage: R.<x> = PolynomialRing(QQ) 

sage: f = x**5 + 1184*x**3 + 1846*x**2 + 956*x + 560 

sage: C = HyperellipticCurve(f) 

sage: Jacobian(C) 

Jacobian of Hyperelliptic Curve over Rational Field defined 

by y^2 = x^5 + 1184*x^3 + 1846*x^2 + 956*x + 560 

 

REFERENCES: 

 

.. [WpJacobianVariety] :wikipedia:`Jacobian_variety` 

""" 

 

############################################################################## 

# Copyright (C) 2013 Volker Braun <vbraun.name@gmail.com> 

# 

# Distributed under the terms of the GNU General Public License (GPL) 

# 

# The full text of the GPL is available at: 

# 

# http://www.gnu.org/licenses/ 

############################################################################## 

 

from sage.rings.all import QQ 

from sage.schemes.elliptic_curves.constructor import EllipticCurve 

 

 

 

 

def Jacobian(X, **kwds): 

""" 

Return the Jacobian. 

 

INPUT: 

 

- ``X`` -- polynomial, algebraic variety, or anything else that 

has a Jacobian elliptic curve. 

 

- ``kwds`` -- optional keyword arguments. 

 

The input ``X`` can be one of the following: 

 

* A polynomial, see :func:`Jacobian_of_equation` for details. 

 

* A curve, see :func:`Jacobian_of_curve` for details. 

 

EXAMPLES:: 

 

sage: R.<u,v,w> = QQ[] 

sage: Jacobian(u^3+v^3+w^3) 

Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

 

sage: C = Curve(u^3+v^3+w^3) 

sage: Jacobian(C) 

Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

 

sage: P2.<u,v,w> = ProjectiveSpace(2, QQ) 

sage: C = P2.subscheme(u^3+v^3+w^3) 

sage: Jacobian(C) 

Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

 

sage: Jacobian(C, morphism=True) 

Scheme morphism: 

From: Closed subscheme of Projective Space of dimension 2 over Rational Field defined by: 

u^3 + v^3 + w^3 

To: Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

Defn: Defined on coordinates by sending (u : v : w) to 

(-u^4*v^4*w - u^4*v*w^4 - u*v^4*w^4 : 

1/2*u^6*v^3 - 1/2*u^3*v^6 - 1/2*u^6*w^3 + 1/2*v^6*w^3 + 1/2*u^3*w^6 - 1/2*v^3*w^6 : 

u^3*v^3*w^3) 

""" 

try: 

return X.jacobian(**kwds) 

except AttributeError: 

pass 

 

morphism = kwds.pop('morphism', False) 

from sage.rings.polynomial.multi_polynomial_element import is_MPolynomial 

if is_MPolynomial(X): 

if morphism: 

from sage.schemes.curves.constructor import Curve 

return Jacobian_of_equation(X, curve=Curve(X), **kwds) 

else: 

return Jacobian_of_equation(X, **kwds) 

 

from sage.schemes.generic.scheme import is_Scheme 

if is_Scheme(X) and X.dimension() == 1: 

return Jacobian_of_curve(X, morphism=morphism, **kwds) 

 

 

def Jacobian_of_curve(curve, morphism=False): 

""" 

Return the Jacobian of a genus-one curve 

 

INPUT: 

 

- ``curve`` -- a one-dimensional algebraic variety of genus one. 

 

OUTPUT: 

 

Its Jacobian elliptic curve. 

 

EXAMPLES:: 

 

sage: R.<u,v,w> = QQ[] 

sage: C = Curve(u^3+v^3+w^3) 

sage: Jacobian(C) 

Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field 

""" 

eqn = None 

try: 

eqn = curve.defining_polynomial() 

except AttributeError: 

pass 

if len(curve.defining_polynomials()) == 1: 

eqn = curve.defining_polynomials()[0] 

if eqn is not None: 

if morphism: 

return Jacobian_of_equation(eqn, curve=curve) 

else: 

return Jacobian_of_equation(eqn) 

raise NotImplementedError('Jacobian for this curve is not implemented') 

 

 

def Jacobian_of_equation(polynomial, variables=None, curve=None): 

r""" 

Construct the Jacobian of a genus-one curve given by a polynomial. 

 

INPUT: 

 

- ``F`` -- a polynomial defining a plane curve of genus one. May 

be homogeneous or inhomogeneous. 

 

- ``variables`` -- list of two or three variables or ``None`` 

(default). The inhomogeneous or homogeneous coordinates. By 

default, all variables in the polynomial are used. 

 

- ``curve`` -- the genus-one curve defined by ``polynomial`` or # 

``None`` (default). If specified, suitable morphism from the 

jacobian elliptic curve to the curve is returned. 

 

OUTPUT: 

 

An elliptic curve in short Weierstrass form isomorphic to the 

curve ``polynomial=0``. If the optional argument ``curve`` is 

specified, a rational multicover from the Jacobian elliptic curve 

to the genus-one curve is returned. 

 

EXAMPLES:: 

 

sage: R.<a,b,c> = QQ[] 

sage: f = a^3+b^3+60*c^3 

sage: Jacobian(f) 

Elliptic Curve defined by y^2 = x^3 - 24300 over Rational Field 

sage: Jacobian(f.subs(c=1)) 

Elliptic Curve defined by y^2 = x^3 - 24300 over Rational Field 

 

If we specify the domain curve the birational covering is returned:: 

 

sage: h = Jacobian(f, curve=Curve(f)); h 

Scheme morphism: 

From: Projective Plane Curve over Rational Field defined by a^3 + b^3 + 60*c^3 

To: Elliptic Curve defined by y^2 = x^3 - 24300 over Rational Field 

Defn: Defined on coordinates by sending (a : b : c) to 

(-216000*a^4*b^4*c - 12960000*a^4*b*c^4 - 12960000*a*b^4*c^4 : 

108000*a^6*b^3 - 108000*a^3*b^6 - 6480000*a^6*c^3 + 6480000*b^6*c^3 + 388800000*a^3*c^6 - 388800000*b^3*c^6 : 

216000*a^3*b^3*c^3) 

 

sage: h([1,-1,0]) 

(0 : 1 : 0) 

 

Plugging in the polynomials defining `h` allows us to verify that 

it is indeed a rational morphism to the elliptic curve:: 

 

sage: E = h.codomain() 

sage: E.defining_polynomial()(h.defining_polynomials()).factor() 

(2519424000000000) * c^3 * b^3 * a^3 * (a^3 + b^3 + 60*c^3) * 

(a^9*b^6 + a^6*b^9 - 120*a^9*b^3*c^3 + 900*a^6*b^6*c^3 - 120*a^3*b^9*c^3 + 

3600*a^9*c^6 + 54000*a^6*b^3*c^6 + 54000*a^3*b^6*c^6 + 3600*b^9*c^6 + 

216000*a^6*c^9 - 432000*a^3*b^3*c^9 + 216000*b^6*c^9) 

 

 

By specifying the variables, we can also construct an elliptic 

curve over a polynomial ring:: 

 

sage: R.<u,v,t> = QQ[] 

sage: Jacobian(u^3+v^3+t, variables=[u,v]) 

Elliptic Curve defined by y^2 = x^3 + (-27/4*t^2) over 

Multivariate Polynomial Ring in u, v, t over Rational Field 

 

TESTS:: 

 

sage: from sage.schemes.elliptic_curves.jacobian import Jacobian_of_equation 

sage: Jacobian_of_equation(f, variables=[a,b,c]) 

Elliptic Curve defined by y^2 = x^3 - 24300 over Rational Field 

""" 

from sage.schemes.toric.weierstrass import WeierstrassForm 

f, g = WeierstrassForm(polynomial, variables=variables) 

try: 

K = polynomial.base_ring() 

f = K(f) 

g = K(g) 

except (TypeError, ValueError): 

pass 

E = EllipticCurve([f, g]) 

if curve is None: 

return E 

X, Y, Z = WeierstrassForm(polynomial, variables=variables, transformation=True) 

from sage.schemes.elliptic_curves.weierstrass_transform import WeierstrassTransformation 

return WeierstrassTransformation(curve, E, [X*Z, Y, Z**3], 1)