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

""" 

Local Generic 

 

Superclass for `p`-adic and power series rings. 

 

AUTHORS: 

 

- David Roe 

""" 

from __future__ import absolute_import 

 

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

# Copyright (C) 2007-2013 David Roe <roed.math@gmail.com> 

# William Stein <wstein@gmail.com> 

# 

# 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 copy import copy 

from sage.rings.ring import CommutativeRing 

from sage.categories.complete_discrete_valuation import CompleteDiscreteValuationRings, CompleteDiscreteValuationFields 

from sage.structure.category_object import check_default_category 

from sage.structure.parent import Parent 

from sage.rings.integer import Integer 

from sage.rings.integer_ring import ZZ 

 

class LocalGeneric(CommutativeRing): 

def __init__(self, base, prec, names, element_class, category=None): 

""" 

Initializes self. 

 

EXAMPLES:: 

 

sage: R = Zp(5) #indirect doctest 

sage: R.precision_cap() 

20 

 

In :trac:`14084`, the category framework has been implemented for p-adic rings:: 

 

sage: TestSuite(R).run() 

sage: K = Qp(7) 

sage: TestSuite(K).run() 

 

TESTS:: 

 

sage: R = Zp(5, 5, 'fixed-mod') 

sage: R._repr_option('element_is_atomic') 

False 

""" 

self._prec = prec 

self.Element = element_class 

default_category = getattr(self, '_default_category', None) 

if self.is_field(): 

category = CompleteDiscreteValuationFields() 

else: 

category = CompleteDiscreteValuationRings() 

category = category.Metric().Complete() 

if default_category is not None: 

category = check_default_category(default_category, category) 

Parent.__init__(self, base, names=(names,), normalize=False, category=category) 

 

def is_capped_relative(self): 

""" 

Returns whether this `p`-adic ring bounds precision in a capped 

relative fashion. 

 

The relative precision of an element is the power of `p` 

modulo which the unit part of that element is defined. In a 

capped relative ring, the relative precision of elements are 

bounded by a constant depending on the ring. 

 

EXAMPLES:: 

 

sage: R = ZpCA(5, 15) 

sage: R.is_capped_relative() 

False 

sage: R(5^7) 

5^7 + O(5^15) 

sage: S = Zp(5, 15) 

sage: S.is_capped_relative() 

True 

sage: S(5^7) 

5^7 + O(5^22) 

""" 

return False 

 

def is_capped_absolute(self): 

""" 

Returns whether this `p`-adic ring bounds precision in a 

capped absolute fashion. 

 

The absolute precision of an element is the power of `p` 

modulo which that element is defined. In a capped absolute 

ring, the absolute precision of elements are bounded by a 

constant depending on the ring. 

 

EXAMPLES:: 

 

sage: R = ZpCA(5, 15) 

sage: R.is_capped_absolute() 

True 

sage: R(5^7) 

5^7 + O(5^15) 

sage: S = Zp(5, 15) 

sage: S.is_capped_absolute() 

False 

sage: S(5^7) 

5^7 + O(5^22) 

""" 

return False 

 

def is_fixed_mod(self): 

""" 

Returns whether this `p`-adic ring bounds precision in a fixed 

modulus fashion. 

 

The absolute precision of an element is the power of `p` 

modulo which that element is defined. In a fixed modulus 

ring, the absolute precision of every element is defined to be 

the precision cap of the parent. This means that some 

operations, such as division by `p`, don't return a well defined 

answer. 

 

EXAMPLES:: 

 

sage: R = ZpFM(5,15) 

sage: R.is_fixed_mod() 

True 

sage: R(5^7,absprec=9) 

5^7 + O(5^15) 

sage: S = ZpCA(5, 15) 

sage: S.is_fixed_mod() 

False 

sage: S(5^7,absprec=9) 

5^7 + O(5^9) 

""" 

return False 

 

def is_floating_point(self): 

""" 

Returns whether this `p`-adic ring bounds precision in a floating 

point fashion. 

 

The relative precision of an element is the power of `p` 

modulo which the unit part of that element is defined. In a 

floating point ring, elements do not store precision, but arithmetic 

operations truncate to a relative precision depending on the ring. 

 

EXAMPLES:: 

 

sage: R = ZpCR(5, 15) 

sage: R.is_floating_point() 

False 

sage: R(5^7) 

5^7 + O(5^22) 

sage: S = ZpFP(5, 15) 

sage: S.is_floating_point() 

True 

sage: S(5^7) 

5^7 

""" 

