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

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

r""" 

Tate's parametrisation of `p`-adic curves with multiplicative reduction 

 

Let `E` be an elliptic curve defined over the `p`-adic numbers `\QQ_p`. 

Suppose that `E` has multiplicative reduction, i.e. that the `j`-invariant 

of `E` has negative valuation, say `n`. Then there exists a parameter 

`q` in `\ZZ_p` of valuation `n` such that the points of `E` defined over 

the algebraic closure `\bar{\QQ}_p` are in bijection with 

`\bar{\QQ}_p^{\times}\,/\, q^{\ZZ}`. More precisely there exists 

the series `s_4(q)` and `s_6(q)` such that the 

`y^2+x y = x^3 + s_4(q) x+s_6(q)` curve is isomorphic to `E` over 

`\bar{\QQ}_p` (or over `\QQ_p` if the reduction is *split* multiplicative). There is `p`-adic analytic map from 

`\bar{\QQ}^{\times}_p` to this curve with kernel `q^{\ZZ}`. 

Points of good reduction correspond to points of valuation 

`0` in `\bar{\QQ}^{\times}_p`. 

 

See chapter V of [Sil1994]_ for more details. 

 

AUTHORS: 

 

- Chris Wuthrich (23/05/2007): first version 

 

- William Stein (2007-05-29): added some examples; editing. 

 

- Chris Wuthrich (04/09): reformatted docstrings. 

 

""" 

 

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

# Copyright (C) 2007 chris wuthrich 

# 

# 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 sage.rings.integer_ring import ZZ 

from sage.rings.padics.factory import Qp 

from sage.structure.sage_object import SageObject 

from sage.structure.richcmp import richcmp, richcmp_method 

from sage.arith.all import LCM 

from sage.modular.modform.constructor import EisensteinForms, CuspForms 

from sage.schemes.elliptic_curves.constructor import EllipticCurve 

from sage.functions.log import log 

from sage.misc.all import denominator, prod 

import sage.matrix.all as matrix 

 

 

@richcmp_method 

class TateCurve(SageObject): 

r""" 

Tate's `p`-adic uniformisation of an elliptic curve with 

multiplicative reduction. 

 

.. note:: 

 

Some of the methods of this Tate curve only work when the 

reduction is split multiplicative over `\QQ_p`. 

 

EXAMPLES:: 

 

sage: e = EllipticCurve('130a1') 

sage: eq = e.tate_curve(5); eq 

5-adic Tate curve associated to the Elliptic Curve defined by y^2 + x*y + y = x^3 - 33*x + 68 over Rational Field 

sage: eq == loads(dumps(eq)) 

True 

 

REFERENCES: [Sil1994]_ 

""" 

def __init__(self, E, p): 

r""" 

INPUT: 

 

- ``E`` - an elliptic curve over the rational numbers 

 

- ``p`` - a prime where `E` has multiplicative reduction, 

i.e., such that `j(E)` has negative valuation. 

 

EXAMPLES:: 

 

sage: e = EllipticCurve('130a1') 

sage: eq = e.tate_curve(2); eq 

2-adic Tate curve associated to the Elliptic Curve defined by y^2 + x*y + y = x^3 - 33*x + 68 over Rational Field 

""" 

if not p.is_prime(): 

raise ValueError("p (=%s) must be a prime" % p) 

if E.j_invariant().valuation(p) >= 0: 

raise ValueError("The elliptic curve must have multiplicative reduction at %s" % p) 

self._p = ZZ(p) 

self._E = E 

self._q = self.parameter() 

 

def __richcmp__(self, other, op): 

r""" 

Compare self and other. 

 

TESTS:: 

 

sage: E = EllipticCurve('35a') 

sage: eq5 = E.tate_curve(5) 

sage: eq7 = E.tate_curve(7) 

sage: eq7 == eq7 

True 

sage: eq7 == eq5 

False 

""" 

if type(self) != type(other): 

return NotImplemented 

 

return richcmp((self._E, self._p), (other._E, other._p), op) 

 

