Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

381

382

383

384

385

386

387

388

389

390

391

392

393

394

395

396

397

398

399

400

401

402

403

404

405

406

407

408

409

410

411

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

455

456

457

458

459

460

461

462

463

464

465

466

467

468

469

470

471

472

473

474

475

476

477

478

479

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

521

522

523

524

525

526

527

528

529

530

531

532

533

534

535

536

537

538

539

540

541

542

543

544

545

546

547

548

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

641

642

643

644

645

646

647

648

649

650

651

652

653

654

655

656

657

658

659

660

661

662

663

664

665

666

667

668

669

670

671

672

673

674

675

676

677

678

679

680

681

682

683

684

685

686

687

688

689

690

691

692

693

694

695

696

697

698

699

700

701

702

703

704

705

706

707

708

709

710

711

712

713

714

715

716

717

718

719

720

721

722

723

724

725

726

727

728

729

730

731

732

733

734

735

736

737

738

739

740

741

742

743

744

745

746

747

748

749

750

751

752

753

754

755

756

757

758

759

760

761

762

763

764

765

766

767

768

769

770

771

772

773

774

775

776

777

778

779

780

781

782

783

784

785

786

787

788

789

790

791

792

793

794

795

796

797

798

799

800

801

802

803

804

805

806

807

808

809

810

811

812

813

814

815

816

817

818

819

820

821

822

823

824

825

826

827

828

829

830

831

832

833

834

835

836

837

838

839

840

841

842

843

844

845

846

847

848

849

850

851

852

853

854

855

856

857

858

859

860

861

862

863

864

865

866

867

868

869

870

871

872

873

874

875

876

877

878

879

880

881

882

883

884

885

886

887

888

889

890

891

892

893

894

895

896

897

898

899

900

901

902

903

904

905

906

907

908

909

910

911

912

913

914

915

916

917

918

919

920

921

922

923

924

925

926

927

928

929

930

931

932

933

934

935

936

937

938

939

940

941

942

943

944

945

946

947

948

949

950

951

952

953

954

955

956

957

958

959

960

961

962

963

964

965

966

967

968

969

970

971

972

973

974

975

976

977

978

979

980

981

982

983

984

985

986

987

988

989

990

991

992

993

994

995

996

997

998

999

1000

1001

1002

1003

1004

1005

1006

1007

1008

1009

1010

1011

1012

1013

1014

1015

1016

1017

1018

1019

1020

1021

1022

1023

1024

1025

1026

1027

1028

1029

1030

1031

1032

1033

1034

1035

1036

1037

1038

1039

1040

1041

1042

1043

1044

1045

1046

1047

1048

1049

1050

1051

1052

1053

1054

1055

1056

1057

1058

1059

1060

1061

1062

1063

1064

1065

1066

1067

1068

1069

1070

1071

1072

1073

1074

1075

1076

1077

1078

1079

1080

1081

1082

1083

1084

1085

1086

1087

1088

1089

1090

1091

1092

1093

1094

1095

1096

1097

1098

1099

1100

1101

1102

1103

1104

1105

1106

1107

1108

1109

1110

1111

1112

1113

1114

1115

1116

1117

1118

1119

1120

1121

1122

1123

1124

1125

1126

1127

1128

1129

1130

1131

1132

1133

1134

1135

1136

1137

1138

1139

1140

1141

1142

1143

1144

1145

1146

1147

1148

1149

1150

1151

1152

1153

1154

1155

1156

1157

1158

1159

1160

1161

1162

1163

1164

1165

1166

1167

1168

1169

1170

1171

1172

1173

1174

1175

1176

1177

1178

1179

1180

1181

1182

1183

1184

1185

1186

1187

1188

1189

1190

1191

1192

1193

1194

1195

1196

1197

1198

1199

1200

1201

1202

1203

1204

1205

1206

1207

1208

1209

1210

1211

1212

1213

1214

1215

1216

1217

1218

1219

1220

1221

1222

1223

1224

1225

1226

1227

1228

1229

1230

1231

1232

1233

1234

1235

1236

1237

1238

1239

1240

1241

1242

1243

1244

1245

1246

1247

1248

1249

1250

1251

1252

1253

1254

1255

1256

1257

1258

1259

1260

1261

1262

1263

1264

1265

1266

1267

1268

1269

1270

1271

1272

1273

1274

1275

1276

1277

1278

1279

1280

1281

1282

1283

1284

1285

1286

1287

1288

1289

1290

1291

1292

1293

1294

1295

1296

1297

1298

1299

1300

1301

1302

1303

1304

1305

1306

1307

1308

1309

1310

1311

1312

1313

1314

1315

1316

1317

1318

1319

1320

1321

1322

1323

1324

1325

1326

1327

1328

1329

1330

1331

1332

1333

1334

1335

1336

1337

1338

1339

1340

1341

1342

""" 

Cython interface to Cremona's ``eclib`` library (also known as ``mwrank``) 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata, _mw 

sage: c = _Curvedata(1,2,3,4,5) 

  

sage: print(c) 

[1,2,3,4,5] 

b2 = 9 b4 = 11 b6 = 29 b8 = 35 

c4 = -183 c6 = -3429 

disc = -10351 (# real components = 1) 

#torsion not yet computed 

  

sage: t= _mw(c) 

sage: t.search(10) 

sage: t 

[[1:2:1]] 

""" 

