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

1622

1623

1624

1625

1626

1627

1628

1629

1630

1631

1632

1633

1634

1635

1636

1637

1638

1639

1640

1641

1642

1643

1644

1645

1646

1647

1648

1649

1650

1651

1652

1653

1654

1655

1656

1657

1658

1659

1660

1661

1662

1663

1664

1665

1666

1667

1668

1669

1670

1671

1672

1673

1674

1675

1676

1677

1678

1679

1680

1681

1682

1683

1684

1685

1686

1687

1688

1689

1690

1691

1692

1693

1694

1695

1696

1697

1698

1699

1700

1701

1702

1703

1704

1705

1706

1707

1708

1709

1710

1711

1712

1713

1714

1715

1716

1717

1718

1719

1720

1721

1722

1723

1724

1725

1726

1727

1728

1729

1730

1731

1732

1733

1734

1735

1736

1737

1738

1739

1740

1741

1742

1743

1744

1745

1746

1747

1748

1749

1750

1751

1752

1753

1754

1755

1756

1757

1758

1759

1760

1761

1762

1763

1764

1765

1766

1767

1768

1769

1770

1771

1772

1773

1774

1775

1776

1777

1778

1779

1780

1781

1782

1783

1784

1785

1786

1787

1788

1789

1790

1791

1792

1793

1794

1795

1796

1797

1798

1799

1800

1801

1802

1803

1804

1805

1806

1807

1808

1809

1810

1811

1812

1813

1814

1815

1816

1817

1818

1819

1820

1821

1822

1823

1824

1825

1826

1827

1828

1829

1830

1831

1832

1833

1834

1835

1836

1837

1838

1839

1840

1841

1842

1843

1844

1845

1846

1847

1848

1849

1850

1851

1852

1853

1854

1855

1856

1857

1858

1859

1860

1861

1862

1863

1864

1865

1866

1867

1868

1869

1870

1871

1872

1873

1874

1875

1876

1877

1878

1879

1880

1881

1882

1883

1884

1885

1886

1887

1888

1889

1890

1891

1892

1893

1894

1895

1896

1897

1898

1899

1900

1901

1902

1903

1904

1905

1906

1907

1908

1909

1910

1911

1912

1913

1914

1915

1916

1917

1918

1919

1920

1921

1922

1923

1924

1925

1926

1927

1928

1929

1930

1931

1932

1933

1934

1935

1936

1937

1938

1939

1940

1941

1942

1943

1944

1945

1946

1947

1948

1949

1950

1951

1952

1953

1954

1955

1956

1957

1958

1959

1960

1961

1962

1963

1964

1965

1966

1967

1968

1969

1970

1971

1972

1973

1974

1975

1976

1977

1978

1979

1980

1981

1982

1983

1984

1985

1986

1987

1988

1989

1990

1991

1992

1993

1994

1995

1996

1997

1998

1999

2000

2001

2002

2003

2004

2005

2006

2007

2008

2009

2010

2011

2012

2013

2014

2015

2016

2017

2018

2019

2020

2021

2022

2023

2024

2025

2026

2027

2028

2029

2030

2031

2032

2033

2034

2035

2036

2037

2038

2039

2040

2041

2042

2043

2044

2045

2046

2047

2048

2049

2050

2051

2052

2053

2054

2055

2056

2057

2058

2059

2060

2061

2062

2063

2064

2065

2066

2067

2068

2069

2070

2071

2072

2073

2074

2075

2076

2077

2078

2079

2080

2081

2082

2083

2084

2085

2086

2087

2088

2089

2090

2091

2092

2093

2094

2095

2096

2097

2098

2099

2100

2101

2102

2103

2104

2105

2106

2107

2108

2109

2110

2111

2112

2113

2114

2115

2116

2117

2118

2119

2120

2121

2122

2123

2124

2125

2126

2127

2128

2129

2130

2131

2132

2133

2134

2135

2136

2137

2138

2139

2140

2141

2142

2143

2144

2145

2146

2147

2148

2149

2150

2151

2152

2153

2154

2155

2156

2157

2158

2159

2160

2161

2162

2163

2164

2165

2166

2167

2168

2169

2170

2171

2172

2173

2174

2175

2176

2177

2178

2179

2180

2181

2182

2183

2184

2185

2186

2187

2188

2189

2190

2191

2192

2193

2194

2195

2196

2197

2198

2199

2200

2201

2202

2203

2204

2205

2206

2207

2208

2209

2210

2211

2212

2213

2214

2215

2216

2217

2218

2219

2220

2221

2222

2223

2224

2225

2226

2227

2228

2229

2230

2231

2232

2233

2234

2235

2236

2237

2238

2239

2240

2241

2242

2243

2244

2245

2246

2247

2248

2249

2250

2251

2252

2253

2254

2255

2256

2257

2258

2259

2260

2261

2262

2263

2264

2265

2266

2267

2268

2269

2270

2271

2272

2273

2274

2275

2276

2277

2278

2279

2280

2281

2282

2283

2284

2285

2286

2287

2288

2289

2290

2291

2292

2293

2294

2295

2296

2297

2298

2299

2300

2301

2302

2303

2304

2305

2306

2307

2308

2309

2310

2311

2312

2313

2314

2315

2316

2317

2318

2319

2320

2321

2322

2323

2324

2325

2326

2327

2328

2329

2330

2331

2332

2333

2334

2335

2336

2337

2338

2339

2340

2341

2342

2343

2344

2345

2346

2347

2348

2349

2350

2351

2352

2353

2354

2355

2356

2357

2358

2359

2360

2361

2362

2363

2364

2365

2366

2367

2368

2369

2370

2371

2372

2373

2374

2375

2376

2377

2378

2379

2380

2381

2382

2383

2384

2385

2386

2387

2388

2389

2390

2391

2392

2393

2394

2395

2396

2397

2398

2399

2400

2401

2402

2403

2404

2405

2406

2407

2408

2409

2410

2411

2412

2413

2414

2415

2416

2417

2418

2419

2420

2421

2422

2423

2424

2425

2426

2427

2428

2429

2430

2431

2432

2433

2434

2435

2436

2437

2438

2439

2440

2441

2442

2443

2444

2445

2446

2447

2448

2449

2450

2451

2452

2453

2454

2455

2456

2457

2458

2459

2460

2461

2462

2463

2464

2465

2466

2467

2468

2469

2470

2471

2472

2473

2474

2475

2476

2477

2478

2479

2480

2481

2482

2483

2484

2485

2486

2487

2488

2489

2490

2491

2492

2493

2494

2495

2496

2497

2498

2499

2500

2501

2502

2503

2504

2505

2506

2507

2508

2509

2510

2511

2512

2513

2514

2515

2516

2517

2518

2519

2520

2521

2522

2523

2524

2525

2526

2527

2528

2529

2530

2531

2532

2533

2534

2535

2536

2537

2538

2539

2540

2541

2542

2543

2544

2545

2546

2547

2548

2549

2550

2551

2552

2553

2554

2555

2556

2557

2558

2559

2560

2561

2562

2563

2564

2565

2566

2567

2568

2569

2570

2571

2572

2573

2574

2575

2576

2577

2578

2579

2580

2581

2582

2583

2584

2585

2586

2587

2588

2589

2590

2591

2592

2593

2594

2595

2596

2597

2598

2599

2600

2601

2602

2603

2604

2605

2606

2607

2608

2609

2610

2611

2612

2613

2614

2615

2616

2617

2618

2619

2620

2621

2622

2623

2624

2625

2626

2627

2628

2629

2630

2631

2632

2633

2634

2635

2636

2637

2638

2639

2640

2641

2642

2643

2644

2645

2646

2647

2648

2649

2650

2651

2652

2653

2654

2655

2656

2657

2658

2659

2660

2661

2662

2663

2664

2665

2666

2667

2668

2669

2670

2671

2672

2673

2674

2675

2676

2677

2678

2679

2680

2681

2682

2683

2684

2685

2686

2687

2688

2689

2690

2691

2692

2693

""" 

Optimized Quadratic Number Field Elements 

  

This file defines a Cython class ``NumberFieldElement_quadratic`` to speed up 

computations in quadratic extensions of `\QQ`. 

  

AUTHORS: 

  

- Robert Bradshaw (2007-09): Initial version 

- David Harvey (2007-10): fix up a few bugs, polish around the edges 

- David Loeffler (2009-05): add more documentation and tests 

- Vincent Delecroix (2012-07): comparisons for quadratic number fields 

(:trac:`13213`), abs, floor and ceil functions (:trac:`13256`) 

  

.. TODO:: 

  

The ``_new()`` method should be overridden in this class to copy the ``D`` 

and ``standard_embedding`` attributes 

""" 

  

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

# Copyright (C) 2007 Robert Bradshaw <robertwb@math.washington.edu> 

# 

# This program is free software: you can redistribute it and/or modify 

# it under the terms of the GNU General Public License 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, division 

  

include "sage/libs/ntl/decl.pxi" 

from cpython.object cimport Py_EQ, Py_NE, Py_LE, Py_GE, Py_LT, Py_GT 

  

from cysignals.signals cimport sig_on, sig_off 

  

from sage.libs.gmp.mpz cimport * 

from sage.libs.gmp.mpq cimport * 

from sage.libs.flint.fmpz cimport * 

from sage.libs.arb.arb cimport * 

from sage.libs.arb.acb cimport * 

from sage.libs.ntl.ntl_ZZ cimport ntl_ZZ 

from sage.libs.ntl.ntl_ZZX cimport ntl_ZZX 

from sage.libs.mpfi cimport * 

  

from sage.structure.parent_base cimport ParentWithBase 

from sage.structure.element cimport Element 

from sage.structure.richcmp cimport rich_to_bool_sgn 

  

from sage.rings.rational cimport Rational 

from sage.rings.integer_ring import ZZ 

from sage.rings.rational_field import QQ 

from sage.rings.real_double import RDF 

from sage.rings.complex_double import CDF 

from sage.categories.morphism cimport Morphism 

from sage.rings.number_field.number_field_element import _inverse_mod_generic 

from sage.rings.real_mpfi cimport RealIntervalField_class 

from sage.rings.complex_interval cimport ComplexIntervalFieldElement 

from sage.rings.real_arb cimport RealBall 

from sage.rings.complex_arb cimport ComplexBall 

  

from sage.libs.gmp.pylong cimport mpz_pythonhash 

  

  

def __make_NumberFieldElement_quadratic0(parent, a, b, denom): 

""" 

Used in unpickling elements of number fields. 

  

TESTS:: 

  

sage: K.<a> = NumberField(x^2-x+13) 

sage: loads(dumps(a)) == a # indirect doctest 

True 

""" 

return NumberFieldElement_quadratic(parent, (a, b, denom)) 

  

def __make_NumberFieldElement_quadratic1(parent, cls, a, b, denom): 

""" 

Used in unpickling elements of number fields. 

  

TESTS:: 

  

sage: K.<a> = NumberField(x^2-x+13) 

sage: loads(dumps(a)) == a # indirect doctest 

True 

  

We test that :trac:`6462` is fixed:: 

  

sage: L = QuadraticField(-11,'a'); OL = L.maximal_order(); w = OL.0 

sage: loads(dumps(w)) == w # indirect doctest 

True 

""" 

return cls(parent, (a, b, denom)) 

  

  

cdef class NumberFieldElement_quadratic(NumberFieldElement_absolute): 

