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

""" 

Arbitrary Precision Complex Numbers using GNU MPC 

  

This is a binding for the MPC arbitrary-precision floating point library. 

It is adaptated from ``real_mpfr.pyx`` and ``complex_number.pyx``. 

  

We define a class :class:`MPComplexField`, where each instance of 

``MPComplexField`` specifies a field of floating-point complex numbers with 

a specified precision shared by the real and imaginary part and a rounding 

mode stating the rounding mode directions specific to real and imaginary 

parts. 

  

Individual floating-point numbers are of class :class:`MPComplexNumber`. 

  

For floating-point representation and rounding mode description see the 

documentation for the :mod:`sage.rings.real_mpfr`. 

  

AUTHORS: 

  

- Philippe Theveny (2008-10-13): initial version. 

  

- Alex Ghitza (2008-11): cache, generators, random element, and many doctests. 

  

- Yann Laigle-Chapuy (2010-01): improves compatibility with CC, updates. 

  

- Jeroen Demeyer (2012-02): reformat documentation, make MPC a standard 

package. 

  

- Travis Scrimshaw (2012-10-18): Added doctests for full coverage. 

  

- Vincent Klein (2017-11-15) : add __mpc__() to class MPComplexNumber. 

MPComplexNumber constructor support gmpy2.mpz, gmpy2.mpq, gmpy2.mpfr 

and gmpy2.mpc parameters. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(42) 

sage: a = MPC(12, '15.64E+32'); a 

12.0000000000 + 1.56400000000e33*I 

sage: a *a *a *a 

5.98338564121e132 - 1.83633318912e101*I 

sage: a + 1 

13.0000000000 + 1.56400000000e33*I 

sage: a / 3 

4.00000000000 + 5.21333333333e32*I 

sage: MPC("infinity + NaN *I") 

+infinity + NaN*I 

""" 

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

# Copyright (C) 2008 Philippe Theveny <thevenyp@loria.fr> 

# 2008 Alex Ghitza 

# 2010 Yann Laigle-Chapuy 

# 2012 Jeroen Demeyer <jdemeyer@cage.ugent.be> 

# 

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

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

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

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

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

  

from __future__ import absolute_import, print_function 

  

  

import sage 

import re 

from . import real_mpfr 

import weakref 

from cpython.object cimport Py_NE 

  

from sage.libs.mpfr cimport * 

from sage.libs.mpc cimport * 

from sage.structure.parent cimport Parent 

from sage.structure.parent_gens cimport ParentWithGens 

from sage.structure.element cimport RingElement, Element, ModuleElement 

from sage.categories.map cimport Map 

from sage.libs.pari.all import pari 

  

from .integer cimport Integer 

from .complex_number cimport ComplexNumber 

from .complex_field import ComplexField_class 

  

from sage.misc.randstate cimport randstate, current_randstate 

from sage.misc.superseded import deprecated_function_alias 

from .real_mpfr cimport RealField_class, RealNumber 

from .real_mpfr import mpfr_prec_min, mpfr_prec_max 

from sage.structure.richcmp cimport rich_to_bool, richcmp 

from sage.categories.fields import Fields 

  

IF HAVE_GMPY2: 

cimport gmpy2 

gmpy2.import_gmpy2() 

  

NumberFieldElement_quadratic = None 

AlgebraicNumber_base = None 

AlgebraicNumber = None 

AlgebraicReal = None 

AA = None 

QQbar = None 

CDF = CLF = RLF = None 

  

def late_import(): 

""" 

Import the objects/modules after build (when needed). 

  

TESTS:: 

  

sage: sage.rings.complex_mpc.late_import() 

""" 

global NumberFieldElement_quadratic 

global AlgebraicNumber_base 

global AlgebraicNumber 

global AlgebraicReal 

global AA, QQbar 

global CLF, RLF, CDF 

if NumberFieldElement_quadratic is None: 

import sage.rings.number_field.number_field_element_quadratic as nfeq 

NumberFieldElement_quadratic = nfeq.NumberFieldElement_quadratic 

import sage.rings.qqbar 

AlgebraicNumber_base = sage.rings.qqbar.AlgebraicNumber_base 

AlgebraicNumber = sage.rings.qqbar.AlgebraicNumber 

AlgebraicReal = sage.rings.qqbar.AlgebraicReal 

AA = sage.rings.qqbar.AA 

QQbar = sage.rings.qqbar.QQbar 

from .real_lazy import CLF, RLF 

from .complex_double import CDF 

  

_mpfr_rounding_modes = ['RNDN', 'RNDZ', 'RNDU', 'RNDD'] 

  

_mpc_rounding_modes = [ 'RNDNN', 'RNDZN', 'RNDUN', 'RNDDN', 

'', '', '', '', '', '', '', '', '', '', '', '', 

'RNDNZ', 'RNDZZ', 'RNDUZ', 'RNDDZ', 

'', '', '', '', '', '', '', '', '', '', '', '', 

'RNDUN', 'RNDZU', 'RNDUU', 'RNDDU', 

'', '', '', '', '', '', '', '', '', '', '', '', 

'RNDDN', 'RNDZD', 'RNDUD', 'RNDDD' ] 

  

cdef inline mpfr_rnd_t rnd_re(mpc_rnd_t rnd): 

""" 

Return the numeric value of the real part rounding mode. This 

is an internal function. 

""" 

return <mpfr_rnd_t>(rnd & 3) 

  

cdef inline mpfr_rnd_t rnd_im(mpc_rnd_t rnd): 

""" 

Return the numeric value of the imaginary part rounding mode. 

This is an internal function. 

""" 

return <mpfr_rnd_t>(rnd >> 4) 

  

sign = '[+-]' 

digit_ten = '[0123456789]' 

exponent_ten = '[e@]' + sign + '?[0123456789]+' 

number_ten = 'inf(?:inity)?|@inf@|nan(?:\([0-9A-Z_]*\))?|@nan@(?:\([0-9A-Z_]*\))?'\ 

'|(?:' + digit_ten + '*\.' + digit_ten + '+|' + digit_ten + '+\.?)(?:' + exponent_ten + ')?' 

imaginary_ten = 'i(?:\s*\*\s*(?:' + number_ten + '))?|(?:' + number_ten + ')\s*\*\s*i' 

complex_ten = '(?P<im_first>(?P<im_first_im_sign>' + sign + ')?\s*(?P<im_first_im_abs>' + imaginary_ten + ')' \ 

'(\s*(?P<im_first_re_sign>' + sign + ')\s*(?P<im_first_re_abs>' + number_ten + '))?)' \ 

'|' \ 

'(?P<re_first>(?P<re_first_re_sign>' + sign + ')?\s*(?P<re_first_re_abs>' + number_ten + ')' \ 

'(\s*(?P<re_first_im_sign>' + sign + ')\s*(?P<re_first_im_abs>' + imaginary_ten + '))?)' 

re_complex_ten = re.compile('^\s*(?:' + complex_ten + ')\s*$', re.I) 

  

cpdef inline split_complex_string(string, int base=10): 

""" 

Split and return in that order the real and imaginary parts 

of a complex in a string. 

  

This is an internal function. 

  

EXAMPLES:: 

  

sage: sage.rings.complex_mpc.split_complex_string('123.456e789') 

('123.456e789', None) 

sage: sage.rings.complex_mpc.split_complex_string('123.456e789*I') 

(None, '123.456e789') 

sage: sage.rings.complex_mpc.split_complex_string('123.+456e789*I') 

('123.', '+456e789') 

sage: sage.rings.complex_mpc.split_complex_string('123.456e789', base=2) 

(None, None) 

""" 

if base == 10: 

number = number_ten 

z = re_complex_ten.match(string) 

else: 

all_digits = "0123456789abcdefghijklmnopqrstuvwxyz" 

digit = '[' + all_digits[0:base] + ']' 

  

# In MPFR, '1e42'-> 10^42, '1p42'->2^42, '1@42'->base^42 

if base == 2: 

exponent = '[e@p]' 

elif base <= 10: 

exponent = '[e@]' 

elif base == 16: 

exponent = '[@p]' 

else: 

exponent = '@' 

exponent += sign + '?' + digit + '+' 

  

# Warning: number, imaginary, and complex should be enclosed in parentheses 

# when used as regexp because of alternatives '|' 

number = '@nan@(?:\([0-9A-Z_]*\))?|@inf@|(?:' + digit + '*\.' + digit + '+|' + digit + '+\.?)(?:' + exponent + ')?' 

if base <= 10: 

number = 'nan(?:\([0-9A-Z_]*\))?|inf(?:inity)?|' + number 

imaginary = 'i(?:\s*\*\s*(?:' + number + '))?|(?:' + number + ')\s*\*\s*i' 

complex = '(?P<im_first>(?P<im_first_im_sign>' + sign + ')?\s*(?P<im_first_im_abs>' + imaginary + ')' \ 

'(\s*(?P<im_first_re_sign>' + sign + ')\s*(?P<im_first_re_abs>' + number + '))?)' \ 

'|' \ 

'(?P<re_first>(?P<re_first_re_sign>' + sign + ')?\s*(?P<re_first_re_abs>' + number + ')' \ 