from __future__ import print_function, absolute_import 

  

import os 

import sys 

  

from cysignals.memory cimport sig_free 

from cysignals.signals cimport sig_on, sig_off 

  

from sage.cpython.string cimport char_to_str, str_to_bytes 

from sage.cpython.string import FS_ENCODING 

from sage.libs.eclib cimport bigint, Curvedata, mw, two_descent 

  

cdef extern from "wrap.cpp": 

### misc functions ### 

long mwrank_get_precision() 

void mwrank_set_precision(long n) 

void mwrank_initprimes(char* pfilename, int verb) 

  

### bigint ### 

bigint* str_to_bigint(char* s) 

char* bigint_to_str(bigint* x) 

  

### Curvedata ### 

char* Curvedata_repr(Curvedata* curve) 

double Curvedata_silverman_bound(Curvedata* curve) 

double Curvedata_cps_bound(Curvedata* curve) 

double Curvedata_height_constant(Curvedata* curve) 

char* Curvedata_getdiscr(Curvedata* curve) 

char* Curvedata_conductor(Curvedata* m) 

char* Curvedata_isogeny_class(Curvedata* E, int verbose) 

  

## mw ## 

int mw_process(Curvedata* curve, mw* m, 

bigint* x, bigint* y, 

bigint* z, int sat) 

char* mw_getbasis(mw* m) 

char* mw_regulator(mw* m) 

int mw_rank(mw* m) 

int mw_saturate(mw* m, bigint* index, char** unsat, 

long sat_bd, int odd_primes_only) 

void mw_search(mw* m, char* h_lim, int moduli_option, int verb) 

  

### two_descent ### 

int two_descent_ok(two_descent* t) 

long two_descent_get_certain(two_descent* t) 

char* two_descent_get_basis(two_descent* t) 

char* two_descent_regulator(two_descent* t) 

long two_descent_get_rank(two_descent* t) 

long two_descent_get_rank_bound(two_descent* t) 

long two_descent_get_selmer_rank(two_descent* t) 

void two_descent_saturate(two_descent* t, long sat_bd) 

  

  

cdef object string_sigoff(char* s): 

sig_off() 

# Makes a python string and deletes what is pointed to by s. 

t = char_to_str(s) 

sig_free(s) 

return t 

  

# set the default 

mwrank_set_precision(50) 

  

def get_precision(): 

""" 

Returns the working floating point precision of mwrank. 

  

OUTPUT: 

  

(int) The current precision in decimal digits. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import get_precision 

sage: get_precision() 

50 

""" 

return mwrank_get_precision() 

  

def set_precision(n): 

""" 

Sets the working floating point precision of mwrank. 

  

INPUT: 

  

- ``n`` (int) -- a positive integer: the number of decimal digits. 

  

OUTPUT: 

  

None. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import set_precision 

sage: set_precision(50) 

  

""" 

mwrank_set_precision(n) 

  

def initprimes(filename, verb=False): 

""" 

Initialises mwrank/eclib's internal prime list. 

  

INPUT: 

  

- ``filename`` (string) -- the name of a file of primes. 

  

- ``verb`` (bool: default ``False``) -- verbose or not? 

  

EXAMPLES:: 

  

sage: file = os.path.join(SAGE_TMP, 'PRIMES') 

sage: with open(file, 'w') as fobj: 

....: _ = fobj.write(' '.join([str(p) for p in prime_range(10^7,10^7+20)])) 

sage: mwrank_initprimes(file, verb=True) 

Computed 78519 primes, largest is 1000253 

reading primes from file ... 

read extra prime 10000019 

finished reading primes from file ... 

Extra primes in list: 10000019 

  

sage: mwrank_initprimes("x" + file, True) 

Traceback (most recent call last): 

... 

IOError: No such file or directory: ... 

""" 

if not os.path.exists(filename): 

raise IOError('No such file or directory: %s' % filename) 

filename = str_to_bytes(filename, FS_ENCODING, 'surrogateescape') 

mwrank_initprimes(filename, verb) 

if verb: 

sys.stdout.flush() 

sys.stderr.flush() 

  

############# bigint ########################################### 

# 

# In mwrank (and eclib) bigint is synonymous with NTL's ZZ class. 

# 

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

  

cdef class _bigint: 

""" 

Cython class wrapping eclib's bigint class. 

""" 

cdef bigint* x 

  

def __init__(self, x="0"): 

""" 

Constructor for bigint class. 

  

INPUT: 

  

- ``x`` -- string or int: a string representing a decimal 

integer, or a Sage integer 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _bigint 

sage: _bigint(123) 

123 

sage: _bigint('123') 

123 

sage: type(_bigint(123)) 

<type 'sage.libs.eclib.mwrank._bigint'> 

""" 

s = str(x) 

if s.isdigit() or s[0] == "-" and s[1:].isdigit(): 

self.x = str_to_bigint(str_to_bytes(s)) 

else: 

