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

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

381

382

383

384

385

386

387

388

389

390

391

392

393

394

395

396

397

398

399

400

401

402

403

404

405

406

407

408

409

410

411

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

455

456

457

458

459

460

461

462

463

464

465

466

467

468

469

470

471

472

473

474

475

476

477

478

479

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

521

522

523

524

525

526

527

528

529

530

531

532

533

534

535

536

537

538

539

540

541

542

543

544

545

546

547

548

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

641

642

643

644

645

646

647

648

649

650

651

652

653

654

655

656

657

658

659

660

661

662

663

664

665

666

667

668

669

670

671

672

673

674

675

676

677

678

679

680

681

682

683

684

685

686

687

688

689

690

691

692

693

694

695

696

697

698

699

700

701

702

703

704

705

706

707

708

709

710

711

712

713

714

715

716

717

718

719

720

721

722

723

724

725

726

727

728

729

730

731

732

733

734

735

736

737

738

739

740

741

742

743

744

745

746

747

748

749

750

751

752

753

754

755

756

757

758

759

760

761

762

763

764

765

766

767

768

769

770

771

772

773

774

775

776

777

778

779

780

781

782

783

784

785

786

787

788

789

790

791

792

793

794

795

796

797

798

799

800

801

802

803

804

805

806

807

808

809

810

811

812

813

814

815

816

817

818

819

820

821

822

823

824

825

826

827

828

829

830

831

832

833

834

835

836

837

838

839

840

841

842

843

844

845

846

847

848

849

850

851

852

853

854

855

856

857

858

859

860

861

862

863

864

865

866

867

868

869

870

871

872

873

874

875

876

877

878

879

880

881

882

883

884

885

886

887

888

889

890

891

892

893

894

895

896

897

898

899

900

901

902

903

904

905

906

907

908

909

910

911

912

913

914

915

916

917

918

919

920

921

922

923

924

925

926

927

928

929

930

931

932

933

934

935

936

937

938

939

940

941

942

943

944

945

946

947

948

949

950

951

952

953

954

955

956

957

958

959

960

961

962

963

964

965

966

967

968

969

970

971

972

973

974

975

976

977

978

979

980

981

982

983

984

985

986

987

988

989

990

991

992

993

994

995

996

997

998

999

1000

1001

1002

1003

1004

1005

1006

1007

1008

1009

1010

1011

1012

1013

1014

1015

1016

1017

1018

1019

1020

1021

1022

1023

1024

1025

1026

1027

1028

1029

1030

1031

1032

1033

1034

1035

1036

1037

1038

1039

1040

1041

1042

1043

1044

1045

1046

1047

1048

1049

1050

1051

1052

1053

1054

1055

1056

1057

1058

1059

1060

1061

1062

1063

1064

1065

1066

1067

1068

1069

1070

1071

1072

1073

1074

1075

1076

1077

1078

1079

1080

1081

1082

1083

1084

1085

1086

1087

1088

1089

1090

1091

1092

1093

1094

1095

1096

1097

1098

1099

1100

1101

1102

1103

1104

1105

1106

1107

1108

1109

1110

1111

1112

1113

1114

1115

1116

1117

1118

1119

1120

1121

1122

1123

1124

1125

1126

1127

1128

1129

1130

1131

1132

1133

1134

1135

1136

1137

1138

1139

1140

1141

1142

1143

1144

1145

1146

1147

1148

1149

1150

1151

1152

1153

1154

1155

1156

1157

1158

1159

1160

1161

1162

1163

1164

1165

1166

1167

1168

1169

1170

1171

1172

1173

1174

1175

1176

1177

1178

1179

1180

1181

1182

1183

1184

1185

1186

1187

1188

1189

1190

1191

1192

1193

1194

1195

1196

1197

1198

1199

1200

1201

1202

1203

1204

1205

1206

1207

1208

1209

1210

1211

1212

1213

1214

1215

1216

1217

1218

1219

1220

1221

1222

1223

1224

1225

1226

1227

1228

1229

1230

1231

1232

1233

1234

1235

1236

1237

1238

1239

1240

1241

1242

1243

1244

1245

1246

1247

1248

1249

1250

1251

1252

1253

1254

1255

1256

1257

1258

1259

1260

1261

1262

1263

1264

1265

1266

1267

1268

1269

1270

1271

1272

1273

1274

1275

1276

1277

1278

1279

1280

1281

1282

1283

1284

1285

1286

1287

1288

1289

1290

1291

1292

1293

1294

1295

1296

1297

1298

1299

1300

1301

1302

1303

1304

1305

1306

1307

1308

1309

1310

1311

1312

1313

1314

1315

1316

1317

1318

1319

1320

1321

1322

1323

1324

1325

1326

1327

1328

1329

1330

1331

1332

1333

1334

1335

1336

1337

1338

1339

1340

1341

1342

1343

1344

1345

1346

1347

1348

1349

1350

1351

1352

1353

1354

1355

1356

1357

1358

1359

1360

1361

1362

1363

1364

1365

1366

1367

1368

1369

1370

1371

1372

1373

1374

1375

1376

1377

1378

1379

1380

1381

1382

1383

1384

1385

1386

1387

1388

1389

1390

1391

1392

1393

1394

1395

1396

1397

1398

1399

1400

1401

1402

1403

1404

1405

1406

1407

1408

1409

1410

1411

1412

1413

1414

1415

1416

1417

1418

1419

1420

1421

1422

1423

1424

1425

1426

1427

1428

1429

1430

1431

1432

1433

1434

1435

1436

1437

1438

1439

1440

1441

1442

1443

1444

1445

1446

1447

1448

1449

1450

1451

1452

1453

1454

1455

1456

1457

1458

1459

1460

1461

1462

1463

1464

1465

1466

1467

1468

1469

1470

1471

1472

1473

1474

1475

1476

1477

1478

1479

1480

1481

1482

1483

1484

1485

""" 

Elliptic curve constructor 

 

AUTHORS: 

 

- William Stein (2005): Initial version 

 

- John Cremona (2008-01): EllipticCurve(j) fixed for all cases 

""" 

from __future__ import absolute_import 

 

#***************************************************************************** 

# Copyright (C) 2005 William Stein <wstein@gmail.com> 

# 

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

# 

# This code is distributed in the hope that it will be useful, 

# but WITHOUT ANY WARRANTY; without even the implied warranty of 

# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 

# General Public License for more details. 

# 

# The full text of the GPL is available at: 

# 

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

#***************************************************************************** 

from six import string_types, integer_types 

 

import sage.rings.all as rings 

 

from sage.rings.finite_rings.integer_mod_ring import is_IntegerModRing 

from sage.rings.polynomial.multi_polynomial_ring import is_MPolynomialRing 

from sage.rings.finite_rings.finite_field_constructor import is_FiniteField 

from sage.rings.number_field.number_field import is_NumberField 

from sage.rings.polynomial.multi_polynomial_element import is_MPolynomial 

from sage.rings.ring import is_Ring 

 

from sage.categories.fields import Fields 

_Fields = Fields() 

 

from sage.structure.sequence import Sequence 

from sage.structure.element import parent 

from sage.structure.factory import UniqueFactory 

from sage.symbolic.ring import SR 

from sage.symbolic.expression import is_SymbolicEquation 

 

 

class EllipticCurveFactory(UniqueFactory): 