r""" 

A NumberFieldElement_quadratic object gives an efficient representation of 

an element of a quadratic extension of `\QQ`. 

  

Elements are represented internally as triples `(a, b, c)` of integers, 

where `{\rm gcd}(a, b, c) = 1` and `c > 0`, representing the element `(a + 

b \sqrt{D}) / c`. Note that if the discriminant `D` is `1 \bmod 4`, 

integral elements do not necessarily have `c = 1`. 

  

TESTS:: 

  

sage: from sage.rings.number_field.number_field_element_quadratic import NumberFieldElement_quadratic 

  

We set up some fields:: 

  

sage: K.<a> = NumberField(x^2+23) 

sage: a.parts() 

(0, 1) 

sage: F.<b> = NumberField(x^2-x+7) 

sage: b.parts() 

(1/2, 3/2) 

  

We construct elements of these fields in various ways - firstly, from 

polynomials:: 

  

sage: NumberFieldElement_quadratic(K, x-1) 

a - 1 

sage: NumberFieldElement_quadratic(F, x-1) 

b - 1 

  

From triples of Integers:: 

  

sage: NumberFieldElement_quadratic(K, (1,2,3)) 

2/3*a + 1/3 

sage: NumberFieldElement_quadratic(F, (1,2,3)) 

4/9*b + 1/9 

sage: NumberFieldElement_quadratic(F, (1,2,3)).parts() 

(1/3, 2/3) 

  

From pairs of Rationals:: 

  

sage: NumberFieldElement_quadratic(K, (1/2,1/3)) 

1/3*a + 1/2 

sage: NumberFieldElement_quadratic(F, (1/2,1/3)) 

2/9*b + 7/18 

sage: NumberFieldElement_quadratic(F, (1/2,1/3)).parts() 

(1/2, 1/3) 

  

Direct from Rationals:: 

  

sage: NumberFieldElement_quadratic(K, 2/3) 

2/3 

sage: NumberFieldElement_quadratic(F, 2/3) 

2/3 

  

This checks a bug when converting from lists:: 

  

sage: w = CyclotomicField(3)([1/2,1]) 

sage: w == w.__invert__().__invert__() 

True 

""" 

  

def __init__(self, parent, f): 

""" 

Standard initialisation function. 

  

EXAMPLES:: 

  

sage: F.<a> = QuadraticField(-7) 

sage: c = a + 7 

sage: type(c) # indirect doctest 

<type 'sage.rings.number_field.number_field_element_quadratic.NumberFieldElement_quadratic'> 

""" 

self.D = parent._D 

cdef Integer a, b, denom 

cdef Rational ad, bd 

  

cdef NumberFieldElement_quadratic gen 

  

if isinstance(f, NumberFieldElement_quadratic): 

self._parent = parent # NOTE: We do *not* call NumberFieldElement_absolute.__init__, for speed reasons. 

mpz_set(self.a, (<NumberFieldElement_quadratic>f).a) 

mpz_set(self.b, (<NumberFieldElement_quadratic>f).b) 

mpz_set(self.denom, (<NumberFieldElement_quadratic>f).denom) 

  

elif type(f) is tuple and len(f) == 2: 

NumberFieldElement_absolute.__init__(self, parent, None) 

ad, bd = f 

mpz_lcm(self.denom, mpq_denref(ad.value), mpq_denref(bd.value)) 

mpz_divexact(self.a, self.denom, mpq_denref(ad.value)) 

mpz_mul(self.a, self.a, mpq_numref(ad.value)) 

mpz_divexact(self.b, self.denom, mpq_denref(bd.value)) 

mpz_mul(self.b, self.b, mpq_numref(bd.value)) 

  

elif type(f) is tuple and len(f) == 3: 

NumberFieldElement_absolute.__init__(self, parent, None) 

a, b, denom = f 

mpz_set(self.a, a.value) 

mpz_set(self.b, b.value) 

mpz_set(self.denom, denom.value) 

self._reduce_c_() 

  

else: 

NumberFieldElement_absolute.__init__(self, parent, f) 

# poly is in gen (which may not be sqrt(d)) 

self._ntl_coeff_as_mpz(self.a, 0) 

self._ntl_coeff_as_mpz(self.b, 1) 

if mpz_cmp_ui(self.a, 0) or mpz_cmp_ui(self.b, 0): 

gen = parent.gen() # should this be cached? 

self._ntl_denom_as_mpz(self.denom) 

if mpz_cmp_ui(self.b, 0): 

mpz_mul(self.a, self.a, gen.denom) 

mpz_addmul(self.a, self.b, gen.a) 

mpz_mul(self.b, self.b, gen.b) 

mpz_mul(self.denom, self.denom, gen.denom) 

else: 

mpz_set_ui(self.denom, 1) 

self._reduce_c_() 

  

# set the attribute standard embedding which is used in the methods 

# __cmp__, sign, real, imag, floor, ceil, ... 

self.standard_embedding = parent._standard_embedding 

  

cdef _new(self): 

""" 

Quickly creates a new initialized NumberFieldElement_quadratic with the 

same parent as self. 

  

EXAMPLES:: 

  

sage: F.<b> = CyclotomicField(3) 

sage: b + b # indirect doctest 

2*b 

""" 

cdef type t = type(self) 

cdef NumberFieldElement_quadratic x = <NumberFieldElement_quadratic>t.__new__(t) 

x._parent = self._parent 

x.standard_embedding = self.standard_embedding 

x.D = self.D 

return x 

  

cdef number_field(self): 

r""" 

Return the number field to which this element belongs. Since this is a 

Cython cdef method, it is not directly accessible by the user, but the 

function "_number_field" calls this one. 

  

EXAMPLES:: 

  

sage: F.<b> = QuadraticField(-7) 

sage: b._number_field() # indirect doctest 

Number Field in b with defining polynomial x^2 + 7 

""" 

return self._parent 

  

def _maxima_init_(self, I=None): 

""" 

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(-1) 

sage: f = 1 + a 

sage: f._maxima_init_() 

'1+%i*1' 

""" 

a = self.parent().gen() 

if a**2 == -1: 

x0, x1 = self 

return str(x0) + "+" + "%i*" + str(x1) 

else: 

NumberFieldElement_absolute._maxima_init_(self, I) 

  

def _polymake_init_(self): 

""" 

EXAMPLES:: 

  

sage: K.<sqrt5> = QuadraticField(5) 

sage: polymake(3+2*sqrt5) # optional - polymake 

3+2r5 

sage: K.<i> = QuadraticField(-1) 

sage: polymake(i) # optional - polymake 

Traceback (most recent call last): 

... 

TypeError: Negative values for the root of the extension ... Bad Thing... 

""" 

x0, x1 = self 

return "new QuadraticExtension({}, {}, {})".format(x0, x1, self.D) 

  

def __copy__(self): 

r""" 

Returns a new copy of self. 

  

TESTS:: 

  

sage: K.<a> = QuadraticField(-3) 

sage: b = a + 3 

sage: c = b.__copy__() 

sage: b 

a + 3 

sage: c 

a + 3 

sage: b is c 

False 

sage: b == c 

True 

""" 

cdef NumberFieldElement_quadratic x = <NumberFieldElement_quadratic>self._new() 

mpz_set(x.a, self.a) 

mpz_set(x.b, self.b) 

mpz_set(x.denom, self.denom) 

return x 

  

  

def __cinit__(self): 

r""" 

Initialisation function. 

  

EXAMPLES:: 

  

sage: QuadraticField(-3, 'a').gen() # indirect doctest 

a 

""" 

mpz_init(self.a) 

mpz_init(self.b) 

mpz_init(self.denom) 

  

  

def __dealloc__(self): 

mpz_clear(self.a) 

mpz_clear(self.b) 

mpz_clear(self.denom) 

  

  

def __reduce__(self): 

""" 

Used for pickling. 

  

TESTS: 

  

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

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

True 

sage: loads(dumps(a/3+5)) == a/3+5 

True 

""" 

cdef Integer a = Integer.__new__(Integer) 

cdef Integer b = Integer.__new__(Integer) 

cdef Integer denom = Integer.__new__(Integer) 

mpz_set(a.value, self.a) 

mpz_set(b.value, self.b) 

mpz_set(denom.value, self.denom) 

return __make_NumberFieldElement_quadratic1, (self._parent, type(self), a, b, denom) 

  

cdef int _randomize(self, num_bound, den_bound, distribution) except -1: 

cdef Integer temp, denom1, denom2 

  

# in theory, we could just generate two random numerators and 

# a random denominator. however, this would mean that we were 

# extraordinarily unlikely to run into results of the form 

# 1/3 + 1/5*sqrt(D), which are often some of the best examples 

# for testing out code. since this is probably the primary use 

# of the random element code, it's worth doing slightly more 

# work to make this possible. 

  

# normalize denominator bound 

if den_bound is None or den_bound < 1: 

den_bound = 1 

  

# generate denominators 

denom1 = <Integer>(ZZ.random_element(x=1, 

y=den_bound+1, 

distribution=distribution)) 

denom2 = <Integer>(ZZ.random_element(x=1, 

y=den_bound+1, 

distribution=distribution)) 

  

# set a, b 

temp = <Integer>(ZZ.random_element(x=num_bound, distribution=distribution)) 

mpz_mul(self.a, temp.value, denom2.value) 

temp = <Integer>(ZZ.random_element(x=num_bound, distribution=distribution)) 

mpz_mul(self.b, temp.value, denom1.value) 

# set denom 

mpz_mul(self.denom, denom1.value, denom2.value) 

  

self._reduce_c_() 

return 0 # No error 

  

  

def _lift_cyclotomic_element(self, new_parent, bint check=True, int rel=0): 

""" 

Creates an element of the passed field from this field. This 

is specific to creating elements in a cyclotomic field from 

elements in another cyclotomic field, in the case that 

self.number_field()._n() divides new_parent()._n(). This 

function aims to make this common coercion extremely fast! 

  

More general coercion (i.e. of zeta6 into CyclotomicField(3)) 

is implemented in the _coerce_from_other_cyclotomic_field 

method of a CyclotomicField. 

  

EXAMPLES:: 

  

sage: C.<zeta4>=CyclotomicField(4) 

sage: CyclotomicField(20)(zeta4+1) # The function _lift_cyclotomic_element does the heavy lifting in the background 

zeta20^5 + 1 

sage: (zeta4+1)._lift_cyclotomic_element(CyclotomicField(40)) # There is rarely a purpose to call this function directly 

zeta40^10 + 1 

  

sage: cf3 = CyclotomicField(3) ; z3 = cf3.0 

sage: cf6 = CyclotomicField(6) ; z6 = cf6.0 

sage: z6._lift_cyclotomic_element(cf3) 

Traceback (most recent call last): 

... 

TypeError: The zeta_order of the new field must be a multiple of the zeta_order of the original. 

sage: cf3(z6) 

zeta3 + 1 

sage: z3._lift_cyclotomic_element(cf6) 

zeta6 - 1 

  

Verify embeddings are respected:: 

  

sage: cf6c = CyclotomicField(6, embedding=CDF(exp(-pi*I/3))) ; z6c = cf6c.0 

sage: cf3(z6c) 

-zeta3 

sage: cf6c(z3) 

-zeta6 

  

AUTHOR: 

  

- Joel B. Mohler (original version) 

  

- Craig Citro (reworked for quadratic elements) 

""" 

if check: 

from .number_field import NumberField_cyclotomic 

if not isinstance(self.number_field(), NumberField_cyclotomic) \ 

or not isinstance(new_parent, NumberField_cyclotomic): 

raise TypeError("The field and the new parent field must both be cyclotomic fields.") 

  

if rel == 0: 

small_order = self.number_field()._n() 

large_order = new_parent._n() 

  

try: 

rel = ZZ(large_order / small_order) 

except TypeError: 

raise TypeError("The zeta_order of the new field must be a multiple of the zeta_order of the original.") 

  

cdef NumberFieldElement_quadratic x2 

cdef int n = self._parent._n() 

  

if new_parent.degree() == 2: 

## since self is a *quadratic* element, we can only get 

## here if self.parent() and new_parent are: 

## - CyclotomicField(3) and CyclotomicField(6) 

## - CyclotomicField(3) and CyclotomicField(3) 

## - CyclotomicField(6) and CyclotomicField(6) 

## - CyclotomicField(4) and CyclotomicField(4) 