'(\s*(?P<re_first_im_sign>' + sign + ')\s*(?P<re_first_im_abs>' + imaginary + '))?)' 

  

z = re.match('^\s*(?:' + complex + ')\s*$', string, re.I) 

  

x, y = None, None 

if z is not None: 

if z.group('im_first') is not None: 

prefix = 'im_first' 

elif z.group('re_first') is not None: 

prefix = 're_first' 

else: 

return None 

  

if z.group(prefix + '_re_abs') is not None: 

x = z.expand('\g<' + prefix + '_re_abs>') 

if z.group(prefix + '_re_sign') is not None: 

x = z.expand('\g<' + prefix + '_re_sign>') + x 

  

if z.group(prefix + '_im_abs') is not None: 

y = re.search('(?P<im_part>' + number + ')', z.expand('\g<' + prefix + '_im_abs>'), re.I) 

if y is None: 

y = '1' 

else: 

y = y.expand('\g<im_part>') 

if z.group(prefix + '_im_sign') is not None: 

y = z.expand('\g<' + prefix + '_im_sign>') + y 

  

return x, y 

  

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

# 

# MPComplex Field 

# 

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

# The complex field is in Cython, so mpc elements will have access to 

# their parent via direct C calls, which will be faster. 

  

cache = {} 

def MPComplexField(prec=53, rnd="RNDNN", names=None): 

""" 

Return the complex field with real and imaginary parts having 

prec *bits* of precision. 

  

EXAMPLES:: 

  

sage: MPComplexField() 

Complex Field with 53 bits of precision 

sage: MPComplexField(100) 

Complex Field with 100 bits of precision 

sage: MPComplexField(100).base_ring() 

Real Field with 100 bits of precision 

sage: i = MPComplexField(200).gen() 

sage: i^2 

-1.0000000000000000000000000000000000000000000000000000000000 

""" 

global cache 

mykey = (prec, rnd) 

if mykey in cache: 

X = cache[mykey] 

C = X() 

if not C is None: 

return C 

C = MPComplexField_class(prec, rnd) 

cache[mykey] = weakref.ref(C) 

return C 

  

  

cdef class MPComplexField_class(sage.rings.ring.Field): 

def __init__(self, int prec=53, rnd="RNDNN"): 

""" 

Initialize ``self``. 

  

INPUT: 

  

- ``prec`` -- (integer) precision; default = 53 

  

prec is the number of bits used to represent the matissa of 

both the real and imaginary part of complex floating-point number. 

  

- ``rnd`` -- (string) the rounding mode; default = ``'RNDNN'`` 

  

Rounding mode is of the form ``'RNDxy'`` where ``x`` and ``y`` are 

the rounding mode for respectively the real and imaginary parts and 

are one of: 

  

- ``'N'`` for rounding to nearest 

- ``'Z'`` for rounding towards zero 

- ``'U'`` for rounding towards plus infinity 

- ``'D'`` for rounding towards minus infinity 

  

For example, ``'RNDZU'`` indicates to round the real part towards 

zero, and the imaginary part towards plus infinity. 

  

EXAMPLES:: 

  

sage: MPComplexField(17) 

Complex Field with 17 bits of precision 

sage: MPComplexField() 

Complex Field with 53 bits of precision 

sage: MPComplexField(1042,'RNDDZ') 

Complex Field with 1042 bits of precision and rounding RNDDZ 

  

ALGORITHMS: Computations are done using the MPC library. 

  

TESTS:: 

  

sage: TestSuite(MPComplexField(17)).run() 

""" 

if prec < mpfr_prec_min() or prec > mpfr_prec_max(): 

raise ValueError("prec (=%s) must be >= %s and <= %s." % ( 

prec, mpfr_prec_min(), mpfr_prec_max())) 

self.__prec = prec 

if not isinstance(rnd, str): 

raise TypeError("rnd must be a string") 

try: 

n = _mpc_rounding_modes.index(rnd) 

except ValueError: 

raise ValueError("rnd (=%s) must be of the form RNDxy"\ 

"where x and y are one of N, Z, U, D" % rnd) 

self.__rnd = n 

self.__rnd_str = rnd 

  

self.__real_field = real_mpfr.RealField(prec, rnd=_mpfr_rounding_modes[rnd_re(n)]) 

self.__imag_field = real_mpfr.RealField(prec, rnd=_mpfr_rounding_modes[rnd_im(n)]) 

  

ParentWithGens.__init__(self, self._real_field(), ('I',), False, category=Fields().Infinite()) 

self._populate_coercion_lists_(coerce_list=[MPFRtoMPC(self._real_field(), self)]) 

  

cdef MPComplexNumber _new(self): 

""" 

Return a new complex number with parent ``self`. 

""" 

cdef MPComplexNumber z 

z = MPComplexNumber.__new__(MPComplexNumber) 

z._parent = self 

mpc_init2(z.value, self.__prec) 

z.init = 1 

return z 

  

def _repr_ (self): 

""" 

Return a string representation of ``self``. 

  

EXAMPLES:: 

  

sage: MPComplexField(200, 'RNDDU') # indirect doctest 

Complex Field with 200 bits of precision and rounding RNDDU 

""" 

s = "Complex Field with %s bits of precision"%self.__prec 

if self.__rnd != MPC_RNDNN: 

s = s + " and rounding %s"%(self.__rnd_str) 

return s 

  

def _latex_(self): 

r""" 

Return a latex representation of ``self``. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(10) 

sage: latex(MPC) # indirect doctest 

\C 

""" 

return "\\C" 

  

def __call__(self, x, im=None): 

""" 

Create a floating-point complex using ``x`` and optionally an imaginary 

part ``im``. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: MPC(2) # indirect doctest 

2.00000000000000 

sage: MPC(0, 1) # indirect doctest 

1.00000000000000*I 

sage: MPC(1, 1) 

1.00000000000000 + 1.00000000000000*I 

sage: MPC(2, 3) 

2.00000000000000 + 3.00000000000000*I 

""" 

if x is None: 

return self.zero() 

# We implement __call__ to gracefully accept the second argument. 

if im is not None: 

x = x, im 

return Parent.__call__(self, x) 

  

def _element_constructor_(self, z): 

""" 

Coerce `z` into this complex field. 

  

EXAMPLES:: 

  

sage: C20 = MPComplexField(20) # indirect doctest 

  

The value can be set with a couple of reals:: 

  

sage: a = C20(1.5625, 17.42); a 

1.5625 + 17.420*I 

sage: a.str(2) 

'1.1001000000000000000 + 10001.011010111000011*I' 

sage: C20(0, 2) 

2.0000*I 

  

Complex number can be coerced into MPComplexNumber:: 

  

sage: C20(14.7+0.35*I) 

14.700 + 0.35000*I 

sage: C20(i*4, 7) 

Traceback (most recent call last): 

... 

TypeError: unable to coerce to a ComplexNumber: <type 'sage.symbolic.expression.Expression'> 

  

Each part can be set with strings (written in base ten):: 

  

sage: C20('1.234', '56.789') 

1.2340 + 56.789*I 

  

The string can represent the whole complex value:: 

  

sage: C20('42 + I * 100') 

42.000 + 100.00*I 

sage: C20('-42 * I') 

- 42.000*I 

  

The imaginary part can be written first:: 

  

sage: C20('100*i+42') 

42.000 + 100.00*I 

  

Use ``'inf'`` for infinity and ``'nan'`` for Not a Number:: 

  

sage: C20('nan+inf*i') 

NaN + +infinity*I 

""" 

cdef MPComplexNumber zz 

zz = self._new() 

zz._set(z) 

return zz 

  

cpdef _coerce_map_from_(self, S): 

""" 

Canonical coercion of `z` to this mpc complex field. 

  

The rings that canonically coerce to this mpc complex field are: 

  

- any mpc complex field with precision that is as large as this one 

- anything that canonically coerces to the mpfr real 

field with this prec and the rounding mode of real part. 

  

EXAMPLES:: 

  

sage: MPComplexField(100)(17, '4.2') + MPComplexField(20)('6.0', -23) # indirect doctest 

23.000 - 18.800*I 

sage: a = MPComplexField(100)(17, '4.2') + MPComplexField(20)('6.0', -23) 

sage: a.parent() 

Complex Field with 20 bits of precision 

""" 

if isinstance(S, RealField_class): 

return MPFRtoMPC(S, self) 

if isinstance(S, sage.rings.integer_ring.IntegerRing_class): 

return INTEGERtoMPC(S, self) 

  

RR = self.__real_field 

if RR.has_coerce_map_from(S): 

return self._coerce_map_via([RR], S) 

  

if isinstance(S, MPComplexField_class) and S.prec() >= self.__prec: 

#FIXME: What map when rounding modes differ but prec is the same ? 

# How to provide commutativity of morphisms ? 

# Change _cmp_ when done 

return MPCtoMPC(S, self) 

  

if isinstance(S, ComplexField_class) and S.prec() >= self.__prec: 

return CCtoMPC(S, self) 

  

late_import() 

if S in [AA, QQbar, CLF, RLF] or (S == CDF and self._prec <= 53): 

return self._generic_coerce_map(S) 

  