return False 

 

def is_lattice_prec(self): 

""" 

Returns whether this `p`-adic ring bounds precision using 

a lattice model. 

 

In lattice precision, relationships between elements 

are stored in a precision object of the parent, which 

allows for optimal precision tracking at the cost of 

increased memory usage and runtime. 

 

EXAMPLES:: 

 

sage: R = ZpCR(5, 15) 

sage: R.is_lattice_prec() 

False 

sage: x = R(25, 8) 

sage: x - x 

O(5^8) 

sage: S = ZpLC(5, 15) 

doctest:...: FutureWarning: This class/method/function is marked as experimental. It, its functionality or its interface might change without a formal deprecation. 

See http://trac.sagemath.org/23505 for details. 

sage: S.is_lattice_prec() 

True 

sage: x = S(25, 8) 

sage: x - x 

O(5^30) 

""" 

return False 

 

def is_lazy(self): 

""" 

Returns whether this `p`-adic ring bounds precision in a lazy 

fashion. 

 

In a lazy ring, elements have mechanisms for computing 

themselves to greater precision. 

 

EXAMPLES:: 

 

sage: R = Zp(5) 

sage: R.is_lazy() 

False 

""" 

return False 

 

def _latex_(self): 

r""" 

Latex. 

 

EXAMPLES:: 

 

sage: latex(Zq(27,names='a')) #indirect doctest 

\mathbf{Z}_{3^{3}} 

""" 

return self._repr_(do_latex = True) 

 

def change(self, **kwds): 

r""" 

Return a new ring with changed attributes. 

 

INPUT: 

 

The following arguments are applied to every ring in the tower: 

 

- ``type`` -- string, the precision type 

- ``p`` -- the prime of the ground ring. Defining polynomials 

will be converted to the new base rings. 

- ``print_mode`` -- string 

- ``print_pos`` -- bool 

- ``print_sep`` -- string 

- ``print_alphabet`` -- dict 

- ``show_prec`` -- bool 

- ``check`` -- bool 

- ``label`` -- string (only for lattice precision) 

 

The following arguments are only applied to the top ring in the tower: 

 

- ``var_name`` -- string 

- ``res_name`` -- string 

- ``unram_name`` -- string 

- ``ram_name`` -- string 

- ``names`` -- string 

- ``modulus`` -- polynomial 

 

The following arguments have special behavior: 

 

- ``prec`` -- integer. If the precision is increased on an extension ring, 

the precision on the base is increased as necessary (respecting ramification). 

If the precision is decreased, the precision of the base is unchanged. 

 

- ``field`` -- bool. If ``True``, switch to a tower of fields via the fraction field. 

If False, switch to a tower of rings of integers. 

 

- ``q`` -- prime power. Replace the initial unramified extension of `\QQ_p` or `\ZZ_p` 

with an unramified extension of residue cardinality `q`. 

If the initial extension is ramified, add in an unramified extension. 

 

- ``base`` -- ring or field. Use a specific base ring instead of recursively 

calling :meth:`change` down the tower. 

 

See the :mod:`constructors <sage.rings.padics.factory>` for more details on the 

meaning of these arguments. 

 

EXAMPLES: 

 

We can use this method to change the precision:: 

 

sage: Zp(5).change(prec=40) 

5-adic Ring with capped relative precision 40 

 

or the precision type:: 

 

sage: Zp(5).change(type="capped-abs") 

5-adic Ring with capped absolute precision 20 

 

or even the prime:: 

 

sage: ZpCA(3).change(p=17) 

17-adic Ring with capped absolute precision 20 

 

You can switch between the ring of integers and its fraction field:: 

 

sage: ZpCA(3).change(field=True) 

3-adic Field with capped relative precision 20 

 

You can also change print modes:: 

 

sage: R = Zp(5).change(prec=5, print_mode='digits') 

sage: repr(~R(17)) 

'...13403' 

 

Changing print mode to 'digits' works for Eisenstein extensions:: 

 

sage: S.<x> = ZZ[] 

sage: W.<w> = Zp(3).extension(x^4 + 9*x^2 + 3*x - 3) 

sage: W.print_mode() 

'series' 

sage: W.change(print_mode='digits').print_mode() 

'digits' 

 

You can change extensions:: 

 

sage: K.<a> = QqFP(125, prec=4) 

sage: K.change(q=64) 

Unramified Extension in a defined by x^6 + x^4 + x^3 + x + 1 with floating precision 4 over 2-adic Field 

sage: R.<x> = QQ[] 

sage: K.change(modulus = x^2 - x + 2, print_pos=False) 

Unramified Extension in a defined by x^2 - x + 2 with floating precision 4 over 5-adic Field 

 

and variable names:: 

 

sage: K.change(names='b') 

Unramified Extension in b defined by x^3 + 3*x + 3 with floating precision 4 over 5-adic Field 

 

and precision:: 

 

sage: Kup = K.change(prec=8); Kup 

Unramified Extension in a defined by x^3 + 3*x + 3 with floating precision 8 over 5-adic Field 

sage: Kup.base_ring() 

5-adic Field with floating precision 8 

 

If you decrease the precision, the precision of the base stays the same:: 

 

sage: Kdown = K.change(prec=2); Kdown 

Unramified Extension in a defined by x^3 + 3*x + 3 with floating precision 2 over 5-adic Field 

sage: Kdown.precision_cap() 

2 

sage: Kdown.base_ring() 

5-adic Field with floating precision 4 

 

Changing the prime works for extensions:: 

 

sage: x = polygen(ZZ) 

sage: R.<a> = Zp(5).extension(x^2 + 2) 

sage: S = R.change(p=7) 

sage: S.defining_polynomial(exact=True) 

x^2 + 2 

sage: A.<y> = Zp(5)[] 

sage: R.<a> = Zp(5).extension(y^2 + 2) 

sage: S = R.change(p=7) 

sage: S.defining_polynomial(exact=True) 

y^2 + 2 

 

:: 

 

sage: R.<a> = Zq(5^3) 

sage: S = R.change(prec=50) 

sage: S.defining_polynomial(exact=True) 

x^3 + 3*x + 3 

 

Changing label for lattice precision (the precision lattice is not copied):: 

 

sage: R = ZpLC(37, (8,11)) 

sage: S = R.change(label = "change"); S 

37-adic Ring with lattice-cap precision (label: change) 

sage: S.change(label = "new") 

37-adic Ring with lattice-cap precision (label: new) 

""" 

