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

1343

1344

1345

1346

1347

1348

1349

1350

1351

1352

1353

1354

1355

1356

1357

1358

1359

1360

1361

1362

1363

1364

1365

1366

1367

1368

1369

1370

1371

1372

1373

1374

1375

1376

1377

1378

1379

1380

1381

1382

1383

1384

1385

1386

1387

1388

1389

1390

1391

1392

1393

1394

1395

1396

1397

1398

1399

1400

1401

1402

1403

1404

1405

1406

1407

1408

1409

1410

1411

1412

1413

1414

1415

1416

1417

1418

1419

1420

1421

1422

1423

1424

1425

1426

1427

1428

1429

1430

1431

1432

1433

1434

1435

1436

1437

1438

1439

1440

1441

1442

1443

1444

1445

1446

1447

1448

1449

1450

1451

1452

1453

1454

1455

1456

1457

1458

1459

1460

1461

1462

1463

1464

1465

1466

1467

1468

1469

1470

1471

1472

1473

1474

1475

1476

1477

1478

1479

1480

1481

1482

1483

1484

1485

1486

1487

1488

1489

1490

1491

1492

1493

1494

1495

1496

1497

1498

1499

1500

1501

1502

1503

1504

1505

1506

1507

1508

1509

1510

1511

1512

1513

1514

1515

1516

1517

1518

1519

1520

1521

1522

1523

1524

1525

1526

1527

1528

1529

1530

1531

1532

1533

1534

1535

1536

1537

1538

1539

1540

1541

1542

1543

1544

1545

1546

1547

1548

1549

1550

1551

1552

1553

1554

1555

1556

1557

1558

1559

1560

1561

1562

1563

1564

1565

1566

1567

1568

1569

1570

1571

1572

1573

1574

1575

1576

1577

1578

1579

1580

1581

1582

1583

1584

1585

1586

1587

1588

1589

1590

1591

1592

1593

1594

1595

1596

1597

1598

1599

1600

1601

1602

1603

1604

1605

1606

1607

1608

1609

1610

1611

1612

1613

1614

1615

1616

1617

1618

1619

1620

1621

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

r""" 

Inductive valuations on polynomial rings 

 

This module provides functionality for inductive valuations, i.e., finite 

chains of :mod:`augmented valuations <sage.rings.valuation.augmented_valuation>` on top of a :mod:`Gauss valuation <sage.rings.valuation.gauss_valuation>`. 

 

AUTHORS: 

 

- Julian Rüth (2016-11-01): initial version 

 

EXAMPLES: 

 

A :mod:`Gauss valuation <sage.rings.valuation.gauss_valuation>` is an example of an inductive valuation:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

 

Generally, an inductive valuation is an augmentation of an inductive valuation, 

i.e., a valuation that was created from a Gauss valuation in a finite number of 

augmentation steps:: 

 

sage: w = v.augmentation(x, 1) 

sage: w = w.augmentation(x^2 + 2*x + 4, 3) 

 

REFERENCES: 

 

Inductive valuations are originally discussed in [Mac1936I]_ and [Mac1936II]_. 

An introduction is also given in Chapter 4 of [Rüt2014]_. 

 

""" 

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

# Copyright (C) 2016-2017 Julian Rüth <julian.rueth@fsfe.org> 

# 

# 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 absolute_import 

 

from .valuation import DiscreteValuation, InfiniteDiscretePseudoValuation 

from .developing_valuation import DevelopingValuation 

 

from sage.misc.cachefunc import cached_method 

from sage.misc.abstract_method import abstract_method 

 

 

class InductiveValuation(DevelopingValuation): 

r""" 

Abstract base class for iterated :mod:`augmented valuations <sage.rings.valuation.augmented_valuation>` on top of a :mod:`Gauss valuation <sage.rings.valuation.gauss_valuation>`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(5)) 

 

TESTS:: 

 

sage: TestSuite(v).run() # long time 

 

""" 

def is_equivalence_unit(self, f, valuations=None): 

r""" 

Return whether the polynomial ``f`` is an equivalence unit, i.e., an 

element of :meth:`~sage.rings.valuation.developing_valuation.DevelopingValuation.effective_degree` 

zero (see [Mac1936II]_ p.497.) 

 

INPUT: 

 

- ``f`` -- a polynomial in the domain of this valuation 

 

EXAMPLES:: 

 

sage: R = Zp(2,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.is_equivalence_unit(x) 

False 

sage: v.is_equivalence_unit(S.zero()) 

False 

sage: v.is_equivalence_unit(2*x + 1) 

True 

 

""" 

f = self.domain().coerce(f) 

 

if f.is_zero(): 

return False 

return self.effective_degree(f, valuations=valuations) == 0 

 

def equivalence_reciprocal(self, f, coefficients=None, valuations=None, check=True): 

r""" 

Return an equivalence reciprocal of ``f``. 

 

An equivalence reciprocal of `f` is a polynomial `h` such that `f\cdot 

h` is equivalent to 1 modulo this valuation (see [Mac1936II]_ p.497.) 

 

INPUT: 

 

- ``f`` -- a polynomial in the domain of this valuation which is an 

:meth:`equivalence_unit` 

 

- ``coefficients`` -- the coefficients of ``f`` in the :meth:`~sage.rings.valuation.developing_valuation.DevelopingValuation.phi`-adic 

expansion if known (default: ``None``) 

 

- ``valuations`` -- the valuations of ``coefficients`` if known 

(default: ``None``) 

 

- ``check`` -- whether or not to check the validity of ``f`` (default: 

``True``) 

 

.. WARNING:: 

 

This method may not work over `p`-adic rings due to problems with 

the xgcd implementation there. 

 

EXAMPLES:: 

 

sage: R = Zp(3,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: f = 3*x + 2 

sage: h = v.equivalence_reciprocal(f); h 

(2 + O(3^5)) 

sage: v.is_equivalent(f*h, 1) 

True 

 

In an extended valuation over an extension field:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v = v.augmentation(x^2 + x + u, 1) 

sage: f = 2*x + u 

sage: h = v.equivalence_reciprocal(f); h 

(u + 1) + O(2^5) 

sage: v.is_equivalent(f*h, 1) 

True 

 

Extending the valuation once more:: 

 

sage: v = v.augmentation((x^2 + x + u)^2 + 2*x*(x^2 + x + u) + 4*x, 3) 

sage: h = v.equivalence_reciprocal(f); h 

(u + 1) + O(2^5) 

sage: v.is_equivalent(f*h, 1) 

True 

 

TESTS: 

 

A case that caused problems at some point:: 

 

sage: K = Qp(2, 4) 

sage: R.<x> = K[] 

sage: L.<a> = K.extension(x^4 + 4*x^3 + 6*x^2 + 4*x + 2) 

sage: R.<t> = L[] 

sage: v = GaussValuation(R) 

sage: w = v.augmentation(t + 1, 5/16) 

sage: w = w.augmentation(t^4 + (a^8 + a^12 + a^14 + a^16 + a^17 + a^19 + a^20 + a^23)*t^3 + (a^6 + a^9 + a^13 + a^15 + a^18 + a^19 + a^21)*t^2 + a^10*t + 1 + a^4 + a^5 + a^8 + a^13 + a^14 + a^15, 17/8) 

sage: f = a^-15*t^2 + (a^-11 + a^-9 + a^-6 + a^-5 + a^-3 + a^-2)*t + a^-15 

sage: f_ = w.equivalence_reciprocal(f) 

sage: w.reduce(f*f_) 

1 

sage: f = f.parent()([f[0], f[1].add_bigoh(1), f[2]]) 

sage: f_ = w.equivalence_reciprocal(f) 

sage: w.reduce(f*f_) 

1 

 

""" 