return self._coerce_map_via([CLF], S) 

  

def __reduce__(self): 

""" 

For pickling. 

  

EXAMPLES:: 

  

sage: C = MPComplexField(prec=200, rnd='RNDDZ') 

sage: loads(dumps(C)) == C 

True 

""" 

return __create__MPComplexField_version0, (self.__prec, self.__rnd_str) 

  

def __richcmp__(left, right, int op): 

""" 

Compare ``self`` and ``other``, ignoring the rounding mode. 

  

EXAMPLES:: 

  

sage: MPComplexField(10) == MPComplexField(11) # indirect doctest 

False 

sage: MPComplexField(10) == MPComplexField(10) 

True 

sage: MPComplexField(10,rnd='RNDZN') == MPComplexField(10,rnd='RNDZU') 

True 

""" 

if left is right: 

return rich_to_bool(op, 0) 

  

if not isinstance(right, MPComplexField_class): 

return op == Py_NE 

  

cdef MPComplexField_class s = <MPComplexField_class>left 

cdef MPComplexField_class o = <MPComplexField_class>right 

return richcmp(s.__prec, o.__prec, op) 

  

def gen(self, n=0): 

""" 

Return the generator of this complex field over its real subfield. 

  

EXAMPLES:: 

  

sage: MPComplexField(34).gen() 

1.00000000*I 

""" 

if n != 0: 

raise IndexError("n must be 0") 

return self(0, 1) 

  

def ngens(self): 

""" 

Return 1, the number of generators of this complex field over its real 

subfield. 

  

EXAMPLES:: 

  

sage: MPComplexField(34).ngens() 

1 

""" 

return 1 

  

cpdef _an_element_(self): 

""" 

Return an element of this complex field. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(20) 

sage: MPC._an_element_() 

1.0000*I 

""" 

return self(0, 1) 

  

def random_element(self, min=0, max=1): 

""" 

Return a random complex number, uniformly distributed with 

real and imaginary parts between min and max (default 0 to 1). 

  

EXAMPLES:: 

  

sage: MPComplexField(100).random_element(-5, 10) # random 

1.9305310520925994224072377281 + 0.94745292506956219710477444855*I 

sage: MPComplexField(10).random_element() # random 

0.12 + 0.23*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

cdef randstate rstate = current_randstate() 

mpc_urandom(z.value, rstate.gmp_state) 

if min == 0 and max == 1: 

return z 

else: 

return (max-min)*z + min*self(1,1) 

  

cpdef bint is_exact(self) except -2: 

""" 

Returns whether or not this field is exact, which is always ``False``. 

  

EXAMPLES:: 

  

sage: MPComplexField(42).is_exact() 

False 

""" 

return False 

  

def is_finite(self): 

""" 

Return ``False``, since the field of complex numbers is not finite. 

  

EXAMPLES:: 

  

sage: MPComplexField(17).is_finite() 

False 

""" 

return False 

  

def characteristic(self): 

""" 

Return 0, since the field of complex numbers has characteristic 0. 

  

EXAMPLES:: 

  

sage: MPComplexField(42).characteristic() 

0 

""" 

return Integer(0) 

  

def name(self): 

""" 

Return the name of the complex field. 

  

EXAMPLES:: 

  

sage: C = MPComplexField(10, 'RNDNZ'); C.name() 

'MPComplexField10_RNDNZ' 

""" 

return "MPComplexField%s_%s"%(self.__prec, self.__rnd_str) 

  

def __hash__(self): 

""" 

Return the hash of ``self``. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: hash(MPC) % 2^32 == hash(MPC.name()) % 2^32 

True 

""" 

return hash(self.name()) 

  

def prec(self): 

""" 

Return the precision of this field of complex numbers. 

  

EXAMPLES:: 

  

sage: MPComplexField().prec() 

53 

sage: MPComplexField(22).prec() 

22 

""" 

return self.__prec 

  

def rounding_mode(self): 

""" 

Return rounding modes used for each part of a complex number. 

  

EXAMPLES:: 

  

sage: MPComplexField().rounding_mode() 

'RNDNN' 

sage: MPComplexField(rnd='RNDZU').rounding_mode() 

'RNDZU' 

""" 

return self.__rnd_str 

  

def rounding_mode_real(self): 

""" 

Return rounding mode used for the real part of complex number. 

  

EXAMPLES:: 

  

sage: MPComplexField(rnd='RNDZU').rounding_mode_real() 

'RNDZ' 

""" 

return _mpfr_rounding_modes[rnd_re(self.__rnd)] 

  

def rounding_mode_imag(self): 

""" 

Return rounding mode used for the imaginary part of complex number. 

  

EXAMPLES:: 

  

sage: MPComplexField(rnd='RNDZU').rounding_mode_imag() 

'RNDU' 

""" 

return _mpfr_rounding_modes[rnd_im(self.__rnd)] 

  

def _real_field(self): 

""" 

Return real field for the real part. 

  

EXAMPLES:: 

  

sage: MPComplexField()._real_field() 

Real Field with 53 bits of precision 

""" 

return self.__real_field 

  

def _imag_field(self): 

""" 

Return real field for the imaginary part. 

  

EXAMPLES:: 

  

sage: MPComplexField(prec=100)._imag_field() 

Real Field with 100 bits of precision 

""" 

return self.__imag_field 

  

  

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

# 

# MPComplex Number -- element of MPComplex Field 

# 

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

  

cdef class MPComplexNumber(sage.structure.element.FieldElement): 

""" 

A floating point approximation to a complex number using any specified 

precision common to both real and imaginary part. 

""" 

cdef MPComplexNumber _new(self): 

""" 

Return a new complex number with same parent as ``self``. 

""" 

cdef MPComplexNumber z 

z = MPComplexNumber.__new__(MPComplexNumber) 

z._parent = self._parent 

mpc_init2(z.value, (<MPComplexField_class>self._parent).__prec) 

z.init = 1 

return z 

  

def __init__(self, MPComplexField_class parent, x, y=None, int base=10): 

""" 

Create a complex number. 

  

INPUT: 

  

- ``x`` -- real part or the complex value in a string 

  

- ``y`` -- imaginary part 

  

- ``base`` -- when ``x`` or ``y`` is a string, base in which the 

number is written 

  

A :class:`MPComplexNumber` should be called by first creating a 

:class:`MPComplexField`, as illustrated in the examples. 

  

EXAMPLES:: 

  

sage: C200 = MPComplexField(200) 

sage: C200(1/3, '0.6789') 

0.33333333333333333333333333333333333333333333333333333333333 + 0.67890000000000000000000000000000000000000000000000000000000*I 

sage: C3 = MPComplexField(3) 

sage: C3('1.2345', '0.6789') 

1.2 + 0.62*I 

sage: C3(3.14159) 

3.0 

  

Rounding modes:: 

  

sage: w = C3(5/2, 7/2); w.str(2) 

'10.1 + 11.1*I' 

sage: MPComplexField(2, rnd="RNDZN")(w).str(2) 

'10. + 100.*I' 

sage: MPComplexField(2, rnd="RNDDU")(w).str(2) 

'10. + 100.*I' 

sage: MPComplexField(2, rnd="RNDUD")(w).str(2) 

'11. + 11.*I' 

sage: MPComplexField(2, rnd="RNDNZ")(w).str(2) 

'10. + 11.*I' 

  

TESTS:: 

  

sage: MPComplexField(42)._repr_option('element_is_atomic') 

False 

""" 

self.init = 0 

if parent is None: 

raise TypeError 

self._parent = parent 

mpc_init2(self.value, parent.__prec) 

self.init = 1 

if x is None: return 

self._set(x, y, base) 

  

def _set(self, z, y=None, base=10): 

""" 

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: r = RealField(100).pi() 

sage: z = MPC(r); z # indirect doctest 

3.1415926535897932384626433833 

sage: MPComplexField(10, rnd='RNDDD')(z) 

3.1 

sage: c = ComplexField(53)(1, r) 

sage: MPC(c) 

1.0000000000000000000000000000 + 3.1415926535897931159979634685*I 

sage: MPC(I) 

1.0000000000000000000000000000*I 

sage: MPC('-0 +i') 

1.0000000000000000000000000000*I 

sage: MPC(1+i) 

1.0000000000000000000000000000 + 1.0000000000000000000000000000*I 

sage: MPC(1/3) 

0.33333333333333333333333333333 

  

:: 

  

sage: MPC(1, r/3) 

1.0000000000000000000000000000 + 1.0471975511965977461542144611*I 

sage: MPC(3, 2) 

3.0000000000000000000000000000 + 2.0000000000000000000000000000*I 

sage: MPC(0, r) 

3.1415926535897932384626433833*I 

sage: MPC('0.625e-26', '0.0000001') 

6.2500000000000000000000000000e-27 + 1.0000000000000000000000000000e-7*I 

  

Conversion from gmpy2 numbers:: 

  

sage: from gmpy2 import * # optional - gmpy2 

sage: MPC(mpc(int(2),int(1))) # optional - gmpy2 

2.0000000000000000000000000000 + 1.0000000000000000000000000000*I 

sage: MPC(mpfr(2.5)) # optional - gmpy2 

