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

""" 

Fixed modulus template for complete discrete valuation rings. 

  

In order to use this template you need to write a linkage file and 

gluing file. For an example see mpz_linkage.pxi (linkage file) and 

padic_fixed_modulus_element.pyx (gluing file). 

  

The linkage file implements a common API that is then used in the 

class FMElement defined here. See sage/libs/linkages/padics/API.pxi 

for the functions needed. 

  

The gluing file does the following: 

  

- ctypedef's celement to be the appropriate type (e.g. mpz_t) 

- includes the linkage file 

- includes this template 

- defines a concrete class inheriting from FMElement, and implements 

any desired extra methods 

  

AUTHORS: 

  

- David Roe (2012-03-01) -- initial version 

""" 

  

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

# Copyright (C) 2007-2012 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/ 

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

  

include "padic_template_element.pxi" 

from cpython.int cimport * 

  

from sage.structure.element cimport Element 

from sage.rings.padics.common_conversion cimport comb_prec, _process_args_and_kwds 

from sage.rings.integer_ring import ZZ 

from sage.rings.rational_field import QQ 

from sage.categories.sets_cat import Sets 

from sage.categories.sets_with_partial_maps import SetsWithPartialMaps 

from sage.categories.homset import Hom 

  

cdef class FMElement(pAdicTemplateElement): 

cdef int _set(self, x, long val, long xprec, absprec, relprec) except -1: 

""" 

Sets the value of this element from given defining data. 

  

This function is intended for use in conversion, and should 

not be called on an element created with :meth:`_new_c`. 

  

INPUT: 

  

- ``x`` -- data defining a `p`-adic element: int, long, 

Integer, Rational, other `p`-adic element... 

  

- ``val`` -- the valuation of the resulting element (unused; 

for compatibility with other `p`-adic precision modes) 

  

- ``xprec -- an inherent precision of ``x`` (unused; for 

compatibility with other `p`-adic precision modes) 

  

- ``absprec`` -- an absolute precision cap for this element 

(unused; for compatibility with other `p`-adic precision 

modes) 

  

- ``relprec`` -- a relative precision cap for this element 

(unused; for compatibility with other `p`-adic precision 

modes) 

  

TESTS:: 

  

sage: R = ZpFM(5) 

sage: a = R(17,5); a #indirect doctest 

2 + 3*5 + O(5^20) 

sage: R = ZpFM(5,5) 

sage: a = R(25/9); a #indirect doctest 

4*5^2 + 2*5^3 + O(5^5) 

""" 

cconstruct(self.value, self.prime_pow) 

if isinstance(x,FMElement) and x.parent() is self.parent(): 

cshift(self.value, (<FMElement>x).value, 0, 0, self.prime_pow, False) 

else: 

cconv(self.value, x, self.prime_pow.prec_cap, 0, self.prime_pow) 

  

cdef FMElement _new_c(self): 

""" 

Creates a new element with the same basic info. 

  

TESTS:: 

  

sage: R = ZpFM(5); R(6) * R(7) #indirect doctest 

2 + 3*5 + 5^2 + O(5^20) 

""" 

cdef type t = type(self) 

cdef FMElement ans = t.__new__(t) 

ans._parent = self._parent 

ans.prime_pow = self.prime_pow 

cconstruct(ans.value, ans.prime_pow) 

return ans 

  

cdef pAdicTemplateElement _new_with_value(self, celement value, long absprec): 

""" 

Creates a new element with a given value and absolute precision. 

  

Used by code that doesn't know the precision type. 

""" 

cdef FMElement ans = self._new_c() 

creduce(ans.value, value, ans.prime_pow.prec_cap, ans.prime_pow) 

return ans 

  

cdef int _get_unit(self, celement value) except -1: 

""" 

Sets ``value`` to the unit of this p-adic element. 

""" 

cremove(value, self.value, 0, self.prime_pow) 

  

cdef int check_preccap(self) except -1: 

""" 

Check that the precision of this element does not exceed the 

precision cap. Does nothing for fixed mod elements. 

  

TESTS:: 

  

sage: ZpFM(5)(1).lift_to_precision(30) # indirect doctest 

1 + O(5^20) 

""" 

pass 

  

def __copy__(self): 

""" 

Return a copy of this element. 

  

EXAMPLES:: 

  

sage: a = ZpFM(5,6)(17); b = copy(a) 

sage: a == b 

True 

sage: a is b 

False 

""" 

cdef FMElement ans = self._new_c() 

ccopy(ans.value, self.value, ans.prime_pow) 

return ans 

  

def __dealloc__(self): 

""" 

Deallocate the underlying data structure. 

  

TESTS:: 

  

sage: R = ZpFM(5) 

sage: a = R(17) 

sage: del(a) 

""" 

cdestruct(self.value, self.prime_pow) 

  

def __reduce__(self): 

""" 

Return a tuple of a function and data that can be used to unpickle this 

element. 

  

EXAMPLES:: 

  

sage: a = ZpFM(5)(-3) 

sage: type(a) 

<type 'sage.rings.padics.padic_fixed_mod_element.pAdicFixedModElement'> 

sage: loads(dumps(a)) == a 

True 

""" 

return unpickle_fme_v2, (self.__class__, self.parent(), cpickle(self.value, self.prime_pow)) 

  

cpdef _neg_(self): 

r""" 

Return the additive inverse of this element. 

  

EXAMPLES:: 

  

sage: R = Zp(7, 4, 'fixed-mod', 'series') 

sage: -R(7) #indirect doctest 

6*7 + 6*7^2 + 6*7^3 + O(7^4) 

""" 

cdef FMElement ans = self._new_c() 