## In all cases, conversion of elements is trivial! 

if n == <int>new_parent._n(): 

conjugate = rel != 1 

else: 

# n = 3, new_n = 6 

conjugate = rel == 4 

x2 = <NumberFieldElement_quadratic>(self._new()) 

x2._parent = new_parent 

mpz_set(x2.a, self.a) 

if conjugate: 

mpz_neg(x2.b, self.b) 

else: 

mpz_set(x2.b, self.b) 

mpz_set(x2.denom, self.denom) 

x2.D = self.D 

return x2 

  

cdef NumberFieldElement x 

cdef ZZX_c elt_num 

cdef ZZ_c elt_den, tmp_coeff 

cdef mpz_t tmp_mpz 

cdef long tmp_const 

  

x = <NumberFieldElement_absolute>NumberFieldElement_absolute.__new__(NumberFieldElement_absolute) 

  

mpz_to_ZZ(&elt_den, self.denom) 

  

mpz_init(tmp_mpz) 

  

## set the two terms in the polynomial 

if n == 4: 

mpz_to_ZZ(&tmp_coeff, self.a) 

ZZX_SetCoeff(elt_num, 0, tmp_coeff) 

mpz_to_ZZ(&tmp_coeff, self.b) 

ZZX_SetCoeff(elt_num, 1, tmp_coeff) 

  

elif n == 3: 

## num[0] = a + b 

mpz_add(tmp_mpz, tmp_mpz, self.a) 

mpz_add(tmp_mpz, tmp_mpz, self.b) 

mpz_to_ZZ(&tmp_coeff, tmp_mpz) 

ZZX_SetCoeff(elt_num, 0, tmp_coeff) 

  

## num[1] = 2*b 

mpz_sub(tmp_mpz, tmp_mpz, self.a) 

tmp_const = 2 

mpz_mul_si(tmp_mpz, tmp_mpz, tmp_const) 

mpz_to_ZZ(&tmp_coeff, tmp_mpz) 

ZZX_SetCoeff(elt_num, 1, tmp_coeff) 

  

elif n == 6: 

## num[0] = a - b 

mpz_add(tmp_mpz, tmp_mpz, self.a) 

mpz_sub(tmp_mpz, tmp_mpz, self.b) 

mpz_to_ZZ(&tmp_coeff, tmp_mpz) 

ZZX_SetCoeff(elt_num, 0, tmp_coeff) 

  

## num[1] = 2*b 

mpz_sub(tmp_mpz, tmp_mpz, self.a) 

tmp_const = -2 

mpz_mul_si(tmp_mpz, tmp_mpz, tmp_const) 

mpz_to_ZZ(&tmp_coeff, tmp_mpz) 

ZZX_SetCoeff(elt_num, 1, tmp_coeff) 

  

mpz_clear(tmp_mpz) 

  

x._parent = <ParentWithBase>new_parent 

x.__fld_numerator, x.__fld_denominator = new_parent.polynomial_ntl() 

x.__denominator = elt_den 

cdef ZZX_c result 

cdef ZZ_c tmp 

cdef int i 

cdef ntl_ZZX _num 

cdef ntl_ZZ _den 

_num, _den = new_parent.polynomial_ntl() 

for i from 0 <= i <= ZZX_deg(elt_num): 

tmp = ZZX_coeff(elt_num, i) 

ZZX_SetCoeff(result, i*rel, tmp) 

ZZX_rem(x.__numerator, result, _num.x) 

(<NumberFieldElement_absolute>x)._reduce_c_() 

return x 

  

def _real_mpfi_(self, R): 

r""" 

Conversion to a real interval field 

  

TESTS:: 

  

sage: K1 = QuadraticField(2) 

sage: RIF(K1.gen()) 

1.414213562373095? 

sage: RIF(K1(1/5)) 

0.2000000000000000? 

sage: RIF(3/5*K1.gen() + 1/5) 

1.048528137423857? 

  

sage: K2 = QuadraticField(2, embedding=-RLF(2).sqrt()) 

sage: RIF(K2.gen()) 

-1.414213562373095? 

  

sage: K3 = QuadraticField(-1) 

sage: RIF(K3.gen()) 

Traceback (most recent call last): 

... 

ValueError: unable to convert complex algebraic number a to real interval 

sage: RIF(K3(2)) 

2 

  

Check that :trac:`21979` is fixed:: 

  

sage: a = QuadraticField(5).gen() 

sage: u = -573147844013817084101/2*a + 1281597540372340914251/2 

sage: RealIntervalField(128)(u) 

0.?e-17 

sage: RealIntervalField(128)(u).is_zero() 

False 

  

This was fixed in :trac:`24371`:: 

  

sage: RIF.convert_map_from(QuadraticField(5)) 

Conversion via _real_mpfi_ method map: 

From: Number Field in a with defining polynomial x^2 - 5 

To: Real Interval Field with 53 bits of precision 

""" 

ans = (<RealIntervalField_class?>R)._new() 

  

if mpz_cmp_ui(self.b, 0): 

if mpz_cmp_ui(self.D.value, 0) < 0: 

raise ValueError(f"unable to convert complex algebraic number {self!r} to real interval") 

mpfi_set_z(ans.value, self.D.value) 

mpfi_sqrt(ans.value, ans.value) 

if not self.standard_embedding: 

mpfi_neg(ans.value, ans.value) 

mpfi_mul_z(ans.value, ans.value, self.b) 

mpfi_add_z(ans.value, ans.value, self.a) 

else: 

mpfi_set_z(ans.value, self.a) 

  

mpfi_div_z(ans.value, ans.value, self.denom) 

return ans 

  

def _complex_mpfi_(self, R): 

r""" 

Conversion to a complex interval field 

  

TESTS:: 

  

sage: K.<a> = QuadraticField(2) 

sage: CIF(a) 

1.414213562373095? 

sage: CIF(K(1/5)) 

0.2000000000000000? 

sage: CIF(1/5 + 3/5*a) 

1.048528137423857? 

  

sage: K.<a> = QuadraticField(-2) 

sage: CIF(a) 

1.414213562373095?*I 

sage: CIF(K(1/5)) 

0.2000000000000000? 

sage: CIF(1/5 + 3/5*a) 

0.2000000000000000? + 0.848528137423857?*I 

  

sage: K.<a> = QuadraticField(-2, embedding=-CLF(-2).sqrt()) 

sage: CIF(a) 

-1.414213562373095?*I 

  

This was fixed in :trac:`24371`:: 

  

sage: CIF.convert_map_from(QuadraticField(-5)) 

Conversion via _complex_mpfi_ method map: 

From: Number Field in a with defining polynomial x^2 + 5 

To: Complex Interval Field with 53 bits of precision 

""" 

ans = <ComplexIntervalFieldElement>ComplexIntervalFieldElement.__new__(ComplexIntervalFieldElement, R) 

  

if mpz_cmp_ui(self.b, 0): 

mpfi_set_z(ans.__re, self.D.value) 

if mpfi_is_neg(ans.__re): 

# Imaginary quadratic 

mpfi_neg(ans.__re, ans.__re) 

mpfi_sqrt(ans.__im, ans.__re) 

if not self.standard_embedding: 

mpfi_neg(ans.__im, ans.__im) 

mpfi_set_z(ans.__re, self.a) 

mpfi_mul_z(ans.__im, ans.__im, self.b) 

mpfi_div_z(ans.__im, ans.__im, self.denom) 

else: 

# Real quadratic 

mpfi_sqrt(ans.__re, ans.__re) 

if not self.standard_embedding: 

mpfi_neg(ans.__re, ans.__re) 

mpfi_mul_z(ans.__re, ans.__re, self.b) 

mpfi_add_z(ans.__re, ans.__re, self.a) 

mpfi_set_ui(ans.__im, 0) 

else: 

mpfi_set_z(ans.__re, self.a) 

mpfi_set_ui(ans.__im, 0) 

  

mpfi_div_z(ans.__re, ans.__re, self.denom) 

return ans 

  

cdef int arb_set_real(self, arb_t x, long prec) except -1: 

"Set x to the real part of this element" 

cdef fmpz_t tmpz 

cdef arb_t rootD 

cdef long prec2 

  

fmpz_init(tmpz) 

  

if mpz_sgn(self.D.value) > 0: 

# To mitigate the effect of cancellations between 

# a and b*sqrt(D) we perform a loop with increasing 

# working precision 

arb_init(rootD) 

prec2 = prec 

sig_on() 

while True: 

fmpz_set_mpz(tmpz, self.a) 

arb_set_fmpz(x, tmpz) 

fmpz_set_mpz(tmpz, self.D.value) 

arb_sqrt_fmpz(rootD, tmpz, prec2) 

fmpz_set_mpz(tmpz, self.b) 

if self.standard_embedding: 

arb_addmul_fmpz(x, rootD, tmpz, prec2) 

else: 

arb_submul_fmpz(x, rootD, tmpz, prec2) 

if arb_rel_accuracy_bits(x) < prec - 4: 

prec2 *= 2 

continue 

else: 

break 

sig_off() 

arb_clear(rootD) 

else: 

fmpz_set_mpz(tmpz, self.a) 

arb_set_fmpz(x, tmpz) 

  

fmpz_set_mpz(tmpz, self.denom) 

arb_div_fmpz(x, x, tmpz, prec) 

fmpz_clear(tmpz) 

return 0 

  

cdef void arb_set_imag(self, arb_t x, long prec): 

"Set x to the imaginary part of this element" 

cdef fmpz_t tmpz 

cdef arb_t rootD 

  

if mpz_sgn(self.D.value) < 0 and mpz_sgn(self.b): 

arb_init(rootD) 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, self.D.value) 

fmpz_neg(tmpz, tmpz) 

arb_sqrt_fmpz(rootD, tmpz, prec) 

fmpz_set_mpz(tmpz, self.b) 

if self.standard_embedding: 

arb_addmul_fmpz(x, rootD, tmpz, prec) 

else: 

arb_submul_fmpz(x, rootD, tmpz, prec) 

fmpz_set_mpz(tmpz, self.denom) 

arb_div_fmpz(x, x, tmpz, prec) 

  

fmpz_clear(tmpz) 

arb_clear(rootD) 

  

def _arb_(self, R): 

r""" 

Conversion to a real ball with parent ``R``. 

  

EXAMPLES:: 

  

sage: a = QuadraticField(7).gen() 

sage: a._arb_(RBF) 

[2.645751311064590 +/- 7.17e-16] 

sage: RBF(a) 

[2.645751311064590 +/- 7.17e-16] 

  

sage: a = QuadraticField(7, embedding=-AA(7).sqrt()).gen() 

sage: a._arb_(RBF) 

[-2.645751311064590 +/- 7.17e-16] 

  

sage: a = QuadraticField(-1).gen() 

sage: RBF(a) 

Traceback (most recent call last): 

... 

ValueError: nonzero imaginary part 

  

This conversion takes care of cancellations:: 

  

sage: K.<a> = QuadraticField(3) 

sage: b = (a - 1)**1000 

sage: RBF(b) 

[3.477155118441024e-136 +/- 8.45e-152] 

  

TESTS: 

  

Check that coercions and conversions go throuh this method:: 

  

sage: RBF.convert_map_from(QuadraticField(5)) 

Conversion via _arb_ method map: 

From: Number Field in a with defining polynomial x^2 - 5 

To: Real ball field with 53 bits of precision 

sage: RBF.coerce_map_from(QuadraticField(5, embedding=-AA(5).sqrt())) 

Conversion via _arb_ method map: 

From: Number Field in a with defining polynomial x^2 - 5 

To: Real ball field with 53 bits of precision 

""" 

cdef RealBall res = RealBall.__new__(RealBall) 

res._parent = R 

if mpz_sgn(self.D.value) < 0 and mpz_sgn(self.b): 

raise ValueError('nonzero imaginary part') 

self.arb_set_real(res.value, R._prec) 