2.5000000000000000000000000000 

sage: MPC(mpq('3/2')) # optional - gmpy2 

1.5000000000000000000000000000 

sage: MPC(mpz(int(5))) # optional - gmpy2 

5.0000000000000000000000000000 

sage: re = mpfr('2.5') # optional - gmpy2 

sage: im = mpz(int(2)) # optional - gmpy2 

sage: MPC(re, im) # optional - gmpy2 

2.5000000000000000000000000000 + 2.0000000000000000000000000000*I 

""" 

# This should not be called except when the number is being created. 

# Complex Numbers are supposed to be immutable. 

cdef RealNumber x 

cdef mpc_rnd_t rnd 

rnd =(<MPComplexField_class>self._parent).__rnd 

if y is None: 

if z is None: return 

if isinstance(z, MPComplexNumber): 

mpc_set(self.value, (<MPComplexNumber>z).value, rnd) 

return 

elif isinstance(z, str): 

a, b = split_complex_string(z, base) 

# set real part 

if a is None: 

mpfr_set_ui(self.value.re, 0, MPFR_RNDN) 

else: 

mpfr_set_str(self.value.re, a, base, rnd_re(rnd)) 

# set imag part 

if b is None: 

if a is None: 

raise TypeError("unable to convert {!r} to a MPComplexNumber".format(z)) 

else: 

mpfr_set_str(self.value.im, b, base, rnd_im(rnd)) 

return 

elif isinstance(z, ComplexNumber): 

mpc_set_fr_fr(self.value, (<ComplexNumber>z).__re, (<ComplexNumber>z).__im, rnd) 

return 

elif isinstance(z, sage.libs.pari.all.pari_gen): 

real, imag = z.real(), z.imag() 

elif isinstance(z, list) or isinstance(z, tuple): 

real, imag = z 

elif isinstance(z, complex): 

real, imag = z.real, z.imag 

elif isinstance(z, sage.symbolic.expression.Expression): 

zz = sage.rings.complex_field.ComplexField(self._parent.prec())(z) 

self._set(zz) 

return 

elif HAVE_GMPY2 and type(z) is gmpy2.mpc: 

mpc_set(self.value, (<gmpy2.mpc>z).c, rnd) 

return 

# then, no imaginary part 

elif HAVE_GMPY2 and type(z) is gmpy2.mpfr: 

mpc_set_fr(self.value, (<gmpy2.mpfr>z).f, rnd) 

return 

elif HAVE_GMPY2 and type(z) is gmpy2.mpq: 

mpc_set_q(self.value, (<gmpy2.mpq>z).q, rnd) 

return 

elif HAVE_GMPY2 and type(z) is gmpy2.mpz: 

mpc_set_z(self.value, (<gmpy2.mpz>z).z, rnd) 

return 

elif isinstance(z, RealNumber): 

zz = sage.rings.real_mpfr.RealField(self._parent.prec())(z) 

mpc_set_fr(self.value, (<RealNumber>zz).value, rnd) 

return 

elif isinstance(z, Integer): 

mpc_set_z(self.value, (<Integer>z).value, rnd) 

return 

elif isinstance(z, (int, long)): 

mpc_set_si(self.value, z, rnd) 

return 

else: 

real = z 

imag = 0 

else: 

real = z 

imag = y 

  

cdef RealField_class R = self._parent._real_field() 

try: 

rr = R(real) 

ii = R(imag) 

mpc_set_fr_fr(self.value, (<RealNumber>rr).value, (<RealNumber>ii).value, rnd) 

  

except TypeError: 

raise TypeError("unable to coerce to a ComplexNumber: %s" % type(real)) 

  

def __reduce__(self): 

""" 

For pickling. 

  

EXAMPLES:: 

  

sage: C = MPComplexField(prec=200, rnd='RNDUU') 

sage: b = C(393.39203845902384098234098230948209384028340) 

sage: loads(dumps(b)) == b; 

True 

sage: C(1) 

1.0000000000000000000000000000000000000000000000000000000000 

sage: b = C(1)/C(0); b 

NaN + NaN*I 

sage: loads(dumps(b)) == b 

True 

sage: b = C(-1)/C(0.); b 

NaN + NaN*I 

sage: loads(dumps(b)) == b 

True 

sage: b = C(-1).sqrt(); b 

1.0000000000000000000000000000000000000000000000000000000000*I 

sage: loads(dumps(b)) == b 

True 

""" 

s = self.str(32) 

return (__create_MPComplexNumber_version0, (self._parent, s, 32)) 

  

def __dealloc__(self): 

if self.init: 

mpc_clear(self.value) 

  

def _repr_(self): 

""" 

Return a string representation of ``self``. 

  

EXAMPLES:: 

  

sage: MPComplexField()(2, -3) # indirect doctest 

2.00000000000000 - 3.00000000000000*I 

""" 

return self.str(truncate=True) 

  

def _latex_(self): 

""" 

Return a latex representation of ``self``. 

  

EXAMPLES:: 

  

sage: latex(MPComplexField()(2, -3)) # indirect doctest 

2.00000000000000 - 3.00000000000000i 

""" 

import re 

s = repr(self).replace('*I', 'i') 

return re.sub(r"e(-?\d+)", r" \\times 10^{\1}", s) 

  

def __hash__(self): 

""" 

Returns the hash of ``self``, which coincides with the python 

complex and float (and often int) types. 

  

This has the drawback that two very close high precision 

numbers will have the same hash, but allows them to play 

nicely with other real types. 

  

EXAMPLES:: 

  

sage: hash(MPComplexField()('1.2', 33)) == hash(complex(1.2, 33)) 

True 

""" 

return hash(complex(self)) 

  

def __getitem__(self, i): 

r""" 

Returns either the real or imaginary component of ``self`` 

depending on the choice of ``i``: real (``i``=0), imaginary (``i``=1). 

  

INPUT: 

  

- ``i`` -- 0 or 1 

  

- ``0`` -- will return the real component of ``self`` 

- ``1`` -- will return the imaginary component of ``self`` 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: a.__getitem__(0) 

2.00000000000000 

sage: a.__getitem__(1) 

1.00000000000000 

  

:: 

  

sage: b = MPC(42,0) 

sage: b 

42.0000000000000 

sage: b.__getitem__(1) 

0.000000000000000 

""" 

if i == 0: 

return self.real() 

elif i == 1: 

return self.imag() 

raise IndexError("i must be between 0 and 1.") 

  

def prec(self): 

""" 

Return precision of this complex number. 

  

EXAMPLES:: 

  

sage: i = MPComplexField(2000).0 

sage: i.prec() 

2000 

""" 

return <MPComplexField_class>(self._parent).__prec 

  

def real(self): 

""" 

Return the real part of ``self``. 

  

EXAMPLES:: 

  

sage: C = MPComplexField(100) 

sage: z = C(2, 3) 

sage: x = z.real(); x 

2.0000000000000000000000000000 

sage: x.parent() 

Real Field with 100 bits of precision 

""" 

cdef RealNumber x 

x = RealNumber(self._parent._real_field()) 

mpfr_set (x.value, self.value.re, (<RealField_class>x._parent).rnd) 

return x 

  

def imag(self): 

""" 

Return imaginary part of ``self``. 

  

EXAMPLES:: 

  

sage: C = MPComplexField(100) 

sage: z = C(2, 3) 

sage: x = z.imag(); x 

3.0000000000000000000000000000 

sage: x.parent() 

Real Field with 100 bits of precision 

""" 

cdef RealNumber y 

y = RealNumber(self._parent._imag_field()) 

mpfr_set (y.value, self.value.im, (<RealField_class>y._parent).rnd) 

return y 

  

def str(self, base=10, **kwds): 

""" 

Return a string of ``self``. 

  

INPUT: 

  

- ``base`` -- (default: 10) base for output 

  

- ``**kwds`` -- other arguments to pass to the ``str()`` 

method of the real numbers in the real and imaginary parts. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(64) 

sage: z = MPC(-4, 3)/7 

sage: z.str() 

'-0.571428571428571428564 + 0.428571428571428571436*I' 

sage: z.str(16) 

'-0.92492492492492490 + 0.6db6db6db6db6db70*I' 

sage: z.str(truncate=True) 

'-0.571428571428571429 + 0.428571428571428571*I' 

sage: z.str(2) 

'-0.1001001001001001001001001001001001001001001001001001001001001001 + 0.01101101101101101101101101101101101101101101101101101101101101110*I' 

""" 

s = "" 

if self.real() != 0: 

s = self.real().str(base, **kwds) 

if self.imag() != 0: 

if mpfr_signbit(self.value.im): 

s += ' - ' + (-self.imag()).str(base, **kwds) + '*I' 

else: 

if s: 

s += ' + ' 

s += self.imag().str(base, **kwds) + '*I' 

if not s: 

return "0" 

return s 

  

def __copy__(self): 

""" 

Return copy of ``self``. 

  

Since ``self`` is immutable, we just return ``self`` again. 

  

EXAMPLES:: 

  

sage: a = MPComplexField()(3.5, 3) 

sage: copy(a) is a 

True 

""" 

return self # since object is immutable. 

  

def __int__(self): 

r""" 

