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

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

r""" 

Modular symbols attached to elliptic curves over `\QQ` 

 

To an elliptic curve `E` over the rational numbers with conductor `N`, 

one can associate a space of modular symbols of level `N`, because `E` 

is known to be modular. The space is two-dimensional and contains a 

subspace on which complex conjugation acts as multiplication by `+1` 

and one on which it acts by `-1`. 

 

There are two implementations of modular symbols, one within ``Sage`` 

and the other in Cremona's ``eclib`` library. One can choose here 

which one is used. 

 

Associated to `E` there is a canonical generator in each space. They are maps 

`[.]^+` and `[.]^{-}`, both `\QQ \to\QQ`. They are normalized such that 

 

.. MATH:: 

 

[r]^{+} \Omega^{+} + [r]^{-}\Omega^{-} = \int_{\infty}^r 2\pi i f(z) dz 

 

where `f` is the newform associated to the isogeny class of `E` and 

`\Omega^{+}` is the smallest positive period of the Néron differential 

of `E` and `\Omega^{-}` is the smallest positive purely imaginary 

period. Note that it depends on `E` rather than on its isogeny class. 

 

From ``eclib`` version v20161230, both plus and minus symbols are 

available and are correctly normalized. In the ``Sage`` 

implementation, the computation of the space provides initial 

generators which are not necessarily correctly normalized; here we 

implement two methods that try to find the correct scaling factor. 

 

Modular symbols are used to compute `p`-adic `L`-functions. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve("19a1") 

sage: m = E.modular_symbol() 

sage: m(0) 

1/3 

sage: m(1/17) 

-2/3 

sage: m2 = E.modular_symbol(-1, implementation="sage") 

sage: m2(0) 

0 

sage: m2(1/5) 

1/2 

 

sage: V = E.modular_symbol_space() 

sage: V 

Modular Symbols subspace of dimension 1 of Modular Symbols space of dimension 2 for Gamma_0(19) of weight 2 with sign 1 over Rational Field 

sage: V.q_eigenform(30) 

q - 2*q^3 - 2*q^4 + 3*q^5 - q^7 + q^9 + 3*q^11 + 4*q^12 - 4*q^13 - 6*q^15 + 4*q^16 - 3*q^17 + q^19 - 6*q^20 + 2*q^21 + 4*q^25 + 4*q^27 + 2*q^28 + 6*q^29 + O(q^30) 

 

For more details on modular symbols consult the following 

 

REFERENCES: 

 

.. [MaTaTe] B. Mazur, J. Tate, and J. Teitelbaum, 

On `p`-adic analogues of the conjectures of Birch and 

Swinnerton-Dyer, Inventiones mathematicae 84, (1986), 1-48. 

 

.. [Crem97] John Cremona, Algorithms for modular elliptic curves, 

Cambridge University Press, 1997. 

 

.. [StWu] William Stein and Christian Wuthrich, Algorithms for the 

Arithmetic of Elliptic Curves using Iwasawa Theory, Mathematics 

of Computation 82 (2013), 1757-1792. 

 

AUTHORS: 

 

- William Stein (2007): first version 

 

- Chris Wuthrich (2008): add scaling and reference to eclib 

 

- John Cremona (2016): reworked eclib interface 

 

""" 

 

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

# Copyright (C) 2007 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 __future__ import print_function 

 

from sage.structure.sage_object import SageObject 

from sage.modular.modsym.all import ModularSymbols 

from sage.libs.eclib.newforms import ECModularSymbol 

from sage.databases.cremona import parse_cremona_label 

 

from sage.arith.all import next_prime, kronecker_symbol, prime_divisors, valuation 

from sage.rings.infinity import unsigned_infinity as infinity 

from sage.rings.integer import Integer 

from sage.modular.cusps import Cusps 

from sage.rings.integer_ring import ZZ 

from sage.rings.rational_field import QQ 

from sage.misc.all import verbose 

 

from sage.schemes.elliptic_curves.constructor import EllipticCurve 

 

oo = Cusps(infinity) 

zero = Integer(0) 

 

def modular_symbol_space(E, sign, base_ring, bound=None): 