return res 

  

def _acb_(self, R): 

r""" 

Conversion to complex ball with parent ``R``. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(-3) 

sage: CBF(1/7 + 3/2 * a) 

[0.1428571428571428 +/- 7.70e-17] + [2.59807621135332 +/- 6.17e-15]*I 

sage: CBF(1/7) + CBF(3/2) * CBF(a) 

[0.1428571428571428 +/- 7.70e-17] + [2.59807621135332 +/- 5.21e-15]*I 

  

sage: a._acb_(CBF) 

[1.732050807568877 +/- 4.16e-16]*I 

  

TESTS: 

  

Check that coercions and conversions go through this method:: 

  

sage: CBF.convert_map_from(QuadraticField(-3)) 

Conversion via _acb_ method map: 

From: Number Field in a with defining polynomial x^2 + 3 

To: Complex ball field with 53 bits of precision 

  

sage: CBF.coerce_map_from(QuadraticField(-3, embedding=-QQbar(-3).sqrt())) 

Conversion via _acb_ method map: 

From: Number Field in a with defining polynomial x^2 + 3 

To: Complex ball field with 53 bits of precision 

""" 

cdef ComplexBall res = ComplexBall.__new__(ComplexBall) 

res._parent = R 

self.arb_set_real(acb_realref(res.value), R._prec) 

self.arb_set_imag(acb_imagref(res.value), R._prec) 

return res 

  

def parts(self): 

""" 

This function returns a pair of rationals `a` and `b` such that self `= 

a+b\sqrt{D}`. 

  

This is much closer to the internal storage format of the 

elements than the polynomial representation coefficients (the output of 

``self.list()``), unless the generator with which this number field was 

constructed was equal to `\sqrt{D}`. See the last example below. 

  

EXAMPLES:: 

  

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

sage: K.discriminant() 

13 

sage: a.parts() 

(0, 1) 

sage: (a/2-4).parts() 

(-4, 1/2) 

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

sage: K.discriminant() 

28 

sage: a.parts() 

(0, 1) 

sage: K.<a> = NumberField(x^2-x+7) 

sage: a.parts() 

(1/2, 3/2) 

sage: a._coefficients() 

[0, 1] 

""" 

cdef Rational ad = <Rational>Rational.__new__(Rational) 

cdef Rational bd = <Rational>Rational.__new__(Rational) 

if mpz_cmp_ui(self.a, 0) == 0: 

mpq_set_ui(ad.value, 0, 1) 

else: 

mpz_set(mpq_numref(ad.value), self.a) 

mpz_set(mpq_denref(ad.value), self.denom) 

mpq_canonicalize(ad.value) 

if mpz_cmp_ui(self.b, 0) == 0: 

mpq_set_ui(bd.value, 0, 1) 

else: 

mpz_set(mpq_numref(bd.value), self.b) 

mpz_set(mpq_denref(bd.value), self.denom) 

mpq_canonicalize(bd.value) 

  

return ad, bd 

  

cdef bint is_sqrt_disc(self): 

r""" 

Returns true if self is `\sqrt{D}`. 

  

EXAMPLES:: 

  

sage: F.<b> = NumberField(x^2 - x + 7) 

sage: b.denominator() # indirect doctest 

1 

""" 

return mpz_cmp_ui(self.denom, 1)==0 and mpz_cmp_ui(self.a, 0)==0 and mpz_cmp_ui(self.b, 1)==0 

  

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

# Comparisons 

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

  

def sign(self): 

r""" 

Returns the sign of self (0 if zero, +1 if positive and -1 if negative). 

  

EXAMPLES:: 

  

sage: K.<sqrt2> = QuadraticField(2, name='sqrt2') 

sage: K(0).sign() 

0 

sage: sqrt2.sign() 

1 

sage: (sqrt2+1).sign() 

1 

sage: (sqrt2-1).sign() 

1 

sage: (sqrt2-2).sign() 

-1 

sage: (-sqrt2).sign() 

-1 

sage: (-sqrt2+1).sign() 

-1 

sage: (-sqrt2+2).sign() 

1 

  

sage: K.<a> = QuadraticField(2, embedding=-1.4142) 

sage: K(0).sign() 

0 

sage: a.sign() 

-1 

sage: (a+1).sign() 

-1 

sage: (a+2).sign() 

1 

sage: (a-1).sign() 

-1 

sage: (-a).sign() 

1 

sage: (-a-1).sign() 

1 

sage: (-a-2).sign() 

-1 

  

sage: K.<b> = NumberField(x^2 + 2*x + 7, 'b', embedding=CC(-1,-sqrt(6))) 

sage: b.sign() 

Traceback (most recent call last): 

... 

ValueError: a complex number has no sign! 

sage: K(1).sign() 

1 

sage: K(0).sign() 

0 

sage: K(-2/3).sign() 

-1 

""" 

cdef mpz_t i, j 

cdef int s = 1, test 

  

if mpz_sgn(self.b) == 0: 

return mpz_sgn(self.a) 

  

if mpz_sgn(self.D.value) == -1: 

raise ValueError("a complex number has no sign!") 

  

if not self.standard_embedding: 

s = -1 

  

if mpz_sgn(self.a) == 0: 

return s*mpz_sgn(self.b) 

  

if mpz_sgn(self.a) == 1: 

if mpz_sgn(self.b) == s: 

return 1 

  

elif mpz_sgn(self.b) == -s: 

return -1 

  

mpz_init_set(i,self.a) 

mpz_mul(i,i,i) 

mpz_init_set(j,self.b) 

mpz_mul(j,j,j) 

mpz_mul(j,j,self.D.value) 

test = mpz_cmp(i,j) 

mpz_clear(i) 

mpz_clear(j) 

if mpz_sgn(self.a) == 1 and mpz_sgn(self.b) == -s: 

return test 

return -test 

  

cpdef _richcmp_(left, _right, int op): 

r""" 

Rich comparison of elements. 

  

TESTS:: 

  

sage: K.<i> = QuadraticField(-1) 

sage: sorted([5*i+1, 2, 3*i+1, 2-i]) 

[3*i + 1, 5*i + 1, -i + 2, 2] 

  

Make some random tests to check that the order is compatible with the 

ones of the real field (RR) and complex field (CC):: 

  

sage: K1 = NumberField(x^2 - 2, 'a', embedding=RR(1.4)) 

sage: K2 = NumberField(x^2 - 2, 'a', embedding=RR(-1.4)) 

sage: for _ in range(500): # long time 

....: for K in K1, K2: 

....: a = K.random_element() 

....: b = K.random_element() 

....: assert (a < b) == (RR(a) < RR(b)) 

....: assert (a > b) == (RR(a) > RR(b)) 

....: assert (a == b) == (RR(a) == RR(b)) 

....: assert (a != b) == (RR(a) != RR(b)) 

....: assert (a >= b) == (RR(a) >= RR(b)) 

....: assert (a <= b) == (RR(a) <= RR(b)) 

  

:: 

  

sage: K1 = NumberField(x^2 + 2, 'a', embedding=CC(0,1)) 

sage: K2 = NumberField(x^2 + 2, 'a', embedding=CC(0,-1)) 

sage: for _ in range(500): # long time 

....: for K in K1, K2: 

....: a = K.random_element() 

....: b = K.random_element() 

....: assert (a < b) == (CC(a) < CC(b)) 

....: assert (a > b) == (CC(a) > CC(b)) 

....: assert (a == b) == (CC(a) == CC(b)) 

....: assert (a != b) == (CC(a) != CC(b)) 

....: assert (a >= b) == (CC(a) >= CC(b)) 

....: assert (a <= b) == (CC(a) <= CC(b)) 

  

The following is tested because of the implementation of 

:func:`Q_to_quadratic_field_element` which was the cause of 

some problems with :trac:`13213`:: 

  

sage: K.<sqrt2> = QuadraticField(2) 

sage: 1/2 + sqrt2 > 0 

True 

  

Two examples from the same number field with its two possible real 

embeddings:: 

  

sage: K.<phi> = NumberField(x^2-x-1, 'phi', embedding=1.618) 

sage: phi > 0 

True 

sage: -phi > 0 

False 

sage: phi - 3 == 2*phi + 1 

False 

sage: fibonacci(10)*phi < fibonacci(11) 

True 

sage: RDF(fibonacci(10)*phi) 

88.99186938124421 

sage: fibonacci(11) 

89 

sage: l = [-2, phi+3, 2*phi-1, 2*phi-5, 0, -phi+2, fibonacci(20)*phi - fibonacci(21)] 

sage: l.sort() 

sage: l 

[-2, 2*phi - 5, 6765*phi - 10946, 0, -phi + 2, 2*phi - 1, phi + 3] 

sage: list(map(RDF, l)) 

[-2.0, -1.7639320225002102, -6.610696073039435e-05, 0.0, 0.3819660112501051, 2.23606797749979, 4.618033988749895] 

  

:: 

  

sage: L.<psi> = NumberField(x^2-x-1, 'psi', embedding=-0.618) 

sage: psi < 0 

True 

sage: 2*psi + 3 == 2*psi + 3 

True 

sage: fibonacci(10)*psi < -fibonacci(9) 

False 

sage: RDF(fibonacci(10)*psi) 

-33.99186938124422 

sage: fibonacci(9) 

34 

sage: l = [-1, psi, 0, fibonacci(20)*psi + fibonacci(19), 3*psi+2] 

sage: l.sort() 

sage: l 

[-1, psi, 0, 6765*psi + 4181, 3*psi + 2] 

sage: list(map(RDF, l)) 

[-1.0, -0.6180339887498949, 0.0, 6.610696073039435e-05, 0.1458980337503153] 

  

For a field with no specified embedding the comparison uses the standard 

embedding:: 

  

sage: K.<sqrt2> = NumberField(x^2-2, 'sqrt2') 

sage: sqrt2 > 1 and sqrt2 < 2 

True 

  

The following examples illustrate the same behavior for a complex 

quadratic field:: 

  

sage: K.<i> = QuadraticField(-1) 

sage: l = [-2, i-3, 2*i-2, 2*i+2, 5*i, 1-3*i, -1+i, 1] 

sage: l.sort() 

sage: l 

[i - 3, -2, 2*i - 2, i - 1, 5*i, -3*i + 1, 1, 2*i + 2] 

sage: list(map(CDF, l)) 

[-3.0 + 1.0*I, -2.0, -2.0 + 2.0*I, -1.0 + 1.0*I, 5.0*I, 1.0 - 3.0*I, 1.0, 2.0 + 2.0*I] 

sage: list(map(CDF, l)) == sorted(map(CDF, l)) 

True 

""" 

# When D > 0 and standard embedding, we compare (a + b * sqrt(D)) / d and (aa + 

# bb * sqrt(D)) / dd using the comparison of (dd*a - d * aa)^2 and (d*bb - dd*b)^2 * D 

# mpz_sgn: returns 1 if > 0, 0 if 0 and -1 if < 0 

cdef mpz_t i, j 

cdef NumberFieldElement_quadratic right = <NumberFieldElement_quadratic> _right 

cdef int test 

  

# inequality and equality 

if mpz_cmp(left.a, right.a) or mpz_cmp(left.b, right.b) or mpz_cmp(left.denom, right.denom): 

if op == Py_EQ: 

return False 

elif op == Py_NE: 

return True 

else: # equality 

if op == Py_EQ or op == Py_LE or op == Py_GE: 

return True 

if op == Py_NE or op == Py_LT or op == Py_GT: 

return False 

  

# comparisons are valid only in *real* quadratic number field 

# when no embedding is specified or in the case of complex embeddings we 

# use a lexicographic order. 

if mpz_sgn(left.D.value) == -1: 

mpz_init(i) 

mpz_init(j) 

mpz_mul(i, left.a, right.denom) 

mpz_mul(j, right.a, left.denom) 

test = mpz_cmp(i,j) 

if test: 