raise ValueError("invalid _bigint: %r"%x) 

  

def __dealloc__(self): 

""" 

Destructor for bigint class (releases memory). 

""" 

del self.x 

  

def __repr__(self): 

""" 

String representation of bigint. 

  

OUTPUT: 

  

(string) the bigint as a string (decimal) 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _bigint 

sage: a = _bigint('123') 

sage: a.__repr__() 

'123' 

sage: a = _bigint('-456') 

sage: a.__repr__() 

'-456' 

""" 

sig_on() 

return string_sigoff(bigint_to_str(self.x)) 

  

  

cdef make_bigint(bigint* x): 

cdef _bigint y 

sig_off() 

y = _bigint.__new__(_bigint) 

y.x = x 

return y 

  

############# Curvedata ################# 

  

cdef class _Curvedata: # cython class wrapping eclib's Curvedata class 

cdef Curvedata* x 

  

def __init__(self, a1, a2, a3, 

a4, a6, min_on_init=0): 

""" 

Constructor for Curvedata class. 

  

INPUT: 

  

- ``a1``, ``a2``, ``a3``, ``a4``, ``a6`` (int) -- integer 

coefficients of a Weierstrass equation (must be 

nonsingular). 

  

- ``min_on_init`` (int, default 0) -- flag controlling whether 

the constructed curve is replaced by a global minimal model. 

If nonzero then this minimisation does take place. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: _Curvedata(1,2,3,4,5) 

[1,2,3,4,5] 

b2 = 9 b4 = 11 b6 = 29 b8 = 35 

c4 = -183 c6 = -3429 

disc = -10351 (# real components = 1) 

#torsion not yet computed 

  

A non-minimal example:: 

  

sage: _Curvedata(0,0,0,0,64) 

[0,0,0,0,64] 

b2 = 0 b4 = 0 b6 = 256 b8 = 0 

c4 = 0 c6 = -55296 

disc = -1769472 (# real components = 1) 

#torsion not yet computed 

  

sage: _Curvedata(0,0,0,0,64,min_on_init=1) 

[0,0,0,0,1] (reduced minimal model) 

b2 = 0 b4 = 0 b6 = 4 b8 = 0 

c4 = 0 c6 = -864 

disc = -432 (# real components = 1) 

#torsion not yet computed 

""" 

cdef _bigint _a1, _a2, _a3, _a4, _a6 

_a1 = _bigint(a1) 

_a2 = _bigint(a2) 

_a3 = _bigint(a3) 

_a4 = _bigint(a4) 

_a6 = _bigint(a6) 

self.x = new Curvedata(_a1.x[0], _a2.x[0], _a3.x[0], _a4.x[0], _a6.x[0], min_on_init) 

if self.discriminant() == 0: 

msg = "Invariants (= {},{},{},{},{}) do not describe an elliptic curve." 

raise ArithmeticError(msg.format(a1, a2, a3, a4, a6)) 

  

def __dealloc__(self): 

""" 

Destructor for Curvedata class. 

""" 

del self.x 

  

def __repr__(self): 

r""" 

String representation of Curvedata 

  

OUTPUT: 

  

(string) the Curvedata as a string 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,2,3,4,5) 

sage: E.__repr__() 

'[1,2,3,4,5]\nb2 = 9\t b4 = 11\t b6 = 29\t b8 = 35\nc4 = -183\t\tc6 = -3429\ndisc = -10351\t(# real components = 1)\n#torsion not yet computed' 

sage: E 

[1,2,3,4,5] 

b2 = 9 b4 = 11 b6 = 29 b8 = 35 

c4 = -183 c6 = -3429 

disc = -10351 (# real components = 1) 

#torsion not yet computed 

""" 

sig_on() 

return string_sigoff(Curvedata_repr(self.x))[:-1] 

  

def silverman_bound(self): 

""" 

The Silverman height bound for this elliptic curve. 

  

OUTPUT: 

  

(float) A non-negative real number `B` such that for every 

rational point on this elliptic curve `E`, `h(P)\le\hat{h}(P) 

+ B`, where `h(P)` is the naive height and `\hat{h}(P)` the 

canonical height. 

  

.. TODO:: 

  

Since eclib can compute this to arbitrary precision it would 

make sense to return a Sage real. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,2,3,4,5) 

sage: E.silverman_bound() 

6.52226179519101... 

sage: type(E.silverman_bound()) 

<type 'float'> 

""" 

return Curvedata_silverman_bound(self.x) 

  

def cps_bound(self): 

""" 

The Cremona-Prickett-Siksek height bound for this elliptic curve. 

  

OUTPUT: 

  

(float) A non-negative real number `B` such that for every 

rational point on this elliptic curve `E`, `h(P)\le\hat{h}(P) 

+ B`, where `h(P)` is the naive height and `\hat{h}(P)` the 

canonical height. 

  

.. TODO:: 

  

Since eclib can compute this to arbitrary precision it would 

make sense to return a Sage real. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,2,3,4,5) 

sage: E.cps_bound() 

0.11912451909250982 

  

Note that this is a better bound than Silverman's in this case:: 

  

sage: E.silverman_bound() 

6.52226179519101... 

""" 