Method for converting ``self`` to type ``int``. 

  

Called by the ``int`` function. Note that calling this method returns 

an error since, in general, complex numbers cannot be coerced into 

integers. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: int(a) 

Traceback (most recent call last): 

... 

TypeError: can't convert complex to int; use int(abs(z)) 

sage: a.__int__() 

Traceback (most recent call last): 

... 

TypeError: can't convert complex to int; use int(abs(z)) 

""" 

raise TypeError("can't convert complex to int; use int(abs(z))") 

  

def __long__(self): 

r""" 

Method for converting ``self`` to type ``long``. 

  

Called by the ``long`` function. Note that calling this method 

returns an error since, in general, complex numbers cannot be 

coerced into integers. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: long(a) 

Traceback (most recent call last): 

... 

TypeError: can't convert complex to long; use long(abs(z)) 

sage: a.__long__() 

Traceback (most recent call last): 

... 

TypeError: can't convert complex to long; use long(abs(z)) 

""" 

raise TypeError("can't convert complex to long; use long(abs(z))") 

  

def __float__(self): 

r""" 

Method for converting ``self`` to type ``float``. 

  

Called by the ``float`` function. Note that calling this method returns 

an error since if the imaginary part of the number is not zero. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(1, 0) 

sage: float(a) 

1.0 

sage: a = MPC(2,1) 

sage: float(a) 

Traceback (most recent call last): 

... 

TypeError: can't convert complex to float; use abs(z) 

sage: a.__float__() 

Traceback (most recent call last): 

... 

TypeError: can't convert complex to float; use abs(z) 

""" 

if mpfr_zero_p(self.value.im): 

return mpfr_get_d(self.value.re,\ 

rnd_re((<MPComplexField_class>self._parent).__rnd)) 

else: 

raise TypeError("can't convert complex to float; use abs(z)") 

  

def __complex__(self): 

r""" 

Method for converting ``self`` to type ``complex``. 

  

Called by the ``complex`` function. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: complex(a) 

(2+1j) 

sage: type(complex(a)) 

<... 'complex'> 

sage: a.__complex__() 

(2+1j) 

""" 

# Fixme: is it the right choice for rounding modes ? 

cdef mpc_rnd_t rnd 

rnd = (<MPComplexField_class>self._parent).__rnd 

return complex(mpfr_get_d(self.value.re, rnd_re(rnd)), mpfr_get_d(self.value.im, rnd_im(rnd))) 

  

def __pari__(self): 

r""" 

Convert ``self`` to a PARI object. 

  

OUTPUT: a PARI ``t_COMPLEX`` object if the input is not purely 

real. If the input is real, a ``t_REAL`` is returned. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: a.__pari__() 

2.00000000000000 + 1.00000000000000*I 

sage: pari(a) 

2.00000000000000 + 1.00000000000000*I 

sage: pari(a).type() 

't_COMPLEX' 

sage: a = MPC(pi) 

sage: pari(a) 

3.14159265358979 

sage: pari(a).type() 

't_REAL' 

sage: a = MPC(-2).sqrt() 

sage: pari(a) 

1.41421356237310*I 

  

The precision is preserved, rounded up to the wordsize:: 

  

sage: MPC = MPComplexField(250) 

sage: MPC(1,2).__pari__().bitprecision() 

256 

sage: MPC(pi).__pari__().bitprecision() 

256 

""" 

if mpfr_zero_p(self.value.re): 

re = pari.PARI_ZERO 

else: 

re = self.real().__pari__() 

if mpfr_zero_p(self.value.im): 

return re 

im = self.imag().__pari__() 

return pari.complex(re, im) 

  

def __mpc__(self): 

""" 

Convert Sage ``MPComplexNumber`` to gmpy2 ``mpc``. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: c = MPC(2,1) 

sage: c.__mpc__() # optional - gmpy2 

mpc('2.0+1.0j') 

sage: from gmpy2 import mpc # optional - gmpy2 

sage: mpc(c) # optional - gmpy2 

mpc('2.0+1.0j') 

sage: MPCF = MPComplexField(42) 

sage: mpc(MPCF(12, 12)).precision # optional - gmpy2 

(42, 42) 

sage: MPCF = MPComplexField(236) 

sage: mpc(MPCF(12, 12)).precision # optional - gmpy2 

(236, 236) 

sage: MPCF = MPComplexField(63) 

sage: x = MPCF('15.64E+128', '15.64E+128') 

sage: y = mpc(x) # optional - gmpy2 

sage: y.precision # optional - gmpy2 

(63, 63) 

sage: MPCF(y) == x # optional - gmpy2 

True 

sage: x = mpc('1.324+4e50j', precision=(70,70)) # optional - gmpy2 

sage: y = MPComplexField(70)(x) # optional - gmpy2 

sage: mpc(y) == x # optional - gmpy2 

True 

  

TESTS:: 

  

sage: c.__mpc__(); raise NotImplementedError("gmpy2 is not installed") 

Traceback (most recent call last): 

... 

NotImplementedError: gmpy2 is not installed 

""" 

IF HAVE_GMPY2: 

return gmpy2.GMPy_MPC_From_mpfr(self.value.re, self.value.im) 

ELSE: 

raise NotImplementedError("gmpy2 is not installed") 

  

cpdef int _cmp_(self, other) except -2: 

r""" 

Compare ``self`` to ``other``. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: b = MPC(1,2) 

sage: a < b 

False 

sage: a > b 

True 

""" 

cdef MPComplexNumber z = <MPComplexNumber>other 

# NaN should compare to nothing 

if mpfr_nan_p (self.value.re) or mpfr_nan_p (self.value.im) or mpfr_nan_p (z.value.re) or mpfr_nan_p (z.value.im): 

return False 

cdef int c = mpc_cmp(self.value, z.value) 

cdef int cre = MPC_INEX_RE(c) 

cdef int cim 

if cre: 

if cre>0: return 1 

else: return -1 

else: 

cim = MPC_INEX_IM(c) 

if cim>0: return 1 

elif cim<0: return -1 

else: return 0 

  

def __nonzero__(self): 

""" 

Return ``True`` if ``self`` is not zero. 

  

This is an internal function; use :meth:`is_zero()` instead. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: z = 1 + MPC(I) 

sage: z.is_zero() 

False 

""" 

return not (mpfr_zero_p(self.value.re) and mpfr_zero_p(self.value.im)) 

  

def is_square(self): 

r""" 

This function always returns true as `\CC` is algebraically closed. 

  

EXAMPLES:: 

  

sage: C200 = MPComplexField(200) 

sage: a = C200(2,1) 

sage: a.is_square() 

True 

  

`\CC` is algebraically closed, hence every element is a square:: 

  

sage: b = C200(5) 

sage: b.is_square() 

True 

""" 

return True 

  

def is_real(self): 

""" 

Return ``True`` if ``self`` is real, i.e. has imaginary part zero. 

  

EXAMPLES:: 

  

sage: C200 = MPComplexField(200) 

sage: C200(1.23).is_real() 

True 

sage: C200(1+i).is_real() 

False 

""" 

return (mpfr_zero_p(self.value.im) != 0) 

  

def is_imaginary(self): 

""" 

Return ``True`` if ``self`` is imaginary, i.e. has real part zero. 

  

EXAMPLES:: 

  

sage: C200 = MPComplexField(200) 

sage: C200(1.23*i).is_imaginary() 

True 

sage: C200(1+i).is_imaginary() 

False 

""" 

return (mpfr_zero_p(self.value.re) != 0) 

  

def algebraic_dependency(self, n, **kwds): 

""" 

Return an irreducible polynomial of degree at most `n` which is 

approximately satisfied by this complex number. 

  

ALGORITHM: Uses the PARI C-library algdep command. 

  

INPUT: Type ``algdep?`` at the top level prompt. All additional 

parameters are passed onto the top-level algdep command. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: z = (1/2)*(1 + sqrt(3.0) * MPC.0); z 

0.500000000000000 + 0.866025403784439*I 

sage: p = z.algebraic_dependency(5) 

sage: p 

x^2 - x + 1 

sage: p(z) 

1.11022302462516e-16 

  

TESTS:: 

  

sage: z.algebraic_dependancy(2) 

doctest:...: DeprecationWarning: algebraic_dependancy is deprecated. Please use algebraic_dependency instead. 

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

x^2 - x + 1 

""" 

from sage.arith.all import algdep 

return algdep(self, n, **kwds) 

  

# Former misspelling 

algebraic_dependancy = deprecated_function_alias(22714, algebraic_dependency) 

  

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

# Basic Arithmetic 

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

  

cpdef _add_(self, right): 

""" 

Add two complex numbers with the same parent. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(30) 

sage: MPC(-1.5, 2) + MPC(0.2, 1) # indirect doctest 