mpz_clear(i) 

mpz_clear(j) 

return rich_to_bool_sgn(op, test) 

mpz_mul(i, left.b, right.denom) 

mpz_mul(j, right.b, left.denom) 

test = mpz_cmp(i,j) 

if test: 

if not left.standard_embedding: 

test = -test 

mpz_clear(i) 

mpz_clear(j) 

return rich_to_bool_sgn(op, test) 

test = mpz_cmp(left.denom, right.denom) 

mpz_clear(i) 

mpz_clear(j) 

return rich_to_bool_sgn(op, test) 

  

# comparison in the real case 

mpz_init(i) 

mpz_mul(i, right.denom, left.a) 

mpz_submul(i, left.denom, right.a) 

  

mpz_init(j) 

mpz_mul(j, left.denom, right.b) 

mpz_submul(j, right.denom, left.b) 

  

if not left.standard_embedding: 

mpz_neg(j, j) 

  

if mpz_sgn(i) == 1: 

if mpz_sgn(j) == 1: 

mpz_mul(i, i, i) 

mpz_mul(j, j, j) 

mpz_mul(j, j, left.D.value) 

test = mpz_cmp(i, j) 

else: 

test = 1 

  

else: 

if mpz_sgn(j) == -1: 

mpz_mul(i, i, i) 

mpz_mul(j, j, j) 

mpz_mul(j, j, left.D.value) 

test = mpz_cmp(j, i) 

else: 

test = -1 

  

mpz_clear(i) 

mpz_clear(j) 

return rich_to_bool_sgn(op, test) 

  

def continued_fraction_list(self): 

r""" 

Return the preperiod and the period of the continued fraction expansion 

of ``self``. 

  

EXAMPLES:: 

  

sage: K.<sqrt2> = QuadraticField(2) 

sage: sqrt2.continued_fraction_list() 

((1,), (2,)) 

sage: (1/2+sqrt2/3).continued_fraction_list() 

((0, 1, 33), (1, 32)) 

  

For rational entries a pair of tuples is also returned but the second 

one is empty:: 

  

sage: K(123/567).continued_fraction_list() 

((0, 4, 1, 1, 1, 1, 3, 2), ()) 

""" 

if mpz_sgn(self.b) == 0: 

return tuple(Rational(self).continued_fraction_list()),() 

  

if mpz_sgn(self.D.value) < 0: 

raise ValueError("the method is only available for positive discriminant") 

  

x = self 

orbit = [] 

quots = [] 

while x not in orbit: 

quots.append(x.floor()) 

orbit.append(x) 

x = ~(x - quots[-1]) 

  

i = orbit.index(x) 

  

return tuple(quots[:i]), tuple(quots[i:]) 

  

def continued_fraction(self): 

r""" 

Return the (finite or ultimately periodic) continued fraction of ``self``. 

  

EXAMPLES:: 

  

sage: K.<sqrt2> = QuadraticField(2) 

sage: cf = sqrt2.continued_fraction(); cf 

[1; (2)*] 

sage: cf.n() 

1.41421356237310 

sage: sqrt2.n() 

1.41421356237309 

sage: cf.value() 

sqrt2 

  

sage: (sqrt2/3 + 1/4).continued_fraction() 

[0; 1, (2, 1, 1, 2, 3, 2, 1, 1, 2, 5, 1, 1, 14, 1, 1, 5)*] 

""" 

t1,t2 = self.continued_fraction_list() 

from sage.rings.continued_fraction import ContinuedFraction_periodic 

return ContinuedFraction_periodic(t1,t2) 

  

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

# Arithmetic 

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

  

cdef void _reduce_c_(self): 

r""" 

Reduces into canonical form. 

  

WARNING: this mutates self. 

""" 

cdef mpz_t gcd 

# cancel out common factors 

mpz_init(gcd) 

mpz_gcd(gcd, self.a, self.denom) 

mpz_gcd(gcd, gcd, self.b) 

if mpz_cmp_si(gcd, 1): # != 0 (i.e. it is not 1) 

mpz_divexact(self.a, self.a, gcd) 

mpz_divexact(self.b, self.b, gcd) 

mpz_divexact(self.denom, self.denom, gcd) 

# make denominator positive 

if mpz_sgn(self.denom) < 0: 

mpz_neg(self.denom, self.denom) 

mpz_neg(self.a, self.a) 

mpz_neg(self.b, self.b) 

mpz_clear(gcd) 

  

  

cpdef _add_(self, other_m): 

""" 

EXAMPLES:: 

  

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

sage: K.discriminant() 

5 

sage: a+a # indirect doctest 

2*a 

sage: s = (a+2)/6; s 

1/6*a + 1/3 

sage: s+a 

7/6*a + 1/3 

sage: s+10 

1/6*a + 31/3 

sage: s+(2*a+5)/7 

19/42*a + 22/21 

sage: s+(1+a)/2 

2/3*a + 5/6 

sage: s+(1+a)/8 

7/24*a + 11/24 

sage: s+(a+5)/6 

1/3*a + 7/6 

sage: (a/3+2/3) + (2*a/3+1/3) 

a + 1 

""" 

cdef NumberFieldElement_quadratic other = <NumberFieldElement_quadratic>other_m 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

cdef mpz_t gcd, tmp 

if mpz_cmp(self.denom, other.denom) == 0: 

mpz_add(res.a, self.a, other.a) 

mpz_add(res.b, self.b, other.b) 

mpz_set(res.denom, self.denom) 

else: 

mpz_init(gcd) 

mpz_gcd(gcd, self.denom, other.denom) 

if mpz_cmp_ui(gcd, 1) == 0: 

mpz_mul(res.a, self.a, other.denom) 

mpz_addmul(res.a, self.denom, other.a) 

mpz_mul(res.b, self.b, other.denom) 

mpz_addmul(res.b, self.denom, other.b) 

mpz_mul(res.denom, self.denom, other.denom) 

else: 

mpz_init(tmp) 

mpz_divexact(tmp, other.denom, gcd) 

mpz_mul(res.a, self.a, tmp) 

mpz_mul(res.b, self.b, tmp) 

mpz_divexact(tmp, self.denom, gcd) 

mpz_addmul(res.a, other.a, tmp) 

mpz_addmul(res.b, other.b, tmp) 

mpz_mul(res.denom, other.denom, tmp) 

mpz_clear(tmp) 

mpz_clear(gcd) 

res._reduce_c_() 

return res 

  

  

cpdef _sub_(self, other_m): 

""" 

EXAMPLES:: 

  

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

sage: b = (a-3)/10; b # indirect doctest 

1/10*a - 3/10 

sage: b-1 

1/10*a - 13/10 

sage: b-a 

-9/10*a - 3/10 

sage: b-1/2 

1/10*a - 4/5 

sage: b-a/15 

1/30*a - 3/10 

""" 

cdef NumberFieldElement_quadratic other = <NumberFieldElement_quadratic>other_m 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

cdef mpz_t gcd, tmp 

if mpz_cmp(self.denom, other.denom) == 0: 

mpz_sub(res.a, self.a, other.a) 

mpz_sub(res.b, self.b, other.b) 

mpz_set(res.denom, self.denom) 

else: 

mpz_init(gcd) 

mpz_gcd(gcd, self.denom, other.denom) 

if mpz_cmp_ui(gcd, 1) == 0: 

mpz_mul(res.a, self.a, other.denom) 

mpz_submul(res.a, self.denom, other.a) 

mpz_mul(res.b, self.b, other.denom) 

mpz_submul(res.b, self.denom, other.b) 

mpz_mul(res.denom, self.denom, other.denom) 

else: 

mpz_init(tmp) 

mpz_divexact(tmp, other.denom, gcd) 

mpz_mul(res.a, self.a, tmp) 

mpz_mul(res.b, self.b, tmp) 

mpz_divexact(tmp, self.denom, gcd) 

mpz_submul(res.a, other.a, tmp) 

mpz_submul(res.b, other.b, tmp) 

mpz_mul(res.denom, other.denom, tmp) 

mpz_clear(tmp) 

mpz_clear(gcd) 

res._reduce_c_() 

return res 

  

  

def __neg__(self): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+163) 

sage: -a 

-a 

sage: -(a+4) 

-a - 4 

sage: b = (a-3)/2 

sage: -b 

-1/2*a + 3/2 

""" 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

mpz_neg(res.a, self.a) 

mpz_neg(res.b, self.b) 

mpz_set(res.denom, self.denom) 

return res 

  

cpdef _mul_(self, other_m): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+23) 

sage: a*a # indirect doctest 

-23 

sage: (a+1)*(a-1) 

-24 

sage: (a+1)*(a+2) 

3*a - 21 

sage: (a+1)/2 * (a+2) 

3/2*a - 21/2 

sage: (a+1)/2 * (a+2)/3 

1/2*a - 7/2 

sage: (2*a+4) * (3*a)/2 

6*a - 69 

  

Verify Karatsuba:: 

  

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

sage: (10^1000 * (a+1)) * K(2+3*a) == 10^1000 * ((a+1) * K(2+3*a)) 

True 

""" 

cdef NumberFieldElement_quadratic other = <NumberFieldElement_quadratic>other_m 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

cdef mpz_t tmp 

  

if mpz_size(self.a) + mpz_size(self.b) < 8: # could I use a macro instead? 

# Do it the traditional way 

mpz_mul(res.a, self.b, other.b) 

mpz_mul(res.a, res.a, self.D.value) 

mpz_addmul(res.a, self.a, other.a) 

  

mpz_mul(res.b, self.a, other.b) 

mpz_addmul(res.b, self.b, other.a) 

  

else: 

# Karatsuba 

sig_on() 

mpz_init(tmp) 

mpz_add(res.a, self.a, self.b) # using res.a as tmp 

mpz_add(tmp, other.a, other.b) 

mpz_mul(res.b, res.a, tmp) # res.b = (self.a + self.b)(other.a + other.b) 

  

mpz_mul(res.a, self.a, other.a) 

mpz_sub(res.b, res.b, res.a) 

mpz_mul(tmp, self.b, other.b) 

mpz_sub(res.b, res.b, tmp) 

mpz_mul(tmp, tmp, self.D.value) 

mpz_add(res.a, res.a, tmp) 

mpz_clear(tmp) 

sig_off() 

  

mpz_mul(res.denom, self.denom, other.denom) 

res._reduce_c_() 

return res 

  

cpdef _rmul_(self, Element _c): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+43) 

sage: (1+a)*3 # indirect doctest 

3*a + 3 

""" 

cdef Rational c = <Rational>_c 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

mpz_mul(res.a, self.a, mpq_numref(c.value)) 

mpz_mul(res.b, self.b, mpq_numref(c.value)) 

mpz_mul(res.denom, self.denom, mpq_denref(c.value)) 

res._reduce_c_() 

return res 

  

cpdef _lmul_(self, Element _c): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+43) 

sage: 5*(a-1/5) # indirect doctest 

5*a - 1 

""" 

cdef Rational c = <Rational>_c 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

mpz_mul(res.a, self.a, mpq_numref(c.value)) 

mpz_mul(res.b, self.b, mpq_numref(c.value)) 

mpz_mul(res.denom, self.denom, mpq_denref(c.value)) 

res._reduce_c_() 

return res 

  

def __invert__(self): 

""" 

EXAMPLES:: 

  

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

sage: ~a 

1/5*a 

sage: ~(a+1) 

1/4*a - 1/4 

sage: (a-1)*(a+1) 

4 

sage: b = ~(5*a-3); b 

5/116*a + 3/116 

sage: b*(5*a-3) 

1 

sage: b = ~((3*a-2)/7); b 

21/41*a + 14/41 

sage: (3*a-2)/7 * b 

1 

  

This fixes ticket :trac:`9357`:: 

  

sage: K.<a> = NumberField(x^2+1) 

sage: d = K(0) 

sage: ~d 

Traceback (most recent call last): 

... 

ZeroDivisionError: number field element division by zero 

sage: K.random_element() / d 

Traceback (most recent call last): 

... 

ZeroDivisionError: number field element division by zero 

""" 