cneg(ans.value, self.value, ans.prime_pow.prec_cap, ans.prime_pow) 

creduce_small(ans.value, ans.value, ans.prime_pow.prec_cap, ans.prime_pow) 

return ans 

  

cpdef _add_(self, _right): 

r""" 

Return the sum of this element and ``_right``. 

  

EXAMPLES:: 

  

sage: R = Zp(7, 4, 'fixed-mod', 'series') 

sage: x = R(1721); x 

6 + 5*7^3 + O(7^4) 

sage: y = R(1373); y 

1 + 4*7^3 + O(7^4) 

sage: x + y #indirect doctest 

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

""" 

cdef FMElement right = _right 

cdef FMElement ans = self._new_c() 

cadd(ans.value, self.value, right.value, ans.prime_pow.prec_cap, ans.prime_pow) 

creduce_small(ans.value, ans.value, ans.prime_pow.prec_cap, ans.prime_pow) 

return ans 

  

cpdef _sub_(self, _right): 

r""" 

Return the difference of this element and ``_right``. 

  

EXAMPLES:: 

  

sage: R = Zp(7, 4, 'fixed-mod', 'series') 

sage: x = R(1721); x 

6 + 5*7^3 + O(7^4) 

sage: y = R(1373); y 

1 + 4*7^3 + O(7^4) 

sage: x - y #indirect doctest 

5 + 7^3 + O(7^4) 

""" 

cdef FMElement right = _right 

cdef FMElement ans = self._new_c() 

csub(ans.value, self.value, right.value, ans.prime_pow.prec_cap, ans.prime_pow) 

creduce_small(ans.value, ans.value, ans.prime_pow.prec_cap, ans.prime_pow) 

return ans 

  

def __invert__(self): 

r""" 

Returns multiplicative inverse of this element. The valuation 

of ``self`` must be zero. 

  

EXAMPLES:: 

  

sage: R = Zp(7, 4, 'fixed-mod', 'series') 

sage: ~R(2) 

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

sage: ~R(0) 

Traceback (most recent call last): 

... 

ValueError: cannot invert non-unit 

sage: ~R(7) 

Traceback (most recent call last): 

... 

ValueError: cannot invert non-unit 

""" 

if not cisunit(self.value, self.prime_pow): 

raise ValueError("cannot invert non-unit") 

cdef FMElement ans = self._new_c() 

cinvert(ans.value, self.value, ans.prime_pow.prec_cap, ans.prime_pow) 

return ans 

  

cpdef _mul_(self, _right): 

r""" 

Return the product of this element and ``_right``. 

  

EXAMPLES:: 

  

sage: R = Zp(7, 4, 'fixed-mod', 'series') 

sage: R(3) * R(2) #indirect doctest 

6 + O(7^4) 

sage: R(1/2) * R(2) 

1 + O(7^4) 

""" 

cdef FMElement right = _right 

cdef FMElement ans = self._new_c() 

cmul(ans.value, self.value, right.value, ans.prime_pow.prec_cap, ans.prime_pow) 

creduce(ans.value, ans.value, ans.prime_pow.prec_cap, ans.prime_pow) 

return ans 

  

cpdef _div_(self, _right): 

r""" 

Return the quotient of this element and ``right``. ``right`` must have 

valuation zero. 

  

EXAMPLES:: 

  

sage: R = Zp(7, 4, 'fixed-mod', 'series') 

sage: R(3) / R(2) #indirect doctest 

5 + 3*7 + 3*7^2 + 3*7^3 + O(7^4) 

sage: R(5) / R(0) 

Traceback (most recent call last): 

... 

ValueError: cannot invert non-unit 

sage: R(7) / R(49) 

Traceback (most recent call last): 

... 

ValueError: cannot invert non-unit 

""" 

cdef FMElement right = _right 

cdef FMElement ans = self._new_c() 

if not cisunit(right.value, self.prime_pow): 

raise ValueError("cannot invert non-unit") 

cdivunit(ans.value, self.value, right.value, ans.prime_pow.prec_cap, ans.prime_pow) 

creduce(ans.value, ans.value, ans.prime_pow.prec_cap, ans.prime_pow) 

return ans 

  

def __pow__(FMElement self, _right, dummy): # NOTE: dummy ignored, always use self.prime_pow.prec_cap 

""" 

Exponentiation by an integer 

  

EXAMPLES:: 

  

sage: R = ZpFM(11, 5) 

sage: R(1/2)^5 

10 + 7*11 + 11^2 + 5*11^3 + 4*11^4 + O(11^5) 

sage: R(1/32) 

10 + 7*11 + 11^2 + 5*11^3 + 4*11^4 + O(11^5) 

sage: R(1/2)^5 == R(1/32) 

True 

sage: R(3)^1000 #indirect doctest 

1 + 4*11^2 + 3*11^3 + 7*11^4 + O(11^5) 

  

TESTS: 

  

We check that :trac:`15640` is resolved:: 

  

sage: R(11)^-1 

Traceback (most recent call last): 

... 

ValueError: cannot invert non-unit 

""" 

cdef FMElement ans = self._new_c() 

cdef Integer right = Integer(_right) 

if right < 0: 

self = ~self 

mpz_neg(right.value, right.value) 

cpow(ans.value, self.value, right.value, self.prime_pow.prec_cap, self.prime_pow) 

return ans 

  

cdef pAdicTemplateElement _lshift_c(self, long shift): 