# We support both print_* and * for *=mode, pos, sep, alphabet 

for atr in ('print_mode', 'print_pos', 'print_sep', 'print_alphabet'): 

if atr in kwds: 

kwds[atr[6:]] = kwds.pop(atr) 

def get_unramified_modulus(q, res_name): 

from sage.rings.finite_rings.finite_field_constructor import FiniteField as GF 

return GF(q, res_name).modulus().change_ring(ZZ) 

n = None 

q = None 

from .padic_base_generic import pAdicBaseGeneric 

if 'q' in kwds and isinstance(self.base_ring(), pAdicBaseGeneric): 

q = kwds.pop('q') 

if not isinstance(q, Integer): 

raise TypeError("q must be an integer") 

p, n = q.is_prime_power(get_data=True) 

if n == 0: 

raise ValueError("q must be a prime power") 

if 'p' in kwds and kwds['p'] != p: 

raise ValueError("q does not match p") 

kwds['p'] = p 

functor, ring = self.construction() 

functor = copy(functor) 

if 'mode' in kwds and 'show_prec' not in kwds: 

new_type = kwds.get('type', self._prec_type()) 

cur_type = self._prec_type() 

cur_mode = self._printer._print_mode() 

cur_show_prec = self._printer._show_prec() 

from .factory import _default_show_prec 

if cur_show_prec == _default_show_prec(cur_type, cur_mode): 

kwds['show_prec'] = _default_show_prec(new_type, kwds['mode']) 

else: 

raise RuntimeError 

p = kwds.get('p', functor.p if hasattr(functor, 'p') else self.prime()) 

curpstr = str(self.prime()) 

functor_dict = getattr(functor, "extras", getattr(functor, "kwds", None)) 

# If we are switching to 'digits', or changing p, need to ensure a large enough alphabet. 

if 'alphabet' not in kwds and (kwds.get('mode') == 'digits' or 

(functor_dict['print_mode'].get('mode') == 'digits' and p > getattr(functor, "p", p))): 

from .padic_printing import _printer_defaults 

kwds['alphabet'] = _printer_defaults.alphabet()[:p] 