f = self.domain().coerce(f) 

 

if check: 

if coefficients is None: 

coefficients = list(self.coefficients(f)) 

if valuations is None: 

valuations = list(self.valuations(f, coefficients=coefficients)) 

if not self.is_equivalence_unit(f, valuations=valuations): 

raise ValueError("f must be an equivalence unit but %r is not"%(f,)) 

 

if coefficients is None: 

e0 = next(self.coefficients(f)) 

else: 

e0 = coefficients[0] 

 

# f is an equivalence unit, its valuation is given by the constant coefficient 

if valuations is None: 

vf = self(e0) 

else: 

vf = valuations[0] 

 

e0 = self.simplify(e0, error=vf) 

s_ = self.equivalence_unit(-vf) 

residue = self.reduce(e0 * s_) 

if not isinstance(self, FinalInductiveValuation): 

assert residue.is_constant() 

residue = residue[0] 

h = self.lift(~residue) * s_ 

 

h = self.simplify(h, -vf) 

 

# it might be the case that f*h has non-zero valuation because h has 

# insufficient precision, so we must not assert that here but only 

# until we lifted to higher precision 

 

# We do not actually need g*phi + h*e0 = 1, it is only important that 

# the RHS is 1 in reduction. 

# This allows us to do two things: 

# - we may lift h to arbitrary precision 

# - we can add anything which times e0 has positive valuation, e.g., we 

# may drop coefficients of positive valuation 

h = h.map_coefficients(lambda c:_lift_to_maximal_precision(c)) 

 

return h 

 

@cached_method 

def mu(self): 

r""" 

Return the valuation of :meth:`~sage.rings.valuation.developing_valuation.DevelopingValuation.phi`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: v.mu() 

0 

 

""" 

return self(self.phi()) 

 

@abstract_method 

def equivalence_unit(self, s, reciprocal=False): 

""" 

Return an equivalence unit of valuation ``s``. 

 

INPUT: 

 

- ``s`` -- an element of the :meth:`~sage.rings.valuation.valuation_space.DiscretePseudoValuationSpace.ElementMethods.value_group` 

 

- ``reciprocal`` -- a boolean (default: ``False``); whether or not to 

return the equivalence unit as the :meth:`equivalence_reciprocal` of 

the equivalence unit of valuation ``-s``. 

 

EXAMPLES:: 

 

sage: S.<x> = Qp(3,5)[] 

sage: v = GaussValuation(S) 

sage: v.equivalence_unit(2) 

(3^2 + O(3^7)) 

sage: v.equivalence_unit(-2) 

(3^-2 + O(3^3)) 

 

Note that this might fail for negative ``s`` if the domain is not 

defined over a field:: 

 

sage: v = ZZ.valuation(2) 

sage: R.<x> = ZZ[] 

sage: w = GaussValuation(R, v) 

sage: w.equivalence_unit(1) 

2 

sage: w.equivalence_unit(-1) 

Traceback (most recent call last): 

... 

ValueError: s must be in the value semigroup of this valuation but -1 is not in Additive Abelian Semigroup generated by 1 

 

""" 

 

@abstract_method 

def augmentation_chain(self): 

r""" 

Return a list with the chain of augmentations down to the underlying 

:mod:`Gauss valuation <sage.rings.valuation.gauss_valuation>`. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.augmentation_chain() 

[Gauss valuation induced by 2-adic valuation] 

 

""" 

 

@abstract_method 

def is_gauss_valuation(self): 

r""" 

Return whether this valuation is a Gauss valuation over the domain. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.is_gauss_valuation() 

True 

 

""" 

 

@abstract_method 

def E(self): 

""" 

Return the ramification index of this valuation over its underlying 

Gauss valuation. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.E() 

1 

 

""" 

 

@abstract_method 

def F(self): 

""" 

Return the residual degree of this valuation over its Gauss extension. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.F() 

1 

 

""" 

 

@abstract_method 

def monic_integral_model(self, G): 

r""" 

Return a monic integral irreducible polynomial which defines the same 

extension of the base ring of the domain as the irreducible polynomial 

``G`` together with maps between the old and the new polynomial. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: v.monic_integral_model(5*x^2 + 1/2*x + 1/4) 

(Ring endomorphism of Univariate Polynomial Ring in x over Rational Field 

Defn: x |--> 1/2*x, 

Ring endomorphism of Univariate Polynomial Ring in x over Rational Field 

Defn: x |--> 2*x, 

x^2 + 1/5*x + 1/5) 

 

""" 

 

@abstract_method 

def element_with_valuation(self, s): 

r""" 

Return a polynomial of minimal degree with valuation ``s``. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: v.element_with_valuation(-2) 

1/4 

 

Depending on the base ring, an element of valuation ``s`` might not 

exist:: 

 

sage: R.<x> = ZZ[] 

sage: v = GaussValuation(R, ZZ.valuation(2)) 

sage: v.element_with_valuation(-2) 

Traceback (most recent call last): 

... 

ValueError: s must be in the value semigroup of this valuation but -2 is not in Additive Abelian Semigroup generated by 1 

 

""" 

 

def _test_element_with_valuation_inductive_valuation(self, **options): 

r""" 

Test the correctness of :meth:`element_with_valuation`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: v._test_element_with_valuation_inductive_valuation() 

 

""" 

tester = self._tester(**options) 

chain = self.augmentation_chain() 

for s in tester.some_elements(self.value_group().some_elements()): 

try: 

R = self.element_with_valuation(s) 

except (ValueError, NotImplementedError): 

# this is often not possible unless the underlying ring of 

# constants is a field 

from sage.categories.fields import Fields 

if self.domain().base() not in Fields(): 

continue 

raise 

tester.assertEqual(self(R), s) 

if chain != [self]: 

base = chain[1] 

if s in base.value_group(): 

S = base.element_with_valuation(s) 

tester.assertEqual(self(S), s) 

tester.assertGreaterEqual(S.degree(), R.degree()) 

 

def _test_EF(self, **options): 

r""" 

Test the correctness of :meth:`E` and :meth:`F`. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v._test_EF() 

 

""" 

tester = self._tester(**options) 

chain = self.augmentation_chain() 

for w,v in zip(chain, chain[1:]): 

from sage.rings.all import infinity, ZZ 

if w(w.phi()) is infinity: 

tester.assertEqual(w.E(), v.E()) 

tester.assertIn(w.E(), ZZ) 

tester.assertIn(w.F(), ZZ) 

tester.assertGreaterEqual(w.E(), v.E()) 

tester.assertGreaterEqual(w.F(), v.F()) 

 

def _test_augmentation_chain(self, **options): 

r""" 

Test the correctness of :meth:`augmentation_chain`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: v._test_augmentation_chain() 

 

""" 

tester = self._tester(**options) 

chain = self.augmentation_chain() 

tester.assertIs(chain[0], self) 

tester.assertTrue(chain[-1].is_gauss_valuation()) 

for w,v in zip(chain, chain[1:]): 

tester.assertGreaterEqual(w, v) 

 

def _test_equivalence_unit(self, **options): 