""" 

Multiplies self by `\pi^{shift}`. 

  

If shift < -self.valuation(), digits will be truncated. See 

:meth:`__rshift__` for details. 

  

EXAMPLES: 

  

We create a fixed modulus ring:: 

  

sage: R = ZpFM(5, 20); a = R(1000); a 

3*5^3 + 5^4 + O(5^20) 

  

Shifting to the right is the same as dividing by a power of 

the uniformizer `\pi` of the `p`-adic ring.:: 

  

sage: a >> 1 

3*5^2 + 5^3 + O(5^20) 

  

Shifting to the left is the same as multiplying by a power of 

`\pi`:: 

  

sage: a << 2 

3*5^5 + 5^6 + O(5^20) 

sage: a*5^2 

3*5^5 + 5^6 + O(5^20) 

  

Shifting by a negative integer to the left is the same as 

right shifting by the absolute value:: 

  

sage: a << -3 

3 + 5 + O(5^20) 

sage: a >> 3 

3 + 5 + O(5^20) 

""" 

if shift < 0: 

return self._rshift_c(-shift) 

elif shift == 0: 

return self 

cdef FMElement ans = self._new_c() 

if shift >= self.prime_pow.prec_cap: 

csetzero(ans.value, ans.prime_pow) 

else: 

cshift(ans.value, self.value, shift, ans.prime_pow.prec_cap, ans.prime_pow, False) 

return ans 

  

cdef pAdicTemplateElement _rshift_c(self, long shift): 

""" 

Divides by `\pi^{shift}`, and truncates. 

  

Note that this operation will insert arbitrary digits (in 

practice, currently all zero) in the least significant digits. 

  

EXAMPLES:: 

  

sage: R = ZpFM(997, 7); a = R(123456878908); a 

964*997 + 572*997^2 + 124*997^3 + O(997^7) 

  

Shifting to the right divides by a power of `\pi`, but 

dropping terms with negative valuation:: 

  

sage: a >> 3 

124 + O(997^7) 

  

A negative shift multiplies by that power of `\pi`:: 

  

sage: a >> -3 

964*997^4 + 572*997^5 + 124*997^6 + O(997^7) 

""" 

if shift < 0: 

return self._lshift_c(-shift) 

elif shift == 0: 

return self 

cdef FMElement ans = self._new_c() 

if shift >= self.prime_pow.prec_cap: 

csetzero(ans.value, ans.prime_pow) 

else: 

cshift(ans.value, self.value, -shift, ans.prime_pow.prec_cap, ans.prime_pow, False) 

return ans 

  

def add_bigoh(self, absprec): 

""" 

Returns a new element truncated modulo `\pi^{\mbox{absprec}}`. 

  

INPUT: 

  

- ``absprec`` -- an integer or infinity 

  

OUTPUT: 

  

- a new element truncated modulo `\pi^{\mbox{absprec}}`. 

  

EXAMPLES:: 

  

sage: R = Zp(7,4,'fixed-mod','series'); a = R(8); a.add_bigoh(1) 

1 + O(7^4) 

  

TESTS: 

  

We handle very large and very small values for ``absprec`` correctly:: 

  

sage: a = R(7) 

sage: a.add_bigoh(2^1000) 

7 + O(7^4) 

sage: a.add_bigoh(-2^1000) 

0 

  

""" 

cdef long aprec 

if absprec is infinity: 

return self 

if isinstance(absprec, int): 

aprec = absprec 

else: 

if not isinstance(absprec, Integer): 

absprec = Integer(absprec) 

if mpz_sgn((<Integer>absprec).value) == -1: 

return self.parent().fraction_field()(0) 

elif mpz_fits_slong_p((<Integer>absprec).value) == 0: 

return self 

else: 

aprec = mpz_get_si((<Integer>absprec).value) 

if aprec < 0: 

return self.parent().fraction_field()(self, absprec) 

elif aprec >= self.prime_pow.prec_cap: 

return self 

cdef FMElement ans = self._new_c() 

creduce(ans.value, self.value, aprec, ans.prime_pow) 

return ans 

  

cpdef bint _is_exact_zero(self) except -1: 

""" 

Tests whether this element is an exact zero, which is always 

False for fixed modulus elements. 

  

EXAMPLES:: 

  

sage: R = Zp(7,4,'fixed-mod','series'); a = R(8); a._is_exact_zero() 

False 

""" 

return False 

  

cpdef bint _is_inexact_zero(self) except -1: 

""" 

Returns True if self is indistinguishable from zero. 

  

EXAMPLES:: 

  

sage: R = ZpFM(7, 5) 

sage: R(14)._is_inexact_zero() 

False 

sage: R(0)._is_inexact_zero() 

True 

""" 

return ciszero(self.value, self.prime_pow) 

  

def is_zero(self, absprec = None): 

r""" 

Returns whether self is zero modulo `\pi^{\mbox{absprec}}`. 

  

INPUT: 

  

- ``absprec`` -- an integer 

  

EXAMPLES:: 

  

sage: R = ZpFM(17, 6) 

sage: R(0).is_zero() 

True 

sage: R(17^6).is_zero() 

True 

sage: R(17^2).is_zero(absprec=2) 

True 

""" 

cdef bint iszero = ciszero(self.value, self.prime_pow) 

if absprec is None: 

return iszero 

if not isinstance(absprec, Integer): 

absprec = Integer(absprec) 

if mpz_cmp_si((<Integer>absprec).value, self.prime_pow.prec_cap) >= 0: 

return iszero 

cdef long val = self.valuation_c() 

return mpz_cmp_si((<Integer>absprec).value, val) <= 0 

  

def __nonzero__(self): 

""" 

Returns True if this element is distinguishable from zero. 

  

For most applications, explicitly specifying the power of p 

modulo which the element is supposed to be nonzero is 

preferrable. 

  

EXAMPLES:: 

  

sage: R = ZpFM(5); a = R(0); b = R(75) 

sage: bool(a), bool(b) # indirect doctest 

(False, True) 

""" 