if mpz_cmp_ui(self.a, 0) == 0 and mpz_cmp_ui(self.b, 0) == 0: 

raise ZeroDivisionError("number field element division by zero") 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

cdef mpz_t tmp, gcd 

mpz_init(tmp) 

mpz_init(gcd) 

  

mpz_gcd(gcd, self.a, self.b) 

if mpz_cmp_si(gcd, 1): # != 0 (i.e. it is not 1) 

# cancel out g (g(a'-b'd)) / (g^2 (a'^2-b'^2d^2)) 

mpz_divexact(res.a, self.a, gcd) 

mpz_divexact(res.b, self.b, gcd) 

mpz_neg(res.b, res.b) 

else: 

mpz_set(res.a, self.a) 

mpz_neg(res.b, self.b) 

  

mpz_pow_ui(res.denom, res.a, 2) 

mpz_pow_ui(tmp, res.b, 2) 

mpz_mul(tmp, tmp, self.D.value) 

mpz_sub(res.denom, res.denom, tmp) 

# need to multiply the leftover g back in 

mpz_mul(res.denom, res.denom, gcd) 

  

mpz_mul(res.a, res.a, self.denom) 

mpz_mul(res.b, res.b, self.denom) 

  

mpz_clear(tmp) 

mpz_clear(gcd) 

  

res._reduce_c_() 

return res 

  

cpdef NumberFieldElement galois_conjugate(self): 

""" 

Return the image of this element under action of the nontrivial 

element of the Galois group of this field. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(23) 

sage: a.galois_conjugate() 

-a 

  

sage: K.<a> = NumberField(x^2 - 5*x + 1) 

sage: a.galois_conjugate() 

-a + 5 

sage: b = 5*a + 1/3 

sage: b.galois_conjugate() 

-5*a + 76/3 

sage: b.norm() == b * b.galois_conjugate() 

True 

sage: b.trace() == b + b.galois_conjugate() 

True 

""" 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

mpz_set(res.a, self.a) 

mpz_neg(res.b, self.b) 

mpz_set(res.denom, self.denom) 

return res 

  

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

# We must override everything that makes uses of self.__numerator/__denominator 

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

  

def __hash__(self): 

""" 

Return hash of this number field element. 

  

For elements in `\ZZ` or `\QQ` the hash coincides with the one in the 

native `\ZZ` or `\QQ`. 

  

EXAMPLES:: 

  

sage: L.<a> = QuadraticField(-7) 

sage: hash(a) 

42082631 

sage: hash(L(1)) 

1 

sage: hash(L(-3)) 

-3 

sage: hash(L(-32/118)) == hash(-32/118) 

True 

""" 

# 1. compute the hash of a/denom as if it was a rational 

# (see the corresponding code in sage/rings/rational.pyx) 

cdef Py_hash_t n = mpz_pythonhash(self.a) 

cdef Py_hash_t d = mpz_pythonhash(self.denom) 

cdef Py_hash_t h = n + (d - 1) * <Py_hash_t>(7461864723258187525) 

  

# 2. mix the hash together with b 

h += 42082631 * mpz_pythonhash(self.b) 

return h 

  

def __nonzero__(self): 

""" 

Check whether this element is not zero. 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+163) 

sage: not a 

False 

sage: not (a-a) 

True 

""" 

return mpz_cmp_ui(self.a, 0) != 0 or mpz_cmp_ui(self.b, 0) != 0 

  

def _integer_(self, Z=None): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+163) 

sage: (a+1-a)._integer_() 

1 

sage: (a+1/2-a)._integer_() 

Traceback (most recent call last): 

... 

TypeError: Unable to coerce 1/2 to an integer 

""" 

cdef Integer res 

if mpz_cmp_ui(self.b, 0) != 0 or mpz_cmp_ui(self.denom, 1) != 0: 

raise TypeError("Unable to coerce %s to an integer" % self) 

else: 

res = Integer.__new__(Integer) 

mpz_set(res.value, self.a) 

return res 

  

def _rational_(self): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+163) 

sage: (a+1/2-a)._rational_() 

1/2 

sage: (a+1/2)._rational_() 

Traceback (most recent call last): 

... 

TypeError: Unable to coerce a + 1/2 to a rational 

""" 

cdef Rational res 

if mpz_cmp_ui(self.b, 0)!=0: 

raise TypeError("Unable to coerce %s to a rational" % self) 

else: 

res = <Rational>Rational.__new__(Rational) 

mpz_set(mpq_numref(res.value), self.a) 

mpz_set(mpq_denref(res.value), self.denom) 

mpq_canonicalize(res.value) 

return res 

  

def _algebraic_(self, parent): 

r""" 

Convert this element to an algebraic number, if possible. 

  

EXAMPLES:: 

  

sage: NF.<i> = QuadraticField(-1) 

sage: QQbar(1+i) 

I + 1 

sage: NF.<sqrt3> = QuadraticField(2) 

sage: AA(sqrt3) 

1.414213562373095? 

""" 

import sage.rings.qqbar as qqbar 

if (parent is qqbar.QQbar 

and list(self._parent.polynomial()) == [1, 0, 1]): 

# AlgebraicNumber.__init__ does a better job than 

# NumberFieldElement._algebraic_ in this case, but 

# QQbar._element_constructor_ calls the latter first. 

return qqbar.AlgebraicNumber(self) 

else: 

return NumberFieldElement._algebraic_(self, parent) 

  

cpdef bint is_one(self): 

r""" 

Check whether this number field element is `1`. 

  

EXAMPLES:: 

  

sage: K = QuadraticField(-2) 

sage: K(1).is_one() 

True 

sage: K(-1).is_one() 

False 

sage: K(2).is_one() 

False 

sage: K(0).is_one() 

False 

sage: K(1/2).is_one() 

False 

sage: K.gen().is_one() 

False 

""" 

return mpz_cmp_ui(self.a, 1) == 0 and \ 

mpz_cmp_ui(self.b, 0) == 0 and \ 

mpz_cmp_ui(self.denom, 1) == 0 

  

cpdef bint is_rational(self): 

r""" 

Check whether this number field element is a rational number. 

  

.. SEEALSO:: 

  

- :meth:`is_integer` to test if this element is an integer 

- :meth:`is_integral` to test if this element is an algebraic integer 

  

EXAMPLES:: 

  

sage: K.<sqrt3> = QuadraticField(3) 

sage: sqrt3.is_rational() 

False 

sage: (sqrt3-1/2).is_rational() 

False 

sage: K(0).is_rational() 

True 

sage: K(-12).is_rational() 

True 

sage: K(1/3).is_rational() 

True 

""" 

return mpz_cmp_ui(self.b, 0) == 0 

  

def is_integer(self): 

r""" 

Check whether this number field element is an integer. 

  

.. SEEALSO:: 

  

- :meth:`is_rational` to test if this element is a rational number 

- :meth:`is_integral` to test if this element is an algebraic integer 

  

EXAMPLES:: 

  

sage: K.<sqrt3> = QuadraticField(3) 

sage: sqrt3.is_integer() 

False 

sage: (sqrt3-1/2).is_integer() 

False 

sage: K(0).is_integer() 

True 

sage: K(-12).is_integer() 

True 

sage: K(1/3).is_integer() 

False 

""" 

return mpz_cmp_ui(self.b, 0) == mpz_cmp_ui(self.denom, 1) == 0 

  

def real(self): 

r""" 

Returns the real part of self, which is either self (if self lives 

it a totally real field) or a rational number. 

  

EXAMPLES:: 

  

sage: K.<sqrt2> = QuadraticField(2) 

sage: sqrt2.real() 

sqrt2 

sage: K.<a> = QuadraticField(-3) 

sage: a.real() 

0 

sage: (a + 1/2).real() 

1/2 

sage: K.<a> = NumberField(x^2 + x + 1) 

sage: a.real() 

-1/2 

sage: parent(a.real()) 

Rational Field 

sage: K.<i> = QuadraticField(-1) 

sage: i.real() 

0 

""" 

cdef Rational res 

if mpz_sgn(self.D.value) > 0: 

return self # totally real 

else: 

res = <Rational>Rational.__new__(Rational) 

mpz_set(mpq_numref(res.value), self.a) 

mpz_set(mpq_denref(res.value), self.denom) 

mpq_canonicalize(res.value) 

return res 

  

def imag(self): 

r""" 

Returns the imaginary part of self. 

  

EXAMPLES:: 

  

sage: K.<sqrt2> = QuadraticField(2) 

sage: sqrt2.imag() 

0 

sage: parent(sqrt2.imag()) 

Rational Field 

  

sage: K.<i> = QuadraticField(-1) 

sage: i.imag() 

1 

sage: parent(i.imag()) 

Rational Field 

  

sage: K.<a> = NumberField(x^2 + x + 1, embedding=CDF.0) 

sage: a.imag() 

1/2*sqrt3 

sage: a.real() 

-1/2 

sage: SR(a) 

1/2*I*sqrt(3) - 1/2 

sage: bool(I*a.imag() + a.real() == a) 

True 

  

TESTS:: 

  

sage: K.<a> = QuadraticField(-9, embedding=-CDF.0) 

sage: a.imag() 

-3 

sage: parent(a.imag()) 

Rational Field 

  

Check that :trac:`22095` is fixed:: 

  

sage: K.<a> = NumberField(x^2 + 2*x + 14, embedding=CC(-1,+3)) 

sage: K13.<sqrt13> = QuadraticField(13) 

sage: K13.zero() 

0 

sage: a.imag() 

sqrt13 

sage: K13.zero() 

0 

""" 

if mpz_sgn(self.D.value) > 0: 

return Rational.__new__(Rational) # = 0 

cdef Integer negD = <Integer>Integer.__new__(Integer) 

mpz_neg(negD.value, self.D.value) 

cdef NumberFieldElement_quadratic q = <NumberFieldElement_quadratic>self._new() 

mpz_set_ui(q.b, 1) 

mpz_set_ui(q.denom, 1) 

cdef Rational res 

if mpz_cmp_ui(negD.value, 1) == 0 or mpz_perfect_square_p(negD.value): 

# D = -1 is the most common case we'll see here 

if self._parent._embedding is None: 

raise ValueError("Embedding must be specified.") 

res = <Rational>Rational.__new__(Rational) 

if mpz_cmp_ui(negD.value, 1) == 0: 

mpz_set(mpq_numref(res.value), self.b) 

else: 

mpz_sqrt(mpq_numref(res.value), negD.value) 

mpz_mul(mpq_numref(res.value), mpq_numref(res.value), self.b) 

mpz_set(mpq_denref(res.value), self.denom) 

mpq_canonicalize(res.value) 

if not self.standard_embedding: 

mpq_neg(res.value, res.value) 

return res 

else: 

# avoid circular import 

if self._parent._embedding is None: 

from .number_field import NumberField 

K = NumberField(QQ['x'].gen()**2 - negD, 'sqrt%s' % negD) 

else: 

from .number_field import QuadraticField 

K = QuadraticField(negD, 'sqrt%s' % negD) 

q = (<NumberFieldElement_quadratic> K._zero_element)._new() 

mpz_set(q.denom, self.denom) 

mpz_set_ui(q.a, 0) 

if self.standard_embedding: 

mpz_set(q.b, self.b) 

else: 

mpz_neg(q.b, self.b) 

return q 

  

cpdef list _coefficients(self): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+41) 

sage: a._coefficients() 

[0, 1] 

sage: K.<a> = NumberField(x^2+x+41) 

sage: a._coefficients() 

[0, 1] 

sage: b = 3*a+1/5 

sage: b._coefficients() 

[1/5, 3] 