r""" 

Creates the space of modular symbols of a given sign over a give base_ring, 

attached to the isogeny class of the elliptic curve ``E``. 

 

INPUT: 

 

- ``E`` - an elliptic curve over `\QQ` 

- ``sign`` - integer, -1, 0, or 1 

- ``base_ring`` - ring 

- ``bound`` - (default: None) maximum number of Hecke operators to 

use to cut out modular symbols factor. If None, use 

enough to provably get the correct answer. 

 

OUTPUT: a space of modular symbols 

 

EXAMPLES:: 

 

sage: import sage.schemes.elliptic_curves.ell_modular_symbols 

sage: E=EllipticCurve('11a1') 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.modular_symbol_space(E,-1,GF(37)) 

sage: M 

Modular Symbols space of dimension 1 for Gamma_0(11) of weight 2 with sign -1 over Finite Field of size 37 

 

""" 

if not sign in [-1,0,1]: 

raise TypeError('sign must -1, 0 or 1') 

_sign = int(sign) 

N = E.conductor() 

M = ModularSymbols(N, sign=sign, base_ring=base_ring) 

if bound is None: 

bound = M.hecke_bound() + 10 

V = M 

p = 2 

target_dim = 1 if sign else 2 

while p <= bound and V.dimension() > target_dim: 

t = V.T(p) 

ap = E.ap(p) 

V = (t - ap).kernel() 

p = next_prime(p) 

 

return V 

 

class ModularSymbol(SageObject): 

r""" 

A modular symbol attached to an elliptic curve, which is the map 

`\QQ\to \QQ` obtained by sending `r` to the normalized 

symmetrized (or anti-symmetrized) integral `\infty` to `r`. 

 

This is as defined in [MaTaTe]_, but normalized to depend on the curve 

and not only its isogeny class as in [StWu]_. 

 

See the documentation of ``E.modular_symbol()`` in elliptic curves 

over the rational numbers for help. 

 

""" 

 

def sign(self): 

r""" 

Return the sign of this elliptic curve modular symbol. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol() 

sage: m.sign() 

1 

sage: m = EllipticCurve('11a1').modular_symbol(sign=-1, implementation="sage") 

sage: m.sign() 

-1 

""" 

return self._sign 

 

def elliptic_curve(self): 

r""" 

Return the elliptic curve of this modular symbol. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol() 

sage: m.elliptic_curve() 

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

 

""" 

return self._E 

 

def base_ring(self): 

r""" 

Return the base ring for this modular symbol. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol() 

sage: m.base_ring() 

Rational Field 

""" 

return self._base_ring 

 

def _repr_(self): 

r""" 

String representation of modular symbols. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol() 

sage: m 

Modular symbol with sign 1 over Rational Field attached to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field 

sage: m = EllipticCurve('43a1').modular_symbol(sign=-1, implementation="sage") 

sage: m 

Modular symbol with sign -1 over Rational Field attached to Elliptic Curve defined by y^2 + y = x^3 + x^2 over Rational Field 

""" 

return "Modular symbol with sign %s over %s attached to %s"%( 

self._sign, self._base_ring, self._E) 

 

 

class ModularSymbolECLIB(ModularSymbol): 

def __init__(self, E, sign): 

r""" 

Modular symbols attached to `E` using ``eclib``. 

 

Note that the normalization used within ``eclib`` differs from the 

normalization chosen here by a factor of 2 in the case of elliptic 

curves with negative discriminant (with one real component) since 

the convention there is to write the above integral as 

$[r]^{+}x+[r]^{-}yi$, where the lattice is $\left<2x,x+yi\right>$, 

so that $\Omega^{+}=2x$ and $\Omega^{-}=2yi$. We 

allow for this below. 

 

INPUT: 

 

- ``E`` - an elliptic curve 

- ``sign`` - an integer, -1 or 1 

 

EXAMPLES:: 

 

sage: import sage.schemes.elliptic_curves.ell_modular_symbols 

sage: E=EllipticCurve('11a1') 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolECLIB(E,+1) 

sage: M 

Modular symbol with sign 1 over Rational Field attached to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field 

sage: M(0) 

1/5 

sage: E=EllipticCurve('11a2') 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolECLIB(E,+1) 

sage: M(0) 

1 

 

This is a rank 1 case with vanishing positive twists:: 

 

sage: E=EllipticCurve('121b1') 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolECLIB(E,+1) 

sage: M(0) 

0 

sage: M(1/7) 

1/2 

 

sage: M = EllipticCurve('121d1').modular_symbol(implementation="eclib") 

sage: M(0) 

2 

 

sage: E = EllipticCurve('15a1') 

sage: [C.modular_symbol(implementation="eclib")(0) for C in E.isogeny_class()] 

[1/4, 1/8, 1/4, 1/2, 1/8, 1/16, 1/2, 1] 

 

Since :trac:`10256`, the interface for negative modular symbols in eclib is available:: 

 

sage: E = EllipticCurve('11a1') 

sage: Mplus = E.modular_symbol(+1); Mplus 

Modular symbol with sign 1 over Rational Field attached to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field 

sage: [Mplus(1/i) for i in [1..11]] 

[1/5, -4/5, -3/10, 7/10, 6/5, 6/5, 7/10, -3/10, -4/5, 1/5, 0] 

sage: Mminus = E.modular_symbol(-1); Mminus 

Modular symbol with sign -1 over Rational Field attached to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field 

sage: [Mminus(1/i) for i in [1..11]] 

[0, 0, 1/2, 1/2, 0, 0, -1/2, -1/2, 0, 0, 0] 

 

 

The scaling factor relative to eclib's normalization is 1/2 for curves of negative discriminant:: 

 

sage: [E.discriminant() for E in cremona_curves([14])] 

[-21952, 941192, -1835008, -28, 25088, 98] 

sage: [E.modular_symbol()._scaling for E in cremona_curves([14])] 

[1/2, 1, 1/2, 1/2, 1, 1] 

 

 

TESTS (for :trac:`10236`):: 

 

sage: E = EllipticCurve('11a1') 

sage: m = E.modular_symbol(implementation="eclib") 

sage: m(1/7) 

7/10 

sage: m(0) 

1/5 

""" 