# For fraction fields of fixed-mod rings, we need to explicitly set show_prec = False 

if 'field' in kwds and 'type' not in kwds: 

if self._prec_type() == 'capped-abs': 

kwds['type'] = 'capped-rel' 

elif self._prec_type() == 'fixed-mod': 

kwds['type'] = 'floating-point' 

kwds['show_prec'] = False # This can be removed once printing of fixed mod elements is changed. 

 

# There are two kinds of functors possible: 

# CompletionFunctor and AlgebraicExtensionFunctor 

# We distinguish them by the presence of ``prec``, 

if hasattr(functor, "prec"): 

functor.extras = copy(functor.extras) 

functor.extras['print_mode'] = copy(functor.extras['print_mode']) 

if 'type' in kwds and kwds['type'] not in functor._dvr_types: 

raise ValueError("completion type must be one of %s"%(", ".join(functor._dvr_types[1:]))) 

if 'field' in kwds: 

field = kwds.pop('field') 

if field: 

ring = ring.fraction_field() 

elif ring.is_field(): 

ring = ring.ring_of_integers() 

for atr in ('p', 'prec', 'type'): 

if atr in kwds: 

setattr(functor, atr, kwds.pop(atr)) 

if q is not None: 

if 'names' in kwds: 

names = kwds.pop('names') 

elif 'unram_name' in kwds: 

names = kwds.pop('unram_name') 

else: 

raise TypeError("You must specify the name of the generator") 

res_name = kwds.pop('res_name', names + '0') 

modulus = kwds.pop('modulus', get_unramified_modulus(q, res_name)) 

implementation = kwds.pop('implementation', 'FLINT') 

# We have to change the way p prints in the default case 

if 'names' in kwds: 

functor.extras['names'] = kwds.pop('names') 

elif functor.extras['names'][0] == curpstr: 

functor.extras['names'] = (str(p),) 

# labels for lattice precision 

if 'label' in kwds: 

functor.extras['label'] = kwds.pop('label') 

elif 'label' in functor.extras and functor.type not in ['lattice-cap','lattice-float']: 

del functor.extras['label'] 

for atr in ('ram_name', 'var_name'): 

if atr in kwds: 

functor.extras['print_mode'][atr] = kwds.pop(atr) 

elif functor.extras['print_mode'].get(atr) == curpstr: 

functor.extras['print_mode'][atr] = str(p) 

if 'check' in kwds: 

functor.extras['check'] = kwds.pop('check') 

for atr in ('mode', 'pos', 'unram_name', 'max_ram_terms', 'max_unram_terms', 'max_terse_terms', 'sep', 'alphabet', 'show_prec'): 

if atr in kwds: 

functor.extras['print_mode'][atr] = kwds.pop(atr) 

if kwds: 

raise ValueError("Extra arguments received: %s"%(", ".join(kwds.keys()))) 

if q is not None: 

# Create an unramified extension 

base = functor(ring) 

from .factory import ExtensionFactory 

modulus = modulus.change_ring(base) 

return ExtensionFactory(base=base, premodulus=modulus, names=names, res_name=res_name, unram=True, implementation=implementation) 

else: 

functor.kwds = copy(functor.kwds) 

functor.kwds['print_mode'] = copy(functor.kwds['print_mode']) 

if 'prec' in kwds: 

# This will need to be modified once lattice precision supports extensions 

prec = kwds.pop('prec') 

baseprec = (prec - 1) // self.e() + 1 

if baseprec > self.base_ring().precision_cap(): 

kwds['prec'] = baseprec 

functor.kwds['prec'] = prec 

from sage.rings.padics.padic_base_generic import pAdicBaseGeneric 

if 'names' in kwds: 

functor.names = [kwds.pop('names')] 

modulus = None 

if 'modulus' in kwds: 

modulus = kwds.pop('modulus') 

if n is not None and modulus.degree() != n: 

raise ValueError("modulus must have degree matching q") 

elif q is not None and self.f() != 1: 

# If q is specified, replace the modulus with one from q. 

modulus = get_unramified_modulus(q, functor.kwds.get('res_name', functor.names[0] + '0')) 

for atr in ('var_name', 'res_name', 'unram_name', 'ram_name'): 

if atr in kwds: 

functor.kwds[atr] = kwds.pop(atr) 

if 'check' in kwds: 