def _repr_(self): 

r""" 

Return print representation. 

 

EXAMPLES:: 

 

sage: e = EllipticCurve('130a1') 

sage: eq = e.tate_curve(2) 

sage: eq._repr_() 

'2-adic Tate curve associated to the Elliptic Curve defined by y^2 + x*y + y = x^3 - 33*x + 68 over Rational Field' 

""" 

return "%s-adic Tate curve associated to the %s" % (self._p, self._E) 

 

def original_curve(self): 

r""" 

Return the elliptic curve the Tate curve was constructed from. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.original_curve() 

Elliptic Curve defined by y^2 + x*y + y = x^3 - 33*x + 68 

over Rational Field 

""" 

return self._E 

 

def prime(self): 

r""" 

Return the residual characteristic `p`. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.original_curve() 

Elliptic Curve defined by y^2 + x*y + y = x^3 - 33*x + 68 

over Rational Field 

sage: eq.prime() 

5 

""" 

return self._p 

 

def parameter(self, prec=20): 

r""" 

Return the Tate parameter `q` such that the curve is isomorphic 

over the algebraic closure of `\QQ_p` to the curve 

`\QQ_p^{\times}/q^{\ZZ}`. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.parameter(prec=5) 

3*5^3 + 3*5^4 + 2*5^5 + 2*5^6 + 3*5^7 + O(5^8) 

""" 

try: 

qE = self._q 

if qE.absolute_precision() >= prec: 

return qE 

except AttributeError: 

pass 

 

E4 = EisensteinForms(weight=4).basis()[0] 

Delta = CuspForms(weight=12).basis()[0] 

j = (E4.q_expansion(prec + 3)) ** 3 / Delta.q_expansion(prec + 3) 

jinv = (1 / j).power_series() 

q_in_terms_of_jinv = jinv.reverse() 

R = Qp(self._p, prec=prec) 

qE = q_in_terms_of_jinv(R(1 / self._E.j_invariant())) 

self._q = qE 

return qE 

 

def __sk(e, k, prec): 

return sum([n ** k * e._q ** n / (1 - e._q ** n) 

for n in range(1, prec + 1)]) 

 

def __delta(e, prec): 

return e._q * prod([(1 - e._q ** n) ** 24 

for n in range(1, prec + 1)]) 

 

def curve(self, prec=20): 

r""" 

Return the `p`-adic elliptic curve of the form 

`y^2+x y = x^3 + s_4 x+s_6`. 

 

This curve with split multiplicative reduction is isomorphic 

to the given curve over the algebraic closure of `\QQ_p`. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.curve(prec=5) 

Elliptic Curve defined by y^2 + (1+O(5^5))*x*y = x^3 + 

(2*5^4+5^5+2*5^6+5^7+3*5^8+O(5^9))*x + 

(2*5^3+5^4+2*5^5+5^7+O(5^8)) over 5-adic 

Field with capped relative precision 5 

""" 

try: 

Eq = self.__curve 

if Eq.a6().absolute_precision() >= prec: 

return Eq 

except AttributeError: 

pass 

 

qE = self.parameter(prec=prec) 

n = qE.valuation() 

precp = (prec / n).floor() + 2 

R = qE.parent() 

 

tate_a4 = -5 * self.__sk(3, precp) 

tate_a6 = (tate_a4 - 7 * self.__sk(5, precp)) / 12 

Eq = EllipticCurve([R.one(), R.zero(), R.zero(), tate_a4, tate_a6]) 

self.__curve = Eq 

return Eq 

 

def _Csquare(self, prec=20): 

r""" 

Return the square of the constant `C` such that the canonical 

Neron differential `\omega` and the canonical differential 

`\frac{du}{u}` on `\QQ^{\times}/q^{\ZZ}` are linked by `\omega 

= C \frac{du}{u}`. 

 

This constant is only a square in `\QQ_p` if the curve has split 

multiplicative reduction. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq._Csquare(prec=5) 

4 + 2*5^2 + 2*5^4 + O(5^5) 

""" 