r""" 

Test the correctness of :meth:`lift_to_key`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: v._test_equivalence_unit() 

 

""" 

tester = self._tester(**options) 

 

if self.is_gauss_valuation(): 

value_group = self.value_group() 

else: 

value_group = self.augmentation_chain()[1].value_group() 

 

for s in tester.some_elements(value_group.some_elements()): 

try: 

R = self.equivalence_unit(s) 

except (ValueError, NotImplementedError): 

# this is often not possible unless the underlying ring of 

# constants is a field 

from sage.categories.fields import Fields 

if self.domain().base() not in Fields(): 

continue 

raise 

tester.assertIs(R.parent(), self.domain()) 

tester.assertEqual(self(R), s) 

tester.assertTrue(self.is_equivalence_unit(R)) 

 

def _test_is_equivalence_unit(self, **options): 

r""" 

Test the correctness of :meth:`is_equivalence_unit`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: v._test_is_equivalence_unit() 

 

""" 

tester = self._tester(**options) 

tester.assertFalse(self.is_equivalence_unit(self.phi())) 

 

def _test_equivalence_reciprocal(self, **options): 

r""" 

Test the correctness of :meth:`equivalence_reciprocal`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: v._test_equivalence_reciprocal() 

 

""" 

tester = self._tester(**options) 

S = tester.some_elements(self.domain().some_elements()) 

for f in S: 

if self.is_equivalence_unit(f): 

try: 

g = self.equivalence_reciprocal(f) 

except (ValueError, NotImplementedError): 

# this is often not possible unless the underlying ring of 

# constants is a field 

from sage.categories.fields import Fields 

if self.domain().base() not in Fields(): 

continue 

raise 

tester.assertEqual(self.reduce(f*g), 1) 

 

def _test_inductive_valuation_inheritance(self, **options): 

r""" 

Test that every instance that is a :class:`InductiveValuation` is 

either a :class:`FiniteInductiveValuation` or a 

:class:`InfiniteInductiveValuation`. Same for 

:class:`FinalInductiveValuation` and 

:class:`NonFinalInductiveValuation`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: v._test_inductive_valuation_inheritance() 

 

""" 

tester = self._tester(**options) 

tester.assertTrue(isinstance(self, InfiniteInductiveValuation) != isinstance(self, FiniteInductiveValuation)) 

tester.assertTrue(isinstance(self, FinalInductiveValuation) != isinstance(self, NonFinalInductiveValuation)) 

 

 

class FiniteInductiveValuation(InductiveValuation, DiscreteValuation): 

r""" 

Abstract base class for iterated :mod:`augmented valuations <sage.rings.valuation.augmented_valuation>` 

on top of a :mod:`Gauss valuation <sage.rings.valuation.gauss_valuation>` which is a discrete valuation, 

i.e., the last key polynomial has finite valuation. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

 

""" 

def __init__(self, parent, phi): 

r""" 

TESTS:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: from sage.rings.valuation.inductive_valuation import FiniteInductiveValuation 

sage: isinstance(v, FiniteInductiveValuation) 

True 

 

""" 

InductiveValuation.__init__(self, parent, phi) 

DiscreteValuation.__init__(self, parent) 

 

def extensions(self, other): 

r""" 

Return the extensions of this valuation to ``other``. 

 

EXAMPLES:: 

 

sage: R.<x> = ZZ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(ZZ)) 

sage: K.<x> = FunctionField(QQ) 

sage: v.extensions(K) 

[Trivial valuation on Rational Field] 

 

""" 

from sage.categories.function_fields import FunctionFields 

if other in FunctionFields() and other.ngens() == 1: 

# extend to K[x] and from there to K(x) 

v = self.extension(self.domain().change_ring(self.domain().base().fraction_field())) 

return [other.valuation(v)] 

return super(FiniteInductiveValuation, self).extensions(other) 

 

 

class NonFinalInductiveValuation(FiniteInductiveValuation, DiscreteValuation): 

r""" 

Abstract base class for iterated :mod:`augmented valuations <sage.rings.valuation.augmented_valuation>` 

on top of a :mod:`Gauss valuation <sage.rings.valuation.gauss_valuation>` which can be extended further 

through :meth:`augmentation`. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v = v.augmentation(x^2 + x + u, 1) 

 

""" 

def __init__(self, parent, phi): 

r""" 

TESTS:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v = v.augmentation(x^2 + x + u, 1) 

sage: from sage.rings.valuation.inductive_valuation import NonFinalInductiveValuation 

sage: isinstance(v, NonFinalInductiveValuation) 

True 

 

""" 

FiniteInductiveValuation.__init__(self, parent, phi) 

DiscreteValuation.__init__(self, parent) 

 

def augmentation(self, phi, mu, check=True): 

r""" 

Return the inductive valuation which extends this valuation by mapping 

``phi`` to ``mu``. 

 

INPUT: 

 

- ``phi`` -- a polynomial in the domain of this valuation; this must be 

a key polynomial, see :meth:`is_key` for properties of key 

polynomials. 

 

- ``mu`` -- a rational number or infinity, the valuation of ``phi`` in 

the extended valuation 

 

- ``check`` -- a boolean (default: ``True``), whether or not to check 

the correctness of the parameters 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v = v.augmentation(x^2 + x + u, 1) 

sage: v = v.augmentation((x^2 + x + u)^2 + 2*x*(x^2 + x + u) + 4*x, 3) 

sage: v 

[ Gauss valuation induced by 2-adic valuation, 

v((1 + O(2^5))*x^2 + (1 + O(2^5))*x + u + O(2^5)) = 1, 

v((1 + O(2^5))*x^4 + (2^2 + O(2^6))*x^3 + (1 + (u + 1)*2 + O(2^5))*x^2 + ((u + 1)*2^2 + O(2^6))*x + (u + 1) + (u + 1)*2 + (u + 1)*2^2 + (u + 1)*2^3 + (u + 1)*2^4 + O(2^5)) = 3 ] 

 

TESTS: 

 

Make sure that we do not make the assumption that the degrees of the 

key polynomials are strictly increasing:: 

 

sage: v_K = QQ.valuation(3) 

sage: A.<t> = QQ[] 

sage: v0 = GaussValuation(A,v_K) 

 

sage: v1 = v0.augmentation(t, 1/12) 

sage: v2 = v1.augmentation(t^12 + 3, 7/6) 

sage: v3 = v2.augmentation(t^12 + 3*t^2 + 3, 9/4) 

sage: v4 = v1.augmentation(t^12 + 3*t^2 + 3, 9/4) 

sage: v3 <= v4 and v3 >= v4  

True 

 

.. SEEALSO:: 

 

:mod:`~sage.rings.valuation.augmented_valuation` 

 

""" 

from .augmented_valuation import AugmentedValuation 

return AugmentedValuation(self, phi, mu, check) 

 

def mac_lane_step(self, G, principal_part_bound=None, assume_squarefree=False, assume_equivalence_irreducible=False, report_degree_bounds_and_caches=False, coefficients=None, valuations=None, check=True): 