functor.kwds['check'] = kwds['check'] 

for atr in ('mode', 'pos', 'max_ram_terms', 'max_unram_terms', 'max_terse_terms', 'sep', 'alphabet', 'show_prec'): 

if atr in kwds: 

functor.kwds['print_mode'][atr] = kwds[atr] 

if 'base' in kwds: 

ring = kwds['base'] 

else: 

if q is not None and self.f() == 1: 

kwds['q'] = q 

ring = ring.change(**kwds) 

if modulus is None: 

if len(functor.polys) != 1: 

raise RuntimeError("Unexpected number of defining polynomials") 

modulus = functor.polys[0] 

if isinstance(modulus.base_ring(), pAdicBaseGeneric): 

modulus.change_ring(ring) 

functor.polys = [modulus] 

return functor(ring) 

 

def precision_cap(self): 

r""" 

Returns the precision cap for this ring. 

 

EXAMPLES:: 

 

sage: R = Zp(3, 10,'fixed-mod'); R.precision_cap() 

10 

sage: R = Zp(3, 10,'capped-rel'); R.precision_cap() 

10 

sage: R = Zp(3, 10,'capped-abs'); R.precision_cap() 

10 

 

.. NOTE:: 

 

This will have different meanings depending on the type of 

local ring. For fixed modulus rings, all elements are 

considered modulo ``self.prime()^self.precision_cap()``. 

For rings with an absolute cap (i.e. the class 

``pAdicRingCappedAbsolute``), each element has a precision 

that is tracked and is bounded above by 

``self.precision_cap()``. Rings with relative caps 

(e.g. the class ``pAdicRingCappedRelative``) are the same 

except that the precision is the precision of the unit 

part of each element. 

""" 

return self._prec 

 

def _precision_cap(self): 

r""" 

Returns the precision cap for this ring, in the format 

used by the factory methods to create the ring. 

 

For most `p`-adic types, this is the same as :meth:`precision_cap`, 

but there is a difference for lattice precision. 

 

EXAMPLES:: 

 

sage: Zp(17,34)._precision_cap() 

34 

""" 

return self._prec 

 

def is_exact(self): 

r""" 

Returns whether this p-adic ring is exact, i.e. False. 

 

INPUT: 

self -- a p-adic ring 

 

OUTPUT: 

boolean -- whether self is exact, i.e. False. 

 

EXAMPLES:: 

 

#sage: R = Zp(5, 3, 'lazy'); R.is_exact() 

#False 

sage: R = Zp(5, 3, 'fixed-mod'); R.is_exact() 

False 

""" 

return False 

 

def residue_characteristic(self): 

r""" 

Returns the characteristic of ``self``'s residue field. 

 

INPUT: 

 

- ``self`` -- a p-adic ring. 

 

OUTPUT: 

 

- integer -- the characteristic of the residue field. 

 

EXAMPLES:: 

 

sage: R = Zp(3, 5, 'capped-rel'); R.residue_characteristic() 

3 

""" 

return self.residue_class_field().characteristic() 

 

def defining_polynomial(self, var = 'x'): 

r""" 

Returns the defining polynomial of this local ring, i.e. just ``x``. 

 

INPUT: 

 

- ``self`` -- a local ring 

- ``var`` -- string (default: ``'x'``) the name of the variable 

 

OUTPUT: 

 

- polynomial -- the defining polynomial of this ring as an extension over its ground ring 

 

EXAMPLES:: 

 

sage: R = Zp(3, 3, 'fixed-mod'); R.defining_polynomial('foo') 

(1 + O(3^3))*foo + (O(3^3)) 

""" 

from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing 

return PolynomialRing(self, var).gen() 

 

def ground_ring(self): 

r""" 

Returns ``self``. 

 

Will be overridden by extensions. 

 

INPUT: 

 

- ``self`` -- a local ring 

 

OUTPUT: 

 

- the ground ring of ``self``, i.e., itself 

 

EXAMPLES:: 

 

sage: R = Zp(3, 5, 'fixed-mod') 

sage: S = Zp(3, 4, 'fixed-mod') 

sage: R.ground_ring() is R 

True 

sage: S.ground_ring() is R 

False 

""" 

return self 

 

def ground_ring_of_tower(self): 

