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

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

r""" 

Galois representations for elliptic curves over number fields. 

 

This file contains the code to compute for which primes the Galois 

representation attached to an elliptic curve (over an arbitrary 

number field) is surjective. The functions in this file are called by 

the ``is_surjective`` and ``non_surjective`` methods of an elliptic curve 

over a number field. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 - 29, 'a'); a = K.gen() 

sage: E = EllipticCurve([1, 0, ((5 + a)/2)**2, 0, 0]) 

sage: rho = E.galois_representation() 

sage: rho.is_surjective(29) # Cyclotomic character not surjective. 

False 

sage: rho.is_surjective(31) # See Section 5.10 of [Serre72]. 

True 

sage: rho.non_surjective() # long time (4s on sage.math, 2014) 

[3, 5, 29] 

 

sage: E = EllipticCurve_from_j(1728).change_ring(K) # CM 

sage: E.galois_representation().non_surjective() # long time (2s on sage.math, 2014) 

[0] 

 

AUTHORS: 

 

- Eric Larson (2012-05-28): initial version. 

- Eric Larson (2014-08-13): added isogeny_bound function. 

- John Cremona (2016, 2017): various efficiency improvements to _semistable_reducible_primes 

 

REFERENCES: 

 

.. [Serre72] Serre. Propriétés galoisiennes des points d'ordre fini des 

courbes elliptiques. Inventiones mathematicae, 1972. 

 

.. [Sutherland12] Sutherland. A local-global principle for rational 

isogenies of prime degree. Journal de Théorie des Nombres de Bordeaux, 

2012. 

 

""" 

 

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

# Copyright (C) 2012 Eric Larson <elarson3@gmail.com> 

# 

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

# as published by the Free Software Foundation; either version 2 of 

# the License, or (at your option) any later version. 

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

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

from six.moves import range 

 

from sage.structure.sage_object import SageObject 

from sage.rings.number_field.number_field import NumberField, QuadraticField 

from sage.schemes.elliptic_curves.cm import cm_j_invariants 

from sage.modules.free_module import VectorSpace 

from sage.rings.finite_rings.finite_field_constructor import GF 

from sage.misc.functional import cyclotomic_polynomial 

from sage.arith.all import legendre_symbol, primes 

from sage.sets.set import Set 

from sage.rings.all import PolynomialRing, Integer, ZZ, QQ, Infinity 

 

class GaloisRepresentation(SageObject): 

r""" 

The compatible family of Galois representation 

attached to an elliptic curve over a number field. 

 

Given an elliptic curve `E` over a number field `K` 

and a rational prime number `p`, the `p^n`-torsion 

`E[p^n]` points of `E` is a representation of the 

absolute Galois group `G_K` of `K`. As `n` varies 

we obtain the Tate module `T_p E` which is 

a representation of `G_K` on a free `\ZZ_p`-module 

of rank `2`. As `p` varies the representations 

are compatible. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 + 1, 'a') 

sage: E = EllipticCurve('11a1').change_ring(K) 

sage: rho = E.galois_representation() 

sage: rho 

Compatible family of Galois representations associated to the Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-10)*x + (-20) over Number Field in a with defining polynomial x^2 + 1 

""" 

 

 

def __init__(self, E): 

r""" 

See ``GaloisRepresentation`` for documentation. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 + 1, 'a') 

sage: E = EllipticCurve('11a1').change_ring(K) 

sage: rho = E.galois_representation() 

sage: rho 

Compatible family of Galois representations associated to the Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-10)*x + (-20) over Number Field in a with defining polynomial x^2 + 1 

sage: loads(rho.dumps()) == rho 

True 

""" 

self.E = E 

 

def __repr__(self): 

r""" 

Return a string representation of the class. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 + 1, 'a') 

sage: E = EllipticCurve('11a1').change_ring(K) 

sage: rho = E.galois_representation() 

sage: rho 

Compatible family of Galois representations associated to the Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-10)*x + (-20) over Number Field in a with defining polynomial x^2 + 1 

 

sage: K.<a> = NumberField(x^2-x+1) 

sage: E = EllipticCurve([0,0,0,a,0]) 

sage: E.galois_representation() 

Compatible family of Galois representations associated to the CM Elliptic Curve defined by y^2 = x^3 + a*x over Number Field in a with defining polynomial x^2 - x + 1 

""" 

if self.E.has_cm(): 

return "Compatible family of Galois representations associated to the CM " + repr(self.E) 

else: 

return "Compatible family of Galois representations associated to the " + repr(self.E) 

 

 

def __eq__(self,other): 