if not sign in [-1,1]: 

raise TypeError('sign must -1 or 1') 

self._sign = ZZ(sign) 

self._E = E 

self._scaling = 1 if E.discriminant()>0 else ZZ(1)/2 

self._use_eclib = True 

self._base_ring = QQ 

# The ECModularSymbol class must be initialized with sign=0 to compute minus symbols 

self._modsym = ECModularSymbol(E, int(sign==1)) 

self.cache = {True: {}, False: {}} 

 

def _call_with_caching(self, r, base_at_infinity=True): 

r""" 

Evaluates the modular symbol at {0,`r`} or {oo,`r`}, caching the computed value. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol(implementation="eclib") 

sage: m._call_with_caching(0) 

1/5 

""" 

cache = self.cache[base_at_infinity] 

try: 

return cache[r] 

except KeyError: 

pass 

c = self._modsym(r, sign=self._sign, base_at_infinity=base_at_infinity) * self._scaling 

cache[r] = c 

return c 

 

def __call__(self, r, base_at_infinity=True): 

r""" 

Evaluates the modular symbol at {0,`r`} or {oo,`r`}. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol(implementation="eclib") 

sage: m(0) 

1/5 

 

""" 

from sage.rings.rational import Rational 

if r != oo: 

r = Rational(r) 

r = r.numer() % r.denom() / r.denom() 

return self._modsym(r, sign=self._sign, base_at_infinity=base_at_infinity) * self._scaling 

 

 

class ModularSymbolSage(ModularSymbol): 

def __init__(self, E, sign, normalize="L_ratio"): 

"""Modular symbols attached to `E` using ``sage``. 

 

INPUT: 

 

- ``E`` -- an elliptic curve 

- ``sign`` -- an integer, -1 or 1 

- ``normalize`` -- either 'L_ratio' (default), 'period', or 

'none'; For 'L_ratio', the modular symbol is correctly 

normalized by comparing it to the quotient of `L(E,1)` by 

the least positive period for the curve and some small 

twists. The normalization 'period' uses the 

integral_period_map for modular symbols and is known to be 

equal to the above normalization up to the sign and a 

possible power of 2. For 'none', the modular symbol is 

almost certainly not correctly normalized, i.e. all values 

will be a fixed scalar multiple of what they should be. But 

the initial computation of the modular symbol is much 

faster, though evaluation of it after computing it won't be 

any faster. 

 

EXAMPLES:: 

 

sage: E=EllipticCurve('11a1') 

sage: import sage.schemes.elliptic_curves.ell_modular_symbols 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolSage(E,+1) 

sage: M 

Modular symbol with sign 1 over Rational Field attached to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field 

sage: M(0) 

1/5 

sage: E=EllipticCurve('11a2') 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolSage(E,+1) 

sage: M(0) 

1 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolSage(E,-1) 

sage: M(1/3) 

1/2 

 

This is a rank 1 case with vanishing positive twists. 

The modular symbol is adjusted by -2:: 

 

sage: E=EllipticCurve('121b1') 

sage: M=sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolSage(E,-1,normalize='L_ratio') 

sage: M(1/3) 

1 

sage: M._scaling 

-1 

 

sage: M = EllipticCurve('121d1').modular_symbol(implementation="sage") 

sage: M(0) 

2 

sage: M = EllipticCurve('121d1').modular_symbol(implementation="sage", normalize='none') 

sage: M(0) 

1 

 

sage: E = EllipticCurve('15a1') 

sage: [C.modular_symbol(implementation="sage", normalize='L_ratio')(0) for C in E.isogeny_class()] 

[1/4, 1/8, 1/4, 1/2, 1/8, 1/16, 1/2, 1] 

sage: [C.modular_symbol(implementation="sage", normalize='period')(0) for C in E.isogeny_class()] 

[1/8, 1/16, 1/8, 1/4, 1/16, 1/32, 1/4, 1/2] 

sage: [C.modular_symbol(implementation="sage", normalize='none')(0) for C in E.isogeny_class()] 

[1, 1, 1, 1, 1, 1, 1, 1] 

 

""" 