r""" 

Returns ``self``. 

 

Will be overridden by extensions. 

 

INPUT: 

 

- ``self`` -- a `p`-adic ring 

 

OUTPUT: 

 

- the ground ring of the tower for ``self``, i.e., itself 

 

EXAMPLES:: 

 

sage: R = Zp(5) 

sage: R.ground_ring_of_tower() 

5-adic Ring with capped relative precision 20 

""" 

return self 

 

def degree(self): 

r""" 

Returns the degree of ``self`` over the ground ring, i.e. 1. 

 

INPUT: 

 

- ``self`` -- a local ring 

 

OUTPUT: 

 

- integer -- the degree of this ring, i.e., 1 

 

EXAMPLES:: 

 

sage: R = Zp(3, 10, 'capped-rel'); R.degree() 

1 

""" 

return Integer(1) 

 

def ramification_index(self, K = None): 

r""" 

Returns the ramification index over the ground ring: 1 unless overridden. 

 

INPUT: 

 

- ``self`` -- a local ring 

 

OUTPUT: 

 

- integer -- the ramification index of this ring: 1 unless overridden. 

 

EXAMPLES:: 

 

sage: R = Zp(3, 5, 'capped-rel'); R.ramification_index() 

1 

""" 

if K is None or K is self: 

return Integer(1) 

else: 

raise ValueError("K should be a subring of self") 

 

def e(self, K = None): 

r""" 

Returns the ramification index over the ground ring: 1 unless overridden. 

 

INPUT: 

 

- ``self`` -- a local ring 

- ``K`` -- a subring of ``self`` (default ``None``) 

 

OUTPUT: 

 

- integer -- the ramification index of this ring: 1 unless overridden. 

 

EXAMPLES:: 

 

sage: R = Zp(3, 5, 'capped-rel'); R.e() 

1 

""" 

return self.ramification_index(K) 

 

def inertia_degree(self, K=None): 

r""" 

Returns the inertia degree over ``K`` (defaults to the ground ring): 1 unless overridden. 

 

INPUT: 

 

- ``self`` -- a local ring 

- ``K`` -- a subring of ``self`` (default None) 

 

OUTPUT: 

 

- integer -- the inertia degree of this ring: 1 unless overridden. 

 

EXAMPLES:: 

 

sage: R = Zp(3, 5, 'capped-rel'); R.inertia_degree() 

1 

""" 

return Integer(1) 

 

def residue_class_degree(self, K=None): 

r""" 

Returns the inertia degree over the ground ring: 1 unless overridden. 

 

INPUT: 

 

- ``self`` -- a local ring 

- ``K`` -- a subring (default ``None``) 

 

OUTPUT: 

 

- integer -- the inertia degree of this ring: 1 unless overridden. 

 

EXAMPLES:: 

 

sage: R = Zp(3, 5, 'capped-rel'); R.residue_class_degree() 

1 

""" 

return self.inertia_degree(K) 

 

def f(self, K=None): 

r""" 

Returns the inertia degree over the ground ring: 1 unless overridden. 

 

INPUT: 

 

- ``self`` -- a local ring 

- ``K`` -- a subring (default ``None``) 

 

OUTPUT: 

 

- integer -- the inertia degree of this ring: 1 unless overridden. 

 

EXAMPLES:: 

 

sage: R = Zp(3, 5, 'capped-rel'); R.f() 

1 

""" 

return self.inertia_degree(K) 

 

def inertia_subring(self): 

r""" 

Returns the inertia subring, i.e. ``self``. 

 

INPUT: 

 

- ``self`` -- a local ring 

 

OUTPUT: 

 

- the inertia subring of self, i.e., itself 

 

EXAMPLES:: 

 

sage: R = Zp(5) 

sage: R.inertia_subring() 

5-adic Ring with capped relative precision 20 

""" 

return self 

 

def maximal_unramified_subextension(self): 

r""" 

Returns the maximal unramified subextension. 

 

INPUT: 

 

- ``self`` -- a local ring 

 

OUTPUT: 

 

- the maximal unramified subextension of ``self`` 

 

EXAMPLES:: 

 

sage: R = Zp(5) 

sage: R.maximal_unramified_subextension() 

5-adic Ring with capped relative precision 20 

""" 

return self.inertia_subring() 

 

# def get_extension(self): 

# r""" 

# Returns the trivial extension of self. 

# """ 

# raise NotImplementedError 

 

def uniformiser(self): 