r""" 

Compares two Galois representations. 

 

We define two compatible families of representations attached 

to elliptic curves to be equal if the curves are isomorphic. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 + 1, 'a'); a = K.gen() 

sage: rho1 = EllipticCurve_from_j(1 + a).galois_representation() 

sage: rho2 = EllipticCurve_from_j(2 + a).galois_representation() 

sage: rho1 == rho1 

True 

sage: rho1 == rho2 

False 

sage: rho1 == 42 

False 

""" 

if type(self) is not type(other): 

return False 

return self.E.is_isomorphic(other.E) 

 

def elliptic_curve(self): 

r""" 

Return the elliptic curve associated to this representation. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 + 1, 'a'); a = K.gen() 

sage: E = EllipticCurve_from_j(a) 

sage: rho = E.galois_representation() 

sage: rho.elliptic_curve() == E 

True 

""" 

return self.E 

 

 

def non_surjective(self, A=100): 

r""" 

Return a list of primes `p` including all primes for which the mod-`p` 

representation might not be surjective. 

 

INPUT: 

 

- ``A`` -- int (a bound on the number of traces of Frobenius to use 

while trying to prove surjectivity). 

 

OUTPUT: 

 

- ``list`` -- A list of primes where mod-`p` representation is 

very likely not surjective. At any prime not in this list, 

the representation is definitely surjective. If `E` has CM, 

the list [0] is returned. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 - 29, 'a'); a = K.gen() 

sage: E = EllipticCurve([1, 0, ((5 + a)/2)**2, 0, 0]) 

sage: rho = E.galois_representation() 

sage: rho.non_surjective() # See Section 5.10 of [Serre72]. 

[3, 5, 29] 

sage: K = NumberField(x**2 + 3, 'a'); a = K.gen() 

sage: E = EllipticCurve([0, -1, 1, -10, -20]).change_ring(K) # X_0(11) 

sage: rho = E.galois_representation() 

sage: rho.non_surjective() # long time (4s on sage.math, 2014) 

[3, 5] 

sage: K = NumberField(x**2 + 1, 'a'); a = K.gen() 

sage: E = EllipticCurve_from_j(1728).change_ring(K) # CM 

sage: rho = E.galois_representation() 

sage: rho.non_surjective() 

[0] 

sage: K = NumberField(x**2 - 5, 'a'); a = K.gen() 

sage: E = EllipticCurve_from_j(146329141248*a - 327201914880) # CM 

sage: rho = E.galois_representation() 

sage: rho.non_surjective() # long time (3s on sage.math, 2014) 

[0] 

 

TESTS: 

 

An example which failed until fixed at :trac:`19229`:: 

 

sage: K.<a> = NumberField(x^2-x+1) 

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

sage: rho = E.galois_representation() 

sage: rho.non_surjective() 

[2, 3] 

 

""" 

if self.E.has_cm(): 

return [0] 

return _non_surjective(self.E, A) 

 

def is_surjective(self, p, A=100): 

r""" 

Return ``True`` if the mod-p representation is (provably) 

surjective onto `Aut(E[p]) = GL_2(\mathbb{F}_p)`. Return 

``False`` if it is (probably) not. 

 

INPUT: 

 

* ``p`` - int - a prime number. 

 

* ``A`` - int - a bound on the number of traces of Frobenius to use 

while trying to prove surjectivity. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 - 29, 'a'); a = K.gen() 

sage: E = EllipticCurve([1, 0, ((5 + a)/2)**2, 0, 0]) 

sage: rho = E.galois_representation() 

sage: rho.is_surjective(29) # Cyclotomic character not surjective. 

False 

sage: rho.is_surjective(7) # See Section 5.10 of [Serre72]. 

True 

 

If `E` is defined over `\QQ`, then the exceptional primes for `E_{/K}` 

are the same as the exceptional primes for `E`, except for those primes 

that are ramified in `K/\QQ` or are less than `[K:\QQ]`:: 

 

sage: K = NumberField(x**2 + 11, 'a') 

sage: E = EllipticCurve([2, 14]) 

sage: rhoQQ = E.galois_representation() 

sage: rhoK = E.change_ring(K).galois_representation() 

sage: rhoQQ.is_surjective(2) == rhoK.is_surjective(2) 

False 

sage: rhoQQ.is_surjective(3) == rhoK.is_surjective(3) 

True 

sage: rhoQQ.is_surjective(5) == rhoK.is_surjective(5) 

True 

 

For CM curves, the mod-p representation is never surjective:: 

 

sage: K.<a> = NumberField(x^2-x+1) 

sage: E = EllipticCurve([0,0,0,0,a]) 

sage: E.has_cm() 

True 

sage: rho = E.galois_representation() 

sage: any(rho.is_surjective(p) for p in [2,3,5,7]) 

False 

""" 

if self.E.has_cm(): 

return False 

return (_exceptionals(self.E, [p], A) == []) 

 

def isogeny_bound(self, A=100): 