cdef double x 

# We declare x so there are *no* Python library 

# calls within the sig_on()/sig_off(). 

sig_on() 

x = Curvedata_cps_bound(self.x) 

sig_off() 

return x 

  

def height_constant(self): 

""" 

A height bound for this elliptic curve. 

  

OUTPUT: 

  

(float) A non-negative real number `B` such that for every 

rational point on this elliptic curve `E`, `h(P)\le\hat{h}(P) 

+ B`, where `h(P)` is the naive height and `\hat{h}(P)` the 

canonical height. This is the minimum of the Silverman and 

Cremona_Prickett-Siksek height bounds. 

  

.. TODO:: 

  

Since eclib can compute this to arbitrary precision it would 

make sense to return a Sage real. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,2,3,4,5) 

sage: E.height_constant() 

0.11912451909250982 

""" 

return Curvedata_height_constant(self.x) 

  

def discriminant(self): 

""" 

The discriminant of this elliptic curve. 

  

OUTPUT: 

  

(Integer) The discriminant. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,2,3,4,5) 

sage: E.discriminant() 

-10351 

sage: E = _Curvedata(100,200,300,400,500) 

sage: E.discriminant() 

-1269581104000000 

sage: ZZ(E.discriminant()) 

-1269581104000000 

""" 

sig_on() 

from sage.rings.all import Integer 

return Integer(string_sigoff(Curvedata_getdiscr(self.x))) 

  

def conductor(self): 

""" 

The conductor of this elliptic curve. 

  

OUTPUT: 

  

(Integer) The conductor. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,2,3,4,5) 

sage: E.discriminant() 

-10351 

sage: E = _Curvedata(100,200,300,400,500) 

sage: E.conductor() 

126958110400 

""" 

sig_on() 

from sage.rings.all import Integer 

return Integer(string_sigoff(Curvedata_conductor(self.x))) 

  

def isogeny_class(self, verbose=False): 

""" 

The isogeny class of this elliptic curve. 

  

OUTPUT: 

  

(tuple) A tuple consisting of (1) a list of the curves in the 

isogeny class, each as a list of its Weierstrass coefficients; 

(2) a matrix of the degrees of the isogenies between the 

curves in the class (prime degrees only). 

  

.. warning:: 

  

The list may not be complete, if the precision is too low. 

Use ``mwrank_set_precision()`` to increase the precision. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,0,1,4,-6) 

sage: E.conductor() 

14 

sage: E.isogeny_class() 

([[1, 0, 1, 4, -6], [1, 0, 1, -36, -70], [1, 0, 1, -1, 0], [1, 0, 1, -171, -874], [1, 0, 1, -11, 12], [1, 0, 1, -2731, -55146]], [[0, 2, 3, 3, 0, 0], [2, 0, 0, 0, 3, 3], [3, 0, 0, 0, 2, 0], [3, 0, 0, 0, 0, 2], [0, 3, 2, 0, 0, 0], [0, 3, 0, 2, 0, 0]]) 

  

""" 

sig_on() 

s = string_sigoff(Curvedata_isogeny_class(self.x, verbose)) 

if verbose: 

sys.stdout.flush() 

sys.stderr.flush() 

return eval(s) 

  

  

############# _mw ################# 

  

cdef class _mw: 

""" 

Cython class wrapping eclib's mw class. 

""" 

cdef mw* x 

cdef Curvedata* curve 

cdef int verb 

  

def __init__(self, _Curvedata curve, verb=False, pp=1, maxr=999): 

