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

r""" 

Jacobian 'morphism' as a class in the Picard group 

 

This module implements the group operation in the Picard group of a 

hyperelliptic curve, represented as divisors in Mumford 

representation, using Cantor's algorithm. 

 

A divisor on the hyperelliptic curve `y^2 + y h(x) = f(x)` 

is stored in Mumford representation, that is, as two polynomials 

`u(x)` and `v(x)` such that: 

 

- `u(x)` is monic, 

 

- `u(x)` divides `f(x) - h(x) v(x) - v(x)^2`, 

 

- `deg(v(x)) < deg(u(x)) \le g`. 

 

REFERENCES: 

 

A readable introduction to divisors, the Picard group, Mumford 

representation, and Cantor's algorithm: 

 

- J. Scholten, F. Vercauteren. An Introduction to Elliptic and 

Hyperelliptic Curve Cryptography and the NTRU Cryptosystem. To 

appear in B. Preneel (Ed.) State of the Art in Applied Cryptography 

- COSIC '03, Lecture Notes in Computer Science, Springer 2004. 

 

A standard reference in the field of cryptography: 

 

- R. Avanzi, H. Cohen, C. Doche, G. Frey, T. Lange, K. Nguyen, and F. 

Vercauteren, Handbook of Elliptic and Hyperelliptic Curve 

Cryptography. CRC Press, 2005. 

 

EXAMPLES: The following curve is the reduction of a curve whose 

Jacobian has complex multiplication. 

 

:: 

 

sage: x = GF(37)['x'].gen() 

sage: H = HyperellipticCurve(x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x); H 

Hyperelliptic Curve over Finite Field of size 37 defined 

by y^2 = x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x 

 

At this time, Jacobians of hyperelliptic curves are handled 

differently than elliptic curves:: 

 

sage: J = H.jacobian(); J 

Jacobian of Hyperelliptic Curve over Finite Field of size 37 defined 

by y^2 = x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x 

sage: J = J(J.base_ring()); J 

Set of rational points of Jacobian of Hyperelliptic Curve over Finite Field 

of size 37 defined by y^2 = x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x 

 

Points on the Jacobian are represented by Mumford's polynomials. 

First we find a couple of points on the curve:: 

 

sage: P1 = H.lift_x(2); P1 

(2 : 11 : 1) 

sage: Q1 = H.lift_x(10); Q1 

(10 : 18 : 1) 

 

Observe that 2 and 10 are the roots of the polynomials in x, 

respectively:: 

 

sage: P = J(P1); P 

(x + 35, y + 26) 

sage: Q = J(Q1); Q 

(x + 27, y + 19) 

 

:: 

 

sage: P + Q 

(x^2 + 25*x + 20, y + 13*x) 

sage: (x^2 + 25*x + 20).roots(multiplicities=False) 

[10, 2] 

 

Frobenius satisfies 

 

.. MATH:: 

 

x^4 + 12*x^3 + 78*x^2 + 444*x + 1369 

 

on the Jacobian of this reduction and the order of the Jacobian is 

`N = 1904`. 

 

:: 

 

sage: 1904*P 

(1) 

sage: 34*P == 0 

True 

sage: 35*P == P 

True 

sage: 33*P == -P 

True 

 

:: 

 

sage: Q*1904 

(1) 

sage: Q*238 == 0 

True 

sage: Q*239 == Q 

True 

sage: Q*237 == -Q 

True 

""" 

 

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

# Copyright (C) 2005 David Kohel <kohel@maths.usyd.edu.au> 

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

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

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

from __future__ import print_function 

 

from sage.misc.all import latex 

 

from sage.structure.element import AdditiveGroupElement 

from sage.structure.richcmp import richcmp, op_NE 

from sage.schemes.generic.morphism import SchemeMorphism 

 

 

def cantor_reduction_simple(a, b, f, genus): 

r""" 

Return the unique reduced divisor linearly equivalent to 

`(a, b)` on the curve `y^2 = f(x).` 

 

See the docstring of 

:mod:`sage.schemes.hyperelliptic_curves.jacobian_morphism` for 

information about divisors, linear equivalence, and reduction. 

 

EXAMPLES:: 

 

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

sage: f = x^5 - x 

sage: H = HyperellipticCurve(f); H 

Hyperelliptic Curve over Rational Field defined by y^2 = x^5 - x 

sage: J = H.jacobian()(QQ); J 

Set of rational points of Jacobian of Hyperelliptic Curve over Rational Field 

defined by y^2 = x^5 - x 

 

The following point is 2-torsion:: 

 

sage: P = J(H.lift_x(-1)); P 

(x + 1, y) 

sage: 2 * P # indirect doctest 

(1) 

""" 