if not sign in [-1,1]: 

raise TypeError('sign must -1 or 1') 

self._sign = ZZ(sign) 

self._E = E 

self._use_eclib = False 

self._normalize = normalize 

self._modsym = E.modular_symbol_space(sign=self._sign) 

self._base_ring = self._modsym.base_ring() 

self._ambient_modsym = self._modsym.ambient_module() 

 

if normalize == "L_ratio": 

self._e = self._modsym.dual_eigenvector() 

self._find_scaling_L_ratio() 

if self._failed_to_scale: 

self._find_scaling_period() # will reset _e and _scaling 

else: 

self._e *= self._scaling 

elif normalize == "period" : 

self._find_scaling_period() # this will set _e and _scaling 

elif normalize == "none": 

self._scaling = 1 

self._e = self._modsym.dual_eigenvector() 

else : 

raise ValueError("no normalization %s known for modular symbols"%normalize) 

 

def _find_scaling_L_ratio(self): 

r""" 

This function is use to set ``_scaling``, the factor used to adjust the 

scalar multiple of the modular symbol. 

If `[0]`, the modular symbol evaluated at 0, is non-zero, we can just scale 

it with respect to the approximation of the L-value. It is known that 

the quotient is a rational number with small denominator. 

Otherwise we try to scale using quadratic twists. 

 

``_scaling`` will be set to a rational non-zero multiple if we succeed and to 1 otherwise. 

Even if we fail we scale at least to make up the difference between the periods 

of the `X_0`-optimal curve and our given curve `E` in the isogeny class. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol(implementation="sage") 

sage: m._scaling 

1/5 

sage: m = EllipticCurve('11a2').modular_symbol(implementation="sage") 

sage: m._scaling 

1 

sage: m = EllipticCurve('11a3').modular_symbol(implementation="sage") 

sage: m._scaling 

1/25 

sage: m = EllipticCurve('37a1').modular_symbol(implementation="sage") 

sage: m._scaling 

1 

sage: m = EllipticCurve('37a1').modular_symbol() 

sage: m._scaling 

1 

sage: m = EllipticCurve('389a1').modular_symbol() 

sage: m._scaling 

1 

sage: m = EllipticCurve('389a1').modular_symbol(implementation="sage") 

sage: m._scaling 

2 

sage: m = EllipticCurve('196a1').modular_symbol(implementation="sage") 

sage: m._scaling 

1/2 

 

Some harder cases fail:: 

 

sage: m = EllipticCurve('121b1').modular_symbol(implementation="sage") 

Warning : Could not normalize the modular symbols, maybe all further results will be multiplied by -1 and a power of 2 

sage: m._scaling 

1 

 

TESTS:: 

 

sage: rk0 = ['11a1', '11a2', '15a1', '27a1', '37b1'] 

sage: for la in rk0: # long time (3s on sage.math, 2011) 

....: E = EllipticCurve(la) 

....: me = E.modular_symbol(implementation="eclib") 

....: ms = E.modular_symbol(implementation="sage") 

....: print("{} {} {}".format(E.lseries().L_ratio()*E.real_components(), me(0), ms(0))) 

1/5 1/5 1/5 

1 1 1 

1/4 1/4 1/4 

1/3 1/3 1/3 

2/3 2/3 2/3 

 

sage: rk1 = ['37a1','43a1','53a1', '91b1','91b2','91b3'] 

sage: [EllipticCurve(la).modular_symbol()(0) for la in rk1] # long time (1s on sage.math, 2011) 

[0, 0, 0, 0, 0, 0] 

sage: for la in rk1: # long time (8s on sage.math, 2011) 

....: E = EllipticCurve(la) 

....: m = E.modular_symbol() 

....: lp = E.padic_lseries(5) 

....: for D in [5,17,12,8]: 

....: ED = E.quadratic_twist(D) 

....: md = sum([kronecker(D,u)*m(ZZ(u)/D) for u in range(D)]) 

....: etaD = lp._quotient_of_periods_to_twist(D) 

....: assert ED.lseries().L_ratio()*ED.real_components() * etaD == md 

 

""" 