""" 

Constructor for mw class. 

  

INPUT: 

  

- ``curve`` (_Curvedata) -- an elliptic curve 

  

- ``verb`` (bool, default ``False``) -- verbosity flag (controls 

amount of output produced in point searches) 

  

- ``pp`` (int, default 1) -- process points flag (if nonzero, 

the points found are processed, so that at all times only a 

`\ZZ`-basis for the subgroup generated by the points found 

so far is stored; if zero, no processing is done and all 

points found are stored). 

  

- ``maxr`` (int, default 999) -- maximum rank (quit point 

searching once the points found generate a subgroup of this 

rank; useful if an upper bound for the rank is already 

known). 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _mw 

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(1,0,1,4,-6) 

sage: EQ = _mw(E) 

sage: EQ 

[] 

sage: type(EQ) 

<type 'sage.libs.eclib.mwrank._mw'> 

  

sage: E = _Curvedata(0,0,1,-7,6) 

sage: EQ = _mw(E) 

sage: EQ.search(2) 

sage: EQ 

[[1:-1:1], [-2:3:1], [-14:25:8]] 

  

Example to illustrate the verbose parameter:: 

  

sage: from sage.libs.eclib.mwrank import _mw 

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(0,0,1,-7,6) 

sage: EQ = _mw(E, verb=False) 

sage: EQ.search(1) 

sage: EQ = _mw(E, verb=True) 

sage: EQ.search(1) 

P1 = [0:1:0] is torsion point, order 1 

P1 = [-3:0:1] is generator number 1 

... 

P4 = [12:35:27] = 1*P1 + -1*P2 + -1*P3 (mod torsion) 

  

The previous command produces the following output:: 

  

P1 = [0:1:0] is torsion point, order 1 

P1 = [-3:0:1] is generator number 1 

saturating up to 20...Checking 2-saturation 

Points have successfully been 2-saturated (max q used = 7) 

Checking 3-saturation 

Points have successfully been 3-saturated (max q used = 7) 

Checking 5-saturation 

Points have successfully been 5-saturated (max q used = 23) 

Checking 7-saturation 

Points have successfully been 7-saturated (max q used = 41) 

Checking 11-saturation 

Points have successfully been 11-saturated (max q used = 17) 

Checking 13-saturation 

Points have successfully been 13-saturated (max q used = 43) 

Checking 17-saturation 

Points have successfully been 17-saturated (max q used = 31) 

Checking 19-saturation 

Points have successfully been 19-saturated (max q used = 37) 

done 

P2 = [-2:3:1] is generator number 2 

saturating up to 20...Checking 2-saturation 

possible kernel vector = [1,1] 

This point may be in 2E(Q): [14:-52:1] 

...and it is! 

Replacing old generator #1 with new generator [1:-1:1] 

Points have successfully been 2-saturated (max q used = 7) 

Index gain = 2^1 

Checking 3-saturation 

Points have successfully been 3-saturated (max q used = 13) 

Checking 5-saturation 

Points have successfully been 5-saturated (max q used = 67) 

Checking 7-saturation 

Points have successfully been 7-saturated (max q used = 53) 

Checking 11-saturation 

Points have successfully been 11-saturated (max q used = 73) 

Checking 13-saturation 

Points have successfully been 13-saturated (max q used = 103) 

Checking 17-saturation 

Points have successfully been 17-saturated (max q used = 113) 

Checking 19-saturation 

Points have successfully been 19-saturated (max q used = 47) 

done (index = 2). 

Gained index 2, new generators = [ [1:-1:1] [-2:3:1] ] 

P3 = [-14:25:8] is generator number 3 

saturating up to 20...Checking 2-saturation 

Points have successfully been 2-saturated (max q used = 11) 

Checking 3-saturation 

Points have successfully been 3-saturated (max q used = 13) 

Checking 5-saturation 

Points have successfully been 5-saturated (max q used = 71) 

Checking 7-saturation 

Points have successfully been 7-saturated (max q used = 101) 

Checking 11-saturation 

Points have successfully been 11-saturated (max q used = 127) 

Checking 13-saturation 

Points have successfully been 13-saturated (max q used = 151) 

Checking 17-saturation 

Points have successfully been 17-saturated (max q used = 139) 

Checking 19-saturation 

Points have successfully been 19-saturated (max q used = 179) 

done (index = 1). 

P4 = [-1:3:1] = -1*P1 + -1*P2 + -1*P3 (mod torsion) 

P4 = [0:2:1] = 2*P1 + 0*P2 + 1*P3 (mod torsion) 

P4 = [2:13:8] = -3*P1 + 1*P2 + -1*P3 (mod torsion) 

P4 = [1:0:1] = -1*P1 + 0*P2 + 0*P3 (mod torsion) 

P4 = [2:0:1] = -1*P1 + 1*P2 + 0*P3 (mod torsion) 

P4 = [18:7:8] = -2*P1 + -1*P2 + -1*P3 (mod torsion) 

P4 = [3:3:1] = 1*P1 + 0*P2 + 1*P3 (mod torsion) 

P4 = [4:6:1] = 0*P1 + -1*P2 + -1*P3 (mod torsion) 

P4 = [36:69:64] = 1*P1 + -2*P2 + 0*P3 (mod torsion) 

P4 = [68:-25:64] = -2*P1 + -1*P2 + -2*P3 (mod torsion) 

P4 = [12:35:27] = 1*P1 + -1*P2 + -1*P3 (mod torsion) 

sage: EQ 

[[1:-1:1], [-2:3:1], [-14:25:8]] 

  

Example to illustrate the process points ``pp`` parameter:: 

  

sage: from sage.libs.eclib.mwrank import _mw 

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: E = _Curvedata(0,0,1,-7,6) 

sage: EQ = _mw(E, pp=1) 

sage: EQ.search(1); EQ 

[[1:-1:1], [-2:3:1], [-14:25:8]] 

sage: EQ = _mw(E, pp=0) 

sage: EQ.search(1); EQ 

[[-3:0:1], [-2:3:1], [-14:25:8], [-1:3:1], [0:2:1], [2:13:8], [1:0:1], [2:0:1], [18:7:8], [3:3:1], [4:6:1], [36:69:64], [68:-25:64], [12:35:27]] 

  

""" 

self.curve = curve.x 

self.x = new mw(curve.x, verb, pp, maxr) 

self.verb = verb 

  

def __dealloc__(self): 

""" 

Destructor for mw class. 

""" 

del self.x 

  

def __repr__(self): 

""" 

String representation of the current basis of this mw group. 

  

OUTPUT: 

  

(string) the current basis of this mw as a string 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: from sage.libs.eclib.mwrank import _mw 

sage: E = _Curvedata(0,0,1,-7,6) 

sage: EQ = _mw(E) 

sage: EQ # indirect doctest 

[] 

sage: EQ.search(2) 

sage: EQ 

[[1:-1:1], [-2:3:1], [-14:25:8]] 

""" 

