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

r""" 

Enumeration of Totally Real Fields: Relative Extensions 

 

This module contains functions to enumerate primitive extensions `L / K`, where 

`K` is a given totally real number field, with given degree and small root 

discriminant. This is a relative analogue of the problem described in 

:mod:`sage.rings.number_field.totallyreal`, and we use a similar approach 

based on a relative version of Hunter's theorem. 

 

In this first simple example, we compute the totally real quadratic 

fields of `F = \QQ(\sqrt{2})` of discriminant `\le 2000`. 

 

:: 

 

sage: ZZx = ZZ['x'] 

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

sage: enumerate_totallyreal_fields_rel(F, 2, 2000) 

[[1600, x^4 - 6*x^2 + 4, xF^2 + xF - 1]] 

 

There is indeed only one such extension, given by `F(\sqrt{5})`. 

 

Next, we list all totally real quadratic extensions of `\QQ(\sqrt 5)` 

with root discriminant `\le 10`. 

 

:: 

 

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

sage: ls = enumerate_totallyreal_fields_rel(F, 2, 10^4) 

sage: ls # random (the second factor is platform-dependent) 

[[725, x^4 - x^3 - 3*x^2 + x + 1, xF^2 + (-1/2*t - 7/2)*xF + 1], 

[1125, x^4 - x^3 - 4*x^2 + 4*x + 1, xF^2 + (-1/2*t - 7/2)*xF + 1/2*t + 3/2], 

[1600, x^4 - 6*x^2 + 4, xF^2 - 2], 

[2000, x^4 - 5*x^2 + 5, xF^2 - 1/2*t - 5/2], 

[2225, x^4 - x^3 - 5*x^2 + 2*x + 4, xF^2 + (-1/2*t + 1/2)*xF - 3/2*t - 7/2], 

[2525, x^4 - 2*x^3 - 4*x^2 + 5*x + 5, xF^2 + (-1/2*t - 1/2)*xF - 1/2*t - 5/2], 

[3600, x^4 - 2*x^3 - 7*x^2 + 8*x + 1, xF^2 - 3], 

[4225, x^4 - 9*x^2 + 4, xF^2 + (-1/2*t - 1/2)*xF - 3/2*t - 9/2], 

[4400, x^4 - 7*x^2 + 11, xF^2 - 1/2*t - 7/2], 

[4525, x^4 - x^3 - 7*x^2 + 3*x + 9, xF^2 + (-1/2*t - 1/2)*xF - 3], 

[5125, x^4 - 2*x^3 - 6*x^2 + 7*x + 11, xF^2 + (-1/2*t - 1/2)*xF - t - 4], 

[5225, x^4 - x^3 - 8*x^2 + x + 11, xF^2 + (-1/2*t - 1/2)*xF - 1/2*t - 7/2], 

[5725, x^4 - x^3 - 8*x^2 + 6*x + 11, xF^2 + (-1/2*t + 1/2)*xF - 1/2*t - 7/2], 

[6125, x^4 - x^3 - 9*x^2 + 9*x + 11, xF^2 + (-1/2*t + 1/2)*xF - t - 4], 

[7225, x^4 - 11*x^2 + 9, xF^2 + (-1)*xF - 4], 

[7600, x^4 - 9*x^2 + 19, xF^2 - 1/2*t - 9/2], 

[7625, x^4 - x^3 - 9*x^2 + 4*x + 16, xF^2 + (-1/2*t - 1/2)*xF - 4], 

[8000, x^4 - 10*x^2 + 20, xF^2 - t - 5], 

[8525, x^4 - 2*x^3 - 8*x^2 + 9*x + 19, xF^2 + (-1)*xF - 1/2*t - 9/2], 

[8725, x^4 - x^3 - 10*x^2 + 2*x + 19, xF^2 + (-1/2*t - 1/2)*xF - 1/2*t - 9/2], 

[9225, x^4 - x^3 - 10*x^2 + 7*x + 19, xF^2 + (-1/2*t + 1/2)*xF - 1/2*t - 9/2]] 

sage: [ f[0] for f in ls ] 

[725, 1125, 1600, 2000, 2225, 2525, 3600, 4225, 4400, 4525, 5125, 5225, 5725, 6125, 7225, 7600, 7625, 8000, 8525, 8725, 9225] 

 

sage: [NumberField(ZZx(x[1]), 't').is_galois() for x in ls] 

[False, True, True, True, False, False, True, True, False, False, False, False, False, True, True, False, False, True, False, False, False] 

 

Eight out of 21 such fields are Galois (with Galois group `C_4` 

or `C_2 \times C_2`); the others have Galois closure of degree 8 

(with Galois group `D_8`). 

 

Finally, we compute the cubic extensions of `\QQ(\zeta_7)^+` with 

discriminant `\le 17 \times 10^9`. 

 

:: 

 

sage: F.<t> = NumberField(ZZx([1,-4,3,1])) 

sage: F.disc() 

49 

sage: enumerate_totallyreal_fields_rel(F, 3, 17*10^9) # not tested, too long time (258s on sage.math, 2013) 

[[16240385609L, x^9 - x^8 - 9*x^7 + 4*x^6 + 26*x^5 - 2*x^4 - 25*x^3 - x^2 + 7*x + 1, xF^3 + (-t^2 - 4*t + 1)*xF^2 + (t^2 + 3*t - 5)*xF + 3*t^2 + 11*t - 5]] # 32-bit 

[[16240385609, x^9 - x^8 - 9*x^7 + 4*x^6 + 26*x^5 - 2*x^4 - 25*x^3 - x^2 + 7*x + 1, xF^3 + (-t^2 - 4*t + 1)*xF^2 + (t^2 + 3*t - 5)*xF + 3*t^2 + 11*t - 5]] # 64-bit 

 

AUTHORS: 

 

- John Voight (2007-11-03): Initial version. 

""" 

 

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