return not ciszero(self.value, self.prime_pow) 

  

def is_equal_to(self, _right, absprec=None): 

r""" 

Returns whether this element is equal to ``right`` modulo `p^{\mbox{absprec}}`. 

  

If ``absprec`` is ``None``, returns if ``self == 0``. 

  

INPUT: 

  

- ``right`` -- a p-adic element with the same parent 

- ``absprec`` -- a positive integer or ``None`` (default: ``None``) 

  

EXAMPLES:: 

  

sage: R = ZpFM(2, 6) 

sage: R(13).is_equal_to(R(13)) 

True 

sage: R(13).is_equal_to(R(13+2^10)) 

True 

sage: R(13).is_equal_to(R(17), 2) 

True 

sage: R(13).is_equal_to(R(17), 5) 

False 

""" 

cdef FMElement right 

cdef long aprec, rprec, sval, rval 

if self.parent() is _right.parent(): 

right = _right 

else: 

right = self.parent()(_right) 

if absprec is None: 

# The default absolute precision is given by the precision cap 

aprec = self.prime_pow.prec_cap 

else: 

if not isinstance(absprec, Integer): 

absprec = Integer(absprec) 

# If absprec is not positive, then self and right are always 

# equal. 

if mpz_sgn((<Integer>absprec).value) < 0: 

return True 

# If absprec is bigger than the precision cap, we use it 

# instead. 

if mpz_cmp_si((<Integer>absprec).value, self.prime_pow.prec_cap) >= 0: 

aprec = self.prime_pow.prec_cap 

else: 

aprec = mpz_get_si((<Integer>absprec).value) 

return ccmp(self.value, 

right.value, 

aprec, 

aprec < self.prime_pow.prec_cap, 

aprec < right.prime_pow.prec_cap, 

self.prime_pow) == 0 

  

cdef int _cmp_units(self, pAdicGenericElement _right) except -2: 

""" 

Comparison of units, used in equality testing. 

  

EXAMPLES:: 

  

sage: R = ZpFM(5) 

sage: a = R(17); b = R(0,3); c = R(85,7); d = R(2, 1) 

sage: any([a == b, a == c, b == c, b == d, c == d, a == d]) # indirect doctest 

False 

sage: all([a == a, b == b, c == c, d == d]) 

True 

""" 

cdef FMElement right = _right 

return ccmp(self.value, right.value, self.prime_pow.prec_cap, False, False, self.prime_pow) 

  

cdef pAdicTemplateElement lift_to_precision_c(self, long absprec): 

""" 

Lifts this element to another with precision at least absprec. 

  

Since fixed modulus elements don't track precision, this 

function just returns the same element. 

  

EXAMPLES:: 

  

sage: R = ZpFM(5); 

sage: a = R(77, 2); a 

2 + 3*5^2 + O(5^20) 

sage: a.lift_to_precision(17) # indirect doctest 

2 + 3*5^2 + O(5^20) 

""" 

return self 

  

def _teichmuller_set_unsafe(self): 

""" 

Sets this element to the Teichmuller representative with the 

same residue. 

  

.. WARNING:: 

  

This function modifies the element, which is not safe. 

Elements are supposed to be immutable. 

  

EXAMPLES:: 

  

sage: R = ZpFM(17,5); a = R(11) 

sage: a 

11 + O(17^5) 

sage: a._teichmuller_set_unsafe(); a 

11 + 14*17 + 2*17^2 + 12*17^3 + 15*17^4 + O(17^5) 

sage: E = a.expansion(lift_mode='teichmuller'); E 

17-adic expansion of 11 + 14*17 + 2*17^2 + 12*17^3 + 15*17^4 + O(17^5) (teichmuller) 

sage: list(E) 

[11 + 14*17 + 2*17^2 + 12*17^3 + 15*17^4 + O(17^5), O(17^5), O(17^5), O(17^5), O(17^5)] 

  

Note that if you set an element which is congruent to 0 you 

get 0 to maximum precision:: 

  

sage: b = R(17*5); b 

5*17 + O(17^5) 

sage: b._teichmuller_set_unsafe(); b 

O(17^5) 

""" 

if cisunit(self.value, self.prime_pow): 

cteichmuller(self.value, self.value, self.prime_pow.prec_cap, self.prime_pow) 

else: 

csetzero(self.value, self.prime_pow) 

  

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

""" 

Returns a polynomial over the base ring that yields this element 

when evaluated at the generator of the parent. 

  

INPUT: 

  

- ``var`` -- string, the variable name for the polynomial 

  

EXAMPLES:: 

  

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

sage: a.polynomial() 

(1 + O(5^20))*x + (O(5^20)) 

sage: a.polynomial(var='y') 

(1 + O(5^20))*y + (O(5^20)) 

sage: (5*a^2 + 25).polynomial() 

(5 + O(5^20))*x^2 + (O(5^20))*x + (5^2 + O(5^20)) 

""" 

R = self.base_ring() 

S = R[var] 

prec = self.precision_absolute() 

e = self.parent().e() 

L = ccoefficients(self.value, 0, self.prime_pow.prec_cap, self.prime_pow) 

if e == 1: 

L = [R(c, prec) for c in L] 

else: 