r""" 

Perform an approximation step towards the squarefree monic non-constant 

integral polynomial ``G`` which is not an :meth:`equivalence unit <InductiveValuation.is_equivalence_unit>`. 

 

This performs the individual steps that are used in 

:meth:`~sage.rings.valuation.valuation.DiscreteValuation.mac_lane_approximants`. 

 

INPUT: 

 

- ``G`` -- a sqaurefree monic non-constant integral polynomial ``G`` 

which is not an :meth:`equivalence unit <InductiveValuation.is_equivalence_unit>` 

 

- ``principal_part_bound`` -- an integer or ``None`` (default: 

``None``), a bound on the length of the principal part, i.e., the 

section of negative slope, of the Newton polygon of ``G`` 

 

- ``assume_squarefree`` -- whether or not to assume that ``G`` is 

squarefree (default: ``False``) 

 

- ``assume_equivalence_irreducible`` -- whether or not to assume that 

``G`` is equivalence irreducible (default: ``False``) 

 

- ``report_degree_bounds_and_caches`` -- whether or not to include internal state with the returned value (used by :meth:`~sage.rings.valuation.valuation.DiscreteValuation.mac_lane_approximants` to speed up sequential calls) 

 

- ``coefficients`` -- the coefficients of ``G`` in the :meth:`~sage.rings.valuation.developing_valuation.DevelopingValuation.phi`-adic expansion if known (default: ``None``) 

 

- ``valauations`` -- the valuations of ``coefficients`` if known 

(default: ``None``) 

 

- ``check`` -- whether to check that ``G`` is a squarefree monic 

non-constant integral polynomial and not an :meth:`equivalence unit <InductiveValuation.is_equivalence_unit>` 

(default: ``True``) 

 

TESTS:: 

 

sage: K.<x> = FunctionField(QQ) 

sage: S.<y> = K[] 

sage: F = y^2 - x^2 - x^3 - 3 

sage: v0 = GaussValuation(K._ring, QQ.valuation(3)) 

sage: v1 = v0.augmentation(K._ring.gen(), 1/3) 

sage: mu0 = K.valuation(v1) 

sage: eta0 = GaussValuation(S, mu0) 

sage: eta1 = eta0.mac_lane_step(F)[0] 

sage: eta2 = eta1.mac_lane_step(F)[0] 

sage: eta2 

[ Gauss valuation induced by Valuation on rational function field induced by [ Gauss valuation induced by 3-adic valuation, v(x) = 1/3 ], v(y + x) = 2/3 ] 

 

""" 

G = self.domain().coerce(G) 

 

if G.is_constant(): 

raise ValueError("G must not be constant") 

 

from itertools import islice 

from sage.misc.misc import verbose 

verbose("Augmenting %s towards %s"%(self, G), level=10) 

 

if not G.is_monic(): 

raise ValueError("G must be monic") 

 

if coefficients is None: 

coefficients = self.coefficients(G) 

if principal_part_bound: 

coefficients = islice(coefficients, 0, principal_part_bound + 1, 1) 

coefficients = list(coefficients) 

if valuations is None: 

valuations = self.valuations(G, coefficients=coefficients) 

if principal_part_bound: 

valuations = islice(valuations, 0, principal_part_bound + 1, 1) 

valuations = list(valuations) 

 

if check and min(valuations) < 0: 

raise ValueError("G must be integral") 

 

if check and self.is_equivalence_unit(G, valuations=valuations): 

raise ValueError("G must not be an equivalence-unit") 

 

if check and not assume_squarefree and not G.is_squarefree(): 

raise ValueError("G must be squarefree") 

 

from sage.rings.all import infinity 

assert self(G) is not infinity # this is a valuation and G is non-zero 

 

ret = [] 

 

F = self.equivalence_decomposition(G, assume_not_equivalence_unit=True, coefficients=coefficients, valuations=valuations, compute_unit=False, degree_bound=principal_part_bound) 

assert len(F), "%s equivalence-decomposes as an equivalence-unit %s"%(G, F) 

if len(F) == 1 and F[0][1] == 1 and F[0][0].degree() == G.degree(): 

assert self.is_key(G, assume_equivalence_irreducible=assume_equivalence_irreducible) 

ret.append((self.augmentation(G, infinity, check=False), G.degree(), principal_part_bound, None, None)) 

else: 

for phi,e in F: 

if G == phi: 

# Something strange happened here: 

# G is not a key (we checked that before) but phi==G is; so phi must have less precision than G 

# this can happen if not all coefficients of G have the same precision 

# if we drop some precision of G then it will be a key (but is 

# that really what we should do?) 

assert not G.base_ring().is_exact() 

prec = min([c.precision_absolute() for c in phi.list()]) 

g = G.map_coefficients(lambda c:c.add_bigoh(prec)) 

assert self.is_key(g) 

ret.append((self.augmentation(g, infinity, check=False), g.degree(), principal_part_bound, None, None)) 

assert len(F) == 1 

break 

 

if phi == self.phi(): 

# a factor phi in the equivalence decomposition means that we 

# found an actual factor of G, i.e., we can set 

# v(phi)=infinity 

# However, this should already have happened in the last step 

# (when this polynomial had -infinite slope in the Newton 

# polygon.) 

if self.is_gauss_valuation(): # unless in the first step 

pass 

else: 

continue 

 

verbose("Determining the augmentation of %s for %s"%(self, phi), level=11) 

old_mu = self(phi) 

w = self.augmentation(phi, old_mu, check=False) 

 

# we made some experiments here: instead of computing the 

# coefficients again from scratch, update the coefficients when 

# phi - self.phi() is a constant. 

# It turned out to be slightly slower than just recomputing the 

# coefficients. The main issue with the approach was that we 

# needed to keep track of all the coefficients and not just of 

# the coefficients up to principal_part_bound. 

 

w_coefficients = w.coefficients(G) 

if principal_part_bound: 

w_coefficients = islice(w_coefficients, 0, principal_part_bound + 1, 1) 

w_coefficients = list(w_coefficients) 

 

w_valuations = w.valuations(G, coefficients=w_coefficients) 

if principal_part_bound: 

w_valuations = islice(w_valuations, 0, principal_part_bound + 1, 1) 

w_valuations = list(w_valuations) 

 

from sage.geometry.newton_polygon import NewtonPolygon 

NP = NewtonPolygon(w.newton_polygon(G, valuations=w_valuations).vertices(), last_slope=0) 

 

verbose("Newton-Polygon for v(phi)=%s : %s"%(self(phi), NP), level=11) 

slopes = NP.slopes(repetition=True) 

multiplicities = {slope : len([s for s in slopes if s == slope]) for slope in slopes} 

slopes = multiplicities.keys() 

if NP.vertices()[0][0] != 0: 

slopes = [-infinity] + slopes 

multiplicities[-infinity] = 1 

 

for i, slope in enumerate(slopes): 

verbose("Slope = %s"%slope, level=12) 

new_mu = old_mu - slope 

new_valuations = [val - (j*slope if slope is not -infinity else (0 if j == 0 else -infinity)) 

for j,val in enumerate(w_valuations)] 

base = self 

if phi.degree() == base.phi().degree(): 

# very frequently, the degree of the key polynomials 

# stagnate for a bit while the valuation of the key 

# polynomial is slowly increased. 

# In this case, we can drop previous key polynomials 

# of the same degree. (They have no influence on the 

# phi-adic expansion.) 

assert new_mu > self(phi) 

if not base.is_gauss_valuation(): 

base = base._base_valuation 

w = base.augmentation(phi, new_mu, check=False) 

assert slope is -infinity or 0 in w.newton_polygon(G).slopes(repetition=False) 

 

from sage.rings.all import ZZ 

assert (phi.degree() / self.phi().degree()) in ZZ 

degree_bound = multiplicities[slope] * phi.degree() 