sig_on() 

return string_sigoff(mw_getbasis(self.x)) 

  

def process(self, point, sat=0): 

""" 

Processes the given point, adding it to the mw group. 

  

INPUT: 

  

- ``point`` (tuple or list) -- tuple or list of 3 integers. 

An ``ArithmeticError`` is raised if the point is not on the 

curve. 

  

- ``sat`` (int, default 0) --saturate at primes up to ``sat``. 

No saturation is done if ``sat=0``. (Note that it is more 

efficient to add several points at once and then saturate 

just once at the end). 

  

.. NOTE:: 

  

The eclib function which implements this only carries out 

any saturation if the rank of the points increases upon 

adding the new point. This is because it is assumed that 

one saturates as ones goes along. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: from sage.libs.eclib.mwrank import _mw 

sage: E = _Curvedata(0,1,1,-2,0) 

sage: EQ = _mw(E) 

  

Initially the list of gens is empty:: 

  

sage: EQ 

[] 

  

We process a point of infinite order:: 

  

sage: EQ.process([-1,1,1]) 

sage: EQ 

[[-1:1:1]] 

  

And another independent one:: 

  

sage: EQ.process([0,-1,1]) 

sage: EQ 

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

  

Processing a point dependent on the current basis will not 

change the basis:: 

  

sage: EQ.process([4,8,1]) 

sage: EQ 

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

  

""" 

if not isinstance(point, (tuple, list)) and len(point) == 3: 

raise TypeError("point must be a list or tuple of length 3.") 

cdef _bigint x,y,z 

sig_on() 

x,y,z = _bigint(point[0]), _bigint(point[1]), _bigint(point[2]) 

r = mw_process(self.curve, self.x, x.x, y.x, z.x, sat) 

sig_off() 

if r != 0: 

raise ArithmeticError("point (=%s) not on curve." % point) 

  

def getbasis(self): 

""" 

Returns the current basis of the mw structure. 

  

OUTPUT: 

  

(string) String representation of the points in the basis of 

the mw group. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: from sage.libs.eclib.mwrank import _mw 

sage: E = _Curvedata(0,1,1,-2,0) 

sage: EQ = _mw(E) 

sage: EQ.search(3) 

sage: EQ.getbasis() 

'[[0:-1:1], [-1:1:1]]' 

sage: EQ.rank() 

2 

""" 

sig_on() 

s = string_sigoff(mw_getbasis(self.x)) 

return s 

  

def regulator(self): 

""" 

Returns the regulator of the current basis of the mw group. 

  

OUTPUT: 

  

(float) The current regulator. 

  

.. TODO:: 

  

``eclib`` computes the regulator to arbitrary precision, and 

the full precision value should be returned. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: from sage.libs.eclib.mwrank import _mw 

sage: E = _Curvedata(0,1,1,-2,0) 

sage: EQ = _mw(E) 

sage: EQ.search(3) 

sage: EQ.getbasis() 

'[[0:-1:1], [-1:1:1]]' 

sage: EQ.rank() 

2 

sage: EQ.regulator() 

0.15246017277240753 

""" 

cdef float f 

sig_on() 

f = float(string_sigoff(mw_regulator(self.x))) 

return f 

  

def rank(self): 

""" 

Returns the rank of the current basis of the mw group. 

  

OUTPUT: 

  

(Integer) The current rank. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: from sage.libs.eclib.mwrank import _mw 

sage: E = _Curvedata(0,1,1,-2,0) 

sage: EQ = _mw(E) 

sage: EQ.search(3) 

sage: EQ.getbasis() 

'[[0:-1:1], [-1:1:1]]' 

sage: EQ.rank() 

2 

""" 

sig_on() 

r = mw_rank(self.x) 

sig_off() 

from sage.rings.all import Integer 

return Integer(r) 

  

def saturate(self, int sat_bd=-1, int odd_primes_only=0): 

""" 

Saturates the current subgroup of the mw group. 

  

INPUT: 

  

- ``sat_bnd`` (int, default -1) -- bound on primes at which to 

saturate. If -1 (default), compute a bound for the primes 

which may not be saturated, and use that. 

  

- ``odd_primes_only`` (bool, default ``False``) -- only do 

saturation at odd primes. (If the points have been found 

via 2-descent they should already be 2-saturated.) 

  

OUTPUT: 

  

(tuple) (success flag, index, list) The success flag will be 1 

unless something failed (usually an indication that the points 

were not saturated but the precision is not high enough to 

divide out successfully). The index is the index of the mw 

group before saturation in the mw group after. The list is a 

string representation of the primes at which saturation was 

not proved or achieved. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: from sage.libs.eclib.mwrank import _mw 

sage: E = _Curvedata(0,1,1,-2,0) 

sage: EQ = _mw(E) 

sage: EQ.process([494, -5720, 6859]) # 3 times another point 

sage: EQ 

[[494:-5720:6859]] 

sage: EQ.saturate() 

(1, 3, '[ ]') 

sage: EQ 

[[-1:1:1]] 

  

If we set the saturation bound at 2, then saturation will fail:: 

  

sage: EQ = _mw(E) 

sage: EQ.process([494, -5720, 6859]) # 3 times another point 

sage: EQ.saturate(sat_bd=2) 

Saturation index bound = 10 

WARNING: saturation at primes p > 2 will not be done; 

points may be unsaturated at primes between 2 and index bound 

Failed to saturate MW basis at primes [ ] 

(0, 1, '[ ]') 

sage: EQ 

[[494:-5720:6859]] 

  

The following output is also seen in the preceding example:: 

  

Saturation index bound = 10 

WARNING: saturation at primes p > 2 will not be done; 

points may be unsaturated at primes between 2 and index bound 

Failed to saturate MW basis at primes [ ] 

  

  

""" 