try: 

Csq = self.__Csquare 

if Csq.absolute_precision() >= prec: 

return Csq 

except AttributeError: 

pass 

 

Eq = self.curve(prec=prec) 

tateCsquare = Eq.c6() * self._E.c4() / Eq.c4() / self._E.c6() 

self.__Csquare = tateCsquare 

return tateCsquare 

 

def E2(self, prec=20): 

r""" 

Return the value of the `p`-adic Eisenstein series of weight 2 

evaluated on the elliptic curve having split multiplicative 

reduction. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.E2(prec=10) 

4 + 2*5^2 + 2*5^3 + 5^4 + 2*5^5 + 5^7 + 5^8 + 2*5^9 + O(5^10) 

 

sage: T = EllipticCurve('14').tate_curve(7) 

sage: T.E2(30) 

2 + 4*7 + 7^2 + 3*7^3 + 6*7^4 + 5*7^5 + 2*7^6 + 7^7 + 5*7^8 + 6*7^9 + 5*7^10 + 2*7^11 + 6*7^12 + 4*7^13 + 3*7^15 + 5*7^16 + 4*7^17 + 4*7^18 + 2*7^20 + 7^21 + 5*7^22 + 4*7^23 + 4*7^24 + 3*7^25 + 6*7^26 + 3*7^27 + 6*7^28 + O(7^30) 

""" 

p = self._p 

Csq = self._Csquare(prec=prec) 

qE = self._q 

n = qE.valuation() 

R = Qp(p, prec) 

e2 = Csq*(1 - 24 * sum([qE**i/(1-qE**i)**2 

for i in range(1, (prec / n).floor() + 5)])) 

return R(e2) 

 

def is_split(self): 

r""" 

Returns True if the given elliptic curve has split multiplicative reduction. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.is_split() 

True 

 

sage: eq = EllipticCurve('37a1').tate_curve(37) 

sage: eq.is_split() 

False 

""" 

return self._Csquare().is_square() 

 

def parametrisation_onto_tate_curve(self, u, prec=20): 

r""" 

Given an element `u` in `\QQ_p^{\times}`, this computes its image on the Tate curve 

under the `p`-adic uniformisation of `E`. 

 

INPUT: 

 

- ``u`` - a non-zero `p`-adic number. 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.parametrisation_onto_tate_curve(1+5+5^2+O(5^10)) 

(5^-2 + 4*5^-1 + 1 + 2*5 + 3*5^2 + 2*5^5 + 3*5^6 + O(5^7) : 

4*5^-3 + 2*5^-1 + 4 + 2*5 + 3*5^4 + 2*5^5 + O(5^6) : 1 + O(5^20)) 

""" 

if u == 1: 

return self.curve(prec=prec)(0) 

 

q = self._q 

un = u * q ** (-(u.valuation() / q.valuation()).floor()) 

 

precn = (prec / q.valuation()).floor() + 4 

 

# formulas in Silverman II (Advanced Topics in the Arithmetic 

# of Elliptic curves, p. 425) 

 

xx = un/(1-un)**2 + sum([q**n*un/(1-q**n*un)**2 + 

q**n/un/(1-q**n/un)**2-2*q**n/(1-q**n)**2 

for n in range(1, precn)]) 

 

yy = un**2/(1-un)**3 + sum([q**(2*n)*un**2/(1-q**n*un)**3 - 

q**n/un/(1-q**n/un)**3+q**n/(1-q**n)**2 

for n in range(1, precn)]) 

 

return self.curve(prec=prec)([xx, yy]) 

 

# From here on all functions need that the curve has split 

# multiplicative reduction. 

 

def L_invariant(self, prec=20): 