r""" 

Construct an elliptic curve. 

 

In Sage, an elliptic curve is always specified by 

(the coefficients of) a long Weierstrass equation 

 

.. MATH:: 

 

y^2 + a_1 xy + a_3 y = x^3 + a_2 x^2 + a_4 x + a_6. 

 

INPUT: 

 

There are several ways to construct an elliptic curve: 

 

- ``EllipticCurve([a1,a2,a3,a4,a6])``: Elliptic curve with given 

`a`-invariants. The invariants are coerced into a common parent. 

If all are integers, they are coerced into the rational numbers. 

 

- ``EllipticCurve([a4,a6])``: Same as above, but `a_1=a_2=a_3=0`. 

 

- ``EllipticCurve(label)``: Returns the elliptic curve over `\QQ` 

from the Cremona database with the given label. The label is a 

string, such as ``"11a"`` or ``"37b2"``. The letters in the 

label *must* be lower case (Cremona's new labeling). 

 

- ``EllipticCurve(R, [a1,a2,a3,a4,a6])``: Create the elliptic 

curve over `R` with given `a`-invariants. Here `R` can be an 

arbitrary commutative ring, although most functionality is only 

implemented over fields. 

 

- ``EllipticCurve(j=j0)`` or ``EllipticCurve_from_j(j0)``: Return 

an elliptic curve with `j`-invariant ``j0``. 

 

- ``EllipticCurve(polynomial)``: Read off the `a`-invariants from 

the polynomial coefficients, see 

:func:`EllipticCurve_from_Weierstrass_polynomial`. 

 

Instead of giving the coefficients as a *list* of length 2 or 5, 

one can also give a *tuple*. 

 

EXAMPLES: 

 

We illustrate creating elliptic curves:: 

 

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

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

 

We create a curve from a Cremona label:: 

 

sage: EllipticCurve('37b2') 

Elliptic Curve defined by y^2 + y = x^3 + x^2 - 1873*x - 31833 over Rational Field 

sage: EllipticCurve('5077a') 

Elliptic Curve defined by y^2 + y = x^3 - 7*x + 6 over Rational Field 

sage: EllipticCurve('389a') 

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

 

Old Cremona labels are allowed:: 

 

sage: EllipticCurve('2400FF') 

Elliptic Curve defined by y^2 = x^3 + x^2 + 2*x + 8 over Rational Field 

 

Unicode labels are allowed:: 

 

sage: EllipticCurve(u'389a') 

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

 

We create curves over a finite field as follows:: 

 

sage: EllipticCurve([GF(5)(0),0,1,-1,0]) 

Elliptic Curve defined by y^2 + y = x^3 + 4*x over Finite Field of size 5 

sage: EllipticCurve(GF(5), [0, 0,1,-1,0]) 

Elliptic Curve defined by y^2 + y = x^3 + 4*x over Finite Field of size 5 

 

Elliptic curves over `\ZZ/N\ZZ` with `N` prime are of type 

"elliptic curve over a finite field":: 

 

sage: F = Zmod(101) 

sage: EllipticCurve(F, [2, 3]) 

Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 101 

sage: E = EllipticCurve([F(2), F(3)]) 

sage: type(E) 

<class 'sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field_with_category'> 

sage: E.category() 

Category of schemes over Ring of integers modulo 101 

 

In contrast, elliptic curves over `\ZZ/N\ZZ` with `N` composite 

are of type "generic elliptic curve":: 

 

sage: F = Zmod(95) 

sage: EllipticCurve(F, [2, 3]) 

Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 95 

sage: E = EllipticCurve([F(2), F(3)]) 

sage: type(E) 

<class 'sage.schemes.elliptic_curves.ell_generic.EllipticCurve_generic_with_category'> 

sage: E.category() 

Category of schemes over Ring of integers modulo 95 

 

The following is a curve over the complex numbers:: 

 

sage: E = EllipticCurve(CC, [0,0,1,-1,0]) 

sage: E 

Elliptic Curve defined by y^2 + 1.00000000000000*y = x^3 + (-1.00000000000000)*x over Complex Field with 53 bits of precision 

sage: E.j_invariant() 

2988.97297297297 

 

We can also create elliptic curves by giving the Weierstrass equation:: 

 

sage: x, y = var('x,y') 

sage: EllipticCurve(y^2 + y == x^3 + x - 9) 

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

 

sage: R.<x,y> = GF(5)[] 

sage: EllipticCurve(x^3 + x^2 + 2 - y^2 - y*x) 

Elliptic Curve defined by y^2 + x*y = x^3 + x^2 + 2 over Finite Field of size 5 

 

We can explicitly specify the `j`-invariant:: 

 

sage: E = EllipticCurve(j=1728); E; E.j_invariant(); E.label() 

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

1728 

'32a2' 

 

sage: E = EllipticCurve(j=GF(5)(2)); E; E.j_invariant() 

Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field of size 5 

2 

 

See :trac:`6657` :: 

 

sage: EllipticCurve(GF(144169),j=1728) 

Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 144169 

 

Elliptic curves over the same ring with the same Weierstrass 

coefficients are identical, even when they are constructed in 

different ways (see :trac:`11474`):: 

 

sage: EllipticCurve('11a3') is EllipticCurve(QQ, [0, -1, 1, 0, 0]) 

True 

 

By default, when a rational value of `j` is given, the constructed 

curve is a minimal twist (minimal conductor for curves with that 

`j`-invariant). This can be changed by setting the optional 

parameter ``minimal_twist``, which is True by default, to False:: 

 

sage: EllipticCurve(j=100) 

Elliptic Curve defined by y^2 = x^3 + x^2 + 3392*x + 307888 over Rational Field 

sage: E =EllipticCurve(j=100); E 

Elliptic Curve defined by y^2 = x^3 + x^2 + 3392*x + 307888 over Rational Field 

sage: E.conductor() 

33129800 

sage: E.j_invariant() 

100 

sage: E =EllipticCurve(j=100, minimal_twist=False); E 

Elliptic Curve defined by y^2 = x^3 + 488400*x - 530076800 over Rational Field 

sage: E.conductor() 

298168200 

sage: E.j_invariant() 

100 

 

Without this option, constructing the curve could take a long time 

since both `j` and `j-1728` have to be factored to compute the 

minimal twist (see :trac:`13100`):: 

 

sage: E = EllipticCurve_from_j(2^256+1,minimal_twist=False) 

sage: E.j_invariant() == 2^256+1 

True 

 

TESTS:: 

 

sage: R = ZZ['u', 'v'] 

sage: EllipticCurve(R, [1,1]) 

Elliptic Curve defined by y^2 = x^3 + x + 1 over Multivariate Polynomial Ring in u, v 

over Integer Ring 

 

We create a curve and a point over ``QQbar`` (see :trac:`6879`):: 

 

sage: E = EllipticCurve(QQbar,[0,1]) 

sage: E(0) 

(0 : 1 : 0) 

sage: E.base_field() 

Algebraic Field 

 

sage: E = EllipticCurve(RR,[1,2]); E; E.base_field() 

Elliptic Curve defined by y^2 = x^3 + 1.00000000000000*x + 2.00000000000000 over Real Field with 53 bits of precision 

Real Field with 53 bits of precision 

sage: EllipticCurve(CC,[3,4]); E; E.base_field() 

Elliptic Curve defined by y^2 = x^3 + 3.00000000000000*x + 4.00000000000000 over Complex Field with 53 bits of precision 

Elliptic Curve defined by y^2 = x^3 + 1.00000000000000*x + 2.00000000000000 over Real Field with 53 bits of precision 

Real Field with 53 bits of precision 

sage: E = EllipticCurve(QQbar,[5,6]); E; E.base_field() 

Elliptic Curve defined by y^2 = x^3 + 5*x + 6 over Algebraic Field 

Algebraic Field 

 

See :trac:`6657` :: 

 

sage: EllipticCurve(3,j=1728) 

Traceback (most recent call last): 

... 

ValueError: First parameter (if present) must be a ring when j is specified 

 

sage: EllipticCurve(GF(5),j=3/5) 

Traceback (most recent call last): 

... 

ValueError: First parameter must be a ring containing 3/5 

 

If the universe of the coefficients is a general field, the object 

constructed has type EllipticCurve_field. Otherwise it is 

EllipticCurve_generic. See :trac:`9816` :: 

 

sage: E = EllipticCurve([QQbar(1),3]); E 

Elliptic Curve defined by y^2 = x^3 + x + 3 over Algebraic Field 

sage: type(E) 

<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'> 

 

sage: E = EllipticCurve([RR(1),3]); E 

Elliptic Curve defined by y^2 = x^3 + 1.00000000000000*x + 3.00000000000000 over Real Field with 53 bits of precision 

sage: type(E) 

<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'> 

 

sage: E = EllipticCurve([i,i]); E 

Elliptic Curve defined by y^2 = x^3 + I*x + I over Symbolic Ring 

sage: type(E) 

<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'> 

sage: E.category() 

Category of schemes over Symbolic Ring 

sage: SR in Fields() 

True 

 

sage: F = FractionField(PolynomialRing(QQ,'t')) 

sage: t = F.gen() 

sage: E = EllipticCurve([t,0]); E 

Elliptic Curve defined by y^2 = x^3 + t*x over Fraction Field of Univariate Polynomial Ring in t over Rational Field 

sage: type(E) 

<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'> 

sage: E.category() 

Category of schemes over Fraction Field of Univariate Polynomial Ring in t over Rational Field 

 

See :trac:`12517`:: 

 

sage: E = EllipticCurve([1..5]) 

sage: EllipticCurve(E.a_invariants()) 

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

 

See :trac:`11773`:: 

 

sage: E = EllipticCurve() 

Traceback (most recent call last): 

... 

TypeError: invalid input to EllipticCurve constructor 

 

""" 