r""" 

Return a list of primes `p` including all primes for which 

the image of the mod-`p` representation is contained in a 

Borel. 

 

.. NOTE:: 

 

For the actual list of primes `p` at which the 

representation is reducible see :meth:`reducible_primes()`. 

 

INPUT: 

 

- ``A`` -- int (a bound on the number of traces of Frobenius to 

use while trying to prove the mod-`p` 

representation is not contained in a Borel). 

 

OUTPUT: 

 

- ``list`` - A list of primes which contains (but may not be 

equal to) all `p` for which the image of the mod-`p` 

representation is contained in a Borel subgroup. At any 

prime not in this list, the image is definitely not 

contained in a Borel. If E has `CM` defined over `K`, the list 

[0] is returned. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 - 29, 'a'); a = K.gen() 

sage: E = EllipticCurve([1, 0, ((5 + a)/2)**2, 0, 0]) 

sage: rho = E.galois_representation() 

sage: rho.isogeny_bound() # See Section 5.10 of [Serre72]. 

[3, 5] 

sage: K = NumberField(x**2 + 1, 'a') 

sage: EllipticCurve_from_j(K(1728)).galois_representation().isogeny_bound() # CM over K 

[0] 

sage: EllipticCurve_from_j(K(0)).galois_representation().isogeny_bound() # CM NOT over K 

[2, 3] 

sage: E = EllipticCurve_from_j(K(2268945/128)) # c.f. [Sutherland12] 

sage: E.galois_representation().isogeny_bound() # No 7-isogeny, but... 

[7] 

 

For curves with rational CM, there are infinitely many primes 

`p` for which the mod-`p` representation is reducible, and [0] 

is returned:: 

 

sage: K.<a> = NumberField(x^2-x+1) 

sage: E = EllipticCurve([0,0,0,0,a]) 

sage: E.has_rational_cm() 

True 

sage: rho = E.galois_representation() 

sage: rho.isogeny_bound() 

[0] 

 

An example (an elliptic curve with everywhere good reduction 

over an imaginary quadratic field with quite large 

discriminant), which failed until fixed at :trac:`21776`:: 

 

sage: K.<a> = NumberField(x^2 - x + 112941801) 

sage: E = EllipticCurve([a+1,a-1,a,-23163076*a + 266044005933275,57560769602038*a - 836483958630700313803]) 

sage: E.conductor().norm() 

1 

sage: GR = E.galois_representation() 

sage: GR.isogeny_bound() 

[] 

 

""" 

if self.E.has_rational_cm(): 

return [0] 

 

E = _over_numberfield(self.E) 

K = E.base_field() 

 

char = lambda P: P.smallest_integer() # cheaper than constructing the residue field 

 

# semistable reducible primes (we are now not in the CM case) 

bad_primes = _semistable_reducible_primes(E) 

 

# primes of additive reduction 

bad_primesK = (K.ideal(E.c4()) + K.ideal(E.discriminant())).prime_factors() 

bad_primes += [char(P) for P in bad_primesK] 

 

# ramified primes 

bad_primes += K.absolute_discriminant().prime_factors() 

 

# remove repeats: 

bad_primes = list(Set(bad_primes)) 

 

return _maybe_borels(E, bad_primes, A) 

 

def reducible_primes(self): 

r""" 

Return a list of primes `p` for which the mod-`p` 

representation is reducible, or [0] for CM curves. 

 

OUTPUT: 

 

- ``list`` - A list of those primes `p` for which the mod-`p` 

representation is contained in a Borel subgroup, i.e. is 

reducible. If E has CM *defined over K*, the list [0] is 

returned (in this case the representation is reducible for 

infinitely many primes). 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 - 29, 'a'); a = K.gen() 

sage: E = EllipticCurve([1, 0, ((5 + a)/2)**2, 0, 0]) 

sage: rho = E.galois_representation() 

sage: rho.isogeny_bound() # See Section 5.10 of [Serre72]. 

[3, 5] 

sage: rho.reducible_primes() 

[3, 5] 

 

sage: K = NumberField(x**2 + 1, 'a') 

sage: EllipticCurve_from_j(K(1728)).galois_representation().isogeny_bound() # CM over K 

[0] 

sage: EllipticCurve_from_j(K(0)).galois_representation().reducible_primes() # CM but NOT over K 

[2, 3] 

sage: E = EllipticCurve_from_j(K(2268945/128)) # c.f. [Sutherland12] 

sage: rho = E.galois_representation() 

sage: rho.isogeny_bound() # ... but there is no 7-isogeny ... 

[7] 

sage: rho.reducible_primes() 

[] 

 

For curves with rational CM, there are infinitely many primes 

`p` for which the mod-`p` representation is reducible, and [0] 

is returned:: 

 

sage: K.<a> = NumberField(x^2-x+1) 

sage: E = EllipticCurve([0,0,0,0,a]) 

sage: E.has_rational_cm() 

True 

sage: rho = E.galois_representation() 

sage: rho.reducible_primes() 

[0] 

""" 