assert degree_bound <= G.degree() 

assert degree_bound >= phi.degree() 

ret.append((w, degree_bound, multiplicities[slope], w_coefficients, new_valuations)) 

 

assert ret 

if not report_degree_bounds_and_caches: 

ret = [v for v,_,_,_,_ in ret] 

return ret 

 

def is_key(self, phi, explain=False, assume_equivalence_irreducible=False): 

r""" 

Return whether ``phi`` is a key polynomial for this valuation, i.e., 

whether it is monic, whether it :meth:`is_equivalence_irreducible`, and 

whether it is :meth:`is_minimal`. 

 

INPUT: 

 

- ``phi`` -- a polynomial in the domain of this valuation 

 

- ``explain`` -- a boolean (default: ``False``), if ``True``, return a 

string explaining why ``phi`` is not a key polynomial 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4, 5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.is_key(x) 

True 

sage: v.is_key(2*x, explain = True) 

(False, 'phi must be monic') 

sage: v.is_key(x^2, explain = True) 

(False, 'phi must be equivalence irreducible') 

 

sage: w = v.augmentation(x, 1) 

sage: w.is_key(x + 1, explain = True) 

(False, 'phi must be minimal') 

 

""" 

phi = self.domain().coerce(phi) 

 

reason = None 

 

if not phi.is_monic(): 

reason = "phi must be monic" 

elif not assume_equivalence_irreducible and not self.is_equivalence_irreducible(phi): 

reason = "phi must be equivalence irreducible" 

elif not self.is_minimal(phi, assume_equivalence_irreducible=True): 

reason = "phi must be minimal" 

 

if explain: 

return reason is None, reason 

else: 

return reason is None 

 

def is_minimal(self, f, assume_equivalence_irreducible=False): 

r""" 

Return whether the polynomial ``f`` is minimal with respect to this 

valuation. 

 

A polynomial `f` is minimal with respect to `v` if it is not a constant 

and any non-zero polynomial `h` which is `v`-divisible by `f` has at 

least the degree of `f`. 

 

A polynomial `h` is `v`-divisible by `f` if there is a polynomial `c` 

such that `fc` :meth:`~sage.rings.valuation.valuation.DiscretePseudoValuation.is_equivalent` to `h`. 

 

ALGORITHM: 

 

Based on Theorem 9.4 of [Mac1936II]_. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4, 5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.is_minimal(x + 1) 

True 

sage: w = v.augmentation(x, 1) 

sage: w.is_minimal(x + 1) 

False 

 

TESTS:: 

 

sage: K = Qp(2, 10) 

sage: R.<x> = K[] 

sage: vp = K.valuation() 

sage: v0 = GaussValuation(R, vp) 

sage: v1 = v0.augmentation(x, 1/4) 

sage: v2 = v1.augmentation(x^4 + 2, 5/4) 

sage: v2.is_minimal(x^5 + x^4 + 2) 

False 

 

Polynomials which are equivalent to the key polynomial are minimal if 

and only if they have the same degree as the key polynomial:: 

 

sage: v2.is_minimal(x^4 + 2) 

True 

sage: v2.is_minimal(x^4 + 4) 

False 

 

""" 

f = self.domain().coerce(f) 

 

if f.is_constant(): 

return False 

 

if not assume_equivalence_irreducible and not self.is_equivalence_irreducible(f): 

# any factor divides f with respect to this valuation 

return False 

 

if not f.is_monic(): 

# divide out the leading factor, it does not change minimality 

v = self 

if not self.domain().base_ring().is_field(): 

domain = self.domain().change_ring(self.domain().base_ring().fraction_field()) 

v = self.extension(domain) 

f = domain(f) 

return v.is_minimal(f / f.leading_coefficient()) 

 

if self.is_gauss_valuation(): 

if self(f) == 0: 

F = self.reduce(f, check=False) 

assert not F.is_constant() 

return F.is_irreducible() 

else: 

assert(self(f) <= 0) # f is monic 

# f is not minimal: 

# Let g be f stripped of its leading term, i.e., g = f - x^n. 

# Then g and f are equivalent with respect to this valuation 

# and in particular g divides f with respect to this valuation 

return False 

 

if self.is_equivalent(self.phi(), f): 

assert f.degree() >= self.phi().degree() 

# If an h divides f with respect to this valuation, then it also divides phi: 

# v(f - c*h) > v(f) = v(c*h) => v(phi - c*h) = v((phi - f) + (f - c*h)) > v(phi) = v(c*h) 

# So if f were not minimal then phi would not be minimal but it is. 

return f.degree() == self.phi().degree() 

 

else: 

tau = self.value_group().index(self._base_valuation.value_group()) 

# see Theorem 9.4 of [Mac1936II] 

return list(self.valuations(f))[-1] == self(f) and \ 

list(self.coefficients(f))[-1].is_constant() and \ 

list(self.valuations(f))[0] == self(f) and \ 

tau.divides(len(list(self.coefficients(f))) - 1) 

 

def _equivalence_reduction(self, f, coefficients=None, valuations=None, degree_bound=None): 

r""" 

Helper method for :meth:`is_equivalence_irreducible` and 

:meth:`equivalence_decomposition` which essentially returns the 

reduction of ``f`` after multiplication with an ``R`` which 

:meth:`is_equivalence_unit`. 

 

This only works when ``f`` is not divisible by :meth:`phi` with respect 

to this valuation. Therefore, we also return the number of times that 

we took out :meth:`phi` of ``f`` before we computed the reduction. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: v._equivalence_reduction(2*x^6 + 4*x^5 + 2*x^4 + 8) 

(1, 4, x^2 + 1) 

 

""" 

f = self.domain().coerce(f) 

 

# base change from R[x] to K[x], so divisions work and sufficient 

# elements of negative valuation exist 

if not self.domain().base_ring().is_field(): 

domain = self.domain().change_ring(self.domain().base_ring().fraction_field()) 

v = self.extension(domain) 

assert self.residue_ring() is v.residue_ring() 

return v._equivalence_reduction(f) 

 

if coefficients is None: 

coefficients = list(self.coefficients(f)) 

if valuations is None: 

valuations = list(self.valuations(f, coefficients=coefficients)) 

valuation = min(valuations) 

for phi_divides in range(len(valuations)): 

# count how many times phi divides f 

if valuations[phi_divides] <= valuation: 

break 

 

if phi_divides: 

from sage.rings.all import PolynomialRing 

R = PolynomialRing(f.parent(), 'phi') 

f = R(coefficients[phi_divides:])(self.phi()) 

valuations = [v-self.mu()*phi_divides for v in valuations[phi_divides:]] 

coefficients = coefficients[phi_divides:] 

valuation = min(valuations) 

 

R = self.equivalence_unit(-valuation) 

R = next(self.coefficients(R)) 

fR_valuations = [v-valuation for v in valuations] 

from sage.rings.all import infinity 

fR_coefficients = [next(self.coefficients(c*R)) if v is not infinity and v == 0 else 0 

for c,v in zip(coefficients,fR_valuations)] 

 

return valuation, phi_divides, self.reduce(f*R, check=False, degree_bound=degree_bound, coefficients=fR_coefficients, valuations=fR_valuations) 

 

def is_equivalence_irreducible(self, f, coefficients=None, valuations=None): 