def create_key_and_extra_args(self, x=None, y=None, j=None, minimal_twist=True, **kwds): 

""" 

Return a ``UniqueFactory`` key and possibly extra parameters. 

 

INPUT: 

 

See the documentation for :class:`EllipticCurveFactory`. 

 

OUTPUT: 

 

A pair ``(key, extra_args)``: 

 

- ``key`` has the form `(R, (a_1, a_2, a_3, a_4, a_6))`, 

representing a ring and the Weierstrass coefficients of an 

elliptic curve over that ring; 

 

- ``extra_args`` is a dictionary containing additional data to 

be inserted into the elliptic curve structure. 

 

EXAMPLES:: 

 

sage: EllipticCurve.create_key_and_extra_args(j=8000) 

((Rational Field, (0, -1, 0, -3, -1)), {}) 

 

When constructing a curve over `\\QQ` from a Cremona or LMFDB 

label, the invariants from the database are returned as 

``extra_args``:: 

 

sage: key, data = EllipticCurve.create_key_and_extra_args('389.a1') 

sage: key 

(Rational Field, (0, 1, 1, -2, 0)) 

sage: data['conductor'] 

389 

sage: data['cremona_label'] 

'389a1' 

sage: data['lmfdb_label'] 

'389.a1' 

sage: data['rank'] 

2 

sage: data['torsion_order'] 

1 

 

User-specified keywords are also included in ``extra_args``:: 

 

sage: key, data = EllipticCurve.create_key_and_extra_args((0, 0, 1, -23737, 960366), rank=4) 

sage: data['rank'] 

4 

 

Furthermore, keywords takes precedence over data from the 

database, which can be used to specify an alternative set of 

generators for the Mordell-Weil group:: 

 

sage: key, data = EllipticCurve.create_key_and_extra_args('5077a1', gens=[[1, -1], [-2, 3], [4, -7]]) 

sage: data['gens'] 

[[1, -1], [-2, 3], [4, -7]] 

sage: E = EllipticCurve.create_object(0, key, **data) 

sage: E.gens() 

[(-2 : 3 : 1), (1 : -1 : 1), (4 : -7 : 1)] 

 

Note that elliptic curves are equal if and only they have the 

same base ring and Weierstrass equation; the data in 

``extra_args`` do not influence comparison of elliptic curves. 

A consequence of this is that passing keyword arguments only 

works when constructing an elliptic curve the first time:: 

 

sage: E = EllipticCurve('433a1', gens=[[-1, 1], [3, 4]]) 

sage: E.gens() 

[(-1 : 1 : 1), (3 : 4 : 1)] 

sage: E = EllipticCurve('433a1', gens=[[-1, 0], [0, 1]]) 

sage: E.gens() 

[(-1 : 1 : 1), (3 : 4 : 1)] 

 

.. WARNING:: 

 

Manually specifying extra data is almost never necessary 

and is not guaranteed to have any effect, as the above 

example shows. Almost no checking is done, so specifying 

incorrect data may lead to wrong results of computations 

instead of errors or warnings. 

 

""" 

R = None 

if is_Ring(x): 

(R, x) = (x, y) 

 

if j is not None: 

if R is not None: 

try: 

j = R(j) 

except (ZeroDivisionError, ValueError, TypeError): 

raise ValueError("First parameter must be a ring containing %s" % j) 

elif x is not None: 

raise ValueError("First parameter (if present) must be a ring when j is specified") 

x = coefficients_from_j(j, minimal_twist) 

 

if is_SymbolicEquation(x): 

x = x.lhs() - x.rhs() 

 

if parent(x) is SR: 

x = x._polynomial_(rings.QQ['x', 'y']) 

 

if is_MPolynomial(x): 

if y is None: 

x = coefficients_from_Weierstrass_polynomial(x) 

else: 

# TODO: This function coefficients_from_cubic() is not defined anywhere! 

x = coefficients_from_cubic(x, y, morphism=False) 

 

if isinstance(x, string_types): 

# Interpret x as a Cremona or LMFDB label. 

from sage.databases.cremona import CremonaDatabase 

x, data = CremonaDatabase().coefficients_and_data(x) 

# User-provided keywords may override database entries. 

data.update(kwds) 

kwds = data 

 

if not isinstance(x, (list, tuple)): 

raise TypeError("invalid input to EllipticCurve constructor") 

 

if len(x) == 2: 

x = (0, 0, 0, x[0], x[1]) 

elif len(x) != 5: 

raise ValueError("sequence of coefficients must have length 2 or 5") 

 

if R is None: 

R = Sequence(x).universe() 

if R in (rings.ZZ,) + integer_types: 

R = rings.QQ 

 

return (R, tuple(R(a) for a in x)), kwds 

 

def create_object(self, version, key, **kwds): 

""" 

Create an object from a ``UniqueFactory`` key. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve.create_object(0, (GF(3), (1, 2, 0, 1, 2))) 

sage: type(E) 

<class 'sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field_with_category'> 

 

.. NOTE:: 

 

Keyword arguments are currently only passed to the 

constructor for elliptic curves over `\\QQ`; elliptic 

curves over other fields do not support them. 

 

""" 

R, x = key 

 

if R is rings.QQ: 

from .ell_rational_field import EllipticCurve_rational_field 

return EllipticCurve_rational_field(x, **kwds) 

elif is_NumberField(R): 

from .ell_number_field import EllipticCurve_number_field 

return EllipticCurve_number_field(R, x) 

elif rings.is_pAdicField(R): 

from .ell_padic_field import EllipticCurve_padic_field 

return EllipticCurve_padic_field(R, x) 

elif is_FiniteField(R) or (is_IntegerModRing(R) and R.characteristic().is_prime()): 

from .ell_finite_field import EllipticCurve_finite_field 

return EllipticCurve_finite_field(R, x) 

elif R in _Fields: 

from .ell_field import EllipticCurve_field 

return EllipticCurve_field(R, x) 

from .ell_generic import EllipticCurve_generic 

return EllipticCurve_generic(R, x) 

 

 

EllipticCurve = EllipticCurveFactory('sage.schemes.elliptic_curves.constructor.EllipticCurve') 

 

 

def EllipticCurve_from_Weierstrass_polynomial(f): 

""" 

Return the elliptic curve defined by a cubic in (long) Weierstrass 

form. 

 

INPUT: 

 

- ``f`` -- a inhomogeneous cubic polynomial in long Weierstrass 

form. 

 

OUTPUT: 

 

The elliptic curve defined by it. 

 

EXAMPLES:: 

 

sage: R.<x,y> = QQ[] 

sage: f = y^2 + 1*x*y + 3*y - (x^3 + 2*x^2 + 4*x + 6) 

sage: EllipticCurve(f) 

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

sage: EllipticCurve(f).a_invariants() 

(1, 2, 3, 4, 6) 

 

The polynomial ring may have extra variables as long as they 

do not occur in the polynomial itself:: 

 

sage: R.<x,y,z,w> = QQ[] 

sage: EllipticCurve(-y^2 + x^3 + 1) 

Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field 

sage: EllipticCurve(-x^2 + y^3 + 1) 

Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field 

sage: EllipticCurve(-w^2 + z^3 + 1) 

Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field 

 

TESTS:: 

 

sage: from sage.schemes.elliptic_curves.constructor import EllipticCurve_from_Weierstrass_polynomial 

sage: EllipticCurve_from_Weierstrass_polynomial(-w^2 + z^3 + 1) 

Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field 

""" 