r""" 

Returns the *mysterious* `\mathcal{L}`-invariant associated 

to an elliptic curve with split multiplicative reduction. 

 

One 

instance where this constant appears is in the exceptional 

case of the `p`-adic Birch and Swinnerton-Dyer conjecture as 

formulated in [MTT]_. See [Col]_ for a detailed discussion. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

REFERENCES: 

 

[MTT]_ 

 

.. [Col] Pierre Colmez, Invariant `\mathcal{L}` et derivees de 

valeurs propres de Frobenius, preprint, 2004. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.L_invariant(prec=10) 

5^3 + 4*5^4 + 2*5^5 + 2*5^6 + 2*5^7 + 3*5^8 + 5^9 + O(5^10) 

""" 

if not self.is_split(): 

raise RuntimeError("The curve must have split multiplicative " 

"reduction") 

qE = self.parameter(prec=prec) 

n = qE.valuation() 

u = qE / self._p ** n 

# the p-adic logarithm of Iwasawa normalised by log(p) = 0 

return log(u) / n 

 

def _isomorphism(self, prec=20): 

r""" 

Return the isomorphism between ``self.curve()`` and the given 

curve in the form of a list ``[u,r,s,t]`` of `p`-adic numbers. 

 

For this to exist the given curve has to have split 

multiplicative reduction over `\QQ_p`. 

 

More precisely, if `E` has coordinates `x` and `y` and the Tate 

curve has coordinates `X`, `Y` with `Y^2 + XY = X^3 + s_4 X +s_6` 

then `X = u^2 x +r` and `Y = u^3 y +s u^2 x +t`. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq._isomorphism(prec=5) 

[2 + 3*5^2 + 2*5^3 + 4*5^4 + O(5^5), 

4 + 3*5 + 4*5^2 + 2*5^3 + O(5^5), 

3 + 2*5 + 5^2 + 5^3 + 2*5^4 + O(5^5), 

2 + 5 + 3*5^2 + 5^3 + 5^4 + O(5^5)] 

""" 

if not self.is_split(): 

raise RuntimeError("The curve must have split multiplicative " 

"reduction") 

C = self._Csquare(prec=prec + 4).sqrt() 

R = Qp(self._p, prec) 

C = R(C) 

s = (C * R(self._E.a1()) - R.one()) / R(2) 

r = (C ** 2 * R(self._E.a2()) + s + s ** 2) / R(3) 

t = (C ** 3 * R(self._E.a3()) - r) / R(2) 

return [C, r, s, t] 

 

def _inverse_isomorphism(self, prec=20): 

r""" 

Return the isomorphism between the given curve and 

``self.curve()`` in the form of a list ``[u,r,s,t]`` of 

`p`-adic numbers. 

 

For this to exist the given curve has to have split 

multiplicative reduction over `\QQ_p`. 

 

More precisely, if `E` has coordinates `x` and `y` and the Tate 

curve has coordinates `X`, `Y` with `Y^2 + XY = X^3 + s_4 X +s_6` 

then `x = u^2 X +r` and `y = u^3 Y +s u^2 X +t`. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq._inverse_isomorphism(prec=5) 

[3 + 2*5 + 3*5^3 + O(5^5), 4 + 2*5 + 4*5^3 + 3*5^4 + O(5^5), 

1 + 5 + 4*5^3 + 2*5^4 + O(5^5), 5 + 2*5^2 + 3*5^4 + O(5^5)] 

""" 

if not self.is_split(): 

raise RuntimeError("The curve must have split multiplicative " 

"reduction") 

u, r, s, t = self._isomorphism(prec=prec) 

return [1 / u, -r / u ** 2, -s / u, (r * s - t) / u ** 3] 

 

def lift(self, P, prec=20): 