""" 

# In terms of the generator... 

cdef NumberFieldElement_quadratic gen = self.number_field().gen() # should this be cached? 

cdef Rational const = <Rational>Rational.__new__(Rational) 

cdef Rational lin = <Rational>Rational.__new__(Rational) 

ad, bd = self.parts() 

if not self: 

return [] 

if not bd: 

return [ad] 

if gen.is_sqrt_disc(): 

return [ad,bd] 

else: 

alpha, beta = gen.parts() 

scale = bd/beta 

return [ad - scale*alpha, scale] 

  

def denominator(self): 

""" 

Return the denominator of self. This is the LCM of the denominators of 

the coefficients of self, and thus it may well be `> 1` even when the 

element is an algebraic integer. 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+x+41) 

sage: a.denominator() 

1 

sage: b = (2*a+1)/6 

sage: b.denominator() 

6 

sage: K(1).denominator() 

1 

sage: K(1/2).denominator() 

2 

sage: K(0).denominator() 

1 

  

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

sage: b = (a + 1)/2 

sage: b.denominator() 

2 

sage: b.is_integral() 

True 

""" 

# In terms of the generator... 

cdef NumberFieldElement_quadratic gen = self.number_field().gen() # should this be cached? 

cdef Integer denom 

if gen.is_sqrt_disc(): 

denom = Integer.__new__(Integer) 

mpz_set(denom.value, self.denom) 

return denom 

else: 

c = self._coefficients() 

if len(c) == 2: 

const, lin = c 

elif len(c) == 1: 

const = c[0] 

lin = Rational(0) 

else: 

const = lin = Rational(0) 

return const.denominator().lcm(lin.denominator()) 

  

def numerator(self): 

""" 

Return self*self.denominator(). 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+x+41) 

sage: b = (2*a+1)/6 

sage: b.denominator() 

6 

sage: b.numerator() 

2*a + 1 

""" 

return self*self.denominator() 

  

  

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

# Some things are so much easier to compute 

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

  

def trace(self): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+x+41) 

sage: a.trace() 

-1 

sage: a.matrix() 

[ 0 1] 

[-41 -1] 

  

The trace is additive:: 

  

sage: K.<a> = NumberField(x^2+7) 

sage: (a+1).trace() 

2 

sage: K(3).trace() 

6 

sage: (a+4).trace() 

8 

sage: (a/3+1).trace() 

2 

""" 

# trace = 2*self.a / self.denom 

cdef Rational res = <Rational>Rational.__new__(Rational) 

if mpz_odd_p(self.denom): 

mpz_mul_2exp(mpq_numref(res.value), self.a, 1) 

mpz_set(mpq_denref(res.value), self.denom) 

else: 

mpz_set(mpq_numref(res.value), self.a) 

mpz_divexact_ui(mpq_denref(res.value), self.denom, 2) 

mpq_canonicalize(res.value) 

return res 

  

  

def norm(self, K=None): 

""" 

Return the norm of self. If the second argument is None, this is the 

norm down to `\QQ`. Otherwise, return the norm down to K (which had 

better be either `\QQ` or this number field). 

  

EXAMPLES:: 

  

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

sage: a.norm() 

3 

sage: a.matrix() 

[ 0 1] 

[-3 1] 

sage: K.<a> = NumberField(x^2+5) 

sage: (1+a).norm() 

6 

  

The norm is multiplicative:: 

  

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

sage: a.norm() 

-3 

sage: K(3).norm() 

9 

sage: (3*a).norm() 

-27 

  

We test that the optional argument is handled sensibly:: 

  

sage: (3*a).norm(QQ) 

-27 

sage: (3*a).norm(K) 

3*a 

sage: (3*a).norm(CyclotomicField(3)) 

Traceback (most recent call last): 

... 

ValueError: no way to embed L into parent's base ring K 

  

""" 

cdef Rational res = <Rational>Rational.__new__(Rational) 

  

if K is None or K == QQ: 

# norm = (a^2 - d b^2) / self.denom^2 

mpz_pow_ui(mpq_numref(res.value), self.a, 2) 

mpz_pow_ui(mpq_denref(res.value), self.b, 2) # use as temp 

mpz_mul(mpq_denref(res.value), mpq_denref(res.value), self.D.value) 

mpz_sub(mpq_numref(res.value), mpq_numref(res.value), mpq_denref(res.value)) 

mpz_pow_ui(mpq_denref(res.value), self.denom, 2) 

mpq_canonicalize(res.value) 

return res 

else: 

return NumberFieldElement.norm(self, K) 

  

def is_integral(self): 

r""" 

Return whether this element is an algebraic integer. 

  

TESTS:: 

  

sage: K.<a> = QuadraticField(-1) 

sage: a.is_integral() 

True 

sage: K(1).is_integral() 

True 

sage: K(1/2).is_integral() 

False 

sage: (a/2).is_integral() 

False 

sage: ((a+1)/2).is_integral() 

False 

sage: ((a+1)/3).is_integral() 

False 

  

sage: K.<a> = QuadraticField(-3) 

sage: a.is_integral() 

True 

sage: K(1).is_integral() 

True 

sage: K(1/2).is_integral() 

False 

sage: (a/2).is_integral() 

False 

sage: ((a+1)/2).is_integral() 

True 

sage: ((a+1)/3).is_integral() 

False 

  

This works for order elements too, see :trac:`24077`:: 

  

sage: O.<w> = EisensteinIntegers() 

sage: w.is_integral() 

True 

sage: for _ in range(20): 

....: assert O.random_element().is_integral() 

""" 

if mpz_cmp_ui(self.denom, 1) == 0: 

return True 

  

# Check for an element of the form x + y sqrt(D) where x and y 

# are half-integers. 

if mpz_even_p(self.a) or mpz_even_p(self.b): 

return False 

if mpz_cmp_ui(self.denom, 2) != 0: 

return False 

# Numbers with half-integral components are integral only for 

# D = 1 mod 4 

return mpz_fdiv_ui(self.D.value, 4) == 1 

  

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

r""" 

The characteristic polynomial of this element over `\QQ`. 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2-x+13) 

sage: a.charpoly() 

x^2 - x + 13 

sage: b = 3-a/2 

sage: f = b.charpoly(); f 

x^2 - 11/2*x + 43/4 

sage: f(b) 

0 

""" 

R = QQ[var] 

return R([self.norm(), -self.trace(), 1]) 

  

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

r""" 

The minimal polynomial of this element over `\QQ`. 

  

INPUT: 

  

- ``var`` -- the minimal polynomial is defined over a polynomial ring 

in a variable with this name. If not specified this defaults to 

``x``. 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+13) 

sage: a.minpoly() 

x^2 + 13 

sage: a.minpoly('T') 

T^2 + 13 

sage: (a+1/2-a).minpoly() 

x - 1/2 

""" 

if self.is_rational(): 

R = QQ[var] 

return R([-self._rational_(), 1]) 

else: 

return self.charpoly(var) 

  

def __abs__(self): 

""" 

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(2, 'a', embedding=-1.4142) 

sage: abs(a) # indirect test 

-a 

sage: abs(a+1) # indirect test 

-a - 1 

sage: abs(a+2) # indirect test 

a + 2 

  

sage: K.<a> = NumberField(x^2+1, embedding=CDF.gen()) 

sage: abs(a+1) 

sqrt(2) 

""" 

if mpz_sgn(self.D.value) == 1: 

if self.sign() >= 0: 

return self 

return -self 

  

# doing that way the parent is the symbolic ring (or IntegerRing if the 

# norm of self is a square). On the other hand, it is coherent with 

# sage.rings.integer.Integer.sqrt 

return (self.real()**2 + self.imag()**2).sqrt() 

  

def floor(self): 

r""" 

Returns the floor of x. 

  

EXAMPLES:: 

  

sage: K.<sqrt2> = QuadraticField(2,name='sqrt2') 

sage: sqrt2.floor() 

1 

sage: (-sqrt2).floor() 

-2 

sage: (13/197 + 3702/123*sqrt2).floor() 

42 

sage: (13/197-3702/123*sqrt2).floor() 

-43 

  

TESTS:: 

  

sage: K2.<sqrt2> = QuadraticField(2) 

sage: K3.<sqrt3> = QuadraticField(3) 

sage: K5.<sqrt5> = QuadraticField(5) 

sage: for _ in range(100): 

....: a = QQ.random_element(1000,20) 

....: b = QQ.random_element(1000,20) 

....: assert floor(a+b*sqrt(2.)) == floor(a+b*sqrt2) 

....: assert floor(a+b*sqrt(3.)) == floor(a+b*sqrt3) 

....: assert floor(a+b*sqrt(5.)) == floor(a+b*sqrt5) 

  

sage: K = QuadraticField(-2) 

sage: l = [K(52), K(-3), K(43/12), K(-43/12)] 

sage: [x.floor() for x in l] 

[52, -3, 3, -4] 

""" 

cdef mpz_t x 

cdef Integer result 

  

if mpz_sgn(self.b) == 0: 

mpz_init_set(x,self.a) 

mpz_fdiv_q(x,x,self.denom) 

result = Integer.__new__(Integer) 

mpz_set(result.value,x) 

mpz_clear(x) 

return result 

  

if not mpz_sgn(self.D.value) == 1: 

raise ValueError("floor is not defined for complex quadratic number field") 

  

mpz_init(x) 

mpz_mul(x,self.b,self.b) 

mpz_mul(x,x,self.D.value) 

mpz_sqrt(x,x) 

if mpz_sgn(self.b) == -1: 

if self.standard_embedding: 

mpz_neg(x,x) 

mpz_sub_ui(x,x,1) 

elif not self.standard_embedding: 

mpz_neg(x,x) 

mpz_sub_ui(x,x,1) 

  

mpz_add(x,x,self.a) # here x = a + floor(sqrt(b^2 D)) or a + floor(-sqrt(b^2 D)) 

mpz_fdiv_q(x,x,self.denom) 

result = Integer.__new__(Integer) 

mpz_set(result.value,x) 

mpz_clear(x) 

return result 

  

def ceil(self): 

r""" 

Returns the ceil. 

  

EXAMPLES:: 

  

sage: K.<sqrt7> = QuadraticField(7, name='sqrt7') 

sage: sqrt7.ceil() 

3 

sage: (-sqrt7).ceil() 

-2 

sage: (1022/313*sqrt7 - 14/23).ceil() 

9 

  

TESTS:: 

  

sage: K2.<sqrt2> = QuadraticField(2) 

sage: K3.<sqrt3> = QuadraticField(3) 

sage: K5.<sqrt5> = QuadraticField(5) 

sage: for _ in range(100): 

....: a = QQ.random_element(1000,20) 

....: b = QQ.random_element(1000,20) 

....: assert ceil(a+b*sqrt(2.)) == ceil(a+b*sqrt2) 

....: assert ceil(a+b*sqrt(3.)) == ceil(a+b*sqrt3) 

....: assert ceil(a+b*sqrt(5.)) == ceil(a+b*sqrt5) 

  

sage: K = QuadraticField(-2) 

sage: l = [K(52), K(-3), K(43/12), K(-43/12)] 

sage: [x.ceil() for x in l] 

[52, -3, 4, -3] 

""" 

x = self.floor() 

if mpz_sgn(self.b) == 0 and mpz_cmp_ui(self.denom,1) == 0: 

return x 

return x+1 

  

def round(self): 