return EllipticCurve(coefficients_from_Weierstrass_polynomial(f)) 

 

def coefficients_from_Weierstrass_polynomial(f): 

""" 

Return the coefficients `[a_1, a_2, a_3, a_4, a_6]` of a cubic in 

Weierstrass form. 

 

EXAMPLES:: 

 

sage: from sage.schemes.elliptic_curves.constructor import coefficients_from_Weierstrass_polynomial 

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

sage: coefficients_from_Weierstrass_polynomial(-w^2 + z^3 + 1) 

[0, 0, 0, 0, 1] 

""" 

R = f.parent() 

cubic_variables = [ x for x in R.gens() if f.degree(x) == 3 ] 

quadratic_variables = [ y for y in R.gens() if f.degree(y) == 2 ] 

try: 

x = cubic_variables[0] 

y = quadratic_variables[0] 

except IndexError: 

raise ValueError('polynomial is not in long Weierstrass form') 

 

a1 = a2 = a3 = a4 = a6 = 0 

x3 = y2 = None 

for coeff, mon in f: 

if mon == x**3: 

x3 = coeff 

elif mon == x**2: 

a2 = coeff 

elif mon == x: 

a4 = coeff 

elif mon == 1: 

a6 = coeff 

elif mon == y**2: 

y2 = -coeff 

elif mon == x*y: 

a1 = -coeff 

elif mon == y: 

a3 = -coeff 

else: 

raise ValueError('polynomial is not in long Weierstrass form') 

 

if x3 != y2: 

raise ValueError('the coefficient of x^3 and -y^2 must be the same') 

elif x3 != 1: 

a1, a2, a3, a4, a6 = a1/x3, a2/x3, a3/x3, a4/x3, a6/x3 

return [a1, a2, a3, a4, a6] 

 

 

def EllipticCurve_from_c4c6(c4, c6): 

""" 

Return an elliptic curve with given `c_4` and 

`c_6` invariants. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve_from_c4c6(17, -2005) 

sage: E 

Elliptic Curve defined by y^2 = x^3 - 17/48*x + 2005/864 over Rational Field 

sage: E.c_invariants() 

(17, -2005) 

""" 

try: 

K = c4.parent() 

except AttributeError: 

K = rings.RationalField() 

if K not in _Fields: 

K = K.fraction_field() 

return EllipticCurve([-K(c4)/K(48), -K(c6)/K(864)]) 

 

 

def EllipticCurve_from_j(j, minimal_twist=True): 

""" 

Return an elliptic curve with given `j`-invariant. 

 

INPUT: 

 

- ``j`` -- an element of some field. 

 

- ``minimal_twist`` (boolean, default True) -- If True and ``j`` is in `\QQ`, the curve returned is a 

minimal twist, i.e. has minimal conductor. If `j` is not in `\QQ` this parameter is ignored. 

 

OUTPUT: 

 

An elliptic curve with `j`-invariant `j`. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve_from_j(0); E; E.j_invariant(); E.label() 

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

0 

'27a3' 

 

sage: E = EllipticCurve_from_j(1728); E; E.j_invariant(); E.label() 

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

1728 

'32a2' 

 

sage: E = EllipticCurve_from_j(1); E; E.j_invariant() 

Elliptic Curve defined by y^2 + x*y = x^3 + 36*x + 3455 over Rational Field 

1 

 

The ``minimal_twist`` parameter (ignored except over `\QQ` and 

True by default) controls whether or not a minimal twist is 

computed:: 

 

sage: EllipticCurve_from_j(100) 

Elliptic Curve defined by y^2 = x^3 + x^2 + 3392*x + 307888 over Rational Field 

sage: _.conductor() 

33129800 

sage: EllipticCurve_from_j(100, minimal_twist=False) 

Elliptic Curve defined by y^2 = x^3 + 488400*x - 530076800 over Rational Field 

sage: _.conductor() 

298168200 

 

Since computing the minimal twist requires factoring both `j` and 

`j-1728` the following example would take a long time without 

setting ``minimal_twist`` to False:: 

 

sage: E = EllipticCurve_from_j(2^256+1,minimal_twist=False) 

sage: E.j_invariant() == 2^256+1 

True 

 

""" 

return EllipticCurve(coefficients_from_j(j, minimal_twist)) 

 

def coefficients_from_j(j, minimal_twist=True): 

""" 

Return Weierstrass coefficients `(a_1, a_2, a_3, a_4, a_6)` for an 

elliptic curve with given `j`-invariant. 

 

INPUT: 

 

See :func:`EllipticCurve_from_j`. 

 

EXAMPLES:: 

 

sage: from sage.schemes.elliptic_curves.constructor import coefficients_from_j 

sage: coefficients_from_j(0) 

[0, 0, 1, 0, 0] 

sage: coefficients_from_j(1728) 

[0, 0, 0, -1, 0] 

sage: coefficients_from_j(1) 

[1, 0, 0, 36, 3455] 

 

The ``minimal_twist`` parameter (ignored except over `\\QQ` and 

True by default) controls whether or not a minimal twist is 

computed:: 

 

sage: coefficients_from_j(100) 

[0, 1, 0, 3392, 307888] 

sage: coefficients_from_j(100, minimal_twist=False) 

[0, 0, 0, 488400, -530076800] 

""" 

try: 

K = j.parent() 

except AttributeError: 

K = rings.RationalField() 

if K not in _Fields: 

K = K.fraction_field() 

 

char=K.characteristic() 

if char==2: 

if j == 0: 

return Sequence([0, 0, 1, 0, 0], universe=K) 

else: 

return Sequence([1, 0, 0, 0, 1/j], universe=K) 

if char == 3: 

if j==0: 

return Sequence([0, 0, 0, 1, 0], universe=K) 

else: 

return Sequence([0, j, 0, 0, -j**2], universe=K) 

 

if K is rings.RationalField(): 

# we construct the minimal twist, i.e. the curve with minimal 

# conductor with this j_invariant: 

if j == 0: 

return Sequence([0, 0, 1, 0, 0], universe=K) # 27a3 

if j == 1728: 

return Sequence([0, 0, 0, -1, 0], universe=K) # 32a2 

 

if not minimal_twist: 

k=j-1728 

return Sequence([0, 0, 0, -3*j*k, -2*j*k**2], universe=K) 

 

n = j.numerator() 

m = n-1728*j.denominator() 

a4 = -3*n*m 

a6 = -2*n*m**2 

 

# Now E=[0,0,0,a4,a6] has j-invariant j=n/d 

from sage.sets.set import Set 

for p in Set(n.prime_divisors()+m.prime_divisors()): 