if self.E.has_rational_cm(): 

return [0] 

 

return [l for l in self.isogeny_bound() if self.E.isogenies_prime_degree(l)] 

 

def _non_surjective(E, patience=100): 

r""" 

Return a list of primes `p` including all primes for which the mod-`p` 

representation might not be surjective. 

 

INPUT: 

 

- ``E`` - EllipticCurve (over a number field). 

 

- ``A`` - int (a bound on the number of traces of Frobenius to use 

while trying to prove surjectivity). 

 

OUTPUT: 

 

- ``list`` - A list of primes where mod-`p` representation is very likely 

not surjective. At any prime not in this list, the representation is 

definitely surjective. If E has CM, a ValueError is raised. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 - 29, 'a'); a = K.gen() 

sage: E = EllipticCurve([1, 0, ((5 + a)/2)**2, 0, 0]) 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._non_surjective(E) # See Section 5.10 of [Serre72]. 

[3, 5, 29] 

sage: E = EllipticCurve_from_j(1728).change_ring(K) # CM 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._non_surjective(E) 

Traceback (most recent call last): 

... 

ValueError: The curve E should not have CM. 

""" 

if E.has_cm(): 

raise ValueError("The curve E should not have CM.") 

 

E = _over_numberfield(E) 

K = E.base_field() 

 

exceptional_primes = [2, 3, 5, 7, 11, 13, 17, 19] 

# The possible primes l unramified in K/QQ for which the image of the mod l 

# Galois representation could be contained in an exceptional subgroup. 

 

# Find the places of additive reduction. 

SA = [] 

for P, n in E.conductor().factor(): 

if n > 1: 

SA.append(P) 

# TODO: All we really require is a list of primes that include all 

# primes at which E has additive reduction. Perhaps we can speed 

# up things by doing something less time-consuming here that produces 

# a list with some extra terms? (Of course, the longer this list is, 

# the slower the rest of the computation is, so it is not clear that 

# this would help...) 

 

char = lambda P: P.smallest_integer() # cheaper than constructing the residue field 

 

bad_primes = exceptional_primes 

bad_primes += [char(P) for P in SA] 

bad_primes += K.discriminant().prime_factors() 

bad_primes += _semistable_reducible_primes(E) 

bad_primes += _possible_normalizers(E, SA) 

 

bad_primes = list(Set(bad_primes)) 

 

return _exceptionals(E, bad_primes, patience) 

 

 

def _maybe_borels(E, L, patience=100): 

r""" 

Determine which primes in L might have an image contained in a 

Borel subgroup, using straight-forward checking of traces of 

Frobenius. 

 

.. NOTE:: 

 

This function will sometimes return primes for which the image 

is not contained in a Borel subgroup. This issue cannot always 

be fixed by increasing patience as it may be a result of a 

failure of a local-global principle for isogenies. 

 

INPUT: 

 

- ``E`` -- EllipticCurve - over a number field. 

 

- ``L`` -- list - a list of prime numbers. 

 

- ``patience`` -- int (a positive integer bounding the number of 

traces of Frobenius to use while trying to 

prove irreducibility). 

 

OUTPUT: 

 

- list -- The list of all primes `\ell` in L for which the 

mod `\ell` image might be contained in a Borel 

subgroup of `GL_2(\mathbf{F}_{\ell})`. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve('11a1') # has a 5-isogeny 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._maybe_borels(E,primes(40)) 

[5] 

 

Example to show that the output may contain primes where the 

representation is in fact reducible. Over `\QQ` the following is 

essentially the unique such example by [Sutherland12]_:: 

 

sage: E = EllipticCurve_from_j(2268945/128) 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._maybe_borels(E, [7, 11]) 

[7] 

 

This curve does possess a 7-isogeny modulo every prime of good 

reduction, but has no rational 7-isogeny:: 

 

sage: E.isogenies_prime_degree(7) 

[] 

 

A number field example: 

 

sage: K.<i> = QuadraticField(-1) 

sage: E = EllipticCurve([1+i, -i, i, -399-240*i, 2627+2869*i]) 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._maybe_borels(E, primes(20)) 

[2, 3] 

 

Here the curve really does possess isogenies of degrees 2 and 3:: 

 

sage: [len(E.isogenies_prime_degree(l)) for l in [2,3]] 

[1, 1] 

 

""" 

E = _over_numberfield(E) 

K = E.base_field() 

 

L = sorted(set(L)) # Remove duplicates from L and makes a copy for output 

 

include_2 = False 

if 2 in L: # c.f. Section 5.3(a) of [Serre72]. 

L.remove(2) 

include_2 = not E.division_polynomial(2).is_irreducible() 

 