r""" 

Returns the round (nearest integer). 

  

EXAMPLES:: 

  

sage: K.<sqrt7> = QuadraticField(7, name='sqrt7') 

sage: sqrt7.round() 

3 

sage: (-sqrt7).round() 

-3 

sage: (12/313*sqrt7 - 1745917/2902921).round() 

0 

sage: (12/313*sqrt7 - 1745918/2902921).round() 

-1 

  

TESTS:: 

  

sage: K2.<sqrt2> = QuadraticField(2) 

sage: K3.<sqrt3> = QuadraticField(3) 

sage: K5.<sqrt5> = QuadraticField(5) 

sage: for _ in range(100): 

....: a = QQ.random_element(1000,20) 

....: b = QQ.random_element(1000,20) 

....: assert round(a) == round(K2(a)), a 

....: assert round(a) == round(K3(a)), a 

....: assert round(a) == round(K5(a)), a 

....: assert round(a+b*sqrt(2.)) == round(a+b*sqrt2), (a, b) 

....: assert round(a+b*sqrt(3.)) == round(a+b*sqrt3), (a, b) 

....: assert round(a+b*sqrt(5.)) == round(a+b*sqrt5), (a, b) 

""" 

n = self.floor() 

test = 2 * (self - n).abs() 

if test < 1: 

return n 

elif test > 1: 

return n + 1 

elif self > 0: 

return n + 1 

else: 

return n 

  

cdef class OrderElement_quadratic(NumberFieldElement_quadratic): 

""" 

Element of an order in a quadratic field. 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2 + 1) 

sage: O2 = K.order(2*a) 

sage: w = O2.1; w 

2*a 

sage: parent(w) 

Order in Number Field in a with defining polynomial x^2 + 1 

""" 

def __init__(self, order, f): 

r""" 

Standard initialisation function. 

  

EXAMPLES:: 

  

sage: OK.<y> = EquationOrder(x^2 + 5) 

sage: v = OK.1 # indirect doctest 

sage: type(v) 

<type 'sage.rings.number_field.number_field_element_quadratic.OrderElement_quadratic'> 

""" 

K = order.number_field() 

NumberFieldElement_quadratic.__init__(self, K, f) 

(<Element>self)._parent = order 

  

def norm(self): 

""" 

The norm of an element of the ring of integers is an Integer. 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2 + 3) 

sage: O2 = K.order(2*a) 

sage: w = O2.gen(1); w 

2*a 

sage: w.norm() 

12 

sage: parent(w.norm()) 

Integer Ring 

""" 

return ZZ(NumberFieldElement_quadratic.norm(self)) 

  

def trace(self): 

""" 

The trace of an element of the ring of integers is an Integer. 

  

EXAMPLES:: 

  

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

sage: R = K.ring_of_integers() 

sage: b = R((1+a)/2) 

sage: b.trace() 

1 

sage: parent(b.trace()) 

Integer Ring 

""" 

return ZZ(NumberFieldElement_quadratic.trace(self)) 

  

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

r""" 

The characteristic polynomial of this element, which is over `\ZZ` 

because this element is an algebraic integer. 

  

EXAMPLES:: 

  

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

sage: R = K.ring_of_integers() 

sage: b = R((5+a)/2) 

sage: f = b.charpoly('x'); f 

x^2 - 5*x + 5 

sage: f.parent() 

Univariate Polynomial Ring in x over Integer Ring 

sage: f(b) 

0 

""" 

R = ZZ[var] 

return R([self.norm(), -self.trace(), 1]) 

  

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

r""" 

The minimal polynomial of this element over `\ZZ`. 

  

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2 + 163) 

sage: R = K.ring_of_integers() 

sage: f = R(a).minpoly('x'); f 

x^2 + 163 

sage: f.parent() 

Univariate Polynomial Ring in x over Integer Ring 

sage: R(5).minpoly() 

x - 5 

""" 

if self.is_rational(): 

R = ZZ[var] 

return R([-self._rational_(), 1]) 

else: 

return self.charpoly() 

  

cdef number_field(self): 

# So few functions actually use self.number_field() for quadratic elements, so 

# it is better *not* to return a cached value (since the call to _parent.number_field()) 

# is expensive. 

return self._parent.number_field() 

  

# We must override these since the basering is now ZZ not QQ. 

cpdef _rmul_(self, Element _c): 

""" 

EXAMPLES:: 

  

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

sage: R = K.ring_of_integers() 

sage: aa = R.gen(1); aa 

1/3*a 

sage: 5 * aa # indirect doctest 

5/3*a 

""" 

cdef Integer c = <Integer>_c 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

mpz_mul(res.a, self.a, c.value) 

mpz_mul(res.b, self.b, c.value) 

mpz_set(res.denom, self.denom) 

res._reduce_c_() 

return res 

  

cpdef _lmul_(self, Element _c): 

""" 

EXAMPLES:: 

  

sage: K.<a> = NumberField(x^2+43) 

sage: R = K.ring_of_integers() 

sage: aa = R.gen(0); aa 

1/2*a + 1/2 

sage: aa*3 # indirect doctest 

3/2*a + 3/2 

""" 

cdef Integer c = <Integer>_c 

cdef NumberFieldElement_quadratic res = <NumberFieldElement_quadratic>self._new() 

mpz_mul(res.a, self.a, c.value) 

mpz_mul(res.b, self.b, c.value) 

mpz_set(res.denom, self.denom) 

res._reduce_c_() 

return res 

  

def __invert__(self): 

r""" 

Implement inversion, checking that the return value has the right parent. 

See trac \#4190. 

  

EXAMPLES:: 

  

sage: K = NumberField(x^2 -x + 2, 'a') 

sage: OK = K.ring_of_integers() 

sage: a = OK(K.gen()) 

sage: (~a).parent() is K 

True 

sage: (~a) in OK 

False 

sage: a**(-1) in OK 

False 

""" 

return self._parent.number_field()(NumberFieldElement_quadratic.__invert__(self)) 

  

def inverse_mod(self, I): 

r""" 

Return an inverse of self modulo the given ideal. 

  

INPUT: 

  

  

- ``I`` - may be an ideal of self.parent(), or an 

element or list of elements of self.parent() generating a nonzero 

ideal. A ValueError is raised if I is non-integral or is zero. 

A ZeroDivisionError is raised if I + (x) != (1). 

  

  

EXAMPLES:: 

  

sage: OE.<w> = EquationOrder(x^2 - x + 2) 

sage: w.inverse_mod(13) == 6*w - 6 

True 

sage: w*(6*w - 6) - 1 

-13 

sage: w.inverse_mod(13).parent() == OE 

True 

sage: w.inverse_mod(2*OE) 

Traceback (most recent call last): 

... 

ZeroDivisionError: w is not invertible modulo Fractional ideal (2) 

""" 

R = self.parent() 

return R(_inverse_mod_generic(self, I)) 

  

  

cdef class Z_to_quadratic_field_element(Morphism): 

""" 

Morphism that coerces from integers to elements of a quadratic number 

field K. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(3) 

sage: phi = K.coerce_map_from(ZZ); phi 

Natural morphism: 

From: Integer Ring 

To: Number Field in a with defining polynomial x^2 - 3 

sage: phi(4) 

4 

sage: phi(5).parent() is K 

True 

""" 

# The zero element of K, lazily initialized 

cdef NumberFieldElement_quadratic zero_element 

  

# TODO: implement __richcmp__ properly so we can have a loads/dumps doctest 

  

def __init__(self, K): 

""" 

``K`` is the target quadratic field 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(3) 

sage: phi = K.coerce_map_from(ZZ) # indirect doctest 

sage: type(phi) 

<type 'sage.rings.number_field.number_field_element_quadratic.Z_to_quadratic_field_element'> 

sage: phi == loads(dumps(phi)) # todo: comparison not implemented 

True 

  

sage: R.<b> = CyclotomicField(6) 

sage: psi = R.coerce_map_from(ZZ) # indirect doctest 

sage: type(psi) 

<type 'sage.rings.number_field.number_field_element_quadratic.Z_to_quadratic_field_element'> 

sage: psi == loads(dumps(psi)) # todo: comparison not implemented 

True 

""" 

import sage.categories.homset 

Morphism.__init__(self, sage.categories.homset.Hom(ZZ, K)) 

  

cpdef Element _call_(self, x): 

r""" 

Evaluate at an integer ``x``. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(3) 

sage: phi = K.coerce_map_from(ZZ) 

sage: a = phi(2); a # indirect doctest 

2 

sage: a.parent() is K 

True 

  

sage: R.<b> = CyclotomicField(6) 

sage: psi = R.coerce_map_from(ZZ) 

sage: b = psi(-42); b # indirect doctest 

-42 

sage: b.parent() is R 

True 

""" 

cdef NumberFieldElement_quadratic y 

if self.zero_element is None: 

self.zero_element = self._codomain.zero() 

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

return self.zero_element 

y = self.zero_element._new() 

mpz_set(y.a, (<Integer> x).value) 

  

# we need to set the denominator to 1 as it is 0 for 

# the zero element of K... (because gcd(0,0) = 0). 

mpz_set_ui(y.denom, 1) 

  

return y 

  

def _repr_type(self): 

r""" 

Return a short name for this morphism. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(3) 

sage: phi = K.coerce_map_from(ZZ) 

sage: phi # indirect doctest 

Natural morphism: 

From: Integer Ring 

To: Number Field in a with defining polynomial x^2 - 3 

  

sage: R.<b> = CyclotomicField(6) 

sage: psi = R.coerce_map_from(ZZ) 

sage: psi # indirect doctest 

Natural morphism: 

From: Integer Ring 

To: Cyclotomic Field of order 6 and degree 2 

""" 

return "Natural" 

  

  

cdef class Q_to_quadratic_field_element(Morphism): 

""" 

Morphism that coerces from rationals to elements of a quadratic number 

field K. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(-3) 

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

Natural morphism: 

From: Rational Field 

To: Number Field in a with defining polynomial x^2 + 3 

sage: f(3/1) 

3 

sage: f(1/2).parent() is K 

True 

""" 

# The zero element of K, lazily initialized 

cdef NumberFieldElement_quadratic zero_element 

  

# TODO: implement __richcmp__ properly so we can have a loads/dumps doctest 

  

def __init__(self, K): 

""" 

``K`` is the target quadratic field 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(3) 

sage: phi = K.coerce_map_from(QQ) # indirect doctest 

sage: type(phi) 

<type 'sage.rings.number_field.number_field_element_quadratic.Q_to_quadratic_field_element'> 

sage: phi == loads(dumps(phi)) # todo: comparison not implemented 

True 

  

sage: R.<b> = CyclotomicField(6) 

sage: psi = R.coerce_map_from(QQ) 

sage: type(psi) 

<type 'sage.rings.number_field.number_field_element_quadratic.Q_to_quadratic_field_element'> 

sage: psi == loads(dumps(psi)) # todo: comparison not implemented 

True 

""" 

import sage.categories.homset 

Morphism.__init__(self, sage.categories.homset.Hom(QQ, K)) 

  

cpdef Element _call_(self, x): 

r""" 

Evaluate at a rational ``x``. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(3) 

sage: phi = K.coerce_map_from(QQ) 

sage: a = phi(2/3); a # indirect doctest 

2/3 

sage: a.parent() is K 

True 

  

sage: R.<b> = CyclotomicField(6) 

sage: psi = R.coerce_map_from(QQ) 

sage: b = psi(-23/15); b # indirect doctest 

-23/15 

sage: b.parent() is R 

True 

""" 

if self.zero_element is None: 

self.zero_element = self._codomain.zero() 

cdef NumberFieldElement_quadratic y = self.zero_element._new() 

mpz_set(y.a, mpq_numref((<Rational>x).value)) 

mpz_set(y.denom, mpq_denref((<Rational>x).value)) 

return y 

  

def _repr_type(self): 

r""" 

Return a short name for this morphism. 

  

EXAMPLES:: 

  

sage: K.<a> = QuadraticField(3) 

sage: phi = K.coerce_map_from(QQ) 

sage: phi # indirect doctest 

Natural morphism: 

From: Rational Field 

To: Number Field in a with defining polynomial x^2 - 3 

  

sage: R.<b> = CyclotomicField(6) 

sage: psi = R.coerce_map_from(QQ) 

sage: psi # indirect doctest 

Natural morphism: 

From: Rational Field 

To: Cyclotomic Field of order 6 and degree 2 

""" 

return "Natural"