e = min(a4.valuation(p)//2,a6.valuation(p)//3) 

if e>0: 

p = p**e 

a4 /= p**2 

a6 /= p**3 

 

# Now E=[0,0,0,a4,a6] is minimal at all p != 2,3 

tw = [-1,2,-2,3,-3,6,-6] 

E1 = EllipticCurve([0,0,0,a4,a6]) 

Elist = [E1] + [E1.quadratic_twist(t) for t in tw] 

Elist.sort(key=lambda E: E.conductor()) 

return Sequence(Elist[0].ainvs()) 

 

# defaults for all other fields: 

if j == 0: 

return Sequence([0, 0, 0, 0, 1], universe=K) 

if j == 1728: 

return Sequence([0, 0, 0, 1, 0], universe=K) 

k=j-1728 

return Sequence([0, 0, 0, -3*j*k, -2*j*k**2], universe=K) 

 

 

def EllipticCurve_from_cubic(F, P=None, morphism=True): 

r""" 

Construct an elliptic curve from a ternary cubic with a rational point. 

 

If you just want the Weierstrass form and are not interested in 

the morphism then it is easier to use the function 

:func:`~sage.schemes.elliptic_curves.jacobian.Jacobian` 

instead. If there is a rational point on the given cubic, this 

function will construct the same elliptic curve but you do not have to 

supply the point ``P``. 

 

INPUT: 

 

- ``F`` -- a homogeneous cubic in three variables with rational 

coefficients, as a polynomial ring element, defining a smooth 

plane cubic curve `C`. 

 

- ``P`` -- a 3-tuple `(x,y,z)` defining a projective point on `C`, 

or ``None``. If ``None`` then a rational flex will be used as a 

base point if one exists, otherwise an error will be raised. 

 

- ``morphism`` -- boolean (default: ``True``). If ``True`` 

returns a birational isomorphism from `C` to a Weierstrass 

elliptic curve `E`, otherwise just returns `E`. 

 

OUTPUT: 

 

Either (when ``morphism``=``False``) an elliptic curve `E` in long 

Weierstrass form isomorphic to the plane cubic curve `C` defined 

by the equation `F=0`. 

 

Or (when ``morphism=True``), a birational isomorphism from `C` to 

the elliptic curve `E`. If the given point is a flex, this is a 

linear isomorphism. 

 

.. NOTE:: 

 

The function 

:func:`~sage.schemes.elliptic_curves.jacobian.Jacobian` may be 

used instead. It constructs the same elliptic curve (which is in 

all cases the Jacobian of `(F=0)`) and needs no base point to be 

provided, but also returns no isomorphism since in general there 

is none: the plane cubic is only isomorphic to its Jacobian when 

it has a rational point. 

 

.. NOTE:: 

 

When ``morphism=True``, a birational isomorphism between the 

curve `F=0` and the Weierstrass curve is returned. If the point 

happens to be a flex, then this is a linear isomorphism. The 

morphism does not necessarily take the given point `P` to the 

point at infinity on `E`, since we always use a rational flex 

on `C` as base-point when one exists. 

 

EXAMPLES: 

 

First we find that the Fermat cubic is isomorphic to the curve 

with Cremona label 27a1:: 

 

sage: R.<x,y,z> = QQ[] 

sage: cubic = x^3+y^3+z^3 

sage: P = [1,-1,0] 

sage: E = EllipticCurve_from_cubic(cubic, P, morphism=False); E 

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

sage: E.cremona_label() 

'27a1' 

sage: EllipticCurve_from_cubic(cubic, [0,1,-1], morphism=False).cremona_label() 

'27a1' 

sage: EllipticCurve_from_cubic(cubic, [1,0,-1], morphism=False).cremona_label() 

'27a1' 

 

Next we find the minimal model and conductor of the Jacobian of the 

Selmer curve:: 

 

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

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

sage: P = [1,-1,0] 

sage: E = EllipticCurve_from_cubic(cubic, P, morphism=False); E 

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

sage: E.minimal_model() 

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

sage: E.conductor() 

24300 

 

We can also get the birational isomorphism to and from the 

Weierstrass form. We start with an example where ``P`` is a flex 

and the equivalence is a linear isomorphism:: 

 

sage: f = EllipticCurve_from_cubic(cubic, P, morphism=True) 

sage: f 

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 - 540*y = x^3 - 97200 over Rational Field 

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

(-c : 3*a : 1/180*a + 1/180*b) 

 

sage: finv = f.inverse(); finv 

Scheme morphism: 

From: Elliptic Curve defined by y^2 - 540*y = x^3 - 97200 over Rational Field 

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

Defn: Defined on coordinates by sending (x : y : z) to 

(1/3*y : -1/3*y + 180*z : -x) 

 

Scheme morphism: 

From: Elliptic Curve defined by y^2 + 2*x*y + 20*y = x^3 - x^2 - 20*x - 400/3 

over Rational Field 

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

a^3 + b^3 + 60*c^3 

Defn: Defined on coordinates by sending (x : y : z) to 

(x + y + 20*z : -x - y : -x) 

 

We verify that `f` maps the chosen point `P=(1,-1,0)` on the cubic 

to the origin of the elliptic curve:: 

 

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

(0 : 1 : 0) 

sage: finv([0,1,0]) 

(-1 : 1 : 0) 

 

To verify the output, we plug in the polynomials to check that 

this indeed transforms the cubic into Weierstrass form:: 

 

sage: cubic(finv.defining_polynomials()) * finv.post_rescaling() 

-x^3 + y^2*z - 540*y*z^2 + 97200*z^3 

 

sage: E.defining_polynomial()(f.defining_polynomials()) * f.post_rescaling() 

a^3 + b^3 + 60*c^3 

 

If the given point is not a flex and the cubic has no rational 

flexes, then the cubic can not be transformed to a Weierstrass 

equation by a linear transformation. The general birational 

transformation is still a birational isomorphism, but is 

quadratic:: 

 

sage: R.<x,y,z> = QQ[] 

sage: cubic = x^2*y + 4*x*y^2 + x^2*z + 8*x*y*z + 4*y^2*z + 9*x*z^2 + 9*y*z^2 

sage: f = EllipticCurve_from_cubic(cubic, [1,-1,1], morphism=True); f 

Scheme morphism: 

From: Projective Plane Curve over Rational Field defined by x^2*y + 4*x*y^2 + x^2*z + 8*x*y*z + 4*y^2*z + 9*x*z^2 + 9*y*z^2 

To: Elliptic Curve defined by y^2 + 7560/19*x*y + 552960000000/2352637*y = x^3 - 3445200/133*x^2 over Rational Field 

Defn: Defined on coordinates by sending (x : y : z) to 

(2527/17280*x^2 + 133/2160*x*y + 133/108000*y^2 + 133/2880*x*z + 931/18000*y*z - 3857/48000*z^2 : -6859/288*x^2 + 323/36*x*y + 359/1800*y^2 + 551/48*x*z + 2813/300*y*z + 24389/800*z^2 : -2352637/99532800000*x^2 - 2352637/124416000000*x*y - 2352637/622080000000*y^2 + 2352637/82944000000*x*z + 2352637/207360000000*y*z - 2352637/276480000000*z^2) 

 

Note that the morphism returned cannot be evaluated directly at 

the given point ``P=(1:-1:1)`` since the polynomials defining it 

all vanish there:: 

 

sage: f([1,-1,1]) 

Traceback (most recent call last): 

... 

ValueError: [0, 0, 0] does not define a valid point since all entries are 0 

 

Using the group law on the codomain elliptic curve, which has rank 

1 and full 2-torsion, and the inverse morphism, we can find many 

points on the cubic. First we find the preimages of multiples of 

the generator:: 

 

sage: E = f.codomain() 

sage: E.label() 

'720e2' 

sage: E.rank() 

1 

sage: R = E.gens()[0]; R 

(-17280000/2527 : 9331200000/6859 : 1) 

sage: finv = f.inverse() 

sage: [finv(k*R) for k in range(1,10)] 

[(-4 : 1 : 0), 

(-1 : 4 : 1), 

(-20 : -55/76 : 1), 

(319/399 : -11339/7539 : 1), 

(159919/14360 : -4078139/1327840 : 1), 

(-27809119/63578639 : 1856146436/3425378659 : 1), 

(-510646582340/56909753439 : 424000923715/30153806197284 : 1), 

(-56686114363679/4050436059492161 : -2433034816977728281/1072927821085503881 : 1), 

(650589589099815846721/72056273157352822480 : -347376189546061993109881/194127383495944026752320 : 1)] 

 

The elliptic curve also has torsion, which we can map back:: 

 

sage: E.torsion_points() 

[(-144000000/17689 : 3533760000000/2352637 : 1), 

(-92160000/17689 : 2162073600000/2352637 : 1), 

(-5760000/17689 : -124070400000/2352637 : 1), 

(0 : 1 : 0)] 

sage: [finv(Q) for Q in E.torsion_points() if Q] 

[(9 : -9/4 : 1), (-9 : 0 : 1), (0 : 1 : 0)] 

 

 

In this example, the given point ``P`` is not a flex but the cubic 

does have a rational flex, ``(-4:0:1)``. We return a linear 

isomorphism which maps this flex to the point at infinity on the 

Weierstrass model:: 

 

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

sage: cubic = a^3+7*b^3+64*c^3 

sage: P = [2,2,-1] 

sage: f = EllipticCurve_from_cubic(cubic, P, morphism=True) 

sage: E = f.codomain(); E 

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

sage: E.minimal_model() 

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

 

sage: f 

Scheme morphism: 

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

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

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

(b : -48*a : -1/5376*a - 1/1344*c) 

 

sage: finv = f.inverse(); finv 

Scheme morphism: 

From: Elliptic Curve defined by y^2 - 258048*y = x^3 - 22196256768 over Rational Field 

To: Projective Plane Curve over Rational Field defined by a^3 + 7*b^3 + 64*c^3 

Defn: Defined on coordinates by sending (x : y : z) to 

(-1/48*y : x : 1/192*y - 1344*z) 

 

sage: cubic(finv.defining_polynomials()) * finv.post_rescaling() 

-x^3 + y^2*z - 258048*y*z^2 + 22196256768*z^3 

 

sage: E.defining_polynomial()(f.defining_polynomials()) * f.post_rescaling() 

a^3 + 7*b^3 + 64*c^3 

 

sage: f(P) 

(5376 : -258048 : 1) 

sage: f([-4,0,1]) 

(0 : 1 : 0) 

 

It is possible to not provide a base point ``P`` provided that the 

cubic has a rational flex. In this case the flexes will be found 

and one will be used as a base point:: 

 

sage: R.<x,y,z> = QQ[] 

sage: cubic = x^3+y^3+z^3 

sage: f = EllipticCurve_from_cubic(cubic, morphism=True) 

sage: f 

Scheme morphism: 

From: Projective Plane Curve over Rational Field defined by x^3 + y^3 + z^3 

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

Defn: Defined on coordinates by sending (x : y : z) to 

(y : -3*x : -1/3*x - 1/3*z) 

 

An error will be raised if no point is given and there are no rational flexes:: 

 

sage: R.<x,y,z> = QQ[] 

sage: cubic = 3*x^3+4*y^3+5*z^3 

sage: EllipticCurve_from_cubic(cubic) 

Traceback (most recent call last): 

... 

ValueError: A point must be given when the cubic has no rational flexes 

 

An example over a finite field, using a flex:: 

 

sage: K = GF(17) 

sage: R.<x,y,z> = K[] 

sage: cubic = 2*x^3+3*y^3+4*z^3 

sage: EllipticCurve_from_cubic(cubic,[0,3,1]) 

Scheme morphism: 

From: Projective Plane Curve over Finite Field of size 17 defined by 2*x^3 + 3*y^3 + 4*z^3 

To: Elliptic Curve defined by y^2 + 16*y = x^3 + 11 over Finite Field of size 17 

Defn: Defined on coordinates by sending (x : y : z) to 

(-x : 4*y : 4*y + 5*z) 

 

An example in characteristic 3:: 

 

sage: K = GF(3) 

sage: R.<x,y,z> = K[] 

sage: cubic = x^3+y^3+z^3+x*y*z 

sage: EllipticCurve_from_cubic(cubic,[0,1,-1]) 

Scheme morphism: 

From: Projective Plane Curve over Finite Field of size 3 defined by x^3 + y^3 + x*y*z + z^3 

To: Elliptic Curve defined by y^2 + x*y = x^3 + 1 over Finite Field of size 3 

Defn: Defined on coordinates by sending (x : y : z) to 

(y + z : -y : x) 

 

An example over a number field, using a non-flex and where there are no rational flexes:: 

 

sage: K.<a> = QuadraticField(-3) 

sage: R.<x,y,z> = K[] 

sage: cubic = 2*x^3+3*y^3+5*z^3 

sage: EllipticCurve_from_cubic(cubic,[1,1,-1]) 

Scheme morphism: 

From: Projective Plane Curve over Number Field in a with defining polynomial x^2 + 3 defined by 2*x^3 + 3*y^3 + 5*z^3 

To: Elliptic Curve defined by y^2 + 1754460/2053*x*y + 5226454388736000/8653002877*y = x^3 + (-652253285700/4214809)*x^2 over Number Field in a with defining polynomial x^2 + 3 

Defn: Defined on coordinates by sending (x : y : z) to 

(-16424/127575*x^2 - 231989/680400*x*y - 14371/64800*y^2 - 26689/81648*x*z - 10265/27216*y*z - 2053/163296*z^2 : 24496/315*x^2 + 119243/840*x*y + 4837/80*y^2 + 67259/504*x*z + 25507/168*y*z + 5135/1008*z^2 : 8653002877/2099914709760000*x^2 + 8653002877/699971569920000*x*y + 8653002877/933295426560000*y^2 + 8653002877/419982941952000*x*z + 8653002877/279988627968000*y*z + 8653002877/335986353561600*z^2) 

 

An example over a function field, using a non-flex:: 

 

sage: K.<t> = FunctionField(QQ) 

sage: R.<x,y,z> = K[] 

sage: cubic = x^3+t*y^3+(1+t)*z^3 

sage: EllipticCurve_from_cubic(cubic,[1,1,-1], morphism=False) 

Elliptic Curve defined by y^2 + ((-236196*t^6-708588*t^5-1180980*t^4-1180980*t^3-708588*t^2-236196*t)/(-1458*t^6-17496*t^5+4374*t^4+29160*t^3+4374*t^2-17496*t-1458))*x*y + ((-459165024*t^14-5969145312*t^13-34207794288*t^12-113872925952*t^11-244304490582*t^10-354331909458*t^9-354331909458*t^8-244304490582*t^7-113872925952*t^6-34207794288*t^5-5969145312*t^4-459165024*t^3)/(-1458*t^14-58320*t^13-841266*t^12-5137992*t^11-11773350*t^10-7709904*t^9+12627738*t^8+25789104*t^7+12627738*t^6-7709904*t^5-11773350*t^4-5137992*t^3-841266*t^2-58320*t-1458))*y = x^3 + ((-118098*t^12-708588*t^11+944784*t^10+11219310*t^9+27871128*t^8+36374184*t^7+27871128*t^6+11219310*t^5+944784*t^4-708588*t^3-118098*t^2)/(-54*t^12-1296*t^11-7452*t^10+6048*t^9+25758*t^8-3888*t^7-38232*t^6-3888*t^5+25758*t^4+6048*t^3-7452*t^2-1296*t-54))*x^2 over Rational function field in t over Rational Field 

 

 

TESTS: 

 

Here is a test for :trac:`21092`:: 

 

sage: R.<x,y,z> = QQ[] 

sage: cubic = -3*x^2*y + 3*x*y^2 + 4*x^2*z + 4*y^2*z - 3*x*z^2 + 3*y*z^2 - 8*z^3 

sage: EllipticCurve_from_cubic(cubic, (-4/5, 4/5, 3/5) ) 

Scheme morphism: 

From: Projective Plane Curve over Rational Field defined by -3*x^2*y + 3*x*y^2 + 4*x^2*z + 4*y^2*z - 3*x*z^2 + 3*y*z^2 - 8*z^3 

To: Elliptic Curve defined by y^2 + 24*x*y + 3024*y = x^3 + 495*x^2 + 36288*x over Rational Field 

Defn: Defined on coordinates by sending (x : y : z) to 

(-1/3*z : 3*x : -1/1008*x + 1/1008*y + 1/378*z) 

""" 

from sage.schemes.curves.constructor import Curve 

from sage.matrix.all import Matrix 

from sage.modules.all import vector 

from sage.schemes.elliptic_curves.weierstrass_transform import \ 

WeierstrassTransformationWithInverse 

 

# check the input 

R = F.parent() 

K = R.base_ring() 

if not is_MPolynomialRing(R): 

raise TypeError('equation must be a polynomial') 

if R.ngens() != 3 or F.nvariables() != 3: 

raise TypeError('equation must be a polynomial in three variables') 

if not F.is_homogeneous(): 

raise TypeError('equation must be a homogeneous polynomial') 

 

C = Curve(F) 

if P: 

try: 

CP = C(P) 

except (TypeError, ValueError): 

raise TypeError('{} does not define a point on a projective curve over {} defined by {}'.format(P,K,F)) 

 

x, y, z = R.gens() 

 

# Test whether P is a flex; if not test whether there are any rational flexes: 

 

hessian = Matrix([[F.derivative(v1, v2) for v1 in R.gens()] for v2 in R.gens()]).det() 

if P and hessian(P)==0: 

flex_point = P 

else: 

flexes = C.intersection(Curve(hessian)).rational_points() 

if flexes: 

flex_point = list(flexes[0]) 

if not P: 

P = flex_point 

CP = C(P) 

else: 

flex_point = None 

 

if flex_point is not None: # first case: base point is a flex 

P = flex_point 

L = tangent_at_smooth_point(C,P) 

dx, dy, dz = [L.coefficient(v) for v in R.gens()] 

 

# find an invertible matrix M such that (0,1,0)M=P and 

# ML'=(0,0,1)' where L=[dx,dy,dx]. Then the linea transform 

# by M takes P to [0,1,0] and L to Z=0: 

 

if P[0]: 

Q1 = [0,-dz,dy] 

Q2 = [0,1,0] if dy else [0,0,1] 

elif P[1]: 

Q1 = [dz,0,-dx] 

Q2 = [1,0,0] if dx else [0,0,1] 

else: 

Q1 = [-dy,dx,0] 

Q2 = [1,0,0] if dx else [0,1,0] 

 

M = Matrix(K,[Q1,P,Q2]) 

# assert M.is_invertible() 

# assert list(vector([0,1,0])*M) == P 

# assert list(M*vector([dx,dy,dz]))[:2] == [0,0] 

 

M = M.transpose() 

F2 = R(M.act_on_polynomial(F)) 

 

# scale and dehomogenise 

a = K(F2.coefficient(x**3)) 

b = K(F2.coefficient(y*y*z)) 

 

F3 = F2([-x, y/b, z*a*b]) / a 

# assert F3.coefficient(x**3) == -1 

# assert F3.coefficient(y*y*z) == 1 

E = EllipticCurve(F3([x,y,1])) 

if not morphism: 

return E 

 

# Construct the (linear) morphism 

M = M * Matrix(K,[[-1,0,0],[0,1/b,0],[0,0,a*b]]) 

inv_defining_poly = [ M[i,0]*x + M[i,1]*y + M[i,2]*z for i in range(3) ] 

inv_post = 1/a 

M = M.inverse() 

fwd_defining_poly = [ M[i,0]*x + M[i,1]*y + M[i,2]*z for i in range(3) ] 

fwd_post = a 

 

else: # Second case: no flexes 

if not P: 

raise ValueError('A point must be given when the cubic has no rational flexes') 

L = tangent_at_smooth_point(C,P) 

Qlist = [Q for Q in C.intersection(Curve(L)).rational_points() if C(Q)!=CP] 

# assert Qlist 

P2 = C(Qlist[0]) 

L2 = tangent_at_smooth_point(C,P2) 

Qlist = [Q for Q in C.intersection(Curve(L2)).rational_points() if C(Q)!=P2] 

# assert Qlist 

P3 = C(Qlist[0]) 

 

# NB This construction of P3 relies on P2 not being a flex. 

# If we want to use a non-flex as P when there are rational 

# flexes this would be a problem. However, the only condition 

# which P3 must satisfy is that it is on the tangent at P2, it 

# need not lie on the cubic. 

 

# send P, P2, P3 to (1:0:0), (0:1:0), (0:0:1) respectively 

M = Matrix(K, [P, list(P2), list(P3)]).transpose() 

F2 = M.act_on_polynomial(F) 

xyzM = [ M[i,0]*x + M[i,1]*y + M[i,2]*z for i in range(3) ] 

# assert F(xyzM)==F2 

 

# substitute x = U^2, y = V*W, z = U*W, and rename (x,y,z)=(U,V,W) 

T1 = [x*x,y*z,x*z] 

S1 = x**2*z 

F3 = F2(T1) // S1 

xyzC = [ t(T1) for t in xyzM ] 

# assert F3 == F(xyzC) // S1 

 

# scale and dehomogenise 

a = K(F3.coefficient(x**3)) 

b = K(F3.coefficient(y*y*z)) 

ab = a*b 

 

T2 = [-x, y/b, ab*z] 

F4 = F3(T2) / a 

# assert F4.coefficient(x**3) == -1 

# assert F4.coefficient(y*y*z) == 1 

xyzW = [ t(T2) for t in xyzC ] 

S2 = a*S1(T2) 

# assert F4 == F(xyzW) // S2 

 

E = EllipticCurve(F4([x,y,1])) 

if not morphism: 

return E 

 

inv_defining_poly = xyzW 

inv_post = 1/S2 

# assert F4==F(inv_defining_poly)*inv_post 

MI = M.inverse() 

xyzI = [ (MI[i,0]*x + MI[i,1]*y + MI[i,2]*z) for i in range(3) ] 

T1I = [x*z,x*y,z*z] # inverse of T1 

xyzIC = [ t(xyzI) for t in T1I ] 

T2I = [-x, b*y, z/ab] # inverse of T2 

xyzIW = [ t(xyzIC) for t in T2I ] 

fwd_defining_poly = xyzIW 

fwd_post = a/(x*z*z)(xyzI) 

# assert F4(fwd_defining_poly)*fwd_post == F 

 

# Construct the morphism 

 

return WeierstrassTransformationWithInverse( 

C, E, fwd_defining_poly, fwd_post, inv_defining_poly, inv_post) 

 

 

def tangent_at_smooth_point(C,P): 

"""Return the tangent at the smooth point `P` of projective curve `C`. 

 

INPUT: 

 

- ``C`` -- a projective plane curve. 

 

- ``P`` -- a 3-tuple `(x,y,z)` defining a projective point on `C`. 

 

OUTPUT: 

 

The linear form defining the tangent at `P` to `C`. 

 

EXAMPLES:: 

 

sage: R.<x,y,z> = QQ[] 

sage: from sage.schemes.elliptic_curves.constructor import tangent_at_smooth_point 

sage: C = Curve(x^3+y^3+60*z^3) 

sage: tangent_at_smooth_point(C, [1,-1,0]) 

x + y 

 

sage: K.<t> = FunctionField(QQ) 

sage: R.<x,y,z> = K[] 

sage: C = Curve(x^3+2*y^3+3*z^3) 

sage: from sage.schemes.elliptic_curves.constructor import tangent_at_smooth_point 

sage: tangent_at_smooth_point(C,[1,1,-1]) 

3*x + 6*y + 9*z 

""" 

# Over function fields such as QQ(t) an error is raised with the 

# default (factor=True). Note that factor=False returns the 

# product of the tangents in case of a multiple point, while here 

# `P` is assumed smooth so factorization is unnecessary, but over 

# QQ (for example) including the factorization gives better 

# results, for example returning x+y instead of 3x+3y in the 

# doctest. 

try: 

return C.tangents(P)[0] 

except NotImplementedError: 

return C.tangents(P,factor=False)[0] 

 

def chord_and_tangent(F, P): 

"""Return the third point of intersection of a cubic with the tangent at one point. 

 

INPUT: 

 

- ``F`` -- a homogeneous cubic in three variables with rational 

coefficients, as a polynomial ring element, defining a smooth 

plane cubic curve. 

 

- ``P`` -- a 3-tuple `(x,y,z)` defining a projective point on the 

curve `F=0`. 

 

OUTPUT: 

 

A point ``Q`` such that ``F(Q)=0``, namely the third point of 

intersection of the tangent at ``P`` with the curve ``F=0``, so 

``Q=P`` if and only if ``P`` is a flex. 

 

EXAMPLES:: 

 

sage: R.<x,y,z> = QQ[] 

sage: from sage.schemes.elliptic_curves.constructor import chord_and_tangent 

sage: F = x^3+y^3+60*z^3 

sage: chord_and_tangent(F, [1,-1,0]) 

(-1 : 1 : 0) 

 

sage: F = x^3+7*y^3+64*z^3 

sage: p0 = [2,2,-1] 

sage: p1 = chord_and_tangent(F, p0); p1 

(5 : -3 : 1) 

sage: p2 = chord_and_tangent(F, p1); p2 

(-1265/314 : 183/314 : 1) 

 

TESTS:: 

 

sage: F(list(p2)) 

0 

sage: list(map(type, p2)) 

[<... 'sage.rings.rational.Rational'>, 

<... 'sage.rings.rational.Rational'>, 

<... 'sage.rings.rational.Rational'>] 

 

See :trac:`16068`:: 

 

sage: F = x**3 - 4*x**2*y - 65*x*y**2 + 3*x*y*z - 76*y*z**2 

sage: chord_and_tangent(F, [0, 1, 0]) 

(0 : 0 : 1) 

 

""" 

from sage.schemes.curves.constructor import Curve 

# check the input 

R = F.parent() 

if not is_MPolynomialRing(R): 

raise TypeError('equation must be a polynomial') 

if R.ngens() != 3: 

raise TypeError('{} is not a polynomial in three variables'.format(F)) 

if not F.is_homogeneous(): 

raise TypeError('{} is not a homogeneous polynomial'.format(F)) 

x, y, z = R.gens() 

if len(P) != 3: 

raise TypeError('{} is not a projective point'.format(P)) 

K = R.base_ring() 

try: 

C = Curve(F) 

P = C(P) 

except (TypeError, ValueError): 

raise TypeError('{} does not define a point on a projective curve over {} defined by {}'.format(P,K,F)) 

 

L = Curve(tangent_at_smooth_point(C,P)) 

Qlist = [Q for Q in C.intersection(L).rational_points() if Q!=P] 

if Qlist: 

return Qlist[0] 

return P 

 

 

def projective_point(p): 

""" 

Return equivalent point with denominators removed 

 

INPUT: 

 

- ``P``, ``Q`` -- list/tuple of projective coordinates. 

 

OUTPUT: 

 

List of projective coordinates. 

 

EXAMPLES:: 

 

sage: from sage.schemes.elliptic_curves.constructor import projective_point 

sage: projective_point([4/5, 6/5, 8/5]) 

[2, 3, 4] 

sage: F = GF(11) 

sage: projective_point([F(4), F(8), F(2)]) 

[4, 8, 2] 

""" 

from sage.rings.integer import GCD_list 

from sage.arith.functions import LCM_list 

try: 

p_gcd = GCD_list([x.numerator() for x in p]) 

p_lcm = LCM_list(x.denominator() for x in p) 

except AttributeError: 

return p 

scale = p_lcm / p_gcd 

return [scale * x for x in p] 

 

 

def are_projectively_equivalent(P, Q, base_ring): 

""" 

Test whether ``P`` and ``Q`` are projectively equivalent. 

 

INPUT: 

 

- ``P``, ``Q`` -- list/tuple of projective coordinates. 

 

- ``base_ring`` -- the base ring. 

 

OUTPUT: 

 

Boolean. 

 

EXAMPLES:: 

 

sage: from sage.schemes.elliptic_curves.constructor import are_projectively_equivalent 

sage: are_projectively_equivalent([0,1,2,3], [0,1,2,2], base_ring=QQ) 

False 

sage: are_projectively_equivalent([0,1,2,3], [0,2,4,6], base_ring=QQ) 

True 

""" 

from sage.matrix.constructor import matrix 

return matrix(base_ring, [P, Q]).rank() < 2 

 

 

def EllipticCurve_from_plane_curve(C, P): 

""" 

Deprecated way to construct an elliptic curve. 

 

Use :meth:`~sage.schemes.elliptic_curves.jacobian.Jacobian` instead. 

 

EXAMPLES:: 

 

sage: R.<x,y,z> = QQ[] 

sage: C = Curve(x^3+y^3+z^3) 

sage: P = C(1,-1,0) 

sage: E = EllipticCurve_from_plane_curve(C,P); E # long time (3s on sage.math, 2013) 

doctest:...: DeprecationWarning: use Jacobian(C) instead 

See http://trac.sagemath.org/3416 for details. 

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

""" 

from sage.misc.superseded import deprecation 

deprecation(3416, 'use Jacobian(C) instead') 

# Note: this function never used the rational point 

from sage.schemes.elliptic_curves.jacobian import Jacobian 

return Jacobian(C) 

 

 

def EllipticCurves_with_good_reduction_outside_S(S=[], proof=None, verbose=False): 

r""" 

Returns a sorted list of all elliptic curves defined over `Q` 

with good reduction outside the set `S` of primes. 

 

INPUT: 

 

- ``S`` - list of primes (default: empty list). 

 

- ``proof`` - True/False (default True): the MW basis for 

auxiliary curves will be computed with this proof flag. 

 

- ``verbose`` - True/False (default False): if True, some details 

of the computation will be output. 

 

.. NOTE:: 

 

Proof flag: The algorithm used requires determining all 

S-integral points on several auxiliary curves, which in turn 

requires the computation of their generators. This is not 

always possible (even in theory) using current knowledge. 

 

The value of this flag is passed to the function which 

computes generators of various auxiliary elliptic curves, in 

order to find their S-integral points. Set to False if the 

default (True) causes warning messages, but note that you can 

then not rely on the set of curves returned being 

complete. 

 

EXAMPLES:: 

 

sage: EllipticCurves_with_good_reduction_outside_S([]) 

[] 

sage: elist = EllipticCurves_with_good_reduction_outside_S([2]) 

sage: elist 

[Elliptic Curve defined by y^2 = x^3 + 4*x over Rational Field, 

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

... 

Elliptic Curve defined by y^2 = x^3 - x^2 - 13*x + 21 over Rational Field] 

sage: len(elist) 

24 

sage: ', '.join([e.label() for e in elist]) 

'32a1, 32a2, 32a3, 32a4, 64a1, 64a2, 64a3, 64a4, 128a1, 128a2, 128b1, 128b2, 128c1, 128c2, 128d1, 128d2, 256a1, 256a2, 256b1, 256b2, 256c1, 256c2, 256d1, 256d2' 

 

Without ``Proof=False``, this example gives two warnings:: 

 

sage: elist = EllipticCurves_with_good_reduction_outside_S([11],proof=False) # long time (14s on sage.math, 2011) 

sage: len(elist) # long time 

12 

sage: ', '.join([e.label() for e in elist]) # long time 

'11a1, 11a2, 11a3, 121a1, 121a2, 121b1, 121b2, 121c1, 121c2, 121d1, 121d2, 121d3' 

 

sage: elist = EllipticCurves_with_good_reduction_outside_S([2,3]) # long time (26s on sage.math, 2011) 

sage: len(elist) # long time 

752 

sage: max([e.conductor() for e in elist]) # long time 

62208 

sage: [N.factor() for N in Set([e.conductor() for e in elist])] # long time 

[2^7, 

2^8, 

2^3 * 3^4, 

2^2 * 3^3, 

2^8 * 3^4, 

2^4 * 3^4, 

2^3 * 3, 

2^7 * 3, 

2^3 * 3^5, 

3^3, 

2^8 * 3, 

2^5 * 3^4, 

2^4 * 3, 

2 * 3^4, 

2^2 * 3^2, 

2^6 * 3^4, 

2^6, 

2^7 * 3^2, 

2^4 * 3^5, 

2^4 * 3^3, 

2 * 3^3, 

2^6 * 3^3, 

2^6 * 3, 

2^5, 

2^2 * 3^4, 

2^3 * 3^2, 

2^5 * 3, 

2^7 * 3^4, 

2^2 * 3^5, 

2^8 * 3^2, 

2^5 * 3^2, 

2^7 * 3^5, 

2^8 * 3^5, 

2^3 * 3^3, 

2^8 * 3^3, 

2^5 * 3^5, 

2^4 * 3^2, 

2 * 3^5, 

2^5 * 3^3, 

2^6 * 3^5, 

2^7 * 3^3, 

3^5, 

2^6 * 3^2] 

""" 

from .ell_egros import (egros_from_jlist, egros_get_j) 

return egros_from_jlist(egros_get_j(S, proof=proof, verbose=verbose), S)