for P in deg_one_primes_iter(K): 

if not (L and patience): # stop if no primes are left, or 

# patience is exhausted 

break 

 

patience -= 1 

 

# Check whether the Frobenius polynomial at P is irreducible 

# modulo each l, dropping l from the list if so. 

 

try: 

trace = E.change_ring(P.residue_field()).trace_of_frobenius() 

except ArithmeticError: # Bad reduction at P. 

continue 

 

determinant = P.norm() 

discriminant = trace**2 - 4 * determinant 

 

for l in L: 

if legendre_symbol(discriminant,l)==-1: 

L.remove(l) 

 

if include_2: 

L = [2] + L 

return L 

 

def _exceptionals(E, L, patience=1000): 

r""" 

Determine which primes in L are exceptional for E, using Proposition 19 

of Section 2.8 of Serre's ``Propriétés Galoisiennes des Points d'Ordre 

Fini des Courbes Elliptiques'' [Serre72]_. 

 

INPUT: 

 

- ``E`` - EllipticCurve - over a number field. 

 

- ``L`` - list - a list of prime numbers. 

 

- ``patience`` - int (a bound on the number of traces of Frobenius to 

use while trying to prove surjectivity). 

 

OUTPUT: 

 

- list -- The list of all primes l in L for which the mod l image 

might fail to be surjective. 

 

EXAMPLES:: 

 

sage: K = NumberField(x**2 - 29, 'a'); a = K.gen() 

sage: E = EllipticCurve([1, 0, ((5 + a)/2)**2, 0, 0]) 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._exceptionals(E, [29, 31]) 

[29] 

 

For CM curves an error is raised:: 

 

sage: E = EllipticCurve_from_j(1728).change_ring(K) # CM 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._exceptionals(E,[2,3,5]) 

Traceback (most recent call last): 

... 

ValueError: The curve E should not have CM. 

""" 

if E.has_cm(): 

raise ValueError("The curve E should not have CM.") 

 

E = _over_numberfield(E) 

K = E.base_field() 

 

output = [] 

 

L = list(set(L)) # Remove duplicates from L. 

 

for l in L: 

if l == 2: # c.f. Section 5.3(a) of [Serre72]. 

if (E.j_invariant() - 1728).is_square(): 

output.append(2) 

elif not E.division_polynomial(2).is_irreducible(): 

output.append(2) 

 

elif l == 3: # c.f. Section 5.3(b) of [Serre72]. 

if K(-3).is_square(): 

output.append(3) 

elif not (K['x'].gen()**3 - E.j_invariant()).is_irreducible(): 

output.append(3) 

elif not E.division_polynomial(3).is_irreducible(): 

output.append(3) 

 

elif (K.discriminant() % l) == 0: 

if not K['x'](cyclotomic_polynomial(l)).is_irreducible(): 

# I.E. if the action on lth roots of unity is not surjective 

# (We want this since as a Galois module, \wedge^2 E[l] 

# is isomorphic to the lth roots of unity.) 

output.append(l) 

 

for l in output: 

L.remove(l) 

if 2 in L: 

L.remove(2) 

if 3 in L: 

L.remove(3) 

 

# If the image is not surjective, then it is contained in one of the 

# maximal subgroups. So, we start by creating a dictionary between primes 

# l in L and possible maximal subgroups in which the mod l image could 

# be contained. This information is stored as a triple whose elements 

# are True/False according to whether the mod l image could be contained 

# in: 

# 0. A Borel or normalizer of split Cartan subgroup. 

# 1. A nonsplit Cartan subgroup or its normalizer. 

# 2. An exceptional subgroup of GL_2. 

 

D = {} 

for l in L: 

D[l] = [True, True, True] 

 

for P in deg_one_primes_iter(K): 

try: 

trace = E.change_ring(P.residue_field()).trace_of_frobenius() 

except ArithmeticError: # Bad reduction at P. 

continue 

 

patience -= 1 

 

determinant = P.norm() 

discriminant = trace**2 - 4 * determinant 

 

unexc = [] # Primes we discover are unexceptional go here. 

 

for l in D: 

tr = GF(l)(trace) 

det = GF(l)(determinant) 

disc = GF(l)(discriminant) 

 

if tr == 0: 

# I.E. if Frob_P could be contained in the normalizer of 

# a Cartan subgroup, but not in the Cartan subgroup. 

continue 

 

if disc == 0: 

# I.E. If the matrix might be non-diagonalizable over F_{p^2}. 

continue 

 

if legendre_symbol(disc, l) == 1: 

# If the matrix is diagonalizable over F_p, it can't be 

# contained in a non-split Cartan subgroup. Since we've 

# gotten rid of the case where it is contained in the 

# of a nonsplit Cartan subgroup but not the Cartan subgroup, 

D[l][1] = False 