""" 

Returns a uniformiser for ``self``, ie a generator for the unique maximal ideal. 

 

EXAMPLES:: 

 

sage: R = Zp(5) 

sage: R.uniformiser() 

5 + O(5^21) 

sage: A = Zp(7,10) 

sage: S.<x> = A[] 

sage: B.<t> = A.ext(x^2+7) 

sage: B.uniformiser() 

t + O(t^21) 

""" 

return self.uniformizer() 

 

def uniformiser_pow(self, n): 

""" 

Returns the `n`th power of the uniformiser of ``self`` (as an element of ``self``). 

 

EXAMPLES:: 

 

sage: R = Zp(5) 

sage: R.uniformiser_pow(5) 

5^5 + O(5^25) 

""" 

return self.uniformizer_pow(n) 

 

def is_finite(self): 

r""" 

Returns whether this ring is finite, i.e. ``False``. 

 

INPUT: 

 

- ``self`` -- a `p`-adic ring 

 

OUTPUT: 

 

- boolean -- whether self is finite, i.e., ``False`` 

 

EXAMPLES:: 

 

sage: R = Zp(3, 10,'fixed-mod'); R.is_finite() 

False 

""" 

return False 

 

def ext(self, *args, **kwds): 

""" 

Constructs an extension of self. See ``extension`` for more details. 

 

EXAMPLES:: 

 

sage: A = Zp(7,10) 

sage: S.<x> = A[] 

sage: B.<t> = A.ext(x^2+7) 

sage: B.uniformiser() 

t + O(t^21) 

""" 

return self.extension(*args, **kwds) 

 

def _test_add_bigoh(self, **options): 

r""" 

Perform tests on ``add_bigoh``. 

 

EXAMPLES:: 

 

sage: K = Qp(3) 

sage: K._test_add_bigoh() 

 

""" 

tester = self._tester(**options) 

for x in tester.some_elements(): 

tester.assertEqual(x.add_bigoh(x.precision_absolute()), x) 

from sage.rings.all import infinity 

tester.assertEqual(x.add_bigoh(infinity), x) 

tester.assertEqual(x.add_bigoh(x.precision_absolute()+1), x) 

 

y = x.add_bigoh(0) 

tester.assertIs(y.parent(), self) 

if self.is_capped_absolute(): 

tester.assertEqual(y.precision_absolute(), 0) 

tester.assertEqual(y, self.zero()) 

elif self.is_capped_relative() or self.is_lattice_prec(): 

tester.assertLessEqual(y.precision_absolute(), 0) 

elif self.is_fixed_mod() or self.is_floating_point(): 

tester.assertGreaterEqual((x-y).valuation(), 0) 

else: 

raise NotImplementedError 

 

# if absprec < 0, then the result is in the fraction field (see #13591) 

y = x.add_bigoh(-1) 

tester.assertIs(y.parent(), self.fraction_field()) 

if not self.is_floating_point() and not self.is_fixed_mod(): 

tester.assertLessEqual(y.precision_absolute(), -1) 

 

# make sure that we handle very large values correctly 

if self._prec_type() != 'lattice-float': # in the lattice-float model, there is no cap 

absprec = Integer(2)**1000 

tester.assertEqual(x.add_bigoh(absprec), x) 

 

def _test_residue(self, **options): 

r""" 

Perform some tests on the residue field of this ring. 

 

EXAMPLES:: 

 

sage: R = Zp(2) 

sage: R._test_residue() 

 

""" 

tester = self._tester(**options) 

tester.assertEqual(self.residue_field().characteristic(), self.residue_characteristic()) 

 

for x in tester.some_elements(): 

errors = [] 

if x.precision_absolute() <= 0: 

from .precision_error import PrecisionError 

errors.append(PrecisionError) 

if x.valuation() < 0: 

errors.append(ValueError) 

if errors: 

with tester.assertRaises(tuple(errors)): 

x.residue() 

continue 

y = x.residue() 

# residue() is in `Z/pZ` which is not identical to the residue field `F_p` 

tester.assertEqual(y.parent().cardinality(), self.residue_field().cardinality()) 

z = self(y) 

tester.assertGreater((x-z).valuation(), 0) 

 

for x in self.residue_field().some_elements(): 

y = self(x) 

if x.is_zero(): 

tester.assertGreater(y.valuation(), 0) 

else: 

tester.assertEqual(y.valuation(), 0) 

z = y.residue() 

tester.assertEqual(x, z)