# Copyright (C) 2007 William Stein and John Voight 

# 

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

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

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

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

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

from __future__ import print_function, absolute_import 

 

from sage.arith.all import binomial, gcd, divisors 

from sage.rings.integer import Integer 

from sage.rings.integer_ring import IntegerRing 

from sage.rings.number_field.totallyreal_data import ZZx, lagrange_degree_3, int_has_small_square_divisor, hermite_constant 

from sage.rings.number_field.number_field import NumberField 

from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing 

from sage.rings.number_field.totallyreal import weed_fields, odlyzko_bound_totallyreal, enumerate_totallyreal_fields_prim 

from sage.libs.pari.all import pari 

from sage.rings.all import ZZ, QQ 

 

import math 

import sys 

 

 

def integral_elements_in_box(K, C): 

r""" 

Return all integral elements of the totally real field `K` whose 

embeddings lie *numerically* within the bounds specified by the 

list `C`. The output is architecture dependent, and one may want 

to expand the bounds that define C by some epsilon. 

 

INPUT: 

 

- `K` -- a totally real number field 

- `C` -- a list [[lower, upper], ...] of lower and upper bounds, 

for each embedding 

 

 

EXAMPLES:: 

 

sage: x = polygen(QQ) 

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

sage: eps = 10e-6 

sage: C = [[0-eps,5+eps],[0-eps,10+eps]] 

sage: ls = sage.rings.number_field.totallyreal_rel.integral_elements_in_box(K, C) 

sage: sorted([ a.trace() for a in ls ]) 

[0, 2, 4, 4, 4, 6, 6, 6, 6, 8, 8, 8, 10, 10, 10, 10, 12, 12, 14] 

sage: len(ls) 

19 

 

sage: v = sage.rings.number_field.totallyreal_rel.integral_elements_in_box(K, C) 

sage: sorted(v) 

[0, -alpha + 2, 1, -alpha + 3, 2, 3, alpha + 2, 4, alpha + 3, 5, alpha + 4, 2*alpha + 3, alpha + 5, 2*alpha + 4, alpha + 6, 2*alpha + 5, 2*alpha + 6, 3*alpha + 5, 2*alpha + 7] 

 

A cubic field:: 

 

sage: x = polygen(QQ) 

sage: K.<a> = NumberField(x^3 - 16*x +16) 

sage: eps = 10e-6 

sage: C = [[0-eps,5+eps]]*3 

sage: v = sage.rings.number_field.totallyreal_rel.integral_elements_in_box(K, C) 

 

Note that the output is platform dependent (sometimes a 5 is listed 

below, and sometimes it isn't):: 

 

sage: sorted(v) 

[-1/2*a + 2, 1/4*a^2 + 1/2*a, 0, 1, 2, 3, 4,...-1/4*a^2 - 1/2*a + 5, 1/2*a + 3, -1/4*a^2 + 5] 

""" 

d = K.degree() 

Z_F = K.maximal_order() 

Foo = K.real_embeddings() 

B = K.reduced_basis() 

 

import numpy 

import numpy.linalg 

L = numpy.array([ [v(b) for b in B] for v in Foo]) 

Linv = numpy.linalg.inv(L) 

Vi = [[C[0][0]],[C[0][1]]] 

for i in range(1,d): 

Vi = sum([ [v + [C[i][0]], v + [C[i][1]]] for v in Vi], []) 

V = numpy.matrix(Linv)*(numpy.matrix(Vi).transpose()) 

j = 0 

while j < 2**d: 

for i in range(d): 

if V[i,j] < V[i,j+1]: 

V[i,j] = math.floor(V[i,j]) 

V[i,j+1] = math.ceil(V[i,j+1]) 

else: 

V[i,j] = math.ceil(V[i,j]) 

V[i,j+1] = math.floor(V[i,j+1]) 

j += 2 

W0 = (Linv*numpy.array([Vi[0]]*d)).transpose() 

W = (Linv*numpy.array([Vi[2**i] for i in range(d)])).transpose() 

for j in range(d): 

for i in range(d): 

if W[i,j] < W0[i,j]: 

W[i,j] = math.floor(W[i,j]) 

W0[i,j] = math.ceil(W0[i,j]) 

else: 

W[i,j] = math.ceil(W[i,j]) 

W0[i,j] = math.floor(W0[i,j]) 

M = [[int(V[i,j]) for i in range(V.shape[0])] for j in range(V.shape[1])] 

M += [[int(W0[i,j]) for j in range(W0.shape[0])] for i in range(W0.shape[0])] 

M += [[int(W[i,j]) for j in range(W.shape[1])] for i in range(W.shape[0])] 

 

from sage.matrix.constructor import matrix 

M = (matrix(IntegerRing(),len(M),len(M[0]), M).transpose()).columns() 

 

i = 0 

while i < len(M): 

j = i+1 

while j < len(M): 

if M[i] == M[j]: 

M.pop(j) 

else: 

j += 1 

i += 1 

 