else: 

# If the matrix is not diagonalizable over F_p, it can't 

# be contained Borel subgroup. 

D[l][0] = False 

 

if det != 0: # c.f. [Serre72], Section 2.8, Prop. 19 

u = trace**2 / det 

if u not in (1, 2, 4) and u**2 - 3 * u + 1 != 0: 

D[l][2] = False 

 

 

if D[l] == [False, False, False]: 

unexc.append(l) 

 

for l in unexc: 

D.pop(l) 

unexc = [] 

 

if (D == {}) or (patience == 0): 

break 

 

for l in D: 

output.append(l) 

 

output.sort() 

return output 

 

 

def _over_numberfield(E): 

r"""Return `E`, defined over a NumberField object. This is necessary 

since if `E` is defined over `\QQ`, then we cannot use SAGE commands 

available for number fields. 

 

INPUT: 

 

- ``E`` - EllipticCurve - over a number field. 

 

OUTPUT: 

 

- If `E` is defined over a NumberField, returns E. 

 

- If `E` is defined over QQ, returns E defined over the NumberField QQ. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve([1, 2]) 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._over_numberfield(E) 

Elliptic Curve defined by y^2 = x^3 + x + 2 over Number Field in a with defining polynomial x 

""" 

 

K = E.base_field() 

 

if K == QQ: 

x = QQ['x'].gen() 

K = NumberField(x, 'a') 

E = E.change_ring(K) 

 

return E 

 

 

def deg_one_primes_iter(K, principal_only=False): 

r""" 

Return an iterator over degree 1 primes of ``K``. 

 

INPUT: 

 

- ``K`` -- a number field 

- ``principal_only`` -- bool; if ``True``, only yield principal primes 

 

OUTPUT: 

 

An iterator over degree 1 primes of `K` up to the given norm, 

optionally yielding only principal primes. 

 

EXAMPLES:: 

 

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

sage: from sage.schemes.elliptic_curves.gal_reps_number_field import deg_one_primes_iter 

sage: it = deg_one_primes_iter(K) 

sage: [next(it) for _ in range(6)] 

[Fractional ideal (2, a + 1), 

Fractional ideal (3, a + 1), 

Fractional ideal (3, a + 2), 

Fractional ideal (-a), 

Fractional ideal (7, a + 3), 

Fractional ideal (7, a + 4)] 

sage: it = deg_one_primes_iter(K, True) 

sage: [next(it) for _ in range(6)] 

[Fractional ideal (-a), 

Fractional ideal (-2*a + 3), 

Fractional ideal (2*a + 3), 

Fractional ideal (a + 6), 

Fractional ideal (a - 6), 

Fractional ideal (-3*a + 4)] 

""" 

# imaginary quadratic fields have no principal primes of norm < disc / 4 

start = K.discriminant().abs() // 4 if principal_only and K.signature() == (0,1) else 2 

 

K_is_Q = (K==QQ) 

from sage.arith.misc import primes 

from sage.rings.infinity import infinity 

for p in primes(start=start, stop=infinity): 

if K_is_Q: 

yield ZZ.ideal(p) 

else: 

for P in K.primes_above(p, degree=1): 

if not principal_only or P.is_principal(): 

yield P 

 

def _semistable_reducible_primes(E, verbose=False): 

r"""Find a list containing all semistable primes l unramified in K/QQ 

for which the Galois image for E could be reducible. 

 

INPUT: 

 

- ``E`` - EllipticCurve - over a number field. 

 

OUTPUT: 

 

A list of primes, which contains all primes `l` unramified in 

`K/\mathbb{QQ}`, such that `E` is semistable at all primes lying 

over `l`, and the Galois image at `l` is reducible. If `E` has CM 

defined over its ground field, a ``ValueError`` is raised. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve([0, -1, 1, -10, -20]) # X_0(11) 

sage: 5 in sage.schemes.elliptic_curves.gal_reps_number_field._semistable_reducible_primes(E) 

True 

 

This example, over a quintic field with Galois group `S_5`, took a 

very long time before :trac:`22343`:: 

 

sage: K.<a> = NumberField(x^5 - 6*x^3 + 8*x - 1) 

sage: E = EllipticCurve(K, [a^3 - 2*a, a^4 - 2*a^3 - 4*a^2 + 6*a + 1, a + 1, -a^3 + a + 1, -a]) 

sage: from sage.schemes.elliptic_curves.gal_reps_number_field import _semistable_reducible_primes 

sage: _semistable_reducible_primes(E) 

[2, 5, 53, 1117] 

""" 

if verbose: print("In _semistable_reducible_primes with E={}".format(E.ainvs())) 

K = E.base_field() 

d = K.degree() 

 

deg_one_primes = deg_one_primes_iter(K, principal_only=True) 

 

bad_primes = set([]) # This will store the output. 

 