E = self._E 

self._scaling = 1 # initial value, may be changed later. 

self._failed_to_scale = False 

 

if self._sign == 1 : 

at0 = self(0) 

if at0 != 0 : 

l1 = self.__lalg__(1) 

if at0 != l1: 

verbose('scale modular symbols by %s'%(l1/at0)) 

self._scaling = l1/at0 

else : 

# if [0] = 0, we can still hope to scale it correctly by considering twists of E 

Dlist = [5,8,12,13,17,21,24,28,29, 33, 37, 40, 41, 44, 53, 56, 57, 60, 61, 65, 69, 73, 76, 77, 85, 88, 89, 92, 93, 97] # a list of positive fundamental discriminants 

j = 0 

at0 = 0 

# computes [0]+ for the twist of E by D until one value is non-zero 

while j < 30 and at0 == 0 : 

D = Dlist[j] 

# the following line checks if the twist of the newform of E by D is a newform 

# this is to avoid that we 'twist back' 

if all( valuation(E.conductor(),ell)<= valuation(D,ell) for ell in prime_divisors(D) ) : 

at0 = sum([kronecker_symbol(D,u) * self(ZZ(u)/D) for u in range(1,abs(D))]) 

j += 1 

if j == 30 and at0 == 0: # curves like "121b1", "225a1", "225e1", "256a1", "256b1", "289a1", "361a1", "400a1", "400c1", "400h1", "441b1", "441c1", "441d1", "441f1 .. will arrive here 

print("Warning : Could not normalize the modular symbols, maybe all further results will be multiplied by -1 and a power of 2") 

self._failed_to_scale = True 

else : 

l1 = self.__lalg__(D) 

if at0 != l1: 

verbose('scale modular symbols by %s found at D=%s '%(l1/at0,D), level=2) 

self._scaling = l1/at0 

 

else : # that is when sign = -1 

Dlist = [-3,-4,-7,-8,-11,-15,-19,-20,-23,-24, -31, -35, -39, -40, -43, -47, -51, -52, -55, -56, -59, -67, -68, -71, -79, -83, -84, -87, -88, -91] # a list of negative fundamental discriminants 

j = 0 

at0 = 0 

while j < 30 and at0 == 0 : 

# computes [0]+ for the twist of E by D until one value is non-zero 

D = Dlist[j] 

if all( valuation(E.conductor(),ell)<= valuation(D,ell) for ell in prime_divisors(D) ) : 

at0 = - sum([kronecker_symbol(D,u) * self(ZZ(u)/D) for u in range(1,abs(D))]) 

j += 1 

if j == 30 and at0 == 0: # no more hope for a normalization 

print("Warning : Could not normalize the modular symbols, maybe all further results will be multiplied by -1 and a power of 2") 

self._failed_to_scale = True 

else : 

l1 = self.__lalg__(D) 

if at0 != l1: 

verbose('scale modular symbols by %s'%(l1/at0)) 

self._scaling = l1/at0 

 

 

def __lalg__(self,D): 

r""" 

For positive `D`, this function evaluates the quotient 

`L(E_D,1)\cdot \sqrt(D)/\Omega_E` where `E_D` is the twist of 

`E` by `D`, `\Omega_E` is the least positive period of `E`. 

For negative `E`, it is the quotient 

`L(E_D,1)\cdot \sqrt(-D)/\Omega^{-}_E` 

where `\Omega^{-}_E` is the least positive imaginary part of a 

non-real period of `E`. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve('11a1') 

sage: m = E.modular_symbol(sign=+1, implementation='sage') 

sage: m.__lalg__(1) 

1/5 

sage: m.__lalg__(3) 

5/2 

 

""" 

from sage.functions.all import sqrt 

# the computation of the L-value could take a lot of time, 

# but then the conductor is so large 