cdef _bigint index 

cdef char* s 

cdef int ok 

sig_on() 

index = _bigint() 

ok = mw_saturate(self.x, index.x, &s, sat_bd, odd_primes_only) 

unsat = string_sigoff(s) 

return ok, index, unsat 

  

def search(self, h_lim, int moduli_option=0, int verb=0): 

""" 

Search for points in the mw group. 

  

INPUT: 

  

- ``h_lim`` (int) -- bound on logarithmic naive height of points 

  

- ``moduli_option`` (int, default 0) -- option for sieving 

strategy. The default (0) uses an adapted version of 

Stoll's ratpoints code and is recommended. 

  

- ``verb`` (int, default 0) -- level of verbosity. If 0, no 

output. If positive, the points are output as found and 

some details of the processing, finding linear relations, 

and partial saturation are output. 

  

.. NOTE:: 

  

The effect of the search is also governed by the class 

options, notably whether the points found are processed: 

meaning that linear relations are found and saturation is 

carried out, with the result that the list of generators 

will always contain a `\ZZ`-span of the saturation of the 

points found, modulo torsion. 

  

OUTPUT: 

  

None. The effect of the search is to update the list of 

generators. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: from sage.libs.eclib.mwrank import _mw 

sage: E = _Curvedata(0,0,1,-19569,-4064513) # 873c1 

sage: EQ = _mw(E) 

sage: EQ = _mw(E) 

sage: for i in [1..11]: 

....: print("{} {} {}".format(i, EQ.search(i), EQ)) 

1 None [] 

2 None [] 

3 None [] 

4 None [] 

5 None [] 

6 None [] 

7 None [] 

8 None [] 

9 None [] 

10 None [] 

11 None [[3639568:106817593:4096]] 

""" 

  

h_lim = str_to_bytes(str(h_lim)) 

  

sig_on() 

mw_search(self.x, h_lim, moduli_option, verb) 

if verb: 

sys.stdout.flush() 

sys.stderr.flush() 

sig_off() 

  

  

############# two_descent ################# 

cdef class _two_descent: 

""" 

Cython class wrapping eclib's two_descent class. 

""" 

cdef two_descent* x 

  

def __init__(self): 

""" 

Constructor for two_descent class. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

""" 

self.x = <two_descent*> 0 

  

def __dealloc__(self): 

""" 

Destructor for two_descent class. 

""" 

del self.x 

  

def do_descent(self, _Curvedata curve, 

int verb = 1, 

int sel = 0, 

int firstlim = 20, 

int secondlim = 8, 

int n_aux = -1, 

int second_descent = 1): 

""" 

Carry out a 2-descent. 

  

INPUT: 

  

- ``curvedata`` (_Curvedata) -- the curve on which to do descent. 

  

- ``verb`` (int, default 1) -- verbosity level. 

  

- ``sel`` (int, default 0) -- Selmer-only flag. If 1, only 

the 2-Selmer group will be computed, with no rational 

points. Useful as a faster way of getting an upper bound on 

the rank. 

  

- ``firstlim`` (int, default 20) -- naive height bound on 

first point search on quartic homogeneous spaces (before 

testing local solubility; very simple search with no 

overheads). 

  

- ``secondlim`` (int, default 8) -- naive height bound on 

second point search on quartic homogeneous spaces (after 

testing local solubility; sieve-assisted search) 

  

- ``n_aux`` (int, default -1) -- If positive, the number of 

auxiliary primes used in sieve-assisted search for quartics. 

If -1 (the default) use a default value (set in the eclib 

code in ``src/qrank/mrank1.cc`` in DEFAULT_NAUX: currently 8). 

Only relevant for curves with no 2-torsion, where full 

2-descent is carried out. Worth increasing for curves 

expected to be of rank > 6 to one or two more than the 

expected rank. 

  

- ``second_descent`` (int, default 1) -- flag specifying 

whether or not a second descent will be carried out (yes if 

1, the default; no if 0). Only relevant for curves with 

2-torsion. Recommended left as the default except for 

experts interested in details of Selmer groups. 

  

OUTPUT: 

  

None 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.getrank() 

3 

sage: D2.getcertain() 

1 

sage: D2.ok() 

1 

""" 

sig_on() 

self.x = new two_descent(curve.x, verb, sel, firstlim, secondlim, n_aux, second_descent) 