L = [R(c, (prec - i - 1) // e + 1) for i, c in enumerate(L)] 

return S(L) 

  

def precision_absolute(self): 

""" 

The absolute precision of this element. 

  

EXAMPLES:: 

  

sage: R = Zp(7,4,'fixed-mod'); a = R(7); a.precision_absolute() 

4 

""" 

cdef Integer ans = Integer.__new__(Integer) 

mpz_set_si(ans.value, self.prime_pow.prec_cap) 

return ans 

  

def precision_relative(self): 

r""" 

The relative precision of this element. 

  

EXAMPLES:: 

  

sage: R = Zp(7,4,'fixed-mod'); a = R(7); a.precision_relative() 

3 

sage: a = R(0); a.precision_relative() 

0 

""" 

cdef Integer ans = Integer.__new__(Integer) 

mpz_set_si(ans.value, self.prime_pow.prec_cap - self.valuation_c()) 

return ans 

  

cpdef pAdicTemplateElement unit_part(FMElement self): 

r""" 

Returns the unit part of self. 

  

If the valuation of self is positive, then the high digits of the 

result will be zero. 

  

EXAMPLES:: 

  

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

sage: R(5).unit_part() 

5 + O(17^4) 

sage: R(18*17).unit_part() 

1 + 17 + O(17^4) 

sage: R(0).unit_part() 

O(17^4) 

sage: type(R(5).unit_part()) 

<type 'sage.rings.padics.padic_fixed_mod_element.pAdicFixedModElement'> 

sage: R = ZpFM(5, 5); a = R(75); a.unit_part() 

3 + O(5^5) 

""" 

cdef FMElement ans = (<FMElement>self)._new_c() 

cremove(ans.value, (<FMElement>self).value, (<FMElement>self).prime_pow.prec_cap, (<FMElement>self).prime_pow) 

return ans 

  

cdef long valuation_c(self): 

""" 

Returns the valuation of this element. 

  

TESTS:: 

  

sage: R = ZpFM(5, 5); R(0).valuation() #indirect doctest 

5 

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

sage: a = R(2*17^2) 

sage: a.valuation() 

2 

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

sage: R(0).valuation() 

4 

sage: R(1).valuation() 

0 

sage: R(2).valuation() 

0 

sage: R(5).valuation() 

1 

sage: R(10).valuation() 

1 

sage: R(25).valuation() 

2 

sage: R(50).valuation() 

2 

""" 

# for backward compatibility 

return cvaluation(self.value, self.prime_pow.prec_cap, self.prime_pow) 

  

cpdef val_unit(self): 

""" 

Returns a 2-tuple, the first element set to the valuation of 

self, and the second to the unit part of self. 

  

If self == 0, then the unit part is O(p^self.parent().precision_cap()). 

  

EXAMPLES:: 

  

sage: R = ZpFM(5,5) 

sage: a = R(75); b = a - a 

sage: a.val_unit() 

(2, 3 + O(5^5)) 

sage: b.val_unit() 

(5, O(5^5)) 

""" 

cdef FMElement unit = self._new_c() 

cdef Integer valuation = Integer.__new__(Integer) 

mpz_set_si(valuation.value, cremove(unit.value, self.value, self.prime_pow.prec_cap, self.prime_pow)) 

return valuation, unit 

  

def __hash__(self): 

""" 

Hashing. 

  

EXAMPLES:: 

  

sage: R = ZpFM(11, 5) 

sage: hash(R(3)) == hash(3) 

True 

""" 

return chash(self.value, 0, self.prime_pow.prec_cap, self.prime_pow) 

  

cdef class pAdicCoercion_ZZ_FM(RingHomomorphism): 

""" 

The canonical inclusion from ZZ to a fixed modulus ring. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ); f 

Ring morphism: 

From: Integer Ring 

To: 5-adic Ring of fixed modulus 5^20 

  

TESTS:: 

  

sage: TestSuite(f).run() 

  

""" 

def __init__(self, R): 

""" 

Initialization. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ); type(f) 

<type 'sage.rings.padics.padic_fixed_mod_element.pAdicCoercion_ZZ_FM'> 

""" 

RingHomomorphism.__init__(self, ZZ.Hom(R)) 

self._zero = R.element_class(R, 0) 

self._section = pAdicConvert_FM_ZZ(R) 

  

cdef dict _extra_slots(self): 

""" 

Helper for copying and pickling. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ) 

sage: g = copy(f) # indirect doctest 

sage: g == f 

True 

sage: g(6) 

1 + 5 + O(5^20) 

sage: g(6) == f(6) 

True 

""" 

_slots = RingHomomorphism._extra_slots(self) 

_slots['_zero'] = self._zero 

_slots['_section'] = self.section() # use method since it copies coercion-internal sections. 

return _slots 

  

cdef _update_slots(self, dict _slots): 

""" 

Helper for copying and pickling. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ) 

sage: g = copy(f) # indirect doctest 

sage: g == f 

True 

sage: g(6) 

1 + 5 + O(5^20) 

sage: g(6) == f(6) 

True 

""" 

self._zero = _slots['_zero'] 

self._section = _slots['_section'] 

RingHomomorphism._update_slots(self, _slots) 

  

cpdef Element _call_(self, x): 

""" 

Evaluation. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ) 

sage: f(0).parent() 

5-adic Ring of fixed modulus 5^20 

sage: f(5) 

5 + O(5^20) 

""" 

if mpz_sgn((<Integer>x).value) == 0: 

return self._zero 

cdef FMElement ans = self._zero._new_c() 

cconv_mpz_t(ans.value, (<Integer>x).value, ans.prime_pow.prec_cap, True, ans.prime_pow) 

return ans 

  

cpdef Element _call_with_args(self, x, args=(), kwds={}): 

""" 

This function is used when some precision cap is passed in (relative or absolute or both). 

  

INPUT: 

  

- ``x`` -- an Integer 

  

- ``absprec``, or the first positional argument -- the maximum 

absolute precision (unused for fixed modulus elements). 

  

- ``relprec``, or the second positional argument -- the 

maximum relative precision (unused for fixed modulus 

elements) 

  

EXAMPLES:: 

  

sage: R = ZpFM(5,4) 

sage: type(R(10,2)) 

<type 'sage.rings.padics.padic_fixed_mod_element.pAdicFixedModElement'> 

sage: R(30,2) 

5 + 5^2 + O(5^4) 

sage: R(30,3,1) 

5 + 5^2 + O(5^4) 

sage: R(30,absprec=2) 

5 + 5^2 + O(5^4) 

sage: R(30,relprec=2) 

5 + 5^2 + O(5^4) 

sage: R(30,absprec=1) 

5 + 5^2 + O(5^4) 

sage: R(30,empty=True) 

5 + 5^2 + O(5^4) 

""" 

if mpz_sgn((<Integer>x).value) == 0: 

return self._zero 

cdef FMElement ans = self._zero._new_c() 

cconv_mpz_t(ans.value, (<Integer>x).value, ans.prime_pow.prec_cap, True, ans.prime_pow) 

return ans 

  

def section(self): 

""" 

Returns a map back to ZZ that approximates an element of this 

`p`-adic ring by an integer. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ).section() 

sage: f(ZpFM(5)(-1)) - 5^20 

-1 

""" 

from sage.misc.constant_function import ConstantFunction 

if not isinstance(self._section.domain, ConstantFunction): 

import copy 

self._section = copy.copy(self._section) 

return self._section 

  

cdef class pAdicConvert_FM_ZZ(RingMap): 

""" 

The map from a fixed modulus ring back to ZZ that returns the smallest 

non-negative integer approximation to its input which is accurate up to the precision. 

  

If the input is not in the closure of the image of ZZ, raises a ValueError. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ).section(); f 

Set-theoretic ring morphism: 

From: 5-adic Ring of fixed modulus 5^20 

To: Integer Ring 

""" 

def __init__(self, R): 

""" 

Initialization. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ).section(); type(f) 

<type 'sage.rings.padics.padic_fixed_mod_element.pAdicConvert_FM_ZZ'> 

sage: f.category() 

Category of homsets of sets 

""" 

if R.degree() > 1 or R.characteristic() != 0 or R.residue_characteristic() == 0: 

RingMap.__init__(self, Hom(R, ZZ, SetsWithPartialMaps())) 

else: 

RingMap.__init__(self, Hom(R, ZZ, Sets())) 

  

cpdef Element _call_(self, _x): 

""" 

Evaluation. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).coerce_map_from(ZZ).section() 

sage: f(ZpFM(5)(-1)) - 5^20 

-1 

sage: f(ZpFM(5)(0)) 

0 

""" 

cdef Integer ans = Integer.__new__(Integer) 

cdef FMElement x = _x 

cconv_mpz_t_out(ans.value, x.value, 0, x.prime_pow.prec_cap, x.prime_pow) 

return ans 

  

cdef class pAdicConvert_QQ_FM(Morphism): 

""" 

The inclusion map from QQ to a fixed modulus ring that is defined 

on all elements with non-negative p-adic valuation. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).convert_map_from(QQ); f 

Generic morphism: 

From: Rational Field 

To: 5-adic Ring of fixed modulus 5^20 

""" 

def __init__(self, R): 

""" 

Initialization. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).convert_map_from(QQ); type(f) 

<type 'sage.rings.padics.padic_fixed_mod_element.pAdicConvert_QQ_FM'> 

""" 

Morphism.__init__(self, Hom(QQ, R, SetsWithPartialMaps())) 

self._zero = R.element_class(R, 0) 

  

cdef dict _extra_slots(self): 

""" 

Helper for copying and pickling. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).convert_map_from(QQ) 

sage: g = copy(f) # indirect doctest 

sage: g == f # todo: comparison not implemented 

True 

sage: g(1/6) 

1 + 4*5 + 4*5^3 + 4*5^5 + 4*5^7 + 4*5^9 + 4*5^11 + 4*5^13 + 4*5^15 + 4*5^17 + 4*5^19 + O(5^20) 

sage: g(1/6) == f(1/6) 

True 

""" 

_slots = Morphism._extra_slots(self) 

_slots['_zero'] = self._zero 

return _slots 

  

cdef _update_slots(self, dict _slots): 

""" 

Helper for copying and pickling. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5).convert_map_from(QQ) 

sage: g = copy(f) # indirect doctest 

sage: g == f # todo: comparison not implemented 

True 

sage: g(1/6) 

1 + 4*5 + 4*5^3 + 4*5^5 + 4*5^7 + 4*5^9 + 4*5^11 + 4*5^13 + 4*5^15 + 4*5^17 + 4*5^19 + O(5^20) 

sage: g(1/6) == f(1/6) 

True 

""" 

self._zero = _slots['_zero'] 

Morphism._update_slots(self, _slots) 

  

cpdef Element _call_(self, x): 

""" 

Evaluation. 

  

EXAMPLES:: 

  

sage: f = ZpFM(5,4).convert_map_from(QQ) 

sage: f(1/7) 

3 + 3*5 + 2*5^3 + O(5^4) 

sage: f(0) 

O(5^4) 

""" 

if mpq_sgn((<Rational>x).value) == 0: 

return self._zero 

cdef FMElement ans = self._zero._new_c() 

cconv_mpq_t(ans.value, (<Rational>x).value, ans.prime_pow.prec_cap, True, ans.prime_pow) 

return ans 

  

cpdef Element _call_with_args(self, x, args=(), kwds={}): 

""" 

This function is used when some precision cap is passed in (relative or absolute or both). 

  

INPUT: 

  

- ``x`` -- a Rational 

  

- ``absprec``, or the first positional argument -- the maximum 

absolute precision (unused for fixed modulus elements). 

  

- ``relprec``, or the second positional argument -- the 

maximum relative precision (unused for fixed modulus 

elements) 

  

EXAMPLES:: 

  

sage: R = ZpFM(5,4) 

sage: type(R(1/7,2)) 

<type 'sage.rings.padics.padic_fixed_mod_element.pAdicFixedModElement'> 

sage: R(1/7,2) 

3 + 3*5 + 2*5^3 + O(5^4) 

sage: R(1/7,3,1) 

3 + 3*5 + 2*5^3 + O(5^4) 

sage: R(1/7,absprec=2) 

3 + 3*5 + 2*5^3 + O(5^4) 

sage: R(1/7,relprec=2) 

3 + 3*5 + 2*5^3 + O(5^4) 

sage: R(1/7,absprec=1) 

3 + 3*5 + 2*5^3 + O(5^4) 

sage: R(1/7,empty=True) 

3 + 3*5 + 2*5^3 + O(5^4) 

""" 

if mpq_sgn((<Rational>x).value) == 0: 

return self._zero 

cdef FMElement ans = self._zero._new_c() 

cconv_mpq_t(ans.value, (<Rational>x).value, ans.prime_pow.prec_cap, True, ans.prime_pow) 

return ans 

  

cdef class pAdicCoercion_FM_frac_field(RingHomomorphism): 

""" 

The canonical inclusion of Zq into its fraction field. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27, implementation='FLINT') 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R); f 

Ring morphism: 

From: Unramified Extension in a defined by x^3 + 2*x + 1 of fixed modulus 3^20 over 3-adic Ring 

To: Unramified Extension in a defined by x^3 + 2*x + 1 with floating precision 20 over 3-adic Field 

  

TESTS:: 

  

sage: TestSuite(f).run() 

  

""" 

def __init__(self, R, K): 

""" 

Initialization. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R); type(f) 

<type 'sage.rings.padics.qadic_flint_FM.pAdicCoercion_FM_frac_field'> 

""" 

RingHomomorphism.__init__(self, R.Hom(K)) 

self._zero = K(0) 

self._section = pAdicConvert_FM_frac_field(K, R) 

  

cpdef Element _call_(self, _x): 

""" 

Evaluation. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R) 

sage: f(a) 

a 

""" 

cdef FMElement x = _x 

if ciszero(x.value, x.prime_pow): 

return self._zero 

cdef FPElement ans = self._zero._new_c() 

ans.ordp = cremove(ans.unit, x.value, x.prime_pow.ram_prec_cap, x.prime_pow) 

return ans 

  

cpdef Element _call_with_args(self, _x, args=(), kwds={}): 

""" 

This function is used when some precision cap is passed in 

(relative or absolute or both). 

  

See the documentation for 

:meth:`pAdicCappedAbsoluteElement.__init__` for more details. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R) 

sage: f(a, 3) 

a 

sage: b = 117*a 

sage: f(b, 3) 

a*3^2 

sage: f(b, 4, 1) 

a*3^2 

sage: f(b, 4, 3) 

a*3^2 + a*3^3 

sage: f(b, absprec=4) 

a*3^2 + a*3^3 

sage: f(b, relprec=3) 

a*3^2 + a*3^3 + a*3^4 

sage: f(b, absprec=1) 

0 

sage: f(R(0)) 

0 

""" 

cdef long aprec, rprec 

cdef FMElement x = _x 

if ciszero(x.value, x.prime_pow): 

return self._zero 

cdef FPElement ans = self._zero._new_c() 

cdef bint reduce = False 

_process_args_and_kwds(&aprec, &rprec, args, kwds, False, x.prime_pow) 

ans.ordp = cremove(ans.unit, x.value, aprec, x.prime_pow) 

if aprec < ans.ordp + rprec: 

rprec = aprec - ans.ordp 

if rprec <= 0: 

return self._zero 

creduce(ans.unit, ans.unit, rprec, x.prime_pow) 

return ans 

  

def section(self): 

""" 

Returns a map back to the ring that converts elements of 

non-negative valuation. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R) 

sage: f.section()(K.gen()) 

a + O(3^20) 

""" 

from sage.misc.constant_function import ConstantFunction 

if not isinstance(self._section.domain, ConstantFunction): 

import copy 

self._section = copy.copy(self._section) 

return self._section 

  

cdef dict _extra_slots(self): 

""" 

Helper for copying and pickling. 

  

TESTS:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R) 

sage: g = copy(f) # indirect doctest 

sage: g 

Ring morphism: 

From: Unramified Extension in a defined by x^3 + 2*x + 1 of fixed modulus 3^20 over 3-adic Ring 

To: Unramified Extension in a defined by x^3 + 2*x + 1 with floating precision 20 over 3-adic Field 

sage: g == f 

True 

sage: g is f 

False 

sage: g(a) 

a 

sage: g(a) == f(a) 

True 

""" 

_slots = RingHomomorphism._extra_slots(self) 

_slots['_zero'] = self._zero 

_slots['_section'] = self.section() # use method since it copies coercion-internal sections. 

return _slots 

  

cdef _update_slots(self, dict _slots): 

""" 

Helper for copying and pickling. 

  

TESTS:: 

  

sage: R.<a> = ZqFM(9) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R) 

sage: g = copy(f) # indirect doctest 

sage: g 

Ring morphism: 

From: Unramified Extension in a defined by x^2 + 2*x + 2 of fixed modulus 3^20 over 3-adic Ring 

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

sage: g == f 

True 

sage: g is f 

False 

sage: g(a) 

a 

sage: g(a) == f(a) 

True 

  

""" 

self._zero = _slots['_zero'] 

self._section = _slots['_section'] 

RingHomomorphism._update_slots(self, _slots) 

  

def is_injective(self): 

r""" 

Return whether this map is injective. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(9) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R) 

sage: f.is_injective() 

True 

  

""" 

return True 

  

def is_surjective(self): 

r""" 

Return whether this map is surjective. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(9) 

sage: K = R.fraction_field() 

sage: f = K.coerce_map_from(R) 

sage: f.is_surjective() 

False 

  

""" 

return False 

  

  

cdef class pAdicConvert_FM_frac_field(Morphism): 

""" 

The section of the inclusion from `\ZZ_q`` to its fraction field. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = R.convert_map_from(K); f 

Generic morphism: 

From: Unramified Extension in a defined by x^3 + 2*x + 1 with floating precision 20 over 3-adic Field 

To: Unramified Extension in a defined by x^3 + 2*x + 1 of fixed modulus 3^20 over 3-adic Ring 

""" 

def __init__(self, K, R): 

""" 

Initialization. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = R.convert_map_from(K); type(f) 

<type 'sage.rings.padics.qadic_flint_FM.pAdicConvert_FM_frac_field'> 

""" 

Morphism.__init__(self, Hom(K, R, SetsWithPartialMaps())) 

self._zero = R(0) 

  

cpdef Element _call_(self, _x): 

""" 

Evaluation. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = R.convert_map_from(K) 

sage: f(K.gen()) 

a + O(3^20) 

""" 

cdef FPElement x = _x 

if x.ordp < 0: raise ValueError("negative valuation") 

if x.ordp >= self._zero.prime_pow.ram_prec_cap: 

return self._zero 

cdef FMElement ans = self._zero._new_c() 

cshift(ans.value, x.unit, x.ordp, ans.prime_pow.ram_prec_cap, ans.prime_pow, x.ordp > 0) 

return ans 

  

cpdef Element _call_with_args(self, _x, args=(), kwds={}): 

""" 

This function is used when some precision cap is passed in 

(relative or absolute or both). 

  

See the documentation for 

:meth:`pAdicCappedAbsoluteElement.__init__` for more details. 

  

EXAMPLES:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = R.convert_map_from(K); a = K(a) 

sage: f(a, 3) 

a + O(3^20) 

sage: b = 117*a 

sage: f(b, 3) 

a*3^2 + O(3^20) 

sage: f(b, 4, 1) 

a*3^2 + O(3^20) 

sage: f(b, 4, 3) 

a*3^2 + a*3^3 + O(3^20) 

sage: f(b, absprec=4) 

a*3^2 + a*3^3 + O(3^20) 

sage: f(b, relprec=3) 

a*3^2 + a*3^3 + a*3^4 + O(3^20) 

sage: f(b, absprec=1) 

O(3^20) 

sage: f(K(0)) 

O(3^20) 

""" 

cdef long aprec, rprec 

cdef FPElement x = _x 

if x.ordp < 0: raise ValueError("negative valuation") 

if x.ordp >= self._zero.prime_pow.ram_prec_cap: 

return self._zero 

cdef FMElement ans = self._zero._new_c() 

_process_args_and_kwds(&aprec, &rprec, args, kwds, True, ans.prime_pow) 

if rprec < aprec - x.ordp: 

aprec = x.ordp + rprec 

sig_on() 

cshift(ans.value, x.unit, x.ordp, aprec, ans.prime_pow, x.ordp > 0) 

sig_off() 

return ans 

  

cdef dict _extra_slots(self): 

""" 

Helper for copying and pickling. 

  

TESTS:: 

  

sage: R.<a> = ZqFM(27) 

sage: K = R.fraction_field() 

sage: f = R.convert_map_from(K) 

sage: a = K(a) 

sage: g = copy(f) # indirect doctest 

sage: g 

Generic morphism: 

From: Unramified Extension in a defined by x^3 + 2*x + 1 with floating precision 20 over 3-adic Field 

To: Unramified Extension in a defined by x^3 + 2*x + 1 of fixed modulus 3^20 over 3-adic Ring 

sage: g == f 

True 

sage: g is f 

False 

sage: g(a) 

a + O(3^20) 

sage: g(a) == f(a) 

True 

""" 

_slots = Morphism._extra_slots(self) 

_slots['_zero'] = self._zero 

return _slots 

  

cdef _update_slots(self, dict _slots): 

""" 

Helper for copying and pickling. 

  

TESTS:: 

  

sage: R.<a> = ZqFM(9) 

sage: K = R.fraction_field() 

sage: f = R.convert_map_from(K) 

sage: a = f(a) 

sage: g = copy(f) # indirect doctest 

sage: g 

Generic morphism: 

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

To: Unramified Extension in a defined by x^2 + 2*x + 2 of fixed modulus 3^20 over 3-adic Ring 

sage: g == f 

True 

sage: g is f 

False 

sage: g(a) 

a + O(3^20) 

sage: g(a) == f(a) 

True 

  

""" 

self._zero = _slots['_zero'] 

Morphism._update_slots(self, _slots) 

  

def unpickle_fme_v2(cls, parent, value): 

""" 

Unpickles a fixed-mod element. 

  

EXAMPLES:: 

  

sage: from sage.rings.padics.padic_fixed_mod_element import pAdicFixedModElement, unpickle_fme_v2 

sage: R = ZpFM(5) 

sage: a = unpickle_fme_v2(pAdicFixedModElement, R, 17*25); a 

2*5^2 + 3*5^3 + O(5^20) 

sage: a.parent() is R 

True 

""" 

cdef FMElement ans = cls.__new__(cls) 

ans._parent = parent 

ans.prime_pow = <PowComputer_?>parent.prime_pow 

cconstruct(ans.value, ans.prime_pow) 

cunpickle(ans.value, value, ans.prime_pow) 

return ans