r""" 

Return whether the polynomial ``f`` is equivalence-irreducible, i.e., 

whether its :meth:`equivalence_decomposition` is trivial. 

 

ALGORITHM: 

 

We use the same algorithm as in :meth:`equivalence_decomposition` we 

just do not lift the result to key polynomials. 

 

INPUT: 

 

- ``f`` -- a non-constant polynomial in the domain of this valuation 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,5) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.is_equivalence_irreducible(x) 

True 

sage: v.is_equivalence_irreducible(x^2) 

False 

sage: v.is_equivalence_irreducible(x^2 + 2) 

False 

 

""" 

f = self.domain().coerce(f) 

 

if not self.domain().base_ring().is_field(): 

domain = self.domain().change_ring(self.domain().base_ring().fraction_field()) 

v = self.extension(domain) 

return v.is_equivalence_irreducible(v.domain()(f)) 

 

if f.is_constant(): 

raise ValueError("f must not be constant") 

 

_, phi_divides, F = self._equivalence_reduction(f, coefficients=coefficients, valuations=valuations) 

if phi_divides == 0: 

return F.is_constant() or F.is_irreducible() 

if phi_divides == 1: 

return F.is_constant() 

if phi_divides > 1: 

return False 

 

def equivalence_decomposition(self, f, assume_not_equivalence_unit=False, coefficients=None, valuations=None, compute_unit=True, degree_bound=None): 

r""" 

Return an equivalence decomposition of ``f``, i.e., a polynomial 

`g(x)=e(x)\prod_i \phi_i(x)` with `e(x)` an :meth:`equivalence unit 

<InductiveValuation.is_equivalence_unit>` and the `\phi_i` :meth:`key 

polynomials <is_key>` such that ``f`` :meth:`~sage.rings.valuation.valuation.DiscretePseudoValuation.is_equivalent` to `g`. 

 

INPUT: 

 

- ``f`` -- a non-zero polynomial in the domain of this valuation 

 

- ``assume_not_equivalence_unit`` -- whether or not to assume that 

``f`` is not an :meth:`equivalence unit <InductiveValuation.is_equivalence_unit>` 

(default: ``False``) 

 

- ``coefficients`` -- the coefficients of ``f`` in the 

:meth:`~sage.rings.valuation.developing_valuation.DevelopingValuation.phi`-adic 

expansion if known (default: ``None``) 

 

- ``valuations`` -- the valuations of ``coefficients`` if known 

(default: ``None``) 

 

- ``compute_unit`` -- whether or not to compute the unit part of the 

decomposition (default: ``True``) 

 

- ``degree_bound`` -- a bound on the degree of the 

:meth:`_equivalence_reduction` of ``f`` (default: ``None``) 

 

ALGORITHM: 

 

We use the algorithm described in Theorem 4.4 of [Mac1936II]_. After 

removing all factors `\phi` from a polynomial `f`, there is an 

equivalence unit `R` such that `Rf` has valuation zero. Now `Rf` can be 

factored as `\prod_i \alpha_i` over the :meth:`~sage.rings.valuation.valuation_space.DiscretePseudoValuationSpace.ElementMethods.residue_field`. Lifting 

all `\alpha_i` to key polynomials `\phi_i` gives `Rf=\prod_i R_i f_i` 

for suitable equivalence units `R_i` (see :meth:`lift_to_key`). Taking 

`R'` an :meth:`~InductiveValuation.equivalence_reciprocal` of `R`, we have `f` equivalent 

to `(R'\prod_i R_i)\prod_i\phi_i`. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,10) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.equivalence_decomposition(S.zero()) 

Traceback (most recent call last): 

... 

ValueError: equivalence decomposition of zero is not defined 

sage: v.equivalence_decomposition(S.one()) 

1 + O(2^10) 

sage: v.equivalence_decomposition(x^2+2) 

((1 + O(2^10))*x)^2 

sage: v.equivalence_decomposition(x^2+1) 

((1 + O(2^10))*x + 1 + O(2^10))^2 

 

A polynomial that is an equivalence unit, is returned as the unit part 

of a :class:`~sage.structure.factorization.Factorization`, leading to a unit 

non-minimal degree:: 

 

sage: w = v.augmentation(x, 1) 

sage: F = w.equivalence_decomposition(x^2+1); F 

(1 + O(2^10))*x^2 + 1 + O(2^10) 

sage: F.unit() 

(1 + O(2^10))*x^2 + 1 + O(2^10) 

 

However, if the polynomial has a non-unit factor, then the unit might 

be replaced by a factor of lower degree:: 

 

sage: f = x * (x^2 + 1) 

sage: F = w.equivalence_decomposition(f); F 

(1 + O(2^10))*x 

sage: F.unit() 

1 + O(2^10) 

 

Examples over an iterated unramified extension:: 

 

sage: v = v.augmentation(x^2 + x + u, 1) 

sage: v = v.augmentation((x^2 + x + u)^2 + 2*x*(x^2 + x + u) + 4*x, 3) 

 

sage: v.equivalence_decomposition(x) 

(1 + O(2^10))*x 

sage: F = v.equivalence_decomposition( v.phi() ) 

sage: len(F) 

1 

sage: F = v.equivalence_decomposition( v.phi() * (x^4 + 4*x^3 + (7 + 2*u)*x^2 + (8 + 4*u)*x + 1023 + 3*u) ) 

sage: len(F) 

2 

 

TESTS:: 

 

sage: R.<x> = QQ[] 

sage: K1.<pi>=NumberField(x^3 - 2) 

sage: K.<alpha>=K1.galois_closure() 

sage: R.<x>=K[] 

sage: vp = QQ.valuation(2) 

sage: vp = vp.extension(K) 

sage: v0 = GaussValuation(R, vp) 

sage: G=x^36 + 36*x^35 + 630*x^34 + 7144*x^33 + 59055*x^32 + 379688*x^31 +1978792*x^30 + 8604440*x^29 + 31895428*x^28 + 102487784*x^27 + 289310720*x^26 + 725361352*x^25 + 1629938380*x^24 + 3307417800*x^23 + 6098786184*x^22+10273444280*x^21 + 15878121214*x^20 + 22596599536*x^19 + 29695703772*x^18 +36117601976*x^17 + 40722105266*x^16 + 42608585080*x^15 + 41395961848*x^14 +37344435656*x^13 + 31267160756*x^12 + 24271543640*x^11 + 17439809008*x^10 + 11571651608*x^9 + 7066815164*x^8 + 3953912472*x^7 + 2013737432*x^6 + 925014888*x^5 + 378067657*x^4 + 134716588*x^3 + 40441790*x^2 + 9532544*x + 1584151 

sage: v1 = v0.mac_lane_step(G)[0] 

sage: V = v1.mac_lane_step(G) 

sage: v2 = V[0] 

sage: F = v2.equivalence_decomposition(G); F 

(x^4 + 2*x^2 + alpha^4 + alpha^3 + 1)^3 * (x^4 + 2*x^2 + 1/2*alpha^4 + alpha^3 + 5*alpha + 1)^3 * (x^4 + 2*x^2 + 3/2*alpha^4 + alpha^3 + 5*alpha + 1)^3 

sage: v2.is_equivalent(F.prod(), G) 

True 

 

""" 

f = self.domain().coerce(f) 

 

if f.is_zero(): 

raise ValueError("equivalence decomposition of zero is not defined") 

 

from sage.structure.factorization import Factorization 

if not assume_not_equivalence_unit and self.is_equivalence_unit(f): 

return Factorization([], unit=f, sort=False) 

 