# We find two primes (of distinct residue characteristics) which are 

# of degree 1, unramified in K/Q, and at which E has good reduction. 

# Each of these primes will give us a nontrivial divisibility constraint 

# on the exceptional primes l. For both of these primes P, we precompute 

# a generator and the characteristic polynomial of Frob_P^12. 

 

precomp = [] 

last_p = 0 # The residue characteristic of the most recent prime. 

 

while len(precomp) < 2: 

P = next(deg_one_primes) 

p = P.norm() 

if p != last_p and (d==1 or P.ramification_index() == 1) and E.has_good_reduction(P): 

precomp.append(P) 

last_p = p 

 

Px, Py = precomp 

x, y = [P.gens_reduced()[0] for P in precomp] 

EmodPx = E.reduction(Px) if d>1 else E.reduction(x) 

EmodPy = E.reduction(Py) if d>1 else E.reduction(y) 

fxpol = EmodPx.frobenius_polynomial() 

fypol = EmodPy.frobenius_polynomial() 

fx12pol = fxpol.adams_operator(12) # roots are 12th powers of those of fxpol 

fy12pol = fypol.adams_operator(12) 

px = x.norm() if d>1 else x 

py = y.norm() if d>1 else x 

Zx = fxpol.parent() 

xpol = x.charpoly() if d>1 else Zx([-x,1]) 

ypol = y.charpoly() if d>1 else Zx([-y,1]) 

 

if verbose: print("Finished precomp, x={} (p={}), y={} (p={})".format(x,px,y,py)) 

 

for w in range(1+d/2): 

if verbose: print("w = {}".format(w)) 

gx = xpol.symmetric_power(w).adams_operator(12).resultant(fx12pol) 

gy = ypol.symmetric_power(w).adams_operator(12).resultant(fy12pol) 

if verbose: print("computed gx and gy") 

 

gxn = Integer(gx.absolute_norm()) if d>1 else gx 

gyn = Integer(gy.absolute_norm()) if d>1 else gy 

gxyn = gxn.gcd(gyn) 

if gxyn: 

xprimes = gxyn.prime_factors() 

if verbose: print("adding prime factors {} of {} to {}".format(xprimes, gxyn, sorted(bad_primes))) 

bad_primes.update(xprimes) 

if verbose: print("...done, bad_primes now {}".format(sorted(bad_primes))) 

continue 

else: 

if verbose: print("gx and gy both 0!") 

 

 

## It is possible that our curve has CM. ## 

 

# Our character must be of the form Nm^K_F for an imaginary 

# quadratic subfield F of K (which is the CM field if E has CM). 

 

# Note that this can only happen when d is even, w=d/2, and K 

# contains (or the Galois closure of K contains?) the 

# imaginary quadratic field F = Q(sqrt(a)) which is the 

# splitting field of both fx12pol and fy12pol. We compute a 

# and relativise K over F: 

 

a = fx12pol.discriminant().squarefree_part() 

 

# Construct a field isomorphic to K but a relative extension over QQ(sqrt(a)). 

 

# See #19229: the names given here, which are not used, should 

# not be the name of the generator of the base field. 

 

rootsa = K(a).sqrt(all=True) # otherwise if a is not a square the 

# returned result is in the symbolic ring! 

try: 

roota = rootsa[0] 

except IndexError: 

raise RuntimeError("error in _semistable_reducible_primes: K={} does not contain sqrt({})".format(K,a)) 

K_rel = K.relativize(roota, ['name1','name2']) 

iso = K_rel.structure()[1] # an isomorphism from K to K_rel 

E_rel = E.change_ring(iso) # same as E but over K_rel 

 

## We try again to find a nontrivial divisibility condition. ## 

 

div = 0 

patience = 5 * K.absolute_degree() 

# Number of Frobenius elements to check before suspecting that E 

# has CM and computing the set of CM j-invariants of K to check. 

# TODO: Is this the best value for this parameter? 

 

while div==0 and patience>0: 

P = next(deg_one_primes) # a prime of K not K_rel 

while E.has_bad_reduction(P): 

P = next(deg_one_primes) 

 

if verbose: print("trying P = {}...".format(P)) 

EmodP = E.reduction(P) 

fpol = EmodP.frobenius_polynomial() 

if verbose: print("...good reduction, frobenius poly = {}".format(fpol)) 

x = iso(P.gens_reduced()[0]).relative_norm() 

xpol = x.charpoly().adams_operator(12) 