a2 = (f - b**2) // a 

a2 = a2.monic() 

b2 = -b % (a2); 

if a2.degree() == a.degree(): 

# XXX 

assert a2.degree() == genus+1 

print("Returning ambiguous form of degree genus+1.") 

return (a2, b2) 

elif a2.degree() > genus: 

return cantor_reduction_simple(a2, b2, f, genus) 

return (a2, b2) 

 

def cantor_reduction(a, b, f, h, genus): 

r""" 

Return the unique reduced divisor linearly equivalent to 

`(a, b)` on the curve `y^2 + y h(x) = f(x)`. 

 

See the docstring of 

:mod:`sage.schemes.hyperelliptic_curves.jacobian_morphism` for 

information about divisors, linear equivalence, and reduction. 

 

EXAMPLES:: 

 

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

sage: f = x^5 - x 

sage: H = HyperellipticCurve(f, x); H 

Hyperelliptic Curve over Rational Field defined by y^2 + x*y = x^5 - x 

sage: J = H.jacobian()(QQ); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Rational Field defined by y^2 + x*y = x^5 - x 

 

The following point is 2-torsion:: 

 

sage: Q = J(H.lift_x(0)); Q 

(x, y) 

sage: 2*Q # indirect doctest 

(1) 

 

The next point is not 2-torsion:: 

 

sage: P = J(H.lift_x(-1)); P 

(x + 1, y - 1) 

sage: 2 * J(H.lift_x(-1)) # indirect doctest 

(x^2 + 2*x + 1, y - 3*x - 4) 

sage: 3 * J(H.lift_x(-1)) # indirect doctest 

(x^2 - 487*x - 324, y - 10754*x - 7146) 

""" 

assert a.degree() < 2*genus+1 

assert b.degree() < a.degree() 

k = f - h*b - b**2 

if 2*a.degree() == k.degree(): 

# must adjust b to include the point at infinity 

g1 = a.degree() 

x = a.parent().gen() 

r = (x**2 + h[g1]*x - f[2*g1]).roots()[0][0] 

b = b + r*(x**g1 - (x**g1) % (a)) 

k = f - h*b - b**2 

assert k % (a) == 0 