if not self.domain().base_ring().is_field(): 

nonfractions = self.domain().base_ring() 

domain = self.domain().change_ring(nonfractions.fraction_field()) 

v = self.extension(domain) 

ret = v.equivalence_decomposition(v.domain()(f)) 

return Factorization([(self._eliminate_denominators(g), e) 

for (g,e) in ret], unit=self._eliminate_denominators(ret.unit())) 

 

valuation, phi_divides, F = self._equivalence_reduction(f, coefficients=coefficients, valuations=valuations, degree_bound=degree_bound) 

F = F.factor() 

from sage.misc.misc import verbose 

verbose("%s factors as %s = %s in reduction"%(f, F.prod(), F), level=20) 

 

unit = self.domain().one() 

if compute_unit: 

R_ = self.equivalence_unit(valuation, reciprocal=True) 

unit = self.lift(self.residue_ring()(F.unit())) * R_ 

F = list(F) 

 

if compute_unit: 

from sage.misc.all import prod 

unit *= self.lift(self.residue_ring()(prod([ psi.leading_coefficient()**e for psi,e in F ]))) 

 

# A potential speedup that we tried to implement here: 

# When F factors as T^n - a, then instead of using any lift of T^n - a 

# we tried to take a lift that approximates well an n-th root of the 

# constant coefficient of f[0]. Doing so saved a few invocations of 

# mac_lane_step but in the end made hardly any difference. 

 

F = [(self.lift_to_key(psi/psi.leading_coefficient()),e) for psi,e in F] 

 

if compute_unit: 

for g,e in F: 

v_g = self(g) 

unit *= self._pow(self.equivalence_unit(-v_g, reciprocal=True), e, error=-v_g*e, effective_degree=0) 

unit = self.simplify(unit, effective_degree=0) 

 

if phi_divides: 

for i,(g,e) in enumerate(F): 

if g == self.phi(): 

F[i] = (self.phi(),e+phi_divides) 

break 

else: 

F.append((self.phi(),phi_divides)) 

 

ret = Factorization(F, unit=unit, sort=False) 

 

if compute_unit: 

assert self.is_equivalent(ret.prod(), f) # this might fail because of leading zeros in inexact rings 

assert self.is_equivalence_unit(ret.unit()) 

 

return ret 

 

def minimal_representative(self, f): 

r""" 

Return a minimal representative for ``f``, i.e., a pair `e, a` such 

that ``f`` :meth:`~sage.rings.valuation.valuation.DiscretePseudoValuation.is_equivalent` to `e a`, `e` is an 

:meth:`equivalence unit <InductiveValuation.is_equivalence_unit>`, and `a` :meth:`is_minimal` and monic. 

 

INPUT: 

 

- ``f`` -- a non-zero polynomial which is not an equivalence unit 

 

OUTPUT: 

 

A factorization which has `e` as its unit and `a` as its unique factor. 

 

ALGORITHM: 

 

We use the algorithm described in the proof of Lemma 4.1 of [Mac1936II]_. 

In the expansion `f=\sum_i f_i\phi^i` take `e=f_i` for the largest `i` 

with `f_i\phi^i` minimal (see :meth:`~sage.rings.valuation.developing_valuation.DevelopingValuation.effective_degree`). 

Let `h` be the :meth:`~InductiveValuation.equivalence_reciprocal` of `e` and take `a` given 

by the terms of minimal valuation in the expansion of `e f`. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,10) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: v.minimal_representative(x + 2) 

(1 + O(2^10))*x 

 

sage: v = v.augmentation(x, 1) 

sage: v.minimal_representative(x + 2) 

(1 + O(2^10))*x + 2 + O(2^11) 

sage: f = x^3 + 6*x + 4 

sage: F = v.minimal_representative(f); F 

(2 + 2^2 + O(2^11)) * ((1 + O(2^10))*x + 2 + O(2^11)) 

sage: v.is_minimal(F[0][0]) 

True 

sage: v.is_equivalent(F.prod(), f) 

True 

 

""" 

f = self.domain().coerce(f) 

 

from sage.categories.fields import Fields 

if not self.domain().base_ring() in Fields(): 

raise NotImplementedError("only implemented for polynomial rings over fields") 

 

if f.is_zero(): 

raise ValueError("zero has no minimal representative") 

 

degree = self.effective_degree(f) 

if degree == 0: 

raise ValueError("equivalence units can not have a minimal representative") 

 

e = list(self.coefficients(f))[degree] 

h = self.equivalence_reciprocal(e).map_coefficients(lambda c:_lift_to_maximal_precision(c)) 

g = h*f 

vg = self(g) 

 

coeffs = [c if v == vg else c.parent().zero() for v,c in zip(self.valuations(g), self.coefficients(g))] 

coeffs[degree] = self.domain().base_ring().one() 

ret = sum([c*self._phi**i for i,c in enumerate(coeffs)]) 

 

assert self.effective_degree(ret) == degree 

assert ret.is_monic() 

assert self.is_minimal(ret) 

 

from sage.structure.factorization import Factorization 

ret = Factorization([(ret, 1)], unit=e, sort=False) 

 

assert self.is_equivalent(ret.prod(), f) # this might fail because of leading zeros 

return ret 

 

@abstract_method 

def lift_to_key(self, F): 

""" 

Lift the irreducible polynomial ``F`` from the 

:meth:`~sage.rings.valuation.valuation_space.DiscretePseudoValuationSpace.ElementMethods.residue_ring` 

to a key polynomial over this valuation. 

 

INPUT: 

 

- ``F`` -- an irreducible non-constant monic polynomial in 

:meth:`~sage.rings.valuation.valuation_space.DiscretePseudoValuationSpace.ElementMethods.residue_ring` 

of this valuation 

 

OUTPUT: 

 

A polynomial `f` in the domain of this valuation which is a key 

polynomial for this valuation and which is such that an 

:meth:`augmentation` with this polynomial adjoins a root of ``F`` to 

the resulting :meth:`~sage.rings.valuation.valuation_space.DiscretePseudoValuationSpace.ElementMethods.residue_ring`. 

 

More specifically, if ``F`` is not the generator of the residue ring, 

then multiplying ``f`` with the :meth:`~InductiveValuation.equivalence_reciprocal` of the 

:meth:`~InductiveValuation.equivalence_unit` of the valuation of ``f``, produces a unit 

which reduces to ``F``. 

 

EXAMPLES:: 

 

sage: R.<u> = Qq(4,10) 

sage: S.<x> = R[] 

sage: v = GaussValuation(S) 

sage: y = v.residue_ring().gen() 

sage: u0 = v.residue_ring().base_ring().gen() 

sage: f = v.lift_to_key(y^2 + y + u0); f 

(1 + O(2^10))*x^2 + (1 + O(2^10))*x + u + O(2^10) 

 

""" 

 

def _eliminate_denominators(self, f): 

r""" 

Return a polynomial in the domain of this valuation that 

:meth:`is_equivalent` to ``f``. 

 

INPUT: 

 

- ``f`` -- a polynomial with coefficients in the fraction field of the 

base ring of the domain of this valuation. 

 

EXAMPLES:: 

 

sage: R.<x> = ZZ[] 

sage: v = GaussValuation(R, ZZ.valuation(2)) 

sage: v._eliminate_denominators(x/3) 

x 

 

In general such a polynomial may not exist:: 

 

sage: w = v.augmentation(x, 1) 

sage: w._eliminate_denominators(x/2) 

Traceback (most recent call last): 

... 

ValueError: element has no approximate inverse in this ring 

 

In general it exists iff the coefficients of minimal valuation in the 

`\phi`-adic expansion of ``f`` do not have denominators of positive 

valuation and if the same is true for these coefficients in their 

expansion; at least if the coefficient ring's residue ring is already a 

field:: 

 

sage: w._eliminate_denominators(x^3/2 + x) 

x 

 

""" 