r""" 

Given a point `P` in the formal group of the elliptic curve `E` with split multiplicative reduction, 

this produces an element `u` in `\QQ_p^{\times}` mapped to the point `P` by the Tate parametrisation. 

The algorithm return the unique such element in `1+p\ZZ_p`. 

 

INPUT: 

 

- ``P`` - a point on the elliptic curve. 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: e = EllipticCurve('130a1') 

sage: eq = e.tate_curve(5) 

sage: P = e([-6,10]) 

sage: l = eq.lift(12*P, prec=10); l 

1 + 4*5 + 5^3 + 5^4 + 4*5^5 + 5^6 + 5^7 + 4*5^8 + 5^9 + O(5^10) 

 

Now we map the lift l back and check that it is indeed right.:: 

 

sage: eq.parametrisation_onto_original_curve(l) 

(4*5^-2 + 2*5^-1 + 4*5 + 3*5^3 + 5^4 + 2*5^5 + 4*5^6 + O(5^7) : 2*5^-3 + 5^-1 + 4 + 4*5 + 5^2 + 3*5^3 + 4*5^4 + O(5^6) : 1 + O(5^20)) 

sage: e5 = e.change_ring(Qp(5,9)) 

sage: e5(12*P) 

(4*5^-2 + 2*5^-1 + 4*5 + 3*5^3 + 5^4 + 2*5^5 + 4*5^6 + O(5^7) : 2*5^-3 + 5^-1 + 4 + 4*5 + 5^2 + 3*5^3 + 4*5^4 + O(5^6) : 1 + O(5^9)) 

""" 

p = self._p 

R = Qp(self._p, prec) 

if not self._E == P.curve(): 

raise ValueError("The point must lie on the original curve.") 

if not self.is_split(): 

raise ValueError("The curve must have split multiplicative reduction.") 

if P.is_zero(): 

return R.one() 

if P[0].valuation(p) >= 0: 

raise ValueError("The point must lie in the formal group.") 

 

Eq = self.curve(prec=prec) 

C, r, s, t = self._isomorphism(prec=prec) 

xx = r + C ** 2 * P[0] 

yy = t + s * C ** 2 * P[0] + C ** 3 * P[1] 

try: 

Eq([xx, yy]) 

except Exception: 

raise RuntimeError("Bug : Point %s does not lie on the curve " % 

(xx, yy)) 

 

tt = -xx / yy 

eqhat = Eq.formal() 

eqlog = eqhat.log(prec + 3) 

z = eqlog(tt) 

u = ZZ.one() 

fac = ZZ.one() 

for i in range(1, 2 * prec + 1): 

fac *= i 

u += z ** i / fac 

return u 

 

def parametrisation_onto_original_curve(self, u, prec=20): 

r""" 

Given an element `u` in `\QQ_p^{\times}`, this computes its image on the original curve 

under the `p`-adic uniformisation of `E`. 

 

INPUT: 

 

- ``u`` - a non-zero `p`-adic number. 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.parametrisation_onto_original_curve(1+5+5^2+O(5^10)) 

(4*5^-2 + 4*5^-1 + 4 + 2*5^3 + 3*5^4 + 2*5^6 + O(5^7) : 

3*5^-3 + 5^-2 + 4*5^-1 + 1 + 4*5 + 5^2 + 3*5^5 + O(5^6) : 

1 + O(5^20)) 

 

Here is how one gets a 4-torsion point on `E` over `\QQ_5`:: 

 

sage: R = Qp(5,10) 

sage: i = R(-1).sqrt() 

sage: T = eq.parametrisation_onto_original_curve(i); T 

(2 + 3*5 + 4*5^2 + 2*5^3 + 5^4 + 4*5^5 + 2*5^7 + 5^8 + 5^9 + O(5^10) : 

3*5 + 5^2 + 5^4 + 3*5^5 + 3*5^7 + 2*5^8 + 4*5^9 + O(5^10) : 1 + O(5^20)) 

sage: 4*T 

(0 : 1 + O(5^20) : 0) 

""" 

if not self.is_split(): 

raise ValueError("The curve must have split multiplicative " 

"reduction.") 

P = self.parametrisation_onto_tate_curve(u, prec=20) 

C, r, s, t = self._inverse_isomorphism(prec=prec) 

xx = r + C ** 2 * P[0] 

yy = t + s * C ** 2 * P[0] + C ** 3 * P[1] 

R = Qp(self._p, prec) 

E_over_Qp = self._E.base_extend(R) 

return E_over_Qp([xx, yy]) 

 