if verb: 

sys.stdout.flush() 

sys.stderr.flush() 

sig_off() 

  

def getrank(self): 

""" 

Returns the rank (after doing a 2-descent). 

  

OUTPUT: 

  

(Integer) the rank (or an upper bound). 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.getrank() 

3 

""" 

cdef int r 

sig_on() 

r = two_descent_get_rank(self.x) 

sig_off() 

from sage.rings.all import Integer 

return Integer(r) 

  

def getrankbound(self): 

""" 

Returns the rank upper bound (after doing a 2-descent). 

  

OUTPUT: 

  

(Integer) an upper bound on the rank. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.getrankbound() 

3 

""" 

cdef int r 

sig_on() 

r = two_descent_get_rank_bound(self.x) 

sig_off() 

from sage.rings.all import Integer 

return Integer(r) 

  

def getselmer(self): 

""" 

Returns the 2-Selmer rank (after doing a 2-descent). 

  

OUTPUT: 

  

(Integer) The 2-Selmer rank. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.getselmer() 

3 

""" 

sig_on() 

r = two_descent_get_selmer_rank(self.x) 

sig_off() 

from sage.rings.all import Integer 

return Integer(r) 

  

def ok(self): 

""" 

Returns the success flag (after doing a 2-descent). 

  

OUTPUT: 

  

(bool) Flag indicating whether or not 2-descent was successful. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.ok() 

1 

""" 

return two_descent_ok(self.x) 

  

def getcertain(self): 

""" 

Returns the certainty flag (after doing a 2-descent). 

  

OUTPUT: 

  

(bool) ``True`` if the rank upper and lower bounds are equal. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.getcertain() 

1 

""" 

return two_descent_get_certain(self.x) 

  

def saturate(self, saturation_bound=0): 

""" 

Carries out saturation of the points found by a 2-descent. 

  

OUTPUT: 

  

None. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.saturate() 

Searching for points (bound = 8)...done: 

found points which generate a subgroup of rank 3 

and regulator 0.417143558758383969817119544618093396749810106098479 

Processing points found during 2-descent...done: 

now regulator = 0.417143558758383969817119544618093396749810106098479 

No saturation being done 

sage: D2.getbasis() 

'[[1:-1:1], [-2:3:1], [-14:25:8]]' 

""" 

sig_on() 

two_descent_saturate(self.x, saturation_bound) 

sig_off() 

  

def getbasis(self): 

""" 

Returns the basis of points found by doing a 2-descent. 

  

If the success and certain flags are 1, this will be a 

`\ZZ/2\ZZ`-basis for `E(\QQ)/2E(\QQ)` (modulo torsion), 

otherwise possibly only for a proper subgroup. 

  

.. NOTE:: 

  

You must call ``saturate()`` first, or a RunTimeError will be raised. 

  

OUTPUT: 

  

(string) String representation of the list of points after 

saturation. 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

sage: D2.saturate() 

Searching for points (bound = 8)...done: 

found points which generate a subgroup of rank 3 

and regulator 0.417143558758383969817119544618093396749810106098479 

Processing points found during 2-descent...done: 

now regulator = 0.417143558758383969817119544618093396749810106098479 

No saturation being done 

sage: D2.getbasis() 

'[[1:-1:1], [-2:3:1], [-14:25:8]]' 

""" 

sig_on() 

return string_sigoff(two_descent_get_basis(self.x)) 

  

def regulator(self): 

""" 

Returns the regulator of the points found by doing a 2-descent. 

  

OUTPUT: 

  

(float) The regulator (of the subgroup found by 2-descent). 

  

EXAMPLES:: 

  

sage: from sage.libs.eclib.mwrank import _Curvedata 

sage: CD = _Curvedata(0,0,1,-7,6) 

sage: from sage.libs.eclib.mwrank import _two_descent 

sage: D2 = _two_descent() 

sage: D2.do_descent(CD) 

Basic pair: I=336, J=-10800 

disc=35092224 

... 

Mordell rank contribution from B=im(eps) = 3 

Selmer rank contribution from B=im(eps) = 3 

Sha rank contribution from B=im(eps) = 0 

Mordell rank contribution from A=ker(eps) = 0 

Selmer rank contribution from A=ker(eps) = 0 

Sha rank contribution from A=ker(eps) = 0 

  

If called before calling ``saturate()``, a bogus value of 1.0 

is returned:: 

  

sage: D2.regulator() 

1.0 

  

After saturation, both ``getbasis()`` and ``regulator()`` 

return the basis and regulator of the subgroup found by 

2-descent:: 

  

sage: D2.saturate() 

Searching for points (bound = 8)...done: 

found points which generate a subgroup of rank 3 

and regulator 0.417143558758383969817119544618093396749810106098479 

Processing points found during 2-descent...done: 

now regulator = 0.417143558758383969817119544618093396749810106098479 

No saturation being done 

sage: D2.getbasis() 

'[[1:-1:1], [-2:3:1], [-14:25:8]]' 

sage: D2.regulator() 

0.417143558758384 

""" 

sig_on() 

return float(string_sigoff(two_descent_regulator(self.x)))