if f in self.domain(): 

return self.domain()(f) 

 

nonfractions = self.domain().base_ring() 

fractions = nonfractions.fraction_field() 

 

extended_domain = self.domain().change_ring(fractions) 

 

g = extended_domain.coerce(f) 

 

w = self.extension(extended_domain) 

# drop coefficients whose valuation is not minimal (recursively) 

valuation = w(g) 

g = w.simplify(g, error=valuation, force=True, phiadic=True) 

 

if g in self.domain(): 

return self.domain()(g) 

 

nonfraction_valuation = self.restriction(nonfractions) 

# if this fails then there is no equivalent polynomial in the domain of this valuation 

ret = g.map_coefficients( 

lambda c: c.numerator()*nonfraction_valuation.inverse(c.denominator(), 

valuation 

+ nonfraction_valuation(c.denominator()) 

- nonfraction_valuation(c.numerator()) 

+ nonfraction_valuation.value_group().gen()), 

nonfractions) 

assert w.is_equivalent(f, ret) 

return ret 

 

 

def _test_eliminate_denominators(self, **options): 

r""" 

Test the correctness of :meth:`_eliminate_denominators`. 

 

EXAMPLES:: 

 

sage: R.<x> = ZZ[] 

sage: v = GaussValuation(R, ZZ.valuation(2)) 

sage: v._test_eliminate_denominators() 

 

""" 

tester = self._tester(**options) 

 

nonfractions = self.domain().base_ring() 

fractions = nonfractions.fraction_field() 

extended_domain = self.domain().change_ring(fractions) 

w = self.extension(extended_domain) 

 

S = tester.some_elements(w.domain().some_elements()) 

for f in S: 

try: 

g = self._eliminate_denominators(f) 

except ValueError: 

continue 

tester.assertTrue(g.parent() is self.domain()) 

tester.assertTrue(w.is_equivalent(f, g)) 

 

def _test_lift_to_key(self, **options): 

r""" 

Test the correctness of :meth:`lift_to_key`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: v._test_lift_to_key() 

 

""" 

tester = self._tester(**options) 

 

try: 

k = self.residue_ring() 

except NotImplementedError: 

from sage.categories.fields import Fields 

if self.domain().base() in Fields(): 

raise 

return 

 

S = tester.some_elements(self.residue_ring().some_elements()) 

for F in S: 

if F.is_monic() and not F.is_constant() and F.is_irreducible(): 

try: 

f = self.lift_to_key(F) 

except NotImplementedError: 

from sage.categories.fields import Fields 

if self.domain().base() in Fields(): 

raise 

continue 

tester.assertIs(f.parent(), self.domain()) 

tester.assertTrue(self.is_key(f)) 

 

# check that augmentation produces a valuation with roots of F 

# in the residue ring 

from sage.rings.all import infinity 

w = self.augmentation(f, infinity) 

F = F.change_ring(w.residue_ring()) 

roots = F.roots(multiplicities=False) 

tester.assertGreaterEqual(len(roots), 1) 

 

# check that f has the right reduction 

if F == F.parent().gen(): 

tester.assertTrue(self.is_equivalent(f, self.phi())) 

else: 

tester.assertEqual(self.reduce(f * self.equivalence_reciprocal(self.equivalence_unit(self(f)))), F) 

 

 

def _test_is_equivalence_irreducible(self, **options): 

r""" 

Test the correctness of :meth:`is_equivalence_irreducible`. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: v._test_is_equivalence_irreducible() 

 

""" 

tester = self._tester(**options) 

S = tester.some_elements(self.domain().some_elements()) 

for f in S: 

if f.is_constant(): continue 

is_equivalence_irreducible = self.is_equivalence_irreducible(f) 

F = self.equivalence_decomposition(f) 

tester.assertEqual(is_equivalence_irreducible, len(F)==0 or (len(F)==1 and F[0][1]==1)) 

if self.is_equivalence_unit(f): 

tester.assertTrue(f.is_constant() or self.is_equivalence_irreducible(f)) 

 

tester.assertTrue(self.is_equivalence_irreducible(self.phi())) 

tester.assertTrue(self.is_equivalence_irreducible(-self.phi())) 

tester.assertFalse(self.is_equivalence_irreducible(self.phi() ** 2)) 

 

 

class FinalInductiveValuation(InductiveValuation): 

r""" 

Abstract base class for an inductive valuation which can not be augmented further. 

 

TESTS:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, valuations.TrivialValuation(QQ)) 

sage: w = v.augmentation(x^2 + x + 1, infinity) 

sage: from sage.rings.valuation.inductive_valuation import FinalInductiveValuation 

sage: isinstance(w, FinalInductiveValuation) 

True 

 

""" 

 

 

class InfiniteInductiveValuation(FinalInductiveValuation, InfiniteDiscretePseudoValuation): 

r""" 

Abstract base class for an inductive valuation which is not discrete, i.e., 

which assigns infinite valuation to its last key polynomial. 

 

EXAMPLES:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: w = v.augmentation(x^2 + x + 1, infinity) 

 

""" 

def __init__(self, parent, base_valuation): 

r""" 

TESTS:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: w = v.augmentation(x^2 + x + 1, infinity) 

sage: from sage.rings.valuation.inductive_valuation import InfiniteInductiveValuation 

sage: isinstance(w, InfiniteInductiveValuation) 

True 

 

""" 

FinalInductiveValuation.__init__(self, parent, base_valuation) 

InfiniteDiscretePseudoValuation.__init__(self, parent) 

 

def change_domain(self, ring): 

r""" 

Return this valuation over ``ring``. 

 

EXAMPLES: 

 

We can turn an infinite valuation into a valuation on the quotient:: 

 

sage: R.<x> = QQ[] 

sage: v = GaussValuation(R, QQ.valuation(2)) 

sage: w = v.augmentation(x^2 + x + 1, infinity) 

sage: w.change_domain(R.quo(x^2 + x + 1)) 

2-adic valuation 

 

""" 

from sage.rings.polynomial.polynomial_quotient_ring import is_PolynomialQuotientRing 

if is_PolynomialQuotientRing(ring) and ring.base() is self.domain() and ring.modulus() == self.phi(): 

return self.restriction(self.domain().base())._extensions_to_quotient(ring, approximants=[self])[0] 

return super(InfiniteInductiveValuation, self).change_domain(ring) 

 

 

def _lift_to_maximal_precision(c): 

r""" 

Lift ``c`` to maximal precision if the parent is not exact. 

 

EXAMPLES:: 

 

sage: R = Zp(2,5) 

sage: x = R(1,2); x 

1 + O(2^2) 

sage: from sage.rings.valuation.inductive_valuation import _lift_to_maximal_precision 

sage: _lift_to_maximal_precision(x) 

1 + O(2^5) 

 

sage: x = 1 

sage: _lift_to_maximal_precision(x) 

1 

 

""" 

return c if c.parent().is_exact() else c.lift_to_precision()