def __padic_sigma_square(e, u, prec): 

return (u - 1) ** 2 / u * prod([((1-e._q**n*u)*(1-e._q**n/u) / 

(1 - e._q ** n) ** 2) ** 2 

for n in range(1, prec + 1)]) 

 

# the following functions are rather functions of the global curve 

# than the local curve 

# we use the same names as for elliptic curves over rationals. 

 

def padic_height(self, prec=20): 

r""" 

Return the canonical `p`-adic height function on the original curve. 

 

INPUT: 

 

- ``prec`` - the `p`-adic precision, default is 20. 

 

OUTPUT: 

 

- A function that can be evaluated on rational points of `E`. 

 

EXAMPLES:: 

 

sage: e = EllipticCurve('130a1') 

sage: eq = e.tate_curve(5) 

sage: h = eq.padic_height(prec=10) 

sage: P=e.gens()[0] 

sage: h(P) 

2*5^-1 + 1 + 2*5 + 2*5^2 + 3*5^3 + 3*5^6 + 5^7 + O(5^8) 

 

Check that it is a quadratic function:: 

 

sage: h(3*P)-3^2*h(P) 

O(5^8) 

""" 

if not self.is_split(): 

raise NotImplementedError("The p-adic height is not implemented for non-split multiplicative reduction.") 

 

p = self._p 

 

# we will have to do it properly with David Harvey's _multiply_point(E, R, Q) 

n = LCM(self._E.tamagawa_numbers()) * (p-1) 

 

# this function is a closure, I don't see how to doctest it (PZ) 

def _height(P, check=True): 

if check: 

assert P.curve() == self._E, "the point P must lie on the curve from which the height function was created" 

Q = n * P 

cQ = denominator(Q[0]) 

uQ = self.lift(Q, prec=prec) 

si = self.__padic_sigma_square(uQ, prec=prec) 

nn = self._q.valuation() 

qEu = self._q / p ** nn 

return -(log(si*self._Csquare()/cQ) + log(uQ)**2/log(qEu)) / n**2 

 

return _height 

 

def padic_regulator(self, prec=20): 

r""" 

Compute the canonical `p`-adic regulator on the extended Mordell-Weil group as in [MTT]_ 

(with the correction of [Wer]_ and sign convention in [SW]_.) 

 

The `p`-adic Birch and Swinnerton-Dyer conjecture predicts 

that this value appears in the formula for the leading term of 

the `p`-adic L-function. 

 

INPUT: 

 

- ``prec`` -- the `p`-adic precision, default is 20. 

 

REFERENCES: 

 

[MTT]_ 

 

.. [Wer] Annette Werner, Local heights on abelian varieties and 

rigid analytic uniformization, Doc. Math. 3 (1998), 301-319. 

 

[SW]_ 

 

EXAMPLES:: 

 

sage: eq = EllipticCurve('130a1').tate_curve(5) 

sage: eq.padic_regulator() 

2*5^-1 + 1 + 2*5 + 2*5^2 + 3*5^3 + 3*5^6 + 5^7 + 3*5^9 + 3*5^10 + 3*5^12 + 4*5^13 + 3*5^15 + 2*5^16 + 3*5^18 + 4*5^19 + O(5^20) 

 

""" 

prec = prec + 4 

 

K = Qp(self._p, prec=prec) 

rank = self._E.rank() 

if rank == 0: 

return K.one() 

 

if not self.is_split(): 

raise NotImplementedError("The p-adic regulator is not implemented for non-split multiplicative reduction.") 

 

basis = self._E.gens() 

M = matrix.matrix(K, rank, rank, 0) 

 

height = self.padic_height(prec=prec) 

point_height = [height(P) for P in basis] 

for i in range(rank): 

for j in range(i + 1, rank): 

M[i, j] = M[j, i] = (- point_height[i] - point_height[j] + height(basis[i] + basis[j])) / 2 

for i in range(rank): 

M[i, i] = point_height[i] 

 

return M.determinant()