# that the computation of modular symbols for E took even longer 

 

E = self._E 

ED = E.quadratic_twist(D) 

lv = ED.lseries().L_ratio() # this is L(ED,1) divided by the Néron period omD of ED 

lv *= ED.real_components() # now it is by the least positive period 

omD = ED.period_lattice().basis()[0] 

if D > 0 : 

om = E.period_lattice().basis()[0] 

q = sqrt(D)*omD/om * 8 

else : 

om = E.period_lattice().basis()[1].imag() 

if E.real_components() == 1: 

om *= 2 

q = sqrt(-D)*omD/om*8 

 

# see padic_lseries.pAdicLeries._quotient_of_periods_to_twist 

# for the explanation of the second factor 

verbose('real approximation is %s'%q) 

return lv/8 * QQ(int(round(q))) 

 

 

def _find_scaling_period(self): 

r""" 

Uses the integral period map of the modular symbol implementation in sage 

in order to determine the scaling. The resulting modular symbol is correct 

only for the `X_0`-optimal curve, at least up to a possible factor +- a 

power of 2. 

 

EXAMPLES:: 

 

sage: E = EllipticCurve('11a1') 

sage: m = sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolSage(E,+1,normalize='period') 

sage: m._e 

(1/5, 1) 

sage: E = EllipticCurve('11a2') 

sage: m = sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolSage(E,+1,normalize='period') 

sage: m._e 

(1, 5) 

sage: E = EllipticCurve('121b2') 

sage: m = sage.schemes.elliptic_curves.ell_modular_symbols.ModularSymbolSage(E,+1,normalize='period') 

sage: m._e 

(0, 11/2, 0, 11/2, 11/2, 0, 0, -3, 2, 1/2, -1, 3/2) 

 

TESTS:: 

 

sage: E = EllipticCurve('19a1') 

sage: m = E.modular_symbol(sign=+1, implementation='sage', normalize='none') 

sage: m._find_scaling_period() 

sage: m._scaling 

1 

 

sage: E = EllipticCurve('19a2') 

sage: m = E.modular_symbol(sign=+1, implementation='sage', normalize='none') 

sage: m._scaling 

1 

sage: m._find_scaling_period() 

sage: m._scaling 

3 

""" 

P = self._modsym.integral_period_mapping() 

self._e = P.matrix().transpose().row(0) 

self._e /= 2 

E = self._E 

try : 

crla = parse_cremona_label(E.label()) 

except RuntimeError: # raised when curve is outside of the table 

print("Warning : Could not normalize the modular symbols, maybe all further results will be multiplied by a rational number.") 

self._scaling = 1 

else : 

cr0 = Integer(crla[0]).str() + crla[1] + '1' 

E0 = EllipticCurve(cr0) 

if self._sign == 1: 

q = E0.period_lattice().basis()[0]/E.period_lattice().basis()[0] 

else: 

q = E0.period_lattice().basis()[1].imag()/E.period_lattice().basis()[1].imag() 

if E0.real_components() == 1: 

q *= 2 

if E.real_components() == 1: 

q /= 2 

q = QQ(int(round(q*200)))/200 

verbose('scale modular symbols by %s'%q) 

self._scaling = q 

c = self(0) # required, to change base point from oo to 0 

if c<0: 

c *= -1 

self._scaling *= -1 

self._at_zero = c 

self._e *= self._scaling 

 

 

def _call_with_caching(self, r): 

r""" 

Evaluates the modular symbol at `r`, caching the computed value. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol(implementation="sage") 

sage: m._call_with_caching(0) 

1/5 

""" 

try: 

return self.__cache[r] 

except AttributeError: 

self.__cache = {} 

except KeyError: 

pass 

w = self._ambient_modsym([oo,r]).element() 

c = (self._e).dot_product(w) 

self.__cache[r] = c 

return c 

 

def __call__(self, r, base_at_infinity=True): 

r""" 

Evaluates the modular symbol at {0,`r`} or {oo,`r`}. 

 

EXAMPLES:: 

 

sage: m = EllipticCurve('11a1').modular_symbol(implementation="sage") 

sage: m(0) 

1/5 

 

""" 

# this next line takes most of the time # zero = weight-2 

w = self._ambient_modsym.modular_symbol([zero, oo, Cusps(r)], check=False) 

c = (self._e).dot_product(w.element()) 

if not base_at_infinity: 

if self._at_zero is None: 

self._at_zero = self(0) 

c -= self._at_zero 

return c