-1.3000000 + 3.0000000*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_add(z.value, self.value, (<MPComplexNumber>right).value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

cpdef _sub_(self, right): 

""" 

Subtract two complex numbers with the same parent. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(30) 

sage: MPC(-1.5, 2) - MPC(0.2, 1) # indirect doctest 

-1.7000000 + 1.0000000*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_sub(z.value, self.value, (<MPComplexNumber>right).value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

cpdef _mul_(self, right): 

""" 

Multiply two complex numbers with the same parent. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(30) 

sage: MPC(-1.5, 2) * MPC(0.2, 1) # indirect doctest 

-2.3000000 - 1.1000000*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_mul(z.value, self.value, (<MPComplexNumber>right).value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

cpdef _div_(self, right): 

""" 

Divide two complex numbers with the same parent. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(30) 

sage: MPC(-1.5, 2) / MPC(0.2, 1) # indirect doctest 

1.6346154 + 1.8269231*I 

sage: MPC(-1, 1) / MPC(0) 

NaN + NaN*I 

""" 

cdef MPComplexNumber z, x 

z = self._new() 

x = <MPComplexNumber>right 

if not x.is_zero(): 

mpc_div(z.value, self.value, x.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

cpdef _neg_(self): 

""" 

Return the negative of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(30) 

sage: - MPC(-1.5, 2) # indirect doctest 

1.5000000 - 2.0000000*I 

sage: - MPC(0) 

0 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_neg(z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def __invert__(self): 

""" 

Return the multiplicative inverse. 

  

EXAMPLES:: 

  

sage: C = MPComplexField() 

sage: a = ~C(5, 1) 

sage: a * C(5, 1) 

1.00000000000000 

""" 

cdef MPComplexNumber z 

z = self._new() 

if not self.is_zero(): 

mpc_ui_div (z.value, 1, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def __neg__(self): 

r""" 

Return the negative of this complex number. 

  

-(a + ib) = -a -i b 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: -a 

-2.00000000000000 - 1.00000000000000*I 

sage: a.__neg__() 

-2.00000000000000 - 1.00000000000000*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_neg(z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def __abs__(self): 

r""" 

Absolute value or modulus of this complex number, 

rounded with the rounding mode of the real part: 

  

.. MATH:: 

  

|a + ib| = \sqrt(a^2 + b^2). 

  

OUTPUT: 

  

A floating-point number in the real field of the real part 

(same precision, same rounding mode). 

  

EXAMPLES: 

  

Note that the absolute value of a complex number with imaginary 

component equal to zero is the absolute value of the real component:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: abs(a) 

2.23606797749979 

sage: a.__abs__() 

2.23606797749979 

sage: float(sqrt(2^2 + 1^1)) 

2.23606797749979 

  

sage: b = MPC(42,0) 

sage: abs(b) 

42.0000000000000 

sage: b.__abs__() 

42.0000000000000 

sage: b 

42.0000000000000 

""" 

cdef RealNumber x 

x = RealNumber(self._parent._real_field()) 

mpc_abs (x.value, self.value, (<RealField_class>x._parent).rnd) 

return x 

  

def norm(self): 

r""" 

Return the norm of a complex number, rounded with the rounding 

mode of the real part. The norm is the square of the absolute 

value: 

  

.. MATH:: 

  

\mathrm{norm}(a + ib) = a^2 + b^2. 

  

OUTPUT: 

  

A floating-point number in the real field of the real part 

(same precision, same rounding mode). 

  

EXAMPLES: 

  

This indeed acts as the square function when the imaginary 

component of self is equal to zero:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2,1) 

sage: a.norm() 

5.00000000000000 

sage: b = MPC(4.2,0) 

sage: b.norm() 

17.6400000000000 

sage: b^2 

17.6400000000000 

""" 

cdef RealNumber x 

x = RealNumber(self._parent._real_field()) 

mpc_norm(x.value, self.value, (<RealField_class>x._parent).rnd) 

return x 

  

def __rdiv__(self, left): 

r""" 

Returns the quotient of ``left`` with ``self``, that is: ``left/self`` 

as a complex number. 

  

INPUT: 

  

- ``left`` -- a complex number 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(2, 2) 

sage: a.__rdiv__(MPC(1)) 

0.250000000000000 - 0.250000000000000*I 

sage: MPC(1)/a 

0.250000000000000 - 0.250000000000000*I 

""" 

return MPComplexNumber(self._parent, left)/self 

  

def __pow__(self, right, modulus): 

""" 

Compute ``self`` raised to the power of exponent, rounded in 

the direction specified by the parent of ``self``. 

  

.. TODO:: 

  

FIXME: Branch cut 

  

EXAMPLES:: 

  

sage: MPC.<i> = MPComplexField(20) 

sage: a = i^2; a 

-1.0000 

sage: a.parent() 

Complex Field with 20 bits of precision 

sage: a = (1+i)^i; a 

0.42883 + 0.15487*I 

sage: (1+i)^(1+i) 

0.27396 + 0.58370*I 

sage: a.parent() 

Complex Field with 20 bits of precision 

sage: i^i 

0.20788 

sage: (2+i)^(0.5) 

1.4553 + 0.34356*I 

""" 

cdef MPComplexNumber z, x, p 

x = <MPComplexNumber>self 

z = x._new() 

  

if isinstance(right, (int,long)): 

mpc_pow_si(z.value, x.value, right, (<MPComplexField_class>x._parent).__rnd) 

elif isinstance(right, Integer): 

mpc_pow_z(z.value, x.value, (<Integer>right).value, (<MPComplexField_class>x._parent).__rnd) 

else: 

try: 

p = (<MPComplexField_class>x._parent)(right) 

except Exception: 

raise ValueError 

mpc_pow(z.value, x.value, p.value, (<MPComplexField_class>x._parent).__rnd) 

  

return z 

  

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

# Trigonometric & hyperbolic functions 

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

  

def cos(self): 

r""" 

Return the cosine of this complex number: 

  

.. MATH:: 

  

\cos(a + ib) = \cos a \cosh b -i \sin a \sinh b. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: cos(u) 

-11.3642347064011 - 24.8146514856342*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_cos (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def sin(self): 

r""" 

Return the sine of this complex number: 

  

.. MATH:: 

  

\sin(a + ib) = \sin a \cosh b + i \cos x \sinh b. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: sin(u) 

24.8313058489464 - 11.3566127112182*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_sin (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def tan(self): 

r""" 

Return the tangent of this complex number: 

  

.. MATH:: 

  

\tan(a + ib) = (\sin 2a + i \sinh 2b)/(\cos 2a + \cosh 2b). 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(-2, 4) 

sage: tan(u) 

0.000507980623470039 + 1.00043851320205*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_tan (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def cosh(self): 

""" 

Return the hyperbolic cosine of this complex number: 

  

.. MATH:: 

  

\cosh(a + ib) = \cosh a \cos b + i \sinh a \sin b. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: cosh(u) 

-2.45913521391738 - 2.74481700679215*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_cosh (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def sinh(self): 

""" 

Return the hyperbolic sine of this complex number: 

  

.. MATH:: 

  

\sinh(a + ib) = \sinh a \cos b + i \cosh a \sin b. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: sinh(u) 

-2.37067416935200 - 2.84723908684883*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_sinh (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def tanh(self): 

r""" 

Return the hyperbolic tangent of this complex number: 

  

.. MATH:: 

  

\tanh(a + ib) = (\sinh 2a + i \sin 2b)/(\cosh 2a + \cos 2b). 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: tanh(u) 

1.00468231219024 + 0.0364233692474037*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_tanh (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def arccos(self): 

""" 

Return the arccosine of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: arccos(u) 

1.11692611683177 - 2.19857302792094*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_acos (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def arcsin(self): 

""" 

Return the arcsine of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: arcsin(u) 

0.453870209963122 + 2.19857302792094*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_asin (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def arctan(self): 

""" 

Return the arctangent of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(-2, 4) 

sage: arctan(u) 

-1.46704821357730 + 0.200586618131234*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_atan (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def arccosh(self): 

""" 

Return the hyperbolic arccos of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: arccosh(u) 

2.19857302792094 + 1.11692611683177*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_acosh (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def arcsinh(self): 

""" 

Return the hyperbolic arcsine of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: arcsinh(u) 

2.18358521656456 + 1.09692154883014*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_asinh (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def arctanh(self): 

""" 

Return the hyperbolic arctangent of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: arctanh(u) 

0.0964156202029962 + 1.37153510396169*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_atanh (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def coth(self): 

""" 

Return the hyperbolic cotangent of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).coth() 

0.86801414289592494863584920892 - 0.21762156185440268136513424361*I 

""" 

return ~(self.tanh()) 

  

def arccoth(self): 

""" 

Return the hyperbolic arccotangent of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).arccoth() 

0.40235947810852509365018983331 - 0.55357435889704525150853273009*I 

""" 

return (~self).arctanh() 

  

def csc(self): 

""" 

Return the cosecant of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).csc() 

0.62151801717042842123490780586 - 0.30393100162842645033448560451*I 

""" 

return ~(self.sin()) 

  

def csch(self): 

""" 

Return the hyperbolic cosecant of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).csch() 

0.30393100162842645033448560451 - 0.62151801717042842123490780586*I 

""" 

return ~(self.sinh()) 

  

def arccsch(self): 

""" 

Return the hyperbolic arcsine of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).arccsch() 

0.53063753095251782601650945811 - 0.45227844715119068206365839783*I 

""" 

return (~self).arcsinh() 

  

def sec(self): 

""" 

Return the secant of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).sec() 

0.49833703055518678521380589177 + 0.59108384172104504805039169297*I 

""" 

return ~(self.cos()) 

  

def sech(self): 

""" 

Return the hyperbolic secant of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).sech() 

0.49833703055518678521380589177 - 0.59108384172104504805039169297*I 

""" 

return ~(self.cosh()) 

  

def arcsech(self): 

""" 

Return the hyperbolic arcsecant of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(100) 

sage: MPC(1,1).arcsech() 

0.53063753095251782601650945811 - 1.1185178796437059371676632938*I 

""" 

return (~self).arccosh() 

  

def cotan(self): 

""" 

Return the cotangent of this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(53) 

sage: (1+MPC(I)).cotan() 

0.217621561854403 - 0.868014142895925*I 

sage: i = MPComplexField(200).0 

sage: (1+i).cotan() 

0.21762156185440268136513424360523807352075436916785404091068 - 0.86801414289592494863584920891627388827343874994609327121115*I 

sage: i = MPComplexField(220).0 

sage: (1+i).cotan() 

0.21762156185440268136513424360523807352075436916785404091068124239 - 0.86801414289592494863584920891627388827343874994609327121115071646*I 

""" 

return ~(self.tan()) 

  

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

# Other functions 

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

  

def argument(self): 

r""" 

The argument (angle) of the complex number, normalized so that 

`-\pi < \theta \leq \pi`. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: i = MPC.0 

sage: (i^2).argument() 

3.14159265358979 

sage: (1+i).argument() 

0.785398163397448 

sage: i.argument() 

1.57079632679490 

sage: (-i).argument() 

-1.57079632679490 

sage: (RR('-0.001') - i).argument() 

-1.57179632646156 

""" 

cdef RealNumber x 

x = RealNumber(self._parent._real_field()) 

mpc_arg(x.value, self.value, (<RealField_class>x._parent).rnd) 

return x 

  

def conjugate(self): 

r""" 

Return the complex conjugate of this complex number: 

  

.. MATH:: 

  

\mathrm{conjugate}(a + ib) = a - ib. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: i = MPC(0, 1) 

sage: (1+i).conjugate() 

1.00000000000000 - 1.00000000000000*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_conj(z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def sqr(self): 

""" 

Return the square of a complex number: 

  

.. MATH:: 

  

(a + ib)^2 = (a^2 - b^2) + 2iab. 

  

EXAMPLES:: 

  

sage: C = MPComplexField() 

sage: a = C(5, 1) 

sage: a.sqr() 

24.0000000000000 + 10.0000000000000*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_sqr (z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def sqrt(self): 

r""" 

Return the square root, taking the branch cut to be the negative real 

axis: 

  

.. MATH:: 

  

\sqrt z = \sqrt{|z|}(\cos(\arg(z)/2) + i \sin(\arg(z)/2)). 

  

EXAMPLES:: 

  

sage: C = MPComplexField() 

sage: a = C(24, 10) 

sage: a.sqrt() 

5.00000000000000 + 1.00000000000000*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_sqrt(z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def exp(self): 

""" 

Return the exponential of this complex number: 

  

.. MATH:: 

  

\exp(a + ib) = \exp(a) (\cos b + i \sin b). 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: exp(u) 

-4.82980938326939 - 5.59205609364098*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_exp(z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def log(self): 

r""" 

Return the logarithm of this complex number with the branch 

cut on the negative real axis: 

  

.. MATH:: 

  

\log(z) = \log |z| + i \arg(z). 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: log(u) 

1.49786613677700 + 1.10714871779409*I 

""" 

cdef MPComplexNumber z 

z = self._new() 

mpc_log(z.value, self.value, (<MPComplexField_class>self._parent).__rnd) 

return z 

  

def __lshift__(self, n): 

""" 

Fast multiplication by `2^n`. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: u<<2 

8.00000000000000 + 16.0000000000000*I 

sage: u<<(-1) 

1.00000000000000 + 2.00000000000000*I 

""" 

cdef MPComplexNumber z, x 

x = <MPComplexNumber>self 

z = x._new() 

mpc_mul_2si(z.value , x.value, n, (<MPComplexField_class>x._parent).__rnd) 

return z 

  

def __rshift__(self, n): 

""" 

Fast division by `2^n`. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(2, 4) 

sage: u>>2 

0.500000000000000 + 1.00000000000000*I 

sage: u>>(-1) 

4.00000000000000 + 8.00000000000000*I 

""" 

cdef MPComplexNumber z, x 

x = <MPComplexNumber>self 

z = x._new() 

mpc_div_2si(z.value , x.value, n, (<MPComplexField_class>x._parent).__rnd) 

return z 

  

def nth_root(self, n, all=False): 

""" 

The `n`-th root function. 

  

INPUT: 

  

- ``all`` - bool (default: ``False``); if ``True``, return a 

list of all `n`-th roots. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(27) 

sage: a.nth_root(3) 

3.00000000000000 

sage: a.nth_root(3, all=True) 

[3.00000000000000, -1.50000000000000 + 2.59807621135332*I, -1.50000000000000 - 2.59807621135332*I] 

""" 

if self.is_zero(): 

return [self] if all else self 

  

cdef RealField_class R = self._parent._real_field() 

cdef mpfr_rnd_t rrnd = R.rnd 

cdef mpc_rnd_t crnd = (<MPComplexField_class>(self._parent)).__rnd 

  

cdef RealNumber a,r 

a = self.argument()/n 

r = self.abs() 

mpfr_rootn_ui(r.value, r.value, n, rrnd) 

  

cdef MPComplexNumber z 

z = self._new() 

mpfr_sin_cos(z.value.im, z.value.re, a.value, rrnd) 

mpc_mul_fr(z.value, z.value, r.value, crnd) 

  

if not all: 

return z 

  

cdef RealNumber t, tt 

t = R.pi()*2/n 

tt= RealNumber(R) 

  

cdef list zlist = [z] 

for i in xrange(1,n): 

z = self._new() 

mpfr_mul_ui(tt.value, t.value, i, rrnd) 

mpfr_add(tt.value, tt.value, a.value, rrnd) 

mpfr_sin_cos(z.value.im, z.value.re, tt.value, rrnd) 

mpc_mul_fr(z.value, z.value, r.value, crnd) 

zlist.append(z) 

  

return zlist 

  

def dilog(self): 

r""" 

Return the complex dilogarithm of ``self``. 

  

The complex dilogarithm, or Spence's function, is defined by 

  

.. MATH:: 

  

Li_2(z) = - \int_0^z \frac{\log|1-\zeta|}{\zeta} d(\zeta) 

= \sum_{k=1}^\infty \frac{z^k}{k^2}. 

  

Note that the series definition can only be used for `|z| < 1`. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: a = MPC(1,0) 

sage: a.dilog() 

1.64493406684823 

sage: float(pi^2/6) 

1.6449340668482262 

  

:: 

  

sage: b = MPC(0,1) 

sage: b.dilog() 

-0.205616758356028 + 0.915965594177219*I 

  

:: 

  

sage: c = MPC(0,0) 

sage: c.dilog() 

0 

""" 

return self._parent(self.__pari__().dilog()) 

  

def eta(self, omit_frac=False): 

r""" 

Return the value of the Dedekind `\eta` function on ``self``, 

intelligently computed using `\mathbb{SL}(2,\ZZ)` transformations. 

  

The `\eta` function is 

  

.. MATH:: 

  

\eta(z) = e^{\pi i z / 12} \prod_{n=1}^{\infty}(1-e^{2\pi inz}) 

  

INPUT: 

  

- ``self`` - element of the upper half plane (if not, 

raises a ``ValueError``). 

  

- ``omit_frac`` - (bool, default: ``False``), if ``True``, 

omit the `e^{\pi i z / 12}` factor. 

  

OUTPUT: a complex number 

  

ALGORITHM: Uses the PARI C library. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: i = MPC.0 

sage: z = 1+i; z.eta() 

0.742048775836565 + 0.198831370229911*I 

""" 

try: 

return self._parent(self.__pari__().eta(not omit_frac)) 

except sage.libs.pari.all.PariError: 

raise ValueError("value must be in the upper half plane") 

  

def gamma(self): 

""" 

Return the Gamma function evaluated at this complex number. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField(30) 

sage: i = MPC.0 

sage: (1+i).gamma() 

0.49801567 - 0.15494983*I 

  

TESTS:: 

  

sage: MPC(0).gamma() 

Infinity 

  

:: 

  

sage: MPC(-1).gamma() 

Infinity 

""" 

try: 

return self._parent(self.__pari__().gamma()) 

except sage.libs.pari.all.PariError: 

from sage.rings.infinity import UnsignedInfinityRing 

return UnsignedInfinityRing.gen() 

  

def gamma_inc(self, t): 

""" 

Return the incomplete Gamma function evaluated at this complex number. 

  

EXAMPLES:: 

  

sage: C, i = MPComplexField(30).objgen() 

sage: (1+i).gamma_inc(2 + 3*i) # abs tol 2e-10 

0.0020969149 - 0.059981914*I 

sage: (1+i).gamma_inc(5) 

-0.0013781309 + 0.0065198200*I 

sage: C(2).gamma_inc(1 + i) 

0.70709210 - 0.42035364*I 

  

""" 

return self._parent(self.__pari__().incgam(t, precision=self.prec())) 

  

def zeta(self): 

""" 

Return the Riemann zeta function evaluated at this complex number. 

  

EXAMPLES:: 

  

sage: i = MPComplexField(30).gen() 

sage: z = 1 + i 

sage: z.zeta() 

0.58215806 - 0.92684856*I 

""" 

return self._parent(self.__pari__().zeta()) 

  

def agm(self, right, algorithm="optimal"): 

""" 

Return the algebro-geometric mean of ``self`` and ``right``. 

  

EXAMPLES:: 

  

sage: MPC = MPComplexField() 

sage: u = MPC(1, 4) 

sage: v = MPC(-2,5) 

sage: u.agm(v, algorithm="pari") 

-0.410522769709397 + 4.60061063922097*I 

sage: u.agm(v, algorithm="principal") 

1.24010691168158 - 0.472193567796433*I 

sage: u.agm(v, algorithm="optimal") 

-0.410522769709397 + 4.60061063922097*I 

""" 

if algorithm=="pari": 

t = self._parent(right).__pari__() 

return self._parent(self.__pari__().agm(t)) 

  

cdef MPComplexNumber a, b, d, s, res 

cdef mpfr_t sn,dn 

cdef mp_exp_t rel_prec 

cdef bint optimal = algorithm == "optimal" 

  

cdef mpc_rnd_t rnd = (<MPComplexField_class>(self._parent)).__rnd 

  

cdef int prec = self._parent.__prec 

  

if optimal or algorithm == "principal": 

if not isinstance(right, MPComplexNumber) or (<MPComplexNumber>right)._parent is not self._parent: 

right = self._parent(right) 

  

res = self._new() 

  

if self.is_zero(): 

return self 

elif (<MPComplexNumber>right).is_zero(): 

return right 

elif (mpfr_cmpabs(self.value.re, (<MPComplexNumber>right).value.re) == 0 and 

mpfr_cmpabs(self.value.im, (<MPComplexNumber>right).value.im) == 0 and 

mpfr_cmp(self.value.re, (<MPComplexNumber>right).value.re) != 0 and 

mpfr_cmp(self.value.im, (<MPComplexNumber>right).value.im) != 0): 

# self = -right 

mpc_set_ui(res.value, 0, rnd) 

return res 

# Do the computations to a bit higher precision so rounding error 

# won't obscure the termination condition. 

a = MPComplexNumber(MPComplexField(prec+5), None) 

b = a._new() 

d = a._new() 

s = a._new() 

  

# Make copies so we don't mutate self or right. 

mpc_set(a.value, self.value, rnd) 

mpc_set(b.value, (<MPComplexNumber>right).value, rnd) 

  

if optimal: 

mpc_add(s.value, a.value, b.value, rnd) 

  

while True: 

# s = a+b 

if not optimal: 

mpc_add(s.value, a.value, b.value, rnd) 

  

# b = sqrt(a*b) 

mpc_mul(d.value, a.value, b.value, rnd) 

mpc_sqrt(b.value, d.value, rnd) 

  

# a = s/2 

mpc_div_2ui(a.value, s.value, 1, rnd) 

  

# d = a - b 

mpc_sub(d.value, a.value, b.value, rnd) 

if d.is_zero(): 

mpc_set(res.value, a.value, rnd) 

return res 

  

if optimal: 

# s = a+b 

mpc_add(s.value, a.value, b.value, rnd) 

if s.is_zero(): 

mpc_set(res.value, a.value, rnd) 

return res 

  

# |s| < |d| 

if s.norm() < d.norm(): 

mpc_swap(d.value, s.value) 

mpc_neg(b.value, b.value, rnd) 

  

rel_prec = min_exp_t(max_exp(a), max_exp(b)) - max_exp(d) 

if rel_prec > prec: 

mpc_set(res.value, a.value, rnd) 

return res 

  

cdef inline mp_exp_t min_exp_t(mp_exp_t a, mp_exp_t b): 

return a if a < b else b 

cdef inline mp_exp_t max_exp_t(mp_exp_t a, mp_exp_t b): 

return a if a > b else b 

  

cdef inline mp_exp_t max_exp(MPComplexNumber z): 

""" 

Quickly return the maximum exponent of the real and complex parts of ``z``, 

which is useful for estimating its magnitude. 

""" 

if mpfr_zero_p(z.value.im): 

return mpfr_get_exp(z.value.re) 

elif mpfr_zero_p(z.value.re): 

return mpfr_get_exp(z.value.im) 

return max_exp_t(mpfr_get_exp(z.value.re), mpfr_get_exp(z.value.im)) 

  

def __create__MPComplexField_version0 (prec, rnd): 

""" 

Create a :class:`MPComplexField`. 

  

EXAMPLES:: 

  

sage: sage.rings.complex_mpc.__create__MPComplexField_version0(200, 'RNDDZ') 

Complex Field with 200 bits of precision and rounding RNDDZ 

""" 

return MPComplexField(prec, rnd) 

  

def __create__MPComplexNumber_version0 (parent, s, base=10): 

""" 

Create a :class:`MPComplexNumber`. 

  

EXAMPLES:: 

  

sage: C = MPComplexField(prec=20, rnd='RNDUU') 

sage: sage.rings.complex_mpc.__create__MPComplexNumber_version0(C, 3.2+2*i) 

3.2001 + 2.0000*I 

""" 

return MPComplexNumber(parent, s, base=base) 

  

# original version of the file had this with only 1 underscore - TCS 

__create_MPComplexNumber_version0 = __create__MPComplexNumber_version0 

  

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

# 

# Morphisms 

# 

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

  

cdef class MPCtoMPC(Map): 

cpdef Element _call_(self, z): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.complex_mpc import * 

sage: C10 = MPComplexField(10) 

sage: C100 = MPComplexField(100) 

sage: f = MPCtoMPC(C100, C10) # indirect doctest 

sage: a = C100(1.2, 24) 

sage: f(a) 

1.2 + 24.*I 

sage: f 

Generic map: 

From: Complex Field with 100 bits of precision 

To: Complex Field with 10 bits of precision 

""" 

cdef MPComplexNumber y 

y = (<MPComplexField_class>self.codomain())._new() 

y._set(z) 

return y 

  

def section(self): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.complex_mpc import * 

sage: C10 = MPComplexField(10) 

sage: C100 = MPComplexField(100) 

sage: f = MPCtoMPC(C100, C10) 

sage: f.section() 

Generic map: 

From: Complex Field with 10 bits of precision 

To: Complex Field with 100 bits of precision 

""" 

return MPCtoMPC(self.codomain(), self.domain()) 

  

cdef class INTEGERtoMPC(Map): 

cpdef Element _call_(self, x): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.complex_mpc import * 

sage: I = IntegerRing() 

sage: C100 = MPComplexField(100) 

sage: f = MPFRtoMPC(I, C100); f # indirect doctest 

Generic map: 

From: Integer Ring 

To: Complex Field with 100 bits of precision 

sage: a = I(625) 

sage: f(a) 

625.00000000000000000000000000 

""" 

cdef MPComplexNumber y 

cdef mpc_rnd_t rnd 

rnd =(<MPComplexField_class>self._parent).__rnd 

y = (<MPComplexField_class>self.codomain())._new() 

mpc_set_z(y.value, (<Integer>x).value, rnd) 

return y 

  

cdef class MPFRtoMPC(Map): 

cpdef Element _call_(self, x): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.complex_mpc import * 

sage: R10 = RealField(10) 

sage: C100 = MPComplexField(100) 

sage: f = MPFRtoMPC(R10, C100); f # indirect doctest 

Generic map: 

From: Real Field with 10 bits of precision 

To: Complex Field with 100 bits of precision 

sage: a = R10(1.625) 

sage: f(a) 

1.6250000000000000000000000000 

""" 

cdef MPComplexNumber y 

# cdef mpc_rnd_t rnd 

# rnd =(<MPComplexField_class>self._parent).__rnd 

y = (<MPComplexField_class>self.codomain())._new() 

# mpc_set_fr(y.value, (<RealNumber>x).value, rnd) 

y._set(x) 

return y 

  

cdef class CCtoMPC(Map): 

cpdef Element _call_(self, z): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.complex_mpc import * 

sage: C10 = ComplexField(10) 

sage: MPC100 = MPComplexField(100) 

sage: f = CCtoMPC(C10, MPC100); f # indirect doctest 

Generic map: 

From: Complex Field with 10 bits of precision 

To: Complex Field with 100 bits of precision 

sage: a = C10(1.625, 42) 

sage: f(a) 

1.6250000000000000000000000000 + 42.000000000000000000000000000*I 

""" 

cdef MPComplexNumber y 

cdef mpc_rnd_t rnd 

rnd =(<MPComplexField_class>self._parent).__rnd 

y = (<MPComplexField_class>self.codomain())._new() 

mpc_set_fr_fr(y.value, (<ComplexNumber>z).__re, (<ComplexNumber>z).__im, rnd) 

return y