a = (k // a).monic() 

b = -(b+h) % (a) 

if a.degree() > genus: 

return cantor_reduction(a, b, f, h, genus) 

return (a, b) 

 

def cantor_composition_simple(D1,D2,f,genus): 

r""" 

Given `D_1` and `D_2` two reduced Mumford 

divisors on the Jacobian of the curve `y^2 = f(x)`, 

computes a representative `D_1 + D_2`. 

 

.. warning:: 

 

The representative computed is NOT reduced! Use 

:func:`cantor_reduction_simple` to reduce it. 

 

EXAMPLES:: 

 

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

sage: f = x^5 + x 

sage: H = HyperellipticCurve(f); H 

Hyperelliptic Curve over Rational Field defined by y^2 = x^5 + x 

 

:: 

 

sage: F.<a> = NumberField(x^2 - 2, 'a') 

sage: J = H.jacobian()(F); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Number Field in a with defining polynomial x^2 - 2 defined 

by y^2 = x^5 + x 

 

:: 

 

sage: P = J(H.lift_x(F(1))); P 

(x - 1, y - a) 

sage: Q = J(H.lift_x(F(0))); Q 

(x, y) 

sage: 2*P + 2*Q # indirect doctest 

(x^2 - 2*x + 1, y - 3/2*a*x + 1/2*a) 

sage: 2*(P + Q) # indirect doctest 

(x^2 - 2*x + 1, y - 3/2*a*x + 1/2*a) 

sage: 3*P # indirect doctest 

(x^2 - 25/32*x + 49/32, y - 45/256*a*x - 315/256*a) 

""" 

a1, b1 = D1 

a2, b2 = D2 

if a1 == a2 and b1 == b2: 

# Duplication law: 

d, h1, h3 = a1.xgcd(2*b1) 

a = (a1 // d)**2 

b = (b1 + h3*((f - b1**2) // d)) % (a) 

else: 

d0, _, h2 = a1.xgcd(a2) 

if d0 == 1: 

a = a1*a2 

b = (b2 + h2*a2*(b1-b2)) % (a) 

else: 

d, l, h3 = d0.xgcd(b1 + b2) 

a = (a1*a2) // (d**2) 

b = ((b2 + l*h2*(b1-b2)*(a2 // d)) + h3*((f - b2**2) // d)) % (a) 

a =a.monic() 

return (a, b) 

 

def cantor_composition(D1,D2,f,h,genus): 

r""" 

EXAMPLES:: 

 

sage: F.<a> = GF(7^2, 'a') 

sage: x = F['x'].gen() 

sage: f = x^7 + x^2 + a 

sage: H = HyperellipticCurve(f, 2*x); H 

Hyperelliptic Curve over Finite Field in a of size 7^2 defined by y^2 + 2*x*y = x^7 + x^2 + a 

sage: J = H.jacobian()(F); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Finite Field in a of size 7^2 defined by y^2 + 2*x*y = x^7 + x^2 + a 

 

:: 

 

sage: Q = J(H.lift_x(F(1))); Q 

(x + 6, y + 2*a + 2) 

sage: 10*Q # indirect doctest 

(x^3 + (3*a + 1)*x^2 + (2*a + 5)*x + a + 5, y + (4*a + 5)*x^2 + (a + 1)*x + 6*a + 3) 

sage: 7*8297*Q 

(1) 

 

:: 

 

sage: Q = J(H.lift_x(F(a+1))); Q 

(x + 6*a + 6, y + 2*a) 

sage: 7*8297*Q # indirect doctest 

(1) 

 

A test over a prime field: 

 

sage: F = GF(next_prime(10^30)) 

sage: x = F['x'].gen() 

sage: f = x^7 + x^2 + 1 

sage: H = HyperellipticCurve(f, 2*x); H 

Hyperelliptic Curve over Finite Field of size 1000000000000000000000000000057 defined by y^2 + 2*x*y = x^7 + x^2 + 1 

sage: J = H.jacobian()(F); J 

verbose 0 (...: multi_polynomial_ideal.py, dimension) Warning: falling back to very slow toy implementation. 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Finite Field of size 1000000000000000000000000000057 defined 

by y^2 + 2*x*y = x^7 + x^2 + 1 

sage: Q = J(H.lift_x(F(1))); Q 

(x + 1000000000000000000000000000056, y + 1000000000000000000000000000056) 

sage: 10*Q # indirect doctest 

(x^3 + 150296037169838934997145567227*x^2 + 377701248971234560956743242408*x + 509456150352486043408603286615, y + 514451014495791237681619598519*x^2 + 875375621665039398768235387900*x + 861429240012590886251910326876) 

sage: 7*8297*Q 

(x^3 + 35410976139548567549919839063*x^2 + 26230404235226464545886889960*x + 681571430588959705539385624700, y + 999722365017286747841221441793*x^2 + 262703715994522725686603955650*x + 626219823403254233972118260890) 

""" 

a1, b1 = D1 

a2, b2 = D2 

if a1 == a2 and b1 == b2: 

# Duplication law: 

d, h1, h3 = a1.xgcd(2*b1 + h) 

a = (a1 // d)**2; 

b = (b1 + h3*((f-h*b1-b1**2) // d)) % (a) 

else: 

d0, _, h2 = a1.xgcd(a2) 

if d0 == 1: 

a = a1*a2; 

b = (b2 + h2*a2*(b1-b2)) % (a) 

else: 

e0 = b1+b2+h 

if e0 == 0: 

a = (a1*a2) // (d0**2) 

b = (b2 + h2*(b1-b2)*(a2 // d0)) % (a) 

else: 

d, l, h3 = d0.xgcd(e0) 

a = (a1*a2) // (d**2) 

b = (b2 + l*h2*(b1-b2)*(a2 // d) + h3*((f-h*b2-b2**2) // d)) % (a) 

a = a.monic() 

return (a, b) 

 

class JacobianMorphism_divisor_class_field(AdditiveGroupElement, SchemeMorphism): 

r""" 

An element of a Jacobian defined over a field, i.e. in 

`J(K) = \mathrm{Pic}^0_K(C)`. 

""" 

def __init__(self, parent, polys, check=True): 

r""" 

Create a new Jacobian element in Mumford representation. 

 

INPUT: parent: the parent Homset polys: Mumford's `u` and 

`v` polynomials check (default: True): if True, ensure that 

polynomials define a divisor on the appropriate curve and are 

reduced 

 

.. warning:: 

 

Not for external use! Use ``J(K)([u, v])`` instead. 

 

EXAMPLES:: 

 

sage: x = GF(37)['x'].gen() 

sage: H = HyperellipticCurve(x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x) 

sage: J = H.jacobian()(GF(37)); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Finite Field of size 37 defined by 

y^2 = x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x 

 

:: 

 

sage: P1 = J(H.lift_x(2)); P1 # indirect doctest 

(x + 35, y + 26) 

sage: P1.parent() 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Finite Field of size 37 defined by 

y^2 = x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x 

sage: type(P1) 

<class 'sage.schemes.hyperelliptic_curves.jacobian_morphism.JacobianMorphism_divisor_class_field'> 

""" 

SchemeMorphism.__init__(self, parent) 

if check: 

C = parent.curve() 

f, h = C.hyperelliptic_polynomials() 

a, b = polys 

if not (b**2 + h*b - f)%a == 0: 

raise ValueError("Argument polys (= %s) must be divisor on curve %s."%( 

polys, C)) 

genus = C.genus() 

if a.degree() > genus: 

polys = cantor_reduction(a, b, f, h, genus) 

self.__polys = polys 

 

def _printing_polys(self): 

r""" 

Internal function formatting Mumford polynomials for printing. 

 

TESTS:: 

 

sage: F.<a> = GF(7^2, 'a') 

sage: x = F['x'].gen() 

sage: f = x^7 + x^2 + a 

sage: H = HyperellipticCurve(f, 2*x) 

sage: J = H.jacobian()(F) 

 

:: 

 

sage: Q = J(H.lift_x(F(1))); Q # indirect doctest 

(x + 6, y + 2*a + 2) 

""" 

a, b = self.__polys 

P = self.parent()._printing_ring 

y = P.gen() 

x = P.base_ring().gen() 

return (a(x), y - b(x)) 

 

def _repr_(self): 

r""" 

Return a string representation of this Mumford divisor. 

 

EXAMPLES:: 

 

sage: F.<a> = GF(7^2, 'a') 

sage: x = F['x'].gen() 

sage: f = x^7 + x^2 + a 

sage: H = HyperellipticCurve(f, 2*x) 

sage: J = H.jacobian()(F) 

 

:: 

 

sage: Q = J(0); Q # indirect doctest 

(1) 

sage: Q = J(H.lift_x(F(1))); Q # indirect doctest 

(x + 6, y + 2*a + 2) 

sage: Q + Q # indirect doctest 

(x^2 + 5*x + 1, y + 3*a*x + 6*a + 2) 

""" 

if self.is_zero(): 

return "(1)" 

a, b = self._printing_polys() 

return "(%s, %s)" % (a, b) 

 

def _latex_(self): 

r""" 

Return a LaTeX string representing this Mumford divisor. 

 

EXAMPLES:: 

 

sage: F.<alpha> = GF(7^2) 

sage: x = F['x'].gen() 

sage: f = x^7 + x^2 + alpha 

sage: H = HyperellipticCurve(f, 2*x) 

sage: J = H.jacobian()(F) 

 

:: 

 

sage: Q = J(0); print(latex(Q)) # indirect doctest 

\left(1\right) 

sage: Q = J(H.lift_x(F(1))); print(latex(Q)) # indirect doctest 

\left(x + 6, y + 2 \alpha + 2\right) 

 

:: 

 

sage: print(latex(Q + Q)) 

\left(x^{2} + 5 x + 1, y + 3 \alpha x + 6 \alpha + 2\right) 

""" 

if self.is_zero(): 

return "\\left(1\\right)" 

a, b = self._printing_polys() 

return "\\left(%s, %s\\right)" % (latex(a), latex(b)) 

 

def scheme(self): 

r""" 

Return the scheme this morphism maps to; or, where this divisor 

lives. 

 

.. warning:: 

 

Although a pointset is defined over a specific field, the 

scheme returned may be over a different (usually smaller) 

field. The example below demonstrates this: the pointset 

is determined over a number field of absolute degree 2 but 

the scheme returned is defined over the rationals. 

 

EXAMPLES:: 

 

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

sage: f = x^5 + x 

sage: H = HyperellipticCurve(f) 

sage: F.<a> = NumberField(x^2 - 2, 'a') 

sage: J = H.jacobian()(F); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Number Field in a with defining polynomial x^2 - 2 defined 

by y^2 = x^5 + x 

 

:: 

 

sage: P = J(H.lift_x(F(1))) 

sage: P.scheme() 

Jacobian of Hyperelliptic Curve over Rational Field defined by y^2 = x^5 + x 

""" 

return self.codomain() 

 

 

def __list__(self): 

r""" 

Return a list `(a(x), b(x))` of the polynomials giving the 

Mumford representation of self. 

 

TESTS:: 

 

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

sage: f = x^5 + x 

sage: H = HyperellipticCurve(f) 

sage: F.<a> = NumberField(x^2 - 2, 'a') 

sage: J = H.jacobian()(F); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Number Field in a with defining polynomial x^2 - 2 defined 

by y^2 = x^5 + x 

 

:: 

 

sage: P = J(H.lift_x(F(1))) 

sage: list(P) # indirect doctest 

[x - 1, a] 

""" 

return list(self.__polys) 

 

def __tuple__(self): 

r""" 

Return a tuple `(a(x), b(x))` of the polynomials giving the 

Mumford representation of self. 

 

TESTS:: 

 

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

sage: f = x^5 + x 

sage: H = HyperellipticCurve(f) 

sage: F.<a> = NumberField(x^2 - 2, 'a') 

sage: J = H.jacobian()(F); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Number Field in a with defining polynomial x^2 - 2 defined 

by y^2 = x^5 + x 

 

:: 

 

sage: P = J(H.lift_x(F(1))) 

sage: tuple(P) # indirect doctest 

(x - 1, a) 

""" 

return tuple(self.__polys) 

 

def __getitem__(self, n): 

r""" 

Return the `n`-th item of the pair `(a(x), b(x))` 

of polynomials giving the Mumford representation of self. 

 

TESTS:: 

 

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

sage: f = x^5 + x 

sage: H = HyperellipticCurve(f) 

sage: F.<a> = NumberField(x^2 - 2, 'a') 

sage: J = H.jacobian()(F); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Number Field in a with defining polynomial x^2 - 2 defined 

by y^2 = x^5 + x 

 

:: 

 

sage: P = J(H.lift_x(F(1))) 

sage: P[0] # indirect doctest 

x - 1 

sage: P[1] # indirect doctest 

a 

sage: P[-1] # indirect doctest 

a 

sage: P[:1] # indirect doctest 

[x - 1] 

""" 

return list(self.__polys)[n] 

 

def _richcmp_(self, other, op): 

r""" 

Compare self and other. 

 

TESTS:: 

 

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

sage: f = x^5 - x 

sage: H = HyperellipticCurve(f); H 

Hyperelliptic Curve over Rational Field defined by y^2 = x^5 - x 

sage: J = H.jacobian()(QQ); J 

Set of rational points of Jacobian of Hyperelliptic Curve over 

Rational Field defined by y^2 = x^5 - x 

 

The following point is 2-torsion:: 

 

sage: P = J(H.lift_x(-1)); P 

(x + 1, y) 

sage: 0 == 2 * P # indirect doctest 

True 

sage: P == P 

True 

 

:: 

 

sage: Q = J(H.lift_x(-1)) 

sage: Q == P 

True 

 

:: 

 

sage: 2 == Q 

False 

sage: P == False 

False 

 

Let's verify the same "points" on different schemes are not equal:: 

 

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

sage: f = x^5 + x 

sage: H2 = HyperellipticCurve(f) 

sage: J2 = H2.jacobian()(QQ) 

 

:: 

 

sage: P1 = J(H.lift_x(0)); P1 

(x, y) 

sage: P2 = J2(H2.lift_x(0)); P2 

(x, y) 

sage: P1 == P2 

False 

""" 

if self.scheme() != other.scheme(): 

return op == op_NE 

# since divisors are internally represented as Mumford divisors, 

# comparing polynomials is well-defined 

return richcmp(self.__polys, other.__polys, op) 

 

def __bool__(self): 

r""" 

Return ``True`` if this divisor is not the additive identity element. 

 

EXAMPLES:: 

 

sage: x = GF(37)['x'].gen() 

sage: H = HyperellipticCurve(x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x) 

sage: J = H.jacobian()(GF(37)) 

 

:: 

 

sage: P1 = J(H.lift_x(2)); P1 

(x + 35, y + 26) 

sage: P1 == 0 # indirect doctest 

False 

sage: P1 - P1 == 0 # indirect doctest 

True 

""" 

return self.__polys[0] != 1 

 

__nonzero__ = __bool__ 

 

def __neg__(self): 

r""" 

Return the additive inverse of this divisor. 

 

EXAMPLES:: 

 

sage: x = GF(37)['x'].gen() 

sage: H = HyperellipticCurve(x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x) 

sage: J = H.jacobian()(GF(37)) 

sage: P1 = J(H.lift_x(2)); P1 

(x + 35, y + 26) 

sage: - P1 # indirect doctest 

(x + 35, y + 11) 

sage: P1 + (-P1) # indirect doctest 

(1) 

 

:: 

 

sage: H2 = HyperellipticCurve(x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x, x) 

sage: J2 = H2.jacobian()(GF(37)) 

sage: P2 = J2(H2.lift_x(2)); P2 

(x + 35, y + 15) 

sage: - P2 # indirect doctest 

(x + 35, y + 24) 

sage: P2 + (-P2) # indirect doctest 

(1) 

 

TESTS: 

 

The following was fixed in :trac:`14264`:: 

 

sage: P.<x> = QQ[] 

sage: f = x^5 - x + 1; h = x 

sage: C = HyperellipticCurve(f,h,'u,v') 

sage: J = C.jacobian() 

sage: K.<t> = NumberField(x^2-2) 

sage: R.<x> = K[] 

sage: Q = J(K)([x^2-t,R(1)]) 

sage: Q 

(u^2 - t, v - 1) 

sage: -Q 

(u^2 - t, v + u + 1) 

sage: Q + (-Q) # indirect doctest 

(1) 

 

""" 

if self.is_zero(): 

return self 

polys = self.__polys 

X = self.parent() 

f, h = X.curve().hyperelliptic_polynomials() 

if h.is_zero(): 

D = (polys[0],-polys[1]) 

else: 

# It is essential that the modulus polys[0] can be converted into 

# the parent of the dividend h. This is not always automatically 

# the case (h can be a rational polynomial and polys[0] can a 

# non-constant polynomial over a number field different from 

# QQ). Hence, we force coercion into a common parent before 

# computing the modulus. See trac #14249 

D = (polys[0],-polys[1]-(h+polys[0]) % (polys[0])) 

return JacobianMorphism_divisor_class_field(X, D, check=False) 

 

def _add_(self,other): 

r""" 

Return a Mumford representative of the divisor self + other. 

 

EXAMPLES:: 

 

sage: x = GF(37)['x'].gen() 

sage: H = HyperellipticCurve(x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x) 

sage: J = H.jacobian()(GF(37)) 

 

:: 

 

sage: P1 = J(H.lift_x(2)); P1 

(x + 35, y + 26) 

sage: P1 + P1 # indirect doctest 

(x^2 + 33*x + 4, y + 13*x) 

""" 

X = self.parent() 

C = X.curve() 

f, h = C.hyperelliptic_polynomials() 

genus = C.genus() 

if h == 0: 

D = cantor_composition_simple(self.__polys, other.__polys, f, genus) 

if D[0].degree() > genus: 

D = cantor_reduction_simple(D[0], D[1], f, genus) 

else: 

D = cantor_composition(self.__polys, other.__polys, f, h, genus) 

if D[0].degree() > genus: 

D = cantor_reduction(D[0], D[1], f, h, genus) 

return JacobianMorphism_divisor_class_field(X, D, check=False) 

 

def _sub_(self, other): 

r""" 

Return a Mumford representative of the divisor self - other. 

 

EXAMPLES:: 

 

sage: x = GF(37)['x'].gen() 

sage: H = HyperellipticCurve(x^5 + 12*x^4 + 13*x^3 + 15*x^2 + 33*x) 

sage: J = H.jacobian()(GF(37)) 

 

:: 

 

sage: P1 = J(H.lift_x(2)); P1 

(x + 35, y + 26) 

sage: P1 - P1 # indirect doctest 

(1) 

 

:: 

 

sage: P2 = J(H.lift_x(4)); P2 

(x + 33, y + 34) 

 

Observe that the `x`-coordinates are the same but the 

`y`-coordinates differ:: 

 

sage: P1 - P2 # indirect doctest 

(x^2 + 31*x + 8, y + 7*x + 12) 

sage: P1 + P2 # indirect doctest 

(x^2 + 31*x + 8, y + 4*x + 18) 

sage: (P1 - P2) - (P1 + P2) + 2*P2 # indirect doctest 

(1) 

""" 

return self + (-other)