div2 = Integer(xpol.resultant(fpol.adams_operator(12)) // x.norm()**12) 

if div2: 

div = div2.isqrt() 

assert div2==div**2 

if verbose: print("...div = {}".format(div)) 

else: 

if verbose: print("...div = 0, continuing") 

patience -= 1 

 

if patience == 0: 

# We suspect that E has CM, so we check: 

if E.has_cm(): 

raise ValueError("In _semistable_reducible_primes, the curve E should not have CM.") 

 

assert div != 0 

# We found our divisibility constraint. 

 

xprimes = div.prime_factors() 

if verbose: print("...adding prime factors {} of {} to {}...".format(xprimes,div, sorted(bad_primes))) 

bad_primes.update(xprimes) 

if verbose: print("...done, bad_primes now {}".format(sorted(bad_primes))) 

 

L = sorted(bad_primes) 

return L 

 

 

def _possible_normalizers(E, SA): 

r"""Find a list containing all primes `l` such that the Galois image at `l` 

is contained in the normalizer of a Cartan subgroup, such that the 

corresponding quadratic character is ramified only at the given primes. 

 

INPUT: 

 

- ``E`` - EllipticCurve - over a number field K. 

 

- ``SA`` - list - a list of primes of K. 

 

OUTPUT: 

 

- list -- A list of primes, which contains all primes `l` such that the 

Galois image at `l` is contained in the normalizer of a Cartan 

subgroup, such that the corresponding quadratic character is 

ramified only at primes in SA. 

 

- If `E` has geometric CM that is not defined over its ground field, a 

ValueError is raised. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve([0,0,0,-56,4848]) 

sage: 5 in sage.schemes.elliptic_curves.gal_reps_number_field._possible_normalizers(E, [ZZ.ideal(2)]) 

True 

 

For CM curves, an error is raised:: 

 

sage: K.<i> = QuadraticField(-1) 

sage: E = EllipticCurve_from_j(1728).change_ring(K) # CM 

sage: sage.schemes.elliptic_curves.gal_reps_number_field._possible_normalizers(E, []) 

Traceback (most recent call last): 

... 

ValueError: The curve E should not have CM. 

 

""" 

if E.has_cm(): 

raise ValueError("The curve E should not have CM.") 

 

E = _over_numberfield(E) 

K = E.base_field() 

SA = [K.ideal(I.gens()) for I in SA] 

 

x = K['x'].gen() 

selmer_group = K.selmer_group(SA, 2) # Generators of the selmer group. 

 

if selmer_group == []: 

return [] 

 

V = VectorSpace(GF(2), len(selmer_group)) 

# We think of this as the character group of the selmer group. 

 

traces_list = [] 

W = V.zero_subspace() 

 

deg_one_primes = deg_one_primes_iter(K) 

 

while W.dimension() < V.dimension() - 1: 

P = next(deg_one_primes) 

 

k = P.residue_field() 

 

defines_valid_character = True 

# A prime P defines a quadratic residue character 

# on the Selmer group. This variable will be set 

# to zero if any elements of the selmer group are 

# zero mod P (i.e. the character is ramified). 

 

splitting_vector = [] # This will be the values of this 

# character on the generators of the Selmer group. 

 

for a in selmer_group: 

abar = k(a) 

if abar == 0: 

# Ramification. 

defines_valid_character = False 

break 

 

if abar.is_square(): 

splitting_vector.append(GF(2)(0)) 

else: 

splitting_vector.append(GF(2)(1)) 

 

if not defines_valid_character: 

continue 

 

if splitting_vector in W: 

continue 

 

try: 

Etilde = E.change_ring(k) 

except ArithmeticError: # Bad reduction. 

continue 

 

tr = Etilde.trace_of_frobenius() 

 

if tr == 0: 

continue 

 

traces_list.append(tr) 

 

W = W + V.span([splitting_vector]) 

 

bad_primes = set([]) 

 

for i in traces_list: 

for p in i.prime_factors(): 

bad_primes.add(p) 

 

# We find the unique vector v in V orthogonal to W: 

v = W.matrix().transpose().kernel().basis()[0] 

 

# We find the element a of the selmer group corresponding to v: 

a = 1 

for i in range(len(selmer_group)): 

if v[i] == 1: 

a *= selmer_group[i] 

 

# Since we've already included the above bad primes, we can assume 

# that the quadratic character corresponding to the exceptional primes 

# we're looking for is given by mapping into Gal(K[\sqrt{a}]/K). 

 

patience = 5 * K.degree() 

# Number of Frobenius elements to check before suspecting that E 

# has CM and computing the set of CM j-invariants of K to check. 

# TODO: Is this the best value for this parameter? 

 

while True: 

P = next(deg_one_primes) 

 

k = P.residue_field() 

 

if not k(a).is_square(): 

try: 

tr = E.change_ring(k).trace_of_frobenius() 

except ArithmeticError: # Bad reduction. 

continue 

 

if tr == 0: 

patience -= 1 

 

else: 

for p in tr.prime_factors(): 

bad_primes.add(p) 

 

bad_primes = sorted(bad_primes) 

return bad_primes