from sage.geometry.lattice_polytope import LatticePolytope 

P = LatticePolytope(M) 

S = [] 

 

try: 

pts = P.points() 

except ValueError: 

return [] 

 

for p in pts: 

theta = sum([ p.list()[i]*B[i] for i in range(d)]) 

inbounds = True 

for i in range(d): 

inbounds = inbounds and Foo[i](theta) >= C[i][0] and Foo[i](theta) <= C[i][1] 

 

if inbounds: 

S.append(theta) 

 

return S 

 

 

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

# Main class 

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

 

eps_global = 10**(-6) 

 

class tr_data_rel: 

r""" 

This class encodes the data used in the enumeration of totally real 

fields for relative extensions. 

 

We do not give a complete description here. For more information, 

see the attached functions; all of these are used internally by the 

functions in totallyreal_rel.py, so see that file for examples and 

further documentation. 

""" 

 

def __init__(self, F, m, B, a=None): 

r""" 

Initialization routine (constructor). 

 

INPUT: 

 

- ``F`` -- number field, the base field 

- ``m`` -- integer, the relative degree 

- ``B`` -- integer, the discriminant bound 

- ``a`` -- list (default: []), the coefficient list to begin with, 

corresponding to ``a[len(a)]*x^n + ... + a[0]x^(n-len(a))``. 

 

OUTPUT: 

 

the data initialized to begin enumeration of totally real fields 

with base field F, degree n, discriminant bounded by B, and starting 

with coefficients a. 

 

EXAMPLES:: 

 

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

sage: T = sage.rings.number_field.totallyreal_rel.tr_data_rel(F, 2, 2000) 

""" 

if a is None: # don't make the stupid noob mistake of putting a=[] 

a = [] # in the function signature above. 

 

# Initialize constants. 

self.m = m 

d = F.degree() 

self.d = d 

self.n = m*d 

self.B = B 

self.gamma = hermite_constant(self.n-self.d) 

 

self.F = F 

self.Z_F = F.maximal_order() 

self.Foo = F.real_embeddings() 

self.dF = abs(F.disc()) 

self.Fx = PolynomialRing(F, 'xF') 

 

self.beta = [[]]*m 

self.gnk = [[]]*m 

 

self.trace_elts = [] 

 

Z_Fbasis = self.Z_F.basis() 

 

# Initialize variables. 

if a == []: 

# No starting input, all polynomials will be found; initialize to zero. 

self.a = [0]*m + [1] 

self.amaxvals = [[]]*m 

anm1s = [[i] for i in range(0,m//2+1)] 

for i in range(1,self.d): 

for j in range(len(anm1s)): 

anm1s[j] = [ anm1s[j] + [i] for i in range(m)] 

anm1s = sum(anm1s, []) 

anm1s = [sum([Z_Fbasis[i]*a[i] for i in range(self.d)]) for a in anm1s] 

# Minimize trace in class. 

import numpy 

for i in range(len(anm1s)): 

Q = [ [ v(m*x) for v in self.Foo] + [0] for x in Z_Fbasis] + [[v(anm1s[i]) for v in self.Foo] + [10**6]] 

pari_string = '['+';'.join([','.join(["%s"%ii for ii in row]) for row in zip(*Q)])+']' 

adj = pari(pari_string).qflll()[self.d] 

anm1s[i] += sum([m*Z_Fbasis[ii]*int(adj[ii])//int(adj[self.d]) for ii in range(self.d)]) 

 

self.amaxvals[m-1] = anm1s 

self.a[m-1] = self.amaxvals[m-1].pop() 

self.k = m-2 

 

bl = math.ceil(1.7719*self.n) 

br = max([1./m*(am1**2).trace() + \ 

self.gamma*(1./(m**d)*self.B/self.dF)**(1./(self.n-d)) for am1 in anm1s]) 

br = math.floor(br) 

T2s = self.F._positive_integral_elements_with_trace([bl,br]) 

self.trace_elts.append([bl,br,T2s]) 

 

elif len(a) <= m+1: 

# First few coefficients have been specified. 

# The value of k is the largest index of the coefficients of a which is 

# currently unknown; e.g., if k == -1, then we can iterate 

# over polynomials, and if k == n-1, then we have finished iterating. 

if a[len(a)-1] != 1: 

raise ValueError("a[len(a)-1](=%s) must be 1 so polynomial is monic"%a[len(a)-1]) 

 

raise NotImplementedError("These have not been checked.") 

 

k = m-len(a) 

self.k = k 

a = [0]*(k+1) + a 

self.amaxvals = [[]]*m 

for i in range(0,n+1): 

self.a[i] = a[i] 

 

# Bounds come from an application of Lagrange multipliers in degrees 2,3. 

self.b_lower = [-1./m*(v(self.a[m-1]) + 

(m-1.)*math.sqrt(v(self.a[m-1])**2 - 2.*(1+1./(m-1))*v(self.a[m-2]))) for v in self.Foo] 

self.b_upper = [-1./m*(v(self.a[m-1]) - 

(m-1.)*math.sqrt(v(self.a[m-1])**2 - 2.*(1+1./(m-1))*v(self.a[m-2]))) for v in self.Foo] 

if k < m-2: 

bminmax = [lagrange_degree_3(n,v(self.a[m-1]),v(self.a[m-2]),v(self.a[m-3])) for v in self.Foo] 

self.b_lower = bminmax[0] 

self.b_upper = bminmax[1] 

 

# Annoying, but must reverse coefficients for numpy. 

gnk = [binomial(j,k+2)*a[j] for j in range(k+2,n+1)] 

self.beta[k+1] = [[self.b_lower] + numpy.roots([v(gnk[i]) for i in range(len(gnk))].reverse()).tolist().sort() + [self.b_upper] for v in self.Foo] 

 

# Now to really initialize gnk. 

self.gnk[k+1] = [[0] + [binomial(j,k+1)*v(a[j]) for j in range (k+2,m+1)] for v in self.Foo] 

else: 

# Bad input! 

raise ValueError("a has length %s > m+1"%len(a)) 

 

def incr(self, f_out, verbose=False, haltk=0): 

r""" 

This function 'increments' the totally real data to the next 

value which satisfies the bounds essentially given by Rolle's 

theorem, and returns the next polynomial in the sequence 

f_out. 

 

The default or usual case just increments the constant 

coefficient; then inductively, if this is outside of the 

bounds we increment the next higher coefficient, and so on. 

 

If there are no more coefficients to be had, returns the zero 

polynomial. 

 

INPUT: 

 

- ``f_out`` -- an integer sequence, to be written with the 

coefficients of the next polynomial 

- ``verbose`` -- boolean or nonnegative integer (default: False) 

print verbosely computational details. It prints extra 

information if ``verbose`` is set to ``2`` or more 

- ``haltk`` -- integer, the level at which to halt the inductive 

coefficient bounds 

 

OUTPUT: 

 

the successor polynomial as a coefficient list. 

""" 

 

import numpy 

 

m = self.m 

n = self.n 

k = self.k 

d = self.d 

 

# If k == -1, we have a full polynomial, so we add 1 to the constant coefficient. 

if k == -1: 

if len(self.amaxvals[0]) > 0 and self.amaxvals[0]: 

self.a[0] = self.amaxvals[0].pop() 

for i in range(0,m): 

f_out[i] = self.a[i] 

return 

else: 

if verbose: 

print(" finished") 

 

# Already reached maximum, so "carry the 1" to find the next value of k. 

k += 1 

while k < m and len(self.amaxvals[k]) == 0: 

k += 1 

if k < m: 

self.a[k] = self.amaxvals[k].pop() 

k -= 1 

 

# If we are working through an initialization routine, treat that. 

elif haltk and k == haltk-1: 

if len(self.maxvals[k]) == 0: 

k += 1 

while k <= m-1 and len(self.amaxvals[k]) == 0: 

k += 1 

if k < m: 

self.a[k] = self.amaxvals[k].pop() 

k -= 1 

 

# If in the previous step we finished all possible values of 

# the lastmost coefficient, so we must compute bounds on the next coefficient. 

# Recall k == n-1 implies iteration is complete. 

while k < m-1: 

# maxoutflag flags a required abort along the way 

maxoutflag = False 

 

# Recall k == -1 means all coefficients are good to go. 

while k >= 0 and (not haltk or k >= haltk): 

if verbose: 

print(k, ":", end="") 

for i in range(self.m + 1): 

print(self.a[i], end="") 

print("") 

 

if k == m-2: 

# We only know the value of a[n-1], the trace. 

bl = max(math.ceil(1.7719*self.n), ((self.a[m-1]**2).trace()*1./m)) 

br = 1./m*(self.a[m-1]**2).trace() + \ 

self.gamma*(1./(m**d)*self.B/self.dF)**(1./(self.n-d)) 

br = math.floor(br) 

 

# Check for trivially empty. 

if bl > br: 

if verbose: 

print(" ", br, ">", bl) 

maxoutflag = 1 

break 

 

if verbose >= 2: 

print(" bl, br:", bl, br) 

 

# Enumerate all elements of Z_F with T_2 <= br 

T2s = [] 

trace_elts_found = False 

for i in range(len(self.trace_elts)): 

tre = self.trace_elts[i] 

if tre[0] <= bl and tre[1] >= br: 

trace_elts_found = True 

if verbose >= 2: 

print(" found copy!") 

for theta in tre[2]: 

if theta.trace() >= bl and theta.trace() <= br: 

T2s.append(theta) 

break 

if not trace_elts_found: 

T2s = self.F._positive_integral_elements_with_trace([bl,br]) 

self.trace_elts.append([bl,br,T2s]) 

 

# Now ensure that T2 satisfies the correct parity condition 

am2s = [] 

for t2 in T2s: 

am2 = (self.a[m-1]**2-t2)/2 

if am2.is_integral(): 

ispositive = True 

for v in self.Foo: 

ispositive = ispositive and v((m-1)*self.a[m-1]**2-2*m*am2) > 0 

if ispositive: 

am2s.append(am2) 

 

if verbose >= 2: 

print(" am2s:", am2s) 

 

# If none survive, break! 

if len(am2s) == 0: 

if verbose: 

print(" finished") 

maxoutflag = 1 

break 

 

self.amaxvals[m-2] = am2s 

self.a[m-2] = self.amaxvals[m-2].pop() 

 

# Initialize the second derivative. 

self.b_lower = [-1./m*(v(self.a[m-1]) + 

(m-1.)*math.sqrt(v(self.a[m-1])**2 - 2.*(1+1./(m-1))*v(self.a[m-2]))) for v in self.Foo] 

self.b_upper = [-1./m*(v(self.a[m-1]) - 

(m-1.)*math.sqrt(v(self.a[m-1])**2 - 2.*(1+1./(m-1))*v(self.a[m-2]))) for v in self.Foo] 

self.beta[k] = [[self.b_lower[i], -self.Foo[i](self.a[m-1])/m, self.b_upper[i]] for i in range(d)] 

self.gnk[k] = [0, (m-1)*self.a[m-1], m*(m-1)/2] 

 

if verbose >= 2: 

print(" betak:", self.beta[k]) 

else: 

# Compute the roots of the derivative. 

self.gnk[k+1][0] = self.a[k+1] 

gnk = self.gnk[k+1] 

self.beta[k] = [numpy.roots([v(gnk[len(gnk)-1-i]) for i in range(len(gnk))]).tolist() for v in self.Foo] 

 

try: 

for i in range(d): 

self.beta[k][i].sort() 

except TypeError: 

if verbose: 

print(" betak:", self.beta[k]) 

maxoutflag = True 

break 

 

# Check for double roots 

for i in range(len(self.beta[k][0])-1): 

if abs(self.beta[k][0][i] - self.beta[k][0][i+1]) < 2*eps_global: 

# This happens reasonably infrequently, so calling 

# the Python routine should be sufficiently fast... 

f = self.Fx(self.gnk[k+1]) 

df = self.Fx(self.gnk[k+2]) 

if gcd(f,df) != 1: 

if verbose: 

print(" gnk has multiple factor!") 

maxoutflag = True 

break 

if maxoutflag: 

break 

 

if k == m-3: 

self.b_lower = [-1./m*(v(self.a[m-1]) + 

(m-1.)*math.sqrt(v(self.a[m-1])**2 - 2.*(1+1./(m-1))*v(self.a[m-2]))) for v in self.Foo] 

self.b_upper = [-1./m*(v(self.a[m-1]) - 

(m-1.)*math.sqrt(v(self.a[m-1])**2 - 2.*(1+1./(m-1))*v(self.a[m-2]))) for v in self.Foo] 

elif k == m-4: 

# New bounds from Lagrange multiplier in degree 3. 

bminmax = [lagrange_degree_3(m,v(self.a[m-1]),v(self.a[m-2]),v(self.a[m-3])) for v in self.Foo] 

self.b_lower = [bminmax[i][0] for i in range(len(bminmax))] 

self.b_upper = [bminmax[i][1] for i in range(len(bminmax))] 

 

self.beta[k] = [[self.b_lower[i]] + self.beta[k][i] + [self.b_upper[i]] for i in range(len(self.beta[k]))] 

 

if verbose >= 2: 

print(" betak:", self.beta[k]) 

 

# Compute next g_(m-(k+1)), k times the formal integral of g_(m-k). 

self.gnk[k] = [self.F.primitive_element()*0] + [self.gnk[k+1][i-1]*(k+1)/i for i in range(1,m-k+1)] 

gnk = self.gnk[k] 

gnks = [ [v(gnk[len(gnk)-1-i]) for i in range(len(gnk))] for v in self.Foo ] 

gnkm1 = self.gnk[k+1] 

gnkm1s = [ [v(gnkm1[len(gnkm1)-1-i]) for i in range(len(gnkm1))] for v in self.Foo ] 

mk = m-(k+1) 

 

if verbose >= 2: 

print(" gnk:", self.gnk[k]) 

print(" gnks:", gnks) 

 

# Compute upper and lower bounds which guarantee one retains 

# a polynomial with all real roots. 

betak = self.beta[k] 

akmin = [-numpy.polyval(gnks[j], betak[j][mk+1]) - \ 

abs(numpy.polyval(gnkm1s[j], betak[j][mk+1]))*eps_global for j in range(self.d)] 

for i in range(1,(mk+1)//2+1): 

# Use the fact that f(z) <= f(x)+|f'(x)|eps if |x-z| < eps 

# for sufficiently small eps, f(z) = 0, and f''(z) < 0. 

akmin = [max(akmin[j], 

-numpy.polyval(gnks[j], betak[j][mk+1-2*i]) - \ 

abs(numpy.polyval(gnkm1s[j], betak[j][mk+1-2*i])*eps_global)) for j in range(self.d)] 

 

akmax = [-numpy.polyval(gnks[j], betak[j][mk]) + \ 

abs(numpy.polyval(gnkm1s[j], betak[j][mk]))*eps_global for j in range(self.d)] 

for i in range(1,mk//2+1): 

akmax = [min(akmax[j], 

-numpy.polyval(gnks[j], betak[j][mk-2*i]) + \ 

abs(numpy.polyval(gnkm1s[j], betak[j][mk-2*i])*eps_global)) for j in range(self.d)] 

 

if verbose >= 2: 

print(" akmin:", akmin) 

print(" akmax:", akmax) 

 

for i in range(self.d): 

if akmin[i] > akmax[i]: 

if verbose: 

print(" ", akmin[i], ">", akmax[i]) 

maxoutflag = 1 

break 

if maxoutflag: 

break 

 

self.amaxvals[k] = integral_elements_in_box(self.F, [[akmin[i],akmax[i]] for i in range(d)]) 

if k == 0: 

a0s = [0, -sum([self.a[i] for i in range(1,m+1)]), 

-sum([self.a[i]*(-1)**i for i in range(1,m+1)]), 

-sum([self.a[i]*2**i for i in range(1,m+1)]), 

-sum([self.a[i]*(-2)**i for i in range(1,m+1)])] 

for a0 in a0s: 

try: 

self.amaxvals[0].remove(a0) 

except Exception: 

pass 

 

if verbose: 

print(" amaxvals[k]:", self.amaxvals[k]) 

if len(self.amaxvals[k]) == 0: 

if verbose: 

print(" finished") 

maxoutflag = True 

break 

self.a[k] = self.amaxvals[k].pop() 

 

self.k -= 1 

k -= 1 

 

if not maxoutflag: 

self.k = k 

for i in range(m): 

f_out[i] = self.a[i] 

return 

else: 

k += 1 

while k < m and len(self.amaxvals[k]) == 0: 

k += 1 

if k < m: 

self.a[k] = self.amaxvals[k].pop() 

k -= 1 

 

# k == n-1, so iteration is complete; return the zero polynomial (of degree n+1). 

self.k = k 

f_out[m] = 0 

return 

 

 

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

# Main routine 

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

 

def enumerate_totallyreal_fields_rel(F, m, B, a = [], verbose=0, 

return_seqs=False, 

return_pari_objects=True): 

r""" 

This function enumerates (primitive) totally real field extensions of 

degree `m>1` of the totally real field F with discriminant `d \leq B`; 

optionally one can specify the first few coefficients, where the sequence ``a`` 

corresponds to a polynomial by 

 

:: 

 

a[d]*x^n + ... + a[0]*x^(n-d) 

 

if ``length(a) = d+1``, so in particular always ``a[d] = 1``. 

 

.. note:: 

 

This is guaranteed to give all primitive such fields, and 

seems in practice to give many imprimitive ones. 

 

INPUT: 

 

- ``F`` -- number field, the base field 

- ``m`` -- integer, the degree 

- ``B`` -- integer, the discriminant bound 

- ``a`` -- list (default: []), the coefficient list to begin with 

- ``verbose`` -- boolean or nonnegative integer or string (default: 0) 

give a verbose description of the computations being performed. If 

``verbose`` is set to ``2`` or more then it outputs some extra 

information. If ``verbose`` is a string then it outputs to a file 

specified by ``verbose`` 

- ``return_seqs`` -- (boolean, default False) If ``True``, then return 

the polynomials as sequences (for easier exporting to a file). This 

also returns a list of four numbers, as explained in the OUTPUT 

section below. 

- ``return_pari_objects`` -- (boolean, default: True) if 

both ``return_seqs`` and ``return_pari_objects`` are ``False`` then 

it returns the elements as Sage objects; otherwise it returns pari 

objects. 

 

OUTPUT: 

 

- the list of fields with entries ``[d,fabs,f]``, where ``d`` is the 

discriminant, ``fabs`` is an absolute defining polynomial, and ``f`` 

is a defining polynomial relative to ``F``, sorted by discriminant. 

 

- if ``return_seqs`` is ``True``, then the first field of the list is 

a list containing the count of four items as explained below 

 

- the first entry gives the number of polynomials tested 

- the second entry gives the number of polynomials with its 

discriminant having a large enough square divisor 

- the third entry is the number of irreducible polynomials 

- the fourth entry is the number of irreducible polynomials with 

discriminant at most ``B`` 

 

EXAMPLES:: 

 

sage: ZZx = ZZ['x'] 

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

sage: enumerate_totallyreal_fields_rel(F, 1, 2000) 

[[1, [-2, 0, 1], xF - 1]] 

sage: enumerate_totallyreal_fields_rel(F, 2, 2000) 

[[1600, x^4 - 6*x^2 + 4, xF^2 + xF - 1]] 

sage: enumerate_totallyreal_fields_rel(F, 2, 2000, return_seqs=True) 

[[9, 6, 5, 0], [[1600, [4, 0, -6, 0, 1], [-1, 1, 1]]]] 

 

TESTS: 

 

Each of the outputs must be elements of Sage if ``return_pari_objects`` 

is set to ``False``:: 

 

sage: type(enumerate_totallyreal_fields_rel(F, 2, 2000)[0][1]) 

<type 'cypari2.gen.Gen'> 

sage: enumerate_totallyreal_fields_rel(F, 2, 2000, return_pari_objects=False)[0][0].parent() 

Integer Ring 

sage: enumerate_totallyreal_fields_rel(F, 2, 2000, return_pari_objects=False)[0][1].parent() 

Univariate Polynomial Ring in x over Rational Field 

sage: enumerate_totallyreal_fields_rel(F, 2, 2000, return_pari_objects=False)[0][2].parent() 

Univariate Polynomial Ring in xF over Number Field in t with defining polynomial x^2 - 2 

sage: enumerate_totallyreal_fields_rel(F, 2, 2000, return_seqs=True)[1][0][1][0].parent() 

Rational Field 

 

AUTHORS: 

 

- John Voight (2007-11-01) 

""" 

 

if not isinstance(m, Integer): 

try: 

m = Integer(m) 

except TypeError: 

raise TypeError("cannot coerce m (= %s) to an integer" % m) 

if (m < 1): 

raise ValueError("m must be at least 1.") 

 

n = F.degree()*m 

 

# Initialize 

S = {} # dictionary of the form {(d, fabs): f, ...} 

dB_odlyzko = odlyzko_bound_totallyreal(n) 

dB = math.ceil(40000*dB_odlyzko**n) 

counts = [0,0,0,0] 

 

# Trivial case 

if m == 1: 

g = pari(F.defining_polynomial()).polrecip().Vec() 

if return_seqs: 

return [[0,0,0,0], [1, [-1, 1], g]] 

elif return_pari_objects: 

return [[1, g, pari('xF-1')]] 

else: 

Px = PolynomialRing(QQ, 'xF') 

return [[ZZ(1), [QQ(_) for _ in g], Px.gen()-1]] 

 

if verbose: 

saveout = sys.stdout 

if isinstance(verbose, str): 

fsock = open(verbose, 'w') 

sys.stdout = fsock 

# Else, print to screen 

f_out = [0]*m + [1] 

T = tr_data_rel(F,m,B,a) 

if verbose == 2: 

T.incr(f_out,verbose) 

else: 

T.incr(f_out) 

 

Fx = PolynomialRing(F, 'xF') 

 

nfF = pari(str(F.defining_polynomial()).replace('x', str(F.primitive_element()) ) ) 

parit = pari(str(F.primitive_element())) 

 

while f_out[m] != 0: 

counts[0] += 1 

if verbose: 

print("==>", f_out, end="") 

 

f_str = '' 

for i in range(len(f_out)): 

f_str += '(' + str(f_out[i]) + ')*x^' + str(i) 

if i < len(f_out)-1: 

f_str += '+' 

nf = pari(f_str) 

if nf.poldegree('t') == 0: 

nf = nf.subst('x', 'x-t') 

nf = nf.polresultant(nfF, parit) 

d = nf.poldisc() 

#counts[0] += 1 

if d > 0 and nf.polsturm() == n: 

da = int_has_small_square_divisor(Integer(d)) 

if d > dB or d <= B*da: 

counts[1] += 1 

if nf.polisirreducible(): 

counts[2] += 1 

[zk,d] = nf.nfbasis_d() 

 

if d <= B: 

if verbose: 

print("has discriminant", d, end="") 

 

# Find a minimal lattice element 

counts[3] += 1 

ng = pari([nf,zk]).polredabs() 

 

# Check if K is contained in the list. 

if (d, ng) in S: 

if verbose: 

print("but is not new") 

else: 

if verbose: 

print("and is new!") 

S[(d, ng)] = Fx(f_out) 

else: 

if verbose: 

print("has discriminant", abs(d), "> B") 

else: 

if verbose: 

print("is not absolutely irreducible") 

else: 

if verbose: 

print("has discriminant", abs(d), "with no large enough square divisor") 

else: 

if verbose: 

if d == 0: 

print("is not squarefree") 

else: 

print("is not totally real") 

if verbose == 2: 

T.incr(f_out,verbose=verbose) 

else: 

T.incr(f_out) 

 

# In the application of Smyth's theorem above, we exclude finitely 

# many possibilities which we must now throw back in. 

if m == 2: 

if Fx([-1,1,1]).is_irreducible(): 

K = F.extension(Fx([-1,1,1]), 'tK') 

Kabs = K.absolute_field('tKabs') 

Kabs_pari = pari(Kabs.defining_polynomial()) 

d = K.absolute_discriminant() 

if abs(d) <= B: 

ng = Kabs_pari.polredabs() 

S[(d, ng)] = Fx([-1,1,1]) 

elif F.degree() == 2: 

for ff in [[1,-7,13,-7,1],[1,-8,14,-7,1]]: 

f = Fx(ff).factor()[0][0] 

K = F.extension(f, 'tK') 

Kabs = K.absolute_field('tKabs') 

Kabs_pari = pari(Kabs.defining_polynomial()) 

d = K.absolute_discriminant() 

if abs(d) <= B: 

ng = Kabs_pari.polredabs() 

S[(d, ng)] = f 

elif m == 3: 

if Fx([-1,6,-5,1]).is_irreducible(): 

K = F.extension(Fx([-1,6,-5,1]), 'tK') 

Kabs = K.absolute_field('tKabs') 

Kabs_pari = pari(Kabs.defining_polynomial()) 

d = K.absolute_discriminant() 

if abs(d) <= B: 

ng = Kabs_pari.polredabs() 

S[(d, ng)] = Fx([-1,6,-5,1]) 

 

# Convert S to a sorted list of triples [d, fabs, f], taking care 

# to use cmp() and not the comparison operators on PARI polynomials. 

S = [[s[0], s[1], t] for s, t in S.items()] 

S.sort(key=lambda x: (x[0], [QQ(x) for x in x[1].polrecip().Vec()])) 

 

# Now check for isomorphic fields 

weed_fields(S) 

 

# Output. 

if verbose: 

print("=" * 80) 

print("Polynomials tested: {}".format(counts[0])) 

print("Polynomials with discriminant with large enough square" 

" divisor: {}".format(counts[1])) 

print("Irreducible polynomials: {}".format(counts[2])) 

print("Polynomials with nfdisc <= B: {}".format(counts[3])) 

for i in range(len(S)): 

print(S[i]) 

if isinstance(verbose, str): 

fsock.close() 

sys.stdout = saveout 

 

# Make sure to return elements that belong to Sage 

if return_seqs: 

return [[ZZ(x) for x in counts], 

[[s[0], [QQ(x) for x in s[1].polrecip().Vec()], s[2].coefficients(sparse=False)] 

for s in S] 

] 

elif return_pari_objects: 

return S 

else: 

Px = PolynomialRing(QQ, 'x') 

return [[s[0], Px([QQ(_) for _ in s[1].list()]), s[2]] for s in S] 

 

def enumerate_totallyreal_fields_all(n, B, verbose=0, return_seqs=False, 

return_pari_objects=True): 

r""" 

Enumerates *all* totally real fields of degree ``n`` with discriminant 

at most ``B``, primitive or otherwise. 

 

INPUT: 

 

- ``n`` -- integer, the degree 

- ``B`` -- integer, the discriminant bound 

- ``verbose`` -- boolean or nonnegative integer or string (default: 0) 

give a verbose description of the computations being performed. If 

``verbose`` is set to ``2`` or more then it outputs some extra 

information. If ``verbose`` is a string then it outputs to a file 

specified by ``verbose`` 

- ``return_seqs`` -- (boolean, default False) If ``True``, then return 

the polynomials as sequences (for easier exporting to a file). This 

also returns a list of four numbers, as explained in the OUTPUT 

section below. 

- ``return_pari_objects`` -- (boolean, default: True) if both 

``return_seqs`` and ``return_pari_objects`` are ``False`` then it 

returns the elements as Sage objects; otherwise it returns pari 

objects. 

 

EXAMPLES:: 

 

sage: enumerate_totallyreal_fields_all(4, 2000) 

[[725, x^4 - x^3 - 3*x^2 + x + 1], 

[1125, x^4 - x^3 - 4*x^2 + 4*x + 1], 

[1600, x^4 - 6*x^2 + 4], 

[1957, x^4 - 4*x^2 - x + 1], 

[2000, x^4 - 5*x^2 + 5]] 

sage: enumerate_totallyreal_fields_all(1, 10) 

[[1, x - 1]] 

 

TESTS: 

 

Each of the outputs must be elements of Sage if ``return_pari_objects`` 

is set to ``False``:: 

 

sage: enumerate_totallyreal_fields_all(2, 10) 

[[5, x^2 - x - 1], [8, x^2 - 2]] 

sage: type(enumerate_totallyreal_fields_all(2, 10)[0][1]) 

<type 'cypari2.gen.Gen'> 

sage: enumerate_totallyreal_fields_all(2, 10, return_pari_objects=False)[0][1].parent() 

Univariate Polynomial Ring in x over Rational Field 

 

 

In practice most of these will be found by 

:func:`~sage.rings.number_field.totallyreal.enumerate_totallyreal_fields_prim`, 

which is guaranteed to return all primitive fields but often returns 

many non-primitive ones as well. For instance, only one of the five 

fields in the example above is primitive, but 

:func:`~sage.rings.number_field.totallyreal.enumerate_totallyreal_fields_prim` 

finds four out of the five (the exception being `x^4 - 6x^2 + 4`). 

 

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

 

sage: enumerate_totallyreal_fields_all(8, 10^6) # long time (about 2 s) 

[] 

""" 

 

S = [] 

counts = [0,0,0,0] 

if len(divisors(n)) > 4: 

raise ValueError("Only implemented for n = p*q with p,q prime") 

for d in divisors(n): 

if d > 1 and d < n: 

Sds = enumerate_totallyreal_fields_prim(d, int(math.floor((1.*B)**(1.*d/n))), verbose=verbose) 

for i in range(len(Sds)): 

if verbose: 

print("=" * 80) 

print("Taking F =", Sds[i][1]) 

F = NumberField(ZZx(Sds[i][1]), 't') 

T = enumerate_totallyreal_fields_rel(F, n/d, B, verbose=verbose, return_seqs=return_seqs) 

if return_seqs: 

for i in range(4): 

counts[i] += T[0][i] 

S += [[t[0],pari(t[1]).Polrev()] for t in T[1]] 

else: 

S += [[t[0],t[1]] for t in T] 

j = i+1 

for E in enumerate_totallyreal_fields_prim(n/d, int(math.floor((1.*B)**(1./d)/(1.*Sds[i][0])**(n*1./d**2)))): 

for EF in F.composite_fields(NumberField(ZZx(E[1]), 'u')): 

if EF.degree() == n and EF.disc() <= B: 

S.append([EF.disc(), pari(EF.absolute_polynomial())]) 

S += enumerate_totallyreal_fields_prim(n, B, verbose=verbose) 

S.sort(key=lambda x: (x[0], [QQ(x) for x in x[1].polrecip().Vec()])) 

weed_fields(S) 

 

# Output. 

if verbose: 

saveout = sys.stdout 

if isinstance(verbose, str): 

fsock = open(verbose, 'w') 

sys.stdout = fsock 

# Else, print to screen 

print("=" * 80) 

print("Polynomials tested: {}".format(counts[0])) 

print("Polynomials with discriminant with large enough square" 

" divisor: {}".format(counts[1])) 

print("Irreducible polynomials: {}".format(counts[2])) 

print("Polynomials with nfdisc <= B: {}".format(counts[3])) 

for i in range(len(S)): 

print(S[i]) 

if isinstance(verbose, str): 

fsock.close() 

sys.stdout = saveout 

 

# Make sure to return elements that belong to Sage 

if return_seqs: 

return [[ZZ(_) for _ in counts], 

[[ZZ(s[0]), [QQ(_) for _ in s[1].polrecip().Vec()]] for s in S]] 

elif return_pari_objects: 

return S 

else: 

Px = PolynomialRing(QQ, 'x') 

return [[ZZ(s[0]), Px([QQ(_) for _ in s[1].list()])] 

for s in S]