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

2694

2695

2696

2697

2698

2699

2700

2701

2702

2703

2704

2705

2706

2707

2708

2709

2710

2711

2712

2713

2714

2715

2716

2717

2718

2719

2720

2721

2722

2723

2724

2725

2726

2727

2728

2729

2730

2731

2732

2733

2734

2735

2736

2737

2738

2739

2740

2741

2742

2743

2744

2745

2746

2747

2748

2749

2750

2751

2752

2753

2754

2755

2756

2757

2758

2759

2760

2761

2762

2763

2764

2765

2766

2767

2768

2769

2770

2771

2772

2773

2774

2775

2776

2777

2778

2779

2780

2781

2782

2783

2784

2785

2786

2787

2788

2789

2790

2791

2792

2793

2794

2795

2796

2797

2798

2799

2800

2801

2802

2803

2804

2805

2806

2807

2808

2809

2810

2811

2812

2813

2814

2815

2816

2817

2818

2819

2820

2821

2822

2823

2824

2825

2826

2827

2828

2829

2830

2831

2832

2833

2834

2835

2836

2837

2838

2839

2840

2841

2842

2843

2844

2845

2846

2847

2848

2849

2850

2851

2852

2853

2854

2855

2856

2857

2858

2859

2860

2861

2862

2863

2864

2865

2866

2867

2868

2869

2870

2871

2872

2873

2874

2875

2876

2877

2878

2879

2880

2881

2882

2883

2884

2885

2886

2887

2888

2889

2890

2891

2892

2893

2894

2895

2896

2897

2898

2899

2900

2901

2902

2903

2904

2905

2906

2907

2908

2909

2910

2911

2912

2913

2914

r""" 

Double Precision Real Numbers 

  

EXAMPLES: 

  

We create the real double vector space of dimension `3`:: 

  

sage: V = RDF^3; V 

Vector space of dimension 3 over Real Double Field 

  

Notice that this space is unique:: 

  

sage: V is RDF^3 

True 

sage: V is FreeModule(RDF, 3) 

True 

sage: V is VectorSpace(RDF, 3) 

True 

  

Also, you can instantly create a space of large dimension:: 

  

sage: V = RDF^10000 

  

TESTS: 

  

Test NumPy conversions:: 

  

sage: RDF(1).__array_interface__ 

{'typestr': '=f8'} 

sage: import numpy 

sage: numpy.array([RDF.pi()]).dtype 

dtype('float64') 

""" 

  

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

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

  

cimport libc.math 

from libc.string cimport memcpy 

from cpython.object cimport * 

from cpython.float cimport * 

  

from cysignals.signals cimport sig_on, sig_off 

  

from sage.ext.stdsage cimport PY_NEW 

from sage.cpython.python_debug cimport if_Py_TRACE_REFS_then_PyObject_INIT 

  

from sage.libs.gsl.all cimport * 

  

gsl_set_error_handler_off() 

  

import math, operator 

  

from cypari2.convert cimport new_gen_from_double 

  

import sage.rings.integer 

import sage.rings.rational 

  

from sage.rings.integer cimport Integer 

from sage.rings.integer_ring import ZZ 

  

from sage.categories.morphism cimport Morphism 

from sage.structure.coerce cimport is_numpy_type 

from sage.misc.randstate cimport randstate, current_randstate 

from sage.structure.richcmp cimport rich_to_bool 

from sage.arith.constants cimport * 

  

IF HAVE_GMPY2: 

cimport gmpy2 

  

  

def is_RealDoubleField(x): 

""" 

Returns ``True`` if ``x`` is the field of real double precision numbers. 

  

EXAMPLES:: 

  

sage: from sage.rings.real_double import is_RealDoubleField 

sage: is_RealDoubleField(RDF) 

True 

sage: is_RealDoubleField(RealField(53)) 

False 

""" 

return isinstance(x, RealDoubleField_class) 

  

cdef class RealDoubleField_class(Field): 

""" 

An approximation to the field of real numbers using double 

precision floating point numbers. Answers derived from calculations 

in this approximation may differ from what they would be if those 

calculations were performed in the true field of real numbers. This 

is due to the rounding errors inherent to finite precision 

calculations. 

  

EXAMPLES:: 

  

sage: RR == RDF 

False 

sage: RDF == RealDoubleField() # RDF is the shorthand 

True 

  

:: 

  

sage: RDF(1) 

1.0 

sage: RDF(2/3) 

0.6666666666666666 

  

A ``TypeError`` is raised if the coercion doesn't make sense:: 

  

sage: RDF(QQ['x'].0) 

Traceback (most recent call last): 

... 

TypeError: cannot coerce nonconstant polynomial to float 

sage: RDF(QQ['x'](3)) 

3.0 

  

One can convert back and forth between double precision real 

numbers and higher-precision ones, though of course there may be 

loss of precision:: 

  

sage: a = RealField(200)(2).sqrt(); a 

1.4142135623730950488016887242096980785696718753769480731767 

sage: b = RDF(a); b 

1.4142135623730951 

sage: a.parent()(b) 

1.4142135623730951454746218587388284504413604736328125000000 

sage: a.parent()(b) == b 

True 

sage: b == RR(a) 

True 

  

""" 

def __init__(self): 

""" 

Initialize ``self``. 

  

TESTS:: 

  

sage: R = RealDoubleField() 

sage: TestSuite(R).run() 

""" 

from sage.categories.fields import Fields 

Field.__init__(self, self, category=Fields().Metric().Complete()) 

self._populate_coercion_lists_(init_no_parent=True, 

convert_method_name='_real_double_') 

  

_element_constructor_ = RealDoubleElement 

  

def __reduce__(self): 

""" 

For pickling. 

  

EXAMPLES:: 

  

sage: loads(dumps(RDF)) is RDF 

True 

""" 

return RealDoubleField, () 

  

cpdef bint is_exact(self) except -2: 

""" 

Returns ``False``, because doubles are not exact. 

  

EXAMPLES:: 

  

sage: RDF.is_exact() 

False 

""" 

return False 

  

def _latex_(self): 

r""" 

Return a latex representation of ``self``. 

  

EXAMPLES:: 

  

sage: latex(RDF) # indirect doctest 

\Bold{R} 

""" 

return "\\Bold{R}" 

  

def _sage_input_(self, sib, coerced): 

r""" 

Produce an expression which will reproduce this value when evaluated. 

  

EXAMPLES:: 

  

sage: sage_input(RDF, verify=True) 

# Verified 

RDF 

sage: from sage.misc.sage_input import SageInputBuilder 

sage: RDF._sage_input_(SageInputBuilder(), False) 

{atomic:RDF} 

""" 

return sib.name('RDF') 

  

def __repr__(self): 

""" 

Return a string representation of ``self``. 

  

EXAMPLES:: 

  

sage: RealDoubleField() # indirect doctest 

Real Double Field 

sage: RDF 

Real Double Field 

""" 

return "Real Double Field" 

  

def _repr_option(self, key): 

""" 

Metadata about the :meth:`_repr_` output. 

  

See :meth:`sage.structure.parent._repr_option` for details. 

  

EXAMPLES:: 

  

sage: RDF._repr_option('element_is_atomic') 

True 

""" 

if key == 'element_is_atomic': 

return True 

return super(RealDoubleField_class, self)._repr_option(key) 

  

def __richcmp__(self, x, op): 

""" 

Compare ``self`` to ``x``. 

  

EXAMPLES:: 

  

sage: RDF == 5 

False 

sage: loads(dumps(RDF)) == RDF 

True 

""" 

if isinstance(x, RealDoubleField_class): 

return rich_to_bool(op, 0) 

if op == Py_NE: 

return True 

return NotImplemented 

  

def construction(self): 

r""" 

Returns the functorial construction of ``self``, namely, completion of 

the rational numbers with respect to the prime at `\infty`. 

  

Also preserves other information that makes this field unique (i.e. 

the Real Double Field). 

  

EXAMPLES:: 

  

sage: c, S = RDF.construction(); S 

Rational Field 

sage: RDF == c(S) 

True 

""" 

from sage.categories.pushout import CompletionFunctor 

return (CompletionFunctor(sage.rings.infinity.Infinity, 

53, 

{'type': 'RDF'}), 

sage.rings.rational_field.QQ) 

  

def complex_field(self): 

""" 

Return the complex field with the same precision as ``self``, i.e., 

the complex double field. 

  

EXAMPLES:: 

  

sage: RDF.complex_field() 

Complex Double Field 

""" 

from sage.rings.complex_double import CDF 

return CDF 

  

def algebraic_closure(self): 

""" 

Return the algebraic closure of ``self``, i.e., the complex double 

field. 

  

EXAMPLES:: 

  

sage: RDF.algebraic_closure() 

Complex Double Field 

""" 

from sage.rings.complex_double import CDF 

return CDF 

  

cpdef _coerce_map_from_(self, S): 

""" 

Canonical coercion of ``S`` to the real double field. 

  

The rings that canonically coerce to the real double field are: 

  

- the real double field itself 

- int, long, integer, and rational rings 

- numpy integers and floatings 

- the real lazy field 

- the MPFR real field with at least 53 bits of precision 

  

EXAMPLES:: 

  

sage: RDF.coerce(5) # indirect doctest 

5.0 

sage: RDF.coerce(9499294r) 

9499294.0 

sage: RDF.coerce(61/3) 

20.333333333333332 

sage: parent(RDF(3) + CDF(5)) 

Complex Double Field 

sage: parent(CDF(5) + RDF(3)) 

Complex Double Field 

sage: CDF.gen(0) + 5.0 

5.0 + 1.0*I 

sage: RLF(2/3) + RDF(1) 

1.6666666666666665 

  

sage: import numpy 

sage: RDF.coerce(numpy.int8('1')) 

1.0 

sage: RDF.coerce(numpy.float64('1')) 

1.0 

  

sage: RDF.coerce(pi) 

Traceback (most recent call last): 

... 

TypeError: no canonical coercion from Symbolic Ring to Real Double Field 

  

Test that :trac:`15695` is fixed (see also :trac:`18076`):: 

  

sage: 1j + numpy.float64(2) 

2.00000000000000 + 1.00000000000000*I 

sage: parent(_) 

Complex Field with 53 bits of precision 

""" 

if S is int or S is float: 

return ToRDF(S) 

  

from .rational_field import QQ 

from .real_lazy import RLF 

if S is ZZ or S is QQ or S is RLF: 

return ToRDF(S) 

  

from .real_mpfr import RR, RealField_class 

if isinstance(S, RealField_class): 

if S.prec() >= 53: 

return ToRDF(S) 

else: 

return None 

elif is_numpy_type(S): 

import numpy 

if issubclass(S, numpy.integer) or issubclass(S, numpy.floating): 

return ToRDF(S) 

else: 

return None 

  

connecting = RR._internal_coerce_map_from(S) 

if connecting is not None: 

return ToRDF(RR) * connecting 

  

def _magma_init_(self, magma): 

r""" 

Return a string representation of ``self`` in the Magma language. 

  

EXAMPLES: 

  

Magma handles precision in decimal digits, so we lose a bit:: 

  

sage: magma(RDF) # optional - magma # indirect doctest 

Real field of precision 15 

sage: 10^15 < 2^53 < 10^16 

True 

  

When we convert back from Magma, we convert to a generic real field 

that has 53 bits of precision:: 

  

sage: magma(RDF).sage() # optional - magma 

Real Field with 53 bits of precision 

""" 

return "RealField(%s : Bits := true)" % self.prec() 

  

def _polymake_init_(self): 

r""" 

Return the polymake representation of the real double field. 

  

EXAMPLES:: 

  

sage: polymake(RDF) #optional - polymake # indirect doctest 

Float 

""" 

return '"Float"' 

  

def precision(self): 

""" 

Return the precision of this real double field in bits. 

  

Always returns 53. 

  

EXAMPLES:: 

  

sage: RDF.precision() 

53 

""" 

return 53 

  

prec = precision 

  

def to_prec(self, prec): 

""" 

Return the real field to the specified precision. As doubles have 

fixed precision, this will only return a real double field if ``prec`` 

is exactly 53. 

  

EXAMPLES:: 

  

sage: RDF.to_prec(52) 

Real Field with 52 bits of precision 

sage: RDF.to_prec(53) 

Real Double Field 

""" 

if prec == 53: 

return self 

else: 

from .real_mpfr import RealField 

return RealField(prec) 

  

  

def gen(self, n=0): 

""" 

Return the generator of the real double field. 

  

EXAMPLES:: 

  

sage: RDF.0 

1.0 

sage: RDF.gens() 

(1.0,) 

""" 

if n != 0: 

raise ValueError("only 1 generator") 

return RealDoubleElement(1) 

  

def ngens(self): 

""" 

Return the number of generators which is always 1. 

  

EXAMPLES:: 

  

sage: RDF.ngens() 

1 

""" 

return 1 

  

def is_finite(self): 

""" 

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

  

Technical note: There exists an upper bound on the double 

representation. 

  

EXAMPLES:: 

  

sage: RDF.is_finite() 

False 

""" 

return False 

  

def characteristic(self): 

""" 

Returns 0, since the field of real numbers has characteristic 0. 

  

EXAMPLES:: 

  

sage: RDF.characteristic() 

0 

""" 

return Integer(0) 

  

cdef _new_c(self, double value): 

cdef RealDoubleElement x 

x = PY_NEW(RealDoubleElement) 

x._value = value 

return x 

  

def random_element(self, double min=-1, double max=1): 

""" 

Return a random element of this real double field in the interval 

``[min, max]``. 

  

EXAMPLES:: 

  

sage: RDF.random_element() 

0.7369454235661859 

sage: RDF.random_element(min=100, max=110) 

102.8159473516245 

""" 

cdef randstate rstate = current_randstate() 

  

return self._new_c((max-min)*rstate.c_rand_double() + min) 

  

def name(self): 

""" 

The name of ``self``. 

  

EXAMPLES:: 

  

sage: RDF.name() 

'RealDoubleField' 

""" 

return "RealDoubleField" 

  

def __hash__(self): 

""" 

Return the hash value of ``self``. 

  

TESTS:: 

  

sage: hash(RDF) % 2^32 == hash(str(RDF)) % 2^32 

True 

""" 

return 1157042230 #return hash(str(self)) 

  

def pi(self): 

r""" 

Returns `\pi` to double-precision. 

  

EXAMPLES:: 

  

sage: RDF.pi() 

3.141592653589793 

sage: RDF.pi().sqrt()/2 

0.8862269254527579 

""" 

return self(M_PI) 

  

def euler_constant(self): 

""" 

Return Euler's gamma constant to double precision. 

  

EXAMPLES:: 

  

sage: RDF.euler_constant() 

0.5772156649015329 

""" 

return self(M_EULER) 

  

def log2(self): 

r""" 

Return `\log(2)` to the precision of this field. 

  

EXAMPLES:: 

  

sage: RDF.log2() 

0.6931471805599453 

sage: RDF(2).log() 

0.6931471805599453 

""" 

return self(M_LN2) 

  

def factorial(self, int n): 

""" 

Return the factorial of the integer `n` as a real number. 

  

EXAMPLES:: 

  

sage: RDF.factorial(100) 

9.332621544394415e+157 

""" 

if n < 0: 

raise ArithmeticError("n must be nonnegative") 

return self(gsl_sf_fact(n)) 

  

def zeta(self, n=2): 

""" 

Return an `n`-th root of unity in the real field, if one 

exists, or raise a ``ValueError`` otherwise. 

  

EXAMPLES:: 

  

sage: RDF.zeta() 

-1.0 

sage: RDF.zeta(1) 

1.0 

sage: RDF.zeta(5) 

Traceback (most recent call last): 

... 

ValueError: No 5th root of unity in self 

""" 

if n == 1: 

return self(1) 

elif n == 2: 

return self(-1) 

raise ValueError("No %sth root of unity in self" % n) 

  

def NaN(self): 

""" 

Return Not-a-Number ``NaN``. 

  

EXAMPLES:: 

  

sage: RDF.NaN() 

NaN 

""" 

return self(0)/self(0) 

  

nan = NaN 

  

def _factor_univariate_polynomial(self, f): 

""" 

Factor the univariate polynomial ``f``. 

  

INPUT: 

  

- ``f`` -- a univariate polynomial defined over the double precision 

real numbers 

  

OUTPUT: 

  

- A factorization of ``f`` over the double precision real numbers 

into a unit and monic irreducible factors 

  

.. NOTE:: 

  

This is a helper method for 

:meth:`sage.rings.polynomial.polynomial_element.Polynomial.factor`. 

  

TESTS:: 

  

sage: R.<x> = RDF[] 

sage: RDF._factor_univariate_polynomial(x) 

x 

sage: RDF._factor_univariate_polynomial(2*x) 

(2.0) * x 

sage: RDF._factor_univariate_polynomial(x^2) 

x^2 

sage: RDF._factor_univariate_polynomial(x^2 + 1) 

x^2 + 1.0 

sage: RDF._factor_univariate_polynomial(x^2 - 1) 

(x - 1.0) * (x + 1.0) 

  

The implementation relies on the ``roots()`` method which often reports 

roots not to be real even though they are:: 

  

sage: f = (x-1)^3 

sage: f.roots(ring=CDF) # abs tol 2e-5 

[(1.0000065719436413, 1), 

(0.9999967140281792 - 5.691454546815028e-06*I, 1), 

(0.9999967140281792 + 5.691454546815028e-06*I, 1)] 

  

This leads to the following incorrect factorization:: 

  

sage: f.factor() # abs tol 2e-5 

(x - 1.0000065719436413) * (x^2 - 1.9999934280563585*x + 0.9999934280995487) 

""" 

roots = f.roots(sage.rings.complex_double.CDF) 

  

# collect real roots and conjugate pairs of non-real roots 

real_roots = [(r, e) for r, e in roots if r.imag().is_zero()] 

non_real_roots = {r: e for r, e in roots if not r.imag().is_zero()} 

assert all([non_real_roots[r.conj()] == e for r, e in non_real_roots.items()]), "Bug in root finding code over RDF - roots must always come in conjugate pairs" 

non_real_roots = [(r, e) for r, e in non_real_roots.items() if r.imag() > 0] 

  

# turn the roots into irreducible factors 

x = f.parent().gen() 

real_factors = [(x - r.real(), e) for r, e in real_roots] 

non_real_factors = [(x**2 - (r + r.conj()).real()*x + (r*r.conj()).real(), e) for r, e in non_real_roots] 

  

# make the factors monic 

from sage.structure.factorization import Factorization 

return Factorization([(g.monic(), e) for g, e in real_factors + non_real_factors], f.leading_coefficient()) 

  

  

cdef class RealDoubleElement(FieldElement): 

""" 

An approximation to a real number using double precision floating 

point numbers. Answers derived from calculations with such 

approximations may differ from what they would be if those 

calculations were performed with true real numbers. This is due to 

the rounding errors inherent to finite precision calculations. 

""" 

  

__array_interface__ = {'typestr': '=f8'} 

  

def __cinit__(self): 

""" 

Initialize ``self`` for cython. 

  

EXAMPLES:: 

  

sage: RDF(2.3) # indirect doctest 

2.3 

""" 

(<Element>self)._parent = _RDF 

  

def __init__(self, x): 

""" 

Create a new ``RealDoubleElement`` with value ``x``. 

  

EXAMPLES:: 

  

sage: RDF(10^100) 

1e+100 

  

TESTS:: 

  

sage: from gmpy2 import * # optional - gmpy2 

sage: RDF(mpz(42)) # optional - gmpy2 

42.0 

sage: RDF(mpq(3/4)) # optional - gmpy2 

0.75 

sage: RDF(mpq('4.1')) # optional - gmpy2 

4.1 

""" 

self._value = float(x) 

  

def _magma_init_(self, magma): 

r""" 

Return a string representation of ``self`` in the Magma language. 

  

EXAMPLES:: 

  

sage: RDF(10.5) 

10.5 

sage: magma(RDF(10.5)) # optional - magma # indirect doctest 

10.5000000000000 

""" 

return "%s!%s" % (self.parent()._magma_init_(magma), self) 

  

def __reduce__(self): 

""" 

For pickling. 

  

EXAMPLES:: 

  

sage: a = RDF(-2.7) 

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

True 

""" 

return RealDoubleElement, (self._value, ) 

  

cdef _new_c(self, double value): 

cdef RealDoubleElement x 

x = PY_NEW(RealDoubleElement) 

x._value = value 

return x 

  

def prec(self): 

""" 

Return the precision of this number in bits. 

  

Always returns 53. 

  

EXAMPLES:: 

  

sage: RDF(0).prec() 

53 

""" 

  

return 53 

  

def ulp(self): 

""" 

Returns the unit of least precision of ``self``, which is the 

weight of the least significant bit of ``self``. This is always 

a strictly positive number. It is also the gap between this 

number and the closest number with larger absolute value that 

can be represented. 

  

EXAMPLES:: 

  

sage: a = RDF(pi) 

sage: a.ulp() 

4.440892098500626e-16 

sage: b = a + a.ulp() 

  

Adding or subtracting an ulp always gives a different number:: 

  

sage: a + a.ulp() == a 

False 

sage: a - a.ulp() == a 

False 

sage: b + b.ulp() == b 

False 

sage: b - b.ulp() == b 

False 

  

Since the default rounding mode is round-to-nearest, adding or 

subtracting something less than half an ulp always gives the 

same number, unless the result has a smaller ulp. The latter 

can only happen if the input number is (up to sign) exactly a 

power of 2:: 

  

sage: a - a.ulp()/3 == a 

True 

sage: a + a.ulp()/3 == a 

True 

sage: b - b.ulp()/3 == b 

True 

sage: b + b.ulp()/3 == b 

True 

sage: c = RDF(1) 

sage: c - c.ulp()/3 == c 

False 

sage: c.ulp() 

2.220446049250313e-16 

sage: (c - c.ulp()).ulp() 

1.1102230246251565e-16 

  

The ulp is always positive:: 

  

sage: RDF(-1).ulp() 

2.220446049250313e-16 

  

The ulp of zero is the smallest positive number in RDF:: 

  

sage: RDF(0).ulp() 

5e-324 

sage: RDF(0).ulp()/2 

0.0 

  

Some special values:: 

  

sage: a = RDF(1)/RDF(0); a 

+infinity 

sage: a.ulp() 

+infinity 

sage: (-a).ulp() 

+infinity 

sage: a = RDF('nan') 

sage: a.ulp() is a 

True 

  

The ulp method works correctly with small numbers:: 

  

sage: u = RDF(0).ulp() 

sage: u.ulp() == u 

True 

sage: x = u * (2^52-1) # largest denormal number 

sage: x.ulp() == u 

True 

sage: x = u * 2^52 # smallest normal number 

sage: x.ulp() == u 

True 

  

""" 

# First, check special values 

if self._value == 0: 

return RealDoubleElement(libc.math.ldexp(1.0, -1074)) 

if gsl_isnan(self._value): 

return self 

if gsl_isinf(self._value): 

return self.abs() 

  

# Normal case 

cdef int e 

libc.math.frexp(self._value, &e) 

e -= 53 

# Correction for denormals 

if e < -1074: 

e = -1074 

return RealDoubleElement(libc.math.ldexp(1.0, e)) 

  

def real(self): 

""" 

Return ``self`` - we are already real. 

  

EXAMPLES:: 

  

sage: a = RDF(3) 

sage: a.real() 

3.0 

""" 

return self 

  

def imag(self): 

""" 

Return the imaginary part of this number, which is zero. 

  

EXAMPLES:: 

  

sage: a = RDF(3) 

sage: a.imag() 

0.0 

""" 

return RealDoubleElement(0) 

  

def __complex__(self): 

""" 

Return ``self`` as a python complex number. 

  

EXAMPLES:: 

  

sage: a = 2303 

sage: RDF(a) 

2303.0 

sage: complex(RDF(a)) 

(2303+0j) 

""" 

return complex(self._value,0) 

  

def _integer_(self, ZZ=None): 

""" 

If this floating-point number is actually an integer, return 

that integer. Otherwise, raise an exception. 

  

EXAMPLES:: 

  

sage: ZZ(RDF(237.0)) # indirect doctest 

237 

sage: ZZ(RDF(0.0/0.0)) 

Traceback (most recent call last): 

... 

ValueError: cannot convert float NaN to integer 

sage: ZZ(RDF(1.0/0.0)) 

Traceback (most recent call last): 

... 

OverflowError: cannot convert float infinity to integer 

sage: ZZ(RDF(-123456789.0)) 

-123456789 

sage: ZZ(RDF((2.0))^290) 

1989292945639146568621528992587283360401824603189390869761855907572637988050133502132224 

sage: ZZ(RDF(-2345.67)) 

Traceback (most recent call last): 

... 

TypeError: Cannot convert non-integral float to integer 

""" 

return Integer(self._value) 

  

def __mpfr__(self): 

""" 

Convert Sage ``RealDoubleElement`` to gmpy2 ``mpfr``. 

  

EXAMPLES:: 

  

sage: RDF(42.2).__mpfr__() # optional - gmpy2 

mpfr('42.200000000000003') 

sage: from gmpy2 import mpfr # optional - gmpy2 

sage: mpfr(RDF(5.1)) # optional - gmpy2 

mpfr('5.0999999999999996') 

  

TESTS:: 

  

sage: RDF().__mpfr__(); raise NotImplementedError("gmpy2 is not installed") 

Traceback (most recent call last): 

... 

NotImplementedError: gmpy2 is not installed 

""" 

IF HAVE_GMPY2: 

return gmpy2.mpfr(self._value) 

ELSE: 

raise NotImplementedError("gmpy2 is not installed") 

  

def _interface_init_(self, I=None): 

""" 

Return ``self`` formatted as a string, suitable as input to another 

computer algebra system. (This is the default function used for 

exporting to other computer algebra systems.) 

  

EXAMPLES:: 

  

sage: s1 = RDF(sin(1)); s1 

0.8414709848078965 

sage: s1._interface_init_() 

'0.8414709848078965' 

sage: s1 == RDF(gp(s1)) 

True 

""" 

return repr(self._value) 

  

def _sage_input_(self, sib, coerced): 

r""" 

Produce an expression which will reproduce this value when evaluated. 

  

EXAMPLES:: 

  

sage: sage_input(RDF(NaN)) 

RDF(NaN) 

sage: sage_input(RDF(-infinity), verify=True) 

# Verified 

-RDF(infinity) 

sage: sage_input(RDF(-infinity)*polygen(RDF)) 

R.<x> = RDF[] 

-RDF(infinity)*x + RDF(NaN) 

sage: sage_input(RDF(pi), verify=True) 

# Verified 

RDF(3.1415926535897931) 

sage: sage_input(RDF(-e), verify=True, preparse=False) 

# Verified 

-RDF(2.718281828459045...) 

sage: sage_input(RDF(pi)*polygen(RDF), verify=True, preparse=None) 

# Verified 

R = RDF['x'] 

x = R.gen() 

3.1415926535897931*x 

sage: from sage.misc.sage_input import SageInputBuilder 

sage: sib = SageInputBuilder() 

sage: RDF(22/7)._sage_input_(sib, True) 

{atomic:3.1428571428571428} 

sage: RDF(22/7)._sage_input_(sib, False) 

{call: {atomic:RDF}({atomic:3.1428571428571428})} 

""" 

cdef int isinf = gsl_isinf(self._value) 

cdef bint isnan = gsl_isnan(self._value) 

if isinf or isnan: 

if isnan: 

v = sib.name('NaN') 

else: 

v = sib.name('infinity') 

v = sib(self.parent())(v) 

if isinf < 0: 

v = -v 

return v 

  

from sage.rings.all import ZZ, RR 

  

cdef bint negative = self._value < 0 

if negative: 

self = -self 

  

# There are five possibilities for printing this floating-point 

# number, ordered from prettiest to ugliest (IMHO). 

# 1) An integer: 42 

# 2) A simple literal: 3.14159 

# 3) A coerced integer: RDF(42) 

# 4) A coerced literal: RDF(3.14159) 

# 5) A coerced RR value: RDF(RR('3.14159')) 

  

# str(self) works via libc, which we don't necessarily trust 

# to produce the best possible answer. So this function prints 

# via RR/MPFR. Without the preparser, input works via libc as 

# well, but we don't have a choice about that. 

  

# To use choice 1 or choice 3, this number must be an integer. 

cdef bint can_use_int_literal = \ 

self.abs() < (Integer(1) << self.prec()) and self in ZZ 

  

self_str = RR(self._value).str(truncate=False, skip_zeroes=True) 

  

# To use choice 2 or choice 4, we must be able to read 

# numbers of this precision as a literal. 

cdef bint can_use_float_literal = \ 

(sib.preparse() or float(self_str) == self) 

  

if can_use_int_literal or can_use_float_literal: 

if can_use_int_literal: 

v = sib.int(self._integer_()) 

else: 

v = sib.float_str(self_str) 

else: 

v = sib(RR(self)) 

if not coerced: 

v = sib(self.parent())(v) 

  

if negative: 

v = -v 

  

return v 

  

def _repr_(self): 

""" 

Return the string representation of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(-2/3) 

-0.6666666666666666 

sage: a = RDF(2); a 

2.0 

sage: a^2 

4.0 

sage: RDF("nan") 

NaN 

sage: RDF(1)/RDF(0) 

+infinity 

sage: RDF(-1)/RDF(0) 

-infinity 

""" 

return double_repr(self._value) 

  

def _latex_(self): 

r""" 

Return a latex representation of ``self``. 

  

EXAMPLES:: 

  

sage: latex(RDF(3.4)) # indirect doctest 

3.4 

sage: latex(RDF(2e-100)) # indirect doctest 

2 \times 10^{-100} 

""" 

s = self.str() 

parts = s.split('e') 

if len(parts) > 1: 

# scientific notation 

if parts[1][0] == '+': 

parts[1] = parts[1][1:] 

s = "%s \\times 10^{%s}" % (parts[0], parts[1]) 

return s 

  

def __hash__(self): 

""" 

Return the hash of ``self``, which coincides with the python float 

(and often int) type. 

  

EXAMPLES:: 

  

sage: hash(RDF(1.2)) == hash(1.2r) 

True 

""" 

return hash(self._value) 

  

def _im_gens_(self, codomain, im_gens): 

""" 

Return the image of ``self`` under the homomorphism from the rational 

field to ``codomain``. 

  

This always just returns ``self`` coerced into the ``codomain``. 

  

EXAMPLES:: 

  

sage: RDF(2.1)._im_gens_(RR, [RR(1)]) 

2.10000000000000 

sage: R = RealField(20) 

sage: RDF(2.1)._im_gens_(R, [R(1)]) 

2.1000 

""" 

return codomain(self) # since 1 |--> 1 

  

def __str__(self): 

""" 

Return the string representation of ``self``, see :meth:`str`. 

  

EXAMPLES:: 

  

sage: print(RDF(-2/3)) 

-0.666666666667 

sage: print(RDF(oo)) 

+infinity 

""" 

return double_str(self._value) 

  

def str(self): 

""" 

Return the informal string representation of ``self``. 

  

EXAMPLES:: 

  

sage: a = RDF('4.5'); a.str() 

'4.5' 

sage: a = RDF('49203480923840.2923904823048'); a.str() 

'4.92034809238e+13' 

sage: a = RDF(1)/RDF(0); a.str() 

'+infinity' 

sage: a = -RDF(1)/RDF(0); a.str() 

'-infinity' 

sage: a = RDF(0)/RDF(0); a.str() 

'NaN' 

  

We verify consistency with ``RR`` (mpfr reals):: 

  

sage: str(RR(RDF(1)/RDF(0))) == str(RDF(1)/RDF(0)) 

True 

sage: str(RR(-RDF(1)/RDF(0))) == str(-RDF(1)/RDF(0)) 

True 

sage: str(RR(RDF(0)/RDF(0))) == str(RDF(0)/RDF(0)) 

True 

""" 

return double_str(self._value) 

  

def __copy__(self): 

""" 

Return copy of ``self``, which since ``self`` is immutable, is just 

``self``. 

  

EXAMPLES:: 

  

sage: r = RDF('-1.6') 

sage: r.__copy__() is r 

True 

""" 

return self 

  

def integer_part(self): 

""" 

If in decimal this number is written ``n.defg``, returns ``n``. 

  

EXAMPLES:: 

  

sage: r = RDF('-1.6') 

sage: a = r.integer_part(); a 

-1 

sage: type(a) 

<type 'sage.rings.integer.Integer'> 

sage: r = RDF(0.0/0.0) 

sage: a = r.integer_part() 

Traceback (most recent call last): 

... 

TypeError: Attempt to get integer part of NaN 

""" 

if gsl_isnan(self._value): 

raise TypeError("Attempt to get integer part of NaN") 

else: 

return Integer(int(self._value)) 

  

def sign_mantissa_exponent(self): 

r""" 

Return the sign, mantissa, and exponent of ``self``. 

  

In Sage (as in MPFR), floating-point numbers of precision `p` 

are of the form `s m 2^{e-p}`, where `s \in \{-1, 1\}`, 

`2^{p-1} \leq m < 2^p`, and `-2^{30} + 1 \leq e \leq 2^{30} - 

1`; plus the special values ``+0``, ``-0``, ``+infinity``, 

``-infinity``, and ``NaN`` (which stands for Not-a-Number). 

  

This function returns `s`, `m`, and `e-p`. For the special values: 

  

- ``+0`` returns ``(1, 0, 0)`` 

- ``-0`` returns ``(-1, 0, 0)`` 

- the return values for ``+infinity``, ``-infinity``, and ``NaN`` are 

not specified. 

  

EXAMPLES:: 

  

sage: a = RDF(exp(1.0)); a 

2.718281828459045 

sage: sign,mantissa,exponent = RDF(exp(1.0)).sign_mantissa_exponent() 

sage: sign,mantissa,exponent 

(1, 6121026514868073, -51) 

sage: sign*mantissa*(2**exponent) == a 

True 

  

The mantissa is always a nonnegative number:: 

  

sage: RDF(-1).sign_mantissa_exponent() 

(-1, 4503599627370496, -52) 

  

TESTS:: 

  

sage: RDF('+0').sign_mantissa_exponent() 

(1, 0, 0) 

sage: RDF('-0').sign_mantissa_exponent() 

(-1, 0, 0) 

""" 

from sage.rings.all import RR 

return RR(self._value).sign_mantissa_exponent() 

  

  

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

# Basic Arithmetic 

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

def __invert__(self): 

""" 

Compute the multiplicative inverse of ``self``. 

  

EXAMPLES:: 

  

sage: a = RDF(-1.5)*RDF(2.5) 

sage: a.__invert__() 

-0.26666666666666666 

sage: ~a 

-0.26666666666666666 

""" 

cdef RealDoubleElement x = <RealDoubleElement>PY_NEW(RealDoubleElement) 

x._value = 1.0 / self._value 

return x 

  

cpdef _add_(self, right): 

""" 

Add two real numbers with the same parent. 

  

EXAMPLES:: 

  

sage: RDF('-1.5') + RDF('2.5') # indirect doctest 

1.0 

""" 

cdef RealDoubleElement x = <RealDoubleElement>PY_NEW(RealDoubleElement) 

x._value = self._value + (<RealDoubleElement>right)._value 

return x 

  

cpdef _sub_(self, right): 

""" 

Subtract two real numbers with the same parent. 

  

EXAMPLES:: 

  

sage: RDF('-1.5') - RDF('2.5') # indirect doctest 

-4.0 

""" 

cdef RealDoubleElement x = <RealDoubleElement>PY_NEW(RealDoubleElement) 

x._value = self._value - (<RealDoubleElement>right)._value 

return x 

  

cpdef _mul_(self, right): 

""" 

Multiply two real numbers with the same parent. 

  

EXAMPLES:: 

  

sage: RDF('-1.5') * RDF('2.5') # indirect doctest 

-3.75 

""" 

cdef RealDoubleElement x = <RealDoubleElement>PY_NEW(RealDoubleElement) 

x._value = self._value * (<RealDoubleElement>right)._value 

return x 

  

cpdef _div_(self, right): 

""" 

Divide ``self`` by ``right``. 

  

EXAMPLES:: 

  

sage: RDF('-1.5') / RDF('2.5') # indirect doctest 

-0.6 

sage: RDF(1)/RDF(0) 

+infinity 

""" 

cdef RealDoubleElement x = <RealDoubleElement>PY_NEW(RealDoubleElement) 

x._value = self._value / (<RealDoubleElement>right)._value 

return x 

  

def __neg__(self): 

""" 

Negates ``self``. 

  

EXAMPLES:: 

  

sage: -RDF('-1.5') 

1.5 

""" 

cdef RealDoubleElement x = <RealDoubleElement>PY_NEW(RealDoubleElement) 

x._value = -self._value 

return x 

  

def conjugate(self): 

r""" 

Returns the complex conjugate of this real number, which is 

the real number itself. 

  

EXAMPLES:: 

  

sage: RDF(4).conjugate() 

4.0 

""" 

return self 

  

def __abs__(self): 

""" 

Returns the absolute value of ``self``. 

  

EXAMPLES:: 

  

sage: abs(RDF(1.5)) 

1.5 

sage: abs(RDF(-1.5)) 

1.5 

sage: abs(RDF(0.0)) 

0.0 

sage: abs(RDF(-0.0)) 

0.0 

""" 

# Use signbit instead of >= to handle -0.0 correctly 

if not libc.math.signbit(self._value): 

return self 

else: 

return self._new_c(-self._value) 

  

cpdef RealDoubleElement abs(RealDoubleElement self): 

""" 

Returns the absolute value of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(1e10).abs() 

10000000000.0 

sage: RDF(-1e10).abs() 

10000000000.0 

""" 

if self._value >= 0: 

return self 

else: 

return self._new_c(-self._value) 

  

def __lshift__(x, y): 

""" 

LShifting a double is not supported; nor is lshifting a 

:class:`RealDoubleElement`. 

  

TESTS:: 

  

sage: RDF(2) << 3 

Traceback (most recent call last): 

... 

TypeError: unsupported operand type(s) for << 

""" 

raise TypeError("unsupported operand type(s) for <<") 

  

def __rshift__(x, y): 

""" 

RShifting a double is not supported; nor is rshifting a 

:class:`RealDoubleElement`. 

  

TESTS:: 

  

sage: RDF(2) >> 3 

Traceback (most recent call last): 

... 

TypeError: unsupported operand type(s) for >> 

""" 

raise TypeError("unsupported operand type(s) for >>") 

  

def multiplicative_order(self): 

r""" 

Returns `n` such that ``self^n == 1``. 

  

Only `\pm 1` have finite multiplicative order. 

  

EXAMPLES:: 

  

sage: RDF(1).multiplicative_order() 

1 

sage: RDF(-1).multiplicative_order() 

2 

sage: RDF(3).multiplicative_order() 

+Infinity 

""" 

if self._value == 1: 

return 1 

elif self._value == -1: 

return 2 

return sage.rings.infinity.infinity 

  

def sign(self): 

""" 

Returns -1,0, or 1 if ``self`` is negative, zero, or positive; 

respectively. 

  

EXAMPLES:: 

  

sage: RDF(-1.5).sign() 

-1 

sage: RDF(0).sign() 

0 

sage: RDF(2.5).sign() 

1 

""" 

if not self._value: 

return 0 

if self._value > 0: 

return 1 

return -1 

  

  

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

# Rounding etc 

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

  

def round(self): 

""" 

Given real number `x`, rounds up if fractional part is greater than 

`0.5`, rounds down if fractional part is less than `0.5`. 

  

EXAMPLES:: 

  

sage: RDF(0.49).round() 

0 

sage: a=RDF(0.51).round(); a 

1 

""" 

return Integer(round(self._value)) 

  

def floor(self): 

""" 

Return the floor of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(2.99).floor() 

2 

sage: RDF(2.00).floor() 

2 

sage: RDF(-5/2).floor() 

-3 

""" 

return Integer(math.floor(self._value)) 

  

def ceil(self): 

""" 

Return the ceiling of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(2.99).ceil() 

3 

sage: RDF(2.00).ceil() 

2 

sage: RDF(-5/2).ceil() 

-2 

""" 

return Integer(math.ceil(self._value)) 

  

ceiling = ceil 

  

def trunc(self): 

""" 

Truncates this number (returns integer part). 

  

EXAMPLES:: 

  

sage: RDF(2.99).trunc() 

2 

sage: RDF(-2.00).trunc() 

-2 

sage: RDF(0.00).trunc() 

0 

""" 

return Integer(int(self._value)) 

  

def frac(self): 

""" 

Return a real number in `(-1, 1)`. It satisfies the relation: 

``x = x.trunc() + x.frac()`` 

  

EXAMPLES:: 

  

sage: RDF(2.99).frac() 

0.9900000000000002 

sage: RDF(2.50).frac() 

0.5 

sage: RDF(-2.79).frac() 

-0.79 

""" 

return self._new_c(self._value - int(self._value)) 

  

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

# Conversions 

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

  

def __float__(self): 

""" 

Return ``self`` as a python float. 

  

EXAMPLES:: 

  

sage: float(RDF(1.5)) 

1.5 

sage: type(float(RDF(1.5))) 

<... 'float'> 

""" 

return self._value 

  

def _rpy_(self): 

""" 

Returns ``self.__float__()`` for rpy to convert into the 

appropriate R object. 

  

EXAMPLES:: 

  

sage: n = RDF(2.0) 

sage: n._rpy_() 

2.0 

sage: type(n._rpy_()) 

<... 'float'> 

""" 

return self.__float__() 

  

def __int__(self): 

""" 

Return integer truncation of this real number. 

  

EXAMPLES:: 

  

sage: int(RDF(2.99)) 

2 

sage: int(RDF(-2.99)) 

-2 

""" 

return int(self._value) 

  

def __long__(self): 

""" 

Returns long integer truncation of this real number. 

  

EXAMPLES:: 

  

sage: int(RDF(10e15)) 

10000000000000000L # 32-bit 

10000000000000000 # 64-bit 

sage: long(RDF(2^100)) == 2^100 

True 

""" 

return long(self._value) 

  

def _complex_mpfr_field_(self, CC): 

""" 

EXAMPLES:: 

  

sage: a = RDF(1/3) 

sage: CC(a) 

0.333333333333333 

sage: a._complex_mpfr_field_(CC) 

0.333333333333333 

  

If we coerce to a higher-precision field the extra bits appear 

random; they are actually 0's in base 2. 

  

:: 

  

sage: a._complex_mpfr_field_(ComplexField(100)) 

0.33333333333333331482961625625 

sage: a._complex_mpfr_field_(ComplexField(100)).str(2) 

'0.01010101010101010101010101010101010101010101010101010100000000000000000000000000000000000000000000000' 

""" 

return CC(self._value) 

  

def _complex_double_(self, CDF): 

""" 

Return ``self`` as a complex double. 

  

EXAMPLES:: 

  

sage: CDF(RDF(1/3)) # indirect doctest 

0.3333333333333333 

""" 

return CDF(self._value) 

  

def __pari__(self): 

""" 

Return a PARI representation of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(1.5).__pari__() 

1.50000000000000 

""" 

return new_gen_from_double(self._value) 

  

  

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

# Comparisons: ==, !=, <, <=, >, >= 

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

  

def is_NaN(self): 

""" 

Check if ``self`` is ``NaN``. 

  

EXAMPLES:: 

  

sage: RDF(1).is_NaN() 

False 

sage: a = RDF(0)/RDF(0) 

sage: a.is_NaN() 

True 

""" 

return gsl_isnan(self._value) 

  

def is_positive_infinity(self): 

r""" 

Check if ``self`` is `+\infty`. 

  

EXAMPLES:: 

  

sage: a = RDF(1)/RDF(0) 

sage: a.is_positive_infinity() 

True 

sage: a = RDF(-1)/RDF(0) 

sage: a.is_positive_infinity() 

False 

""" 

return gsl_isinf(self._value) > 0 

  

def is_negative_infinity(self): 

r""" 

Check if ``self`` is `-\infty`. 

  

EXAMPLES:: 

  

sage: a = RDF(2)/RDF(0) 

sage: a.is_negative_infinity() 

False 

sage: a = RDF(-3)/RDF(0) 

sage: a.is_negative_infinity() 

True 

""" 

return gsl_isinf(self._value) < 0 

  

def is_infinity(self): 

r""" 

Check if ``self`` is `\infty`. 

  

EXAMPLES:: 

  

sage: a = RDF(2); b = RDF(0) 

sage: (a/b).is_infinity() 

True 

sage: (b/a).is_infinity() 

False 

""" 

return gsl_isinf(self._value) 

  

cpdef _richcmp_(left, right, int op): 

""" 

Rich comparison of ``left`` and ``right``. 

  

EXAMPLES:: 

  

sage: RDF(2) < RDF(0) 

False 

sage: RDF(2) == RDF(4/2) 

True 

sage: RDF(-2) > RDF(-4) 

True 

  

TESTS: 

  

Check comparisons with ``NaN`` (:trac:`16515`):: 

  

sage: n = RDF('NaN') 

sage: n == n 

False 

sage: n == RDF(1) 

False 

""" 

# We really need to use the correct operators, to deal 

# correctly with NaNs. 

cdef double x = (<RealDoubleElement>left)._value 

cdef double y = (<RealDoubleElement>right)._value 

if op == Py_LT: 

return x < y 

elif op == Py_LE: 

return x <= y 

elif op == Py_EQ: 

return x == y 

elif op == Py_NE: 

return x != y 

elif op == Py_GT: 

return x > y 

else: 

return x >= y 

  

  

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

# Special Functions 

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

  

def NaN(self): 

""" 

Return Not-a-Number ``NaN``. 

  

EXAMPLES:: 

  

sage: RDF.NaN() 

NaN 

""" 

return self(0)/self(0) 

  

nan = NaN 

  

def sqrt(self, extend=True, all=False): 

""" 

The square root function. 

  

INPUT: 

  

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

square root in a complex field if necessary if ``self`` is negative; 

otherwise raise a ``ValueError``. 

  

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

list of all square roots. 

  

EXAMPLES:: 

  

sage: r = RDF(4.0) 

sage: r.sqrt() 

2.0 

sage: r.sqrt()^2 == r 

True 

  

:: 

  

sage: r = RDF(4344) 

sage: r.sqrt() 

65.90902821313632 

sage: r.sqrt()^2 - r 

0.0 

  

:: 

  

sage: r = RDF(-2.0) 

sage: r.sqrt() 

1.4142135623730951*I 

  

:: 

  

sage: RDF(2).sqrt(all=True) 

[1.4142135623730951, -1.4142135623730951] 

sage: RDF(0).sqrt(all=True) 

[0.0] 

sage: RDF(-2).sqrt(all=True) 

[1.4142135623730951*I, -1.4142135623730951*I] 

""" 

if self._value >= 0: 

x = self._new_c(libc.math.sqrt(self._value)) 

if all: 

if x.is_zero(): 

return [x] 

else: 

return [x, -x] 

else: 

return x 

if not extend: 

raise ValueError("negative number %s does not have a square root in the real field" % self) 

import sage.rings.complex_double 

return self._complex_double_(sage.rings.complex_double.CDF).sqrt(all=all) 

  

def is_square(self): 

""" 

Return whether or not this number is a square in this field. For 

the real numbers, this is ``True`` if and only if ``self`` is 

non-negative. 

  

EXAMPLES:: 

  

sage: RDF(3.5).is_square() 

True 

sage: RDF(0).is_square() 

True 

sage: RDF(-4).is_square() 

False 

""" 

return self._value >= 0 

  

def is_integer(self): 

""" 

Return True if this number is a integer 

  

EXAMPLES:: 

  

sage: RDF(3.5).is_integer() 

False 

sage: RDF(3).is_integer() 

True 

""" 

return self._value in ZZ 

  

  

def cube_root(self): 

""" 

Return the cubic root (defined over the real numbers) of ``self``. 

  

EXAMPLES:: 

  

sage: r = RDF(125.0); r.cube_root() 

5.000000000000001 

sage: r = RDF(-119.0) 

sage: r.cube_root()^3 - r # rel tol 1 

-1.4210854715202004e-14 

""" 

return self.nth_root(3) 

  

  

def nth_root(self, int n): 

""" 

Return the `n^{th}` root of ``self``. 

  

INPUT: 

  

- ``n`` -- an integer 

  

OUTPUT: 

  

The output is a complex double if ``self`` is negative and `n` is even, 

otherwise it is a real double. 

  

EXAMPLES:: 

  

sage: r = RDF(-125.0); r.nth_root(3) 

-5.000000000000001 

sage: r.nth_root(5) 

-2.6265278044037674 

sage: RDF(-2).nth_root(5)^5 # rel tol 1e-15 

-2.000000000000001 

sage: RDF(-1).nth_root(5)^5 

-1.0 

sage: RDF(3).nth_root(10)^10 

2.9999999999999982 

sage: RDF(-1).nth_root(2) 

6.123233995736757e-17 + 1.0*I 

sage: RDF(-1).nth_root(4) 

0.7071067811865476 + 0.7071067811865475*I 

""" 

if n == 0: 

return RealDoubleElement(float('nan')) 

if self._value < 0: 

if GSL_IS_EVEN(n): 

return self._complex_double_(sage.rings.complex_double.CDF).nth_root(n) 

else: 

return - ( (-self) ** (float(1)/n) ) 

else: 

return self ** (float(1)/n) 

  

cdef __pow_double(self, double exponent, double sign): 

""" 

If ``sign == 1`` or ``self >= 0``, return ``self ^ exponent``. 

If ``sign == -1`` and ``self < 0``, return ``- abs(self) ^ exponent``. 

""" 

cdef double v = self._value 

if v >= 0: 

if v == 1: 

return self 

elif exponent == 0: 

return self._new_c(1.0) 

elif v == 0: 

if exponent < 0: 

raise ZeroDivisionError("0.0 cannot be raised to a negative power") 

return self 

sign = 1.0 

else: # v < 0 

expmod2 = libc.math.fmod(exponent, 2.0) 

if expmod2 == 0.0: 

pass 

elif expmod2 == 1.0: 

sign = -1.0 

else: 

raise ValueError("negative number cannot be raised to a fractional power") 

v = -v 

return self._new_c(sign * gsl_sf_exp(gsl_sf_log(v) * exponent)) 

  

cpdef _pow_(self, other): 

""" 

Return ``self`` raised to the real double power ``other``. 

  

EXAMPLES:: 

  

sage: a = RDF('1.23456') 

sage: a^a 

1.2971114817819216 

  

TESTS:: 

  

sage: RDF(0) ^ RDF(0.5) 

0.0 

sage: RDF(0) ^ (1/2) 

0.0 

sage: RDF(0) ^ RDF(0) 

1.0 

sage: RDF(0) ^ RDF(-1) 

Traceback (most recent call last): 

... 

ZeroDivisionError: 0.0 cannot be raised to a negative power 

sage: RDF(-1) ^ RDF(0) 

1.0 

sage: RDF(-1) ^ RDF(1) 

-1.0 

sage: RDF(-1) ^ RDF(0.5) 

Traceback (most recent call last): 

... 

ValueError: negative number cannot be raised to a fractional power 

""" 

return self.__pow_double((<RealDoubleElement>other)._value, 1) 

  

cpdef _pow_int(self, n): 

""" 

Return ``self`` raised to the integer power ``n``. 

  

TESTS:: 

  

sage: RDF(1) ^ (2^1000) 

1.0 

sage: RDF(1) ^ (2^1000 + 1) 

1.0 

sage: RDF(1) ^ (-2^1000) 

1.0 

sage: RDF(1) ^ (-2^1000 + 1) 

1.0 

sage: RDF(-1) ^ (2^1000) 

1.0 

sage: RDF(-1) ^ (2^1000 + 1) 

-1.0 

sage: RDF(-1) ^ (-2^1000) 

1.0 

sage: RDF(-1) ^ (-2^1000 + 1) 

-1.0 

  

:: 

  

sage: base = RDF(1.0000000000000002) 

sage: base._pow_int(0) 

1.0 

sage: base._pow_int(1) 

1.0000000000000002 

sage: base._pow_int(2) 

1.0000000000000004 

sage: base._pow_int(3) 

1.0000000000000007 

sage: base._pow_int(2^57) 

78962960182680.42 

sage: base._pow_int(2^57 + 1) 

78962960182680.42 

  

:: 

  

sage: base = RDF(-1.0000000000000002) 

sage: base._pow_int(0) 

1.0 

sage: base._pow_int(1) 

-1.0000000000000002 

sage: base._pow_int(2) 

1.0000000000000004 

sage: base._pow_int(3) 

-1.0000000000000007 

sage: base._pow_int(2^57) 

78962960182680.42 

sage: base._pow_int(2^57 + 1) 

-78962960182680.42 

""" 

return self.__pow_double(n, -1.0 if (n & 1) else 1.0) 

  

cdef _pow_long(self, long n): 

""" 

Compute ``self`` raised to the power ``n``. 

  

EXAMPLES:: 

  

sage: RDF('1.23456') ^ 20 

67.64629770385... 

sage: RDF(3) ^ 32 

1853020188851841.0 

sage: RDF(2)^(-1024) 

5.562684646268003e-309 

  

TESTS:: 

  

sage: base = RDF(1.0000000000000002) 

sage: base ^ RDF(2^31) 

1.000000476837272 

sage: base ^ (2^57) 

78962960182680.42 

sage: base ^ RDF(2^57) 

78962960182680.42 

""" 

if -2048 <= n <= 2048: 

# For small exponents, it is possible that the powering 

# is exact either because the base is a power of 2 

# (e.g. 2.0^1000) or because the exact result has few 

# significant digits (e.g. 3.0^10). Here, we use the 

# square-and-multiply algorithm by GSL. 

return self._new_c(gsl_pow_int(self._value, <int>n)) 

# If the exponent is sufficiently large in absolute value, the 

# result cannot be exact (except if the base is -1.0, 0.0 or 

# 1.0 but those cases are handled by __pow_double too). The 

# log-and-exp algorithm from __pow_double will be more precise 

# than square-and-multiply. 

  

# We do need to take care of the sign since the conversion 

# of n to double might change an odd number to an even number. 

return self.__pow_double(<double>n, -1.0 if (n & 1) else 1.0) 

  

cdef _log_base(self, double log_of_base): 

if self._value == 0: 

return RDF(-1)/RDF(0) 

elif self._value < 0: 

return RDF.NaN() 

sig_on() 

a = self._new_c(gsl_sf_log(self._value) / log_of_base) 

sig_off() 

return a 

  

def log(self, base=None): 

""" 

Return the logarithm. 

  

INPUT: 

  

- ``base`` -- integer or ``None`` (default). The base of the 

logarithm. If ``None`` is specified, the base is `e` (the so-called 

natural logarithm). 

  

OUTPUT: 

  

The logarithm of ``self``. If ``self`` is positive, a double 

floating point number. Infinity if ``self`` is zero. A 

imaginary complex floating point number if ``self`` is 

negative. 

  

EXAMPLES:: 

  

sage: RDF(2).log() 

0.6931471805599453 

sage: RDF(2).log(2) 

1.0 

sage: RDF(2).log(pi) 

0.6055115613982801 

sage: RDF(2).log(10) 

0.30102999566398114 

sage: RDF(2).log(1.5) 

1.7095112913514547 

sage: RDF(0).log() 

-infinity 

sage: RDF(-1).log() 

3.141592653589793*I 

sage: RDF(-1).log(2) # rel tol 1e-15 

4.532360141827194*I 

  

TESTS: 

  

Make sure that we can take the log of small numbers accurately 

and the fix doesn't break preexisting values (:trac:`12557`):: 

  

sage: R = RealField(128) 

sage: def check_error(x): 

....: x = RDF(x) 

....: log_RDF = x.log() 

....: log_RR = R(x).log() 

....: diff = R(log_RDF) - log_RR 

....: if abs(diff) < log_RDF.ulp(): 

....: return True 

....: print("logarithm check failed for %s (diff = %s ulp)"% \ 

....: (x, diff/log_RDF.ulp())) 

....: return False 

sage: all( check_error(2^x) for x in range(-100,100) ) 

True 

sage: all( check_error(x) for x in sxrange(0.01, 2.00, 0.01) ) 

True 

sage: all( check_error(x) for x in sxrange(0.99, 1.01, 0.001) ) 

True 

sage: RDF(1.000000001).log() 

1.000000082240371e-09 

sage: RDF(1e-17).log() 

-39.14394658089878 

sage: RDF(1e-50).log() 

-115.12925464970229 

""" 

if self < 0: 

from sage.rings.complex_double import CDF 

return CDF(self).log(base) 

if base is None: 

return self._log_base(1) 

else: 

if isinstance(base, RealDoubleElement): 

return self._log_base(base._log_base(1)) 

else: 

return self._log_base(gsl_sf_log(float(base))) 

  

def log2(self): 

""" 

Return log to the base 2 of ``self``. 

  

EXAMPLES:: 

  

sage: r = RDF(16.0) 

sage: r.log2() 

4.0 

  

:: 

  

sage: r = RDF(31.9); r.log2() 

4.9954845188775066 

""" 

if self < 0: 

from sage.rings.complex_double import CDF 

return CDF(self).log(2) 

sig_on() 

a = self._new_c(gsl_sf_log(self._value) * M_1_LN2) 

sig_off() 

return a 

  

  

def log10(self): 

""" 

Return log to the base 10 of ``self``. 

  

EXAMPLES:: 

  

sage: r = RDF('16.0'); r.log10() 

1.2041199826559248 

sage: r.log() / RDF(log(10)) 

1.2041199826559246 

sage: r = RDF('39.9'); r.log10() 

1.6009728956867482 

""" 

if self < 0: 

from sage.rings.complex_double import CDF 

return CDF(self).log(10) 

sig_on() 

a = self._new_c(gsl_sf_log(self._value) * M_1_LN10) 

sig_off() 

return a 

  

def logpi(self): 

r""" 

Return log to the base `\pi` of ``self``. 

  

EXAMPLES:: 

  

sage: r = RDF(16); r.logpi() 

2.4220462455931204 

sage: r.log() / RDF(log(pi)) 

2.4220462455931204 

sage: r = RDF('39.9'); r.logpi() 

3.2203023346075152 

""" 

if self < 0: 

from sage.rings.complex_double import CDF 

return CDF(self).log(M_PI) 

sig_on() 

a = self._new_c(gsl_sf_log(self._value) * M_1_LNPI) 

sig_off() 

return a 

  

def exp(self): 

r""" 

Return `e^\mathtt{self}`. 

  

EXAMPLES:: 

  

sage: r = RDF(0.0) 

sage: r.exp() 

1.0 

  

:: 

  

sage: r = RDF('32.3') 

sage: a = r.exp(); a 

106588847274864.47 

sage: a.log() 

32.3 

  

:: 

  

sage: r = RDF('-32.3') 

sage: r.exp() 

9.381844588498685e-15 

  

:: 

  

sage: RDF(1000).exp() 

+infinity 

""" 

sig_on() 

a = self._new_c(gsl_sf_exp(self._value)) 

sig_off() 

return a 

  

def exp2(self): 

""" 

Return `2^\mathtt{self}`. 

  

EXAMPLES:: 

  

sage: r = RDF(0.0) 

sage: r.exp2() 

1.0 

  

:: 

  

sage: r = RDF(32.0) 

sage: r.exp2() 

4294967295.9999967 

  

:: 

  

sage: r = RDF(-32.3) 

sage: r.exp2() 

1.8911724825302065e-10 

""" 

sig_on() 

a = self._new_c(gsl_sf_exp(self._value * M_LN2)) 

sig_off() 

return a 

  

def exp10(self): 

r""" 

Return `10^\mathtt{self}`. 

  

EXAMPLES:: 

  

sage: r = RDF(0.0) 

sage: r.exp10() 

1.0 

  

:: 

  

sage: r = RDF(32.0) 

sage: r.exp10() 

1.0000000000000069e+32 

  

:: 

  

sage: r = RDF(-32.3) 

sage: r.exp10() 

5.011872336272702e-33 

""" 

sig_on() 

a = self._new_c(gsl_sf_exp(self._value * M_LN10)) 

sig_off() 

return a 

  

def cos(self): 

""" 

Return the cosine of ``self``. 

  

EXAMPLES:: 

  

sage: t=RDF.pi()/2 

sage: t.cos() 

6.123233995736757e-17 

""" 

return self._new_c(gsl_sf_cos(self._value)) 

  

def sin(self): 

""" 

Return the sine of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(2).sin() 

0.9092974268256817 

""" 

return self._new_c(gsl_sf_sin(self._value)) 

  

def dilog(self): 

r""" 

Return the dilogarithm of ``self``. 

  

This is defined by the 

series `\sum_n x^n/n^2` for `|x| \le 1`. When the absolute 

value of ``self`` is greater than 1, the returned value is the 

real part of (the analytic continuation to `\CC` of) the 

dilogarithm of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(1).dilog() # rel tol 1.0e-13 

1.6449340668482264 

sage: RDF(2).dilog() # rel tol 1.0e-13 

2.46740110027234 

""" 

return self._new_c(gsl_sf_dilog(self._value)) 

  

def restrict_angle(self): 

r""" 

Return a number congruent to ``self`` mod `2\pi` that lies in 

the interval `(-\pi, \pi]`. 

  

Specifically, it is the unique `x \in (-\pi, \pi]` such 

that ```self`` `= x + 2\pi n` for some `n \in \ZZ`. 

  

EXAMPLES:: 

  

sage: RDF(pi).restrict_angle() 

3.141592653589793 

sage: RDF(pi + 1e-10).restrict_angle() 

-3.1415926534897936 

sage: RDF(1+10^10*pi).restrict_angle() 

0.9999977606... 

""" 

return self._new_c(gsl_sf_angle_restrict_symm(self._value)) 

  

def tan(self): 

""" 

Return the tangent of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/3 

sage: q.tan() 

1.7320508075688767 

sage: q = RDF.pi()/6 

sage: q.tan() 

0.5773502691896256 

""" 

cdef double denom 

cos = gsl_sf_cos(self._value) 

a = self._new_c(gsl_sf_sin(self._value) / cos) 

return a 

  

def sincos(self): 

""" 

Return a pair consisting of the sine and cosine of ``self``. 

  

EXAMPLES:: 

  

sage: t = RDF.pi()/6 

sage: t.sincos() 

(0.49999999999999994, 0.8660254037844387) 

""" 

return self.sin(), self.cos() 

  

def hypot(self, other): 

r""" 

Computes the value `\sqrt{s^2 + o^2}` where `s` is ``self`` and `o` 

is ``other`` in such a way as to avoid overflow. 

  

EXAMPLES:: 

  

sage: x = RDF(4e300); y = RDF(3e300); 

sage: x.hypot(y) 

5e+300 

sage: sqrt(x^2+y^2) # overflow 

+infinity 

""" 

sig_on() 

a = self._new_c(gsl_sf_hypot(self._value, float(other))) 

sig_off() 

return a 

  

def arccos(self): 

""" 

Return the inverse cosine of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/3 

sage: i = q.cos() 

sage: i.arccos() == q 

True 

""" 

return self._new_c(libc.math.acos(self._value)) 

  

def arcsin(self): 

""" 

Return the inverse sine of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/5 

sage: i = q.sin() 

sage: i.arcsin() == q 

True 

""" 

return self._new_c(libc.math.asin(self._value)) 

  

def arctan(self): 

""" 

Return the inverse tangent of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/5 

sage: i = q.tan() 

sage: i.arctan() == q 

True 

""" 

return self._new_c(libc.math.atan(self._value)) 

  

  

def cosh(self): 

""" 

Return the hyperbolic cosine of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/12 

sage: q.cosh() 

1.0344656400955106 

""" 

return self._new_c(gsl_ldexp( gsl_sf_exp(self._value) + gsl_sf_exp(-self._value), -1)) # (e^x + e^-x)/2 

  

def sinh(self): 

""" 

Return the hyperbolic sine of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/12 

sage: q.sinh() 

0.26480022760227073 

""" 

return self._new_c(gsl_ldexp( gsl_sf_expm1(self._value) - gsl_sf_expm1(-self._value), -1)) # (e^x - e^-x)/2 

  

def tanh(self): 

""" 

Return the hyperbolic tangent of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/12 

sage: q.tanh() 

0.25597778924568454 

""" 

return self.sinh() / self.cosh() 

  

def acosh(self): 

""" 

Return the hyperbolic inverse cosine of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/2 

sage: i = q.cosh(); i 

2.5091784786580567 

sage: abs(i.acosh()-q) < 1e-15 

True 

""" 

return self._new_c(gsl_acosh(self._value)) 

  

def arcsinh(self): 

""" 

Return the hyperbolic inverse sine of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/2 

sage: i = q.sinh(); i 

2.3012989023072947 

sage: abs(i.arcsinh()-q) < 1e-15 

True 

""" 

return self._new_c(gsl_asinh(self._value)) 

  

def arctanh(self): 

""" 

Return the hyperbolic inverse tangent of ``self``. 

  

EXAMPLES:: 

  

sage: q = RDF.pi()/2 

sage: i = q.tanh(); i 

0.9171523356672744 

sage: i.arctanh() - q # rel tol 1 

4.440892098500626e-16 

""" 

return self._new_c(gsl_atanh(self._value)) 

  

def sech(self): 

r""" 

Return the hyperbolic secant of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(pi).sech() 

0.08626673833405443 

sage: CDF(pi).sech() 

0.08626673833405443 

""" 

return 1/self.cosh() 

  

def csch(self): 

r""" 

Return the hyperbolic cosecant of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(pi).csch() 

0.08658953753004694 

sage: CDF(pi).csch() # rel tol 1e-15 

0.08658953753004696 

""" 

return 1/self.sinh() 

  

def coth(self): 

r""" 

Return the hyperbolic cotangent of ``self``. 

  

EXAMPLES:: 

  

sage: RDF(pi).coth() 

1.003741873197321 

sage: CDF(pi).coth() 

1.0037418731973213 

""" 

return self.cosh() / self.sinh() 

  

def agm(self, other): 

r""" 

Return the arithmetic-geometric mean of ``self`` and ``other``. The 

arithmetic-geometric mean is the common limit of the sequences 

`u_n` and `v_n`, where `u_0` is ``self``, 

`v_0` is other, `u_{n+1}` is the arithmetic mean 

of `u_n` and `v_n`, and `v_{n+1}` is the 

geometric mean of `u_n` and `v_n`. If any operand is negative, the 

return value is ``NaN``. 

  

EXAMPLES:: 

  

sage: a = RDF(1.5) 

sage: b = RDF(2.3) 

sage: a.agm(b) 

1.8786484558146697 

  

The arithmetic-geometric mean always lies between the geometric and 

arithmetic mean:: 

  

sage: sqrt(a*b) < a.agm(b) < (a+b)/2 

True 

""" 

cdef double a = self._value 

cdef double b = other 

cdef double eps = 2.0**-51 

if a < 0 or b < 0: 

return self._parent.nan() 

while True: 

a1 = (a+b)/2 

b1 = libc.math.sqrt(a*b) 

if abs((b1/a1)-1) < eps: return self._new_c(a1) 

a, b = a1, b1 

  

def erf(self): 

""" 

Return the value of the error function on ``self``. 

  

EXAMPLES:: 

  

sage: RDF(6).erf() 

1.0 

""" 

return self._new_c(gsl_sf_erf(self._value)) 

  

def gamma(self): 

""" 

Return the value of the Euler gamma function on ``self``. 

  

EXAMPLES:: 

  

sage: RDF(6).gamma() 

120.0 

sage: RDF(1.5).gamma() # rel tol 1e-15 

0.8862269254527584 

""" 

sig_on() 

a = self._new_c(gsl_sf_gamma(self._value)) 

sig_off() 

return a 

  

def zeta(self): 

r""" 

Return the Riemann zeta function evaluated at this real number. 

  

.. NOTE:: 

  

PARI is vastly more efficient at computing the Riemann zeta 

function. See the example below for how to use it. 

  

EXAMPLES:: 

  

sage: RDF(2).zeta() # rel tol 1e-15 

1.6449340668482269 

sage: RDF.pi()^2/6 

1.6449340668482264 

sage: RDF(-2).zeta() 

0.0 

sage: RDF(1).zeta() 

+infinity 

""" 

if self._value == 1: 

return self._new_c(1)/self._new_c(0) 

return self._new_c(gsl_sf_zeta(self._value)) 

  

def algebraic_dependency(self, n): 

""" 

Return a polynomial of degree at most `n` which is 

approximately satisfied by this number. 

  

.. NOTE:: 

  

The resulting polynomial need not be irreducible, and indeed 

usually won't be if this number is a good approximation to an 

algebraic number of degree less than `n`. 

  

ALGORITHM: 

  

Uses the PARI C-library ``algdep`` command. 

  

EXAMPLES:: 

  

sage: r = sqrt(RDF(2)); r 

1.4142135623730951 

sage: r.algebraic_dependency(5) 

x^2 - 2 

""" 

return sage.arith.all.algdep(self,n) 

  

algdep = algebraic_dependency 

  

cdef class ToRDF(Morphism): 

def __init__(self, R): 

""" 

Fast morphism from anything with a ``__float__`` method to an ``RDF`` 

element. 

  

EXAMPLES:: 

  

sage: f = RDF.coerce_map_from(ZZ); f 

Native morphism: 

From: Integer Ring 

To: Real Double Field 

sage: f(4) 

4.0 

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

Native morphism: 

From: Rational Field 

To: Real Double Field 

sage: f(1/2) 

0.5 

sage: f = RDF.coerce_map_from(int); f 

Native morphism: 

From: Set of Python objects of class 'int' 

To: Real Double Field 

sage: f(3r) 

3.0 

sage: f = RDF.coerce_map_from(float); f 

Native morphism: 

From: Set of Python objects of class 'float' 

To: Real Double Field 

sage: f(3.5) 

3.5 

""" 

from sage.categories.homset import Hom 

if isinstance(R, type): 

from sage.structure.parent import Set_PythonType 

R = Set_PythonType(R) 

Morphism.__init__(self, Hom(R, RDF)) 

  

cpdef Element _call_(self, x): 

""" 

Send ``x`` to the image under this map. 

  

EXAMPLES:: 

  

sage: f = RDF.coerce_map_from(float) 

sage: f(3.5) # indirect doctest 

3.5 

""" 

cdef RealDoubleElement r = <RealDoubleElement>PY_NEW(RealDoubleElement) 

r._value = PyFloat_AsDouble(x) 

return r 

  

def _repr_type(self): 

""" 

Return the representation type of ``self``. 

  

EXAMPLES:: 

  

sage: RDF.coerce_map_from(float)._repr_type() 

'Native' 

""" 

return "Native" 

  

  

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

# unique objects 

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

cdef RealDoubleField_class _RDF 

_RDF = RealDoubleField_class() 

  

RDF = _RDF # external interface 

  

def RealDoubleField(): 

""" 

Return the unique instance of the 

:class:`real double field<RealDoubleField_class>`. 

  

EXAMPLES:: 

  

sage: RealDoubleField() is RealDoubleField() 

True 

""" 

global _RDF 

return _RDF 

  

def is_RealDoubleElement(x): 

""" 

Check if ``x`` is an element of the real double field. 

  

EXAMPLES:: 

  

sage: from sage.rings.real_double import is_RealDoubleElement 

sage: is_RealDoubleElement(RDF(3)) 

True 

sage: is_RealDoubleElement(RIF(3)) 

False 

""" 

return isinstance(x, RealDoubleElement) 

  

  

################# FAST CREATION CODE ###################### 

########### Based on fast integer creation code ######### 

######## There is nothing to see here, move along ####### 

  

# We use a global element to steal all the references 

# from. DO NOT INITIALIZE IT AGAIN and DO NOT REFERENCE IT! 

cdef RealDoubleElement global_dummy_element 

global_dummy_element = RealDoubleElement(0) 

  

# A global pool for performance when elements are rapidly created and destroyed. 

# It operates on the following principles: 

# 

# - The pool starts out empty. 

# - When an new element is needed, one from the pool is returned 

# if available, otherwise a new RealDoubleElement object is created 

# - When an element is collected, it will add it to the pool 

# if there is room, otherwise it will be deallocated. 

DEF element_pool_size = 50 

  

cdef PyObject* element_pool[element_pool_size] 

cdef int element_pool_count = 0 

  

# used for profiling the pool 

cdef int total_alloc = 0 

cdef int use_pool = 0 

  

  

cdef PyObject* fast_tp_new(type t, args, kwds): 

global element_pool, element_pool_count, total_alloc, use_pool 

  

cdef PyObject* new 

  

# for profiling pool usage 

# total_alloc += 1 

  

# If there is a ready real double in the pool, we will 

# decrement the counter and return that. 

  

if element_pool_count > 0: 

  

# for profiling pool usage 

# use_pool += 1 

  

element_pool_count -= 1 

new = <PyObject *> element_pool[element_pool_count] 

  

# Otherwise, we have to create one. 

  

else: 

  

# allocate enough room for the RealDoubleElement, 

# sizeof_RealDoubleElement is sizeof(RealDoubleElement). 

# The use of PyObject_Malloc directly assumes 

# that RealDoubleElements are not garbage collected, i.e. 

# they do not possess references to other Python 

# objects (As indicated by the Py_TPFLAGS_HAVE_GC flag). 

# See below for a more detailed description. 

  

new = <PyObject*>PyObject_Malloc( sizeof(RealDoubleElement) ) 

  

# Now set every member as set in z, the global dummy RealDoubleElement 

# created before this tp_new started to operate. 

  

memcpy(new, (<void*>global_dummy_element), sizeof(RealDoubleElement) ) 

  

# This line is only needed if Python is compiled in debugging mode 

# './configure --with-pydebug' or SAGE_DEBUG=yes. If that is the 

# case a Python object has a bunch of debugging fields which are 

# initialized with this macro. 

if_Py_TRACE_REFS_then_PyObject_INIT( 

new, Py_TYPE(global_dummy_element)) 

  

# The global_dummy_element may have a reference count larger than 

# one, but it is expected that newly created objects have a 

# reference count of one. This is potentially unneeded if 

# everybody plays nice, because the gobal_dummy_element has only 

# one reference in that case. 

  

# Objects from the pool have reference count zero, so this 

# needs to be set in this case. 

  

new.ob_refcnt = 1 

  

return new 

  

cdef void fast_tp_dealloc(PyObject* o): 

  

# If there is room in the pool for a used integer object, 

# then put it in rather than deallocating it. 

  

global element_pool, element_pool_count 

  

if element_pool_count < element_pool_size: 

  

# And add it to the pool. 

element_pool[element_pool_count] = o 

element_pool_count += 1 

return 

  

# Free the object. This assumes that Py_TPFLAGS_HAVE_GC is not 

# set. If it was set another free function would need to be 

# called. 

  

PyObject_Free(o) 

  

  

from sage.misc.allocator cimport hook_tp_functions 

hook_tp_functions(global_dummy_element, <newfunc>(&fast_tp_new), <destructor>(&fast_tp_dealloc), False) 

  

  

cdef double_repr(double x): 

""" 

Convert a double to a string with maximum precision. 

""" 

if gsl_finite(x): 

return repr(x) 

cdef int v = gsl_isinf(x) 

if v > 0: 

return "+infinity" 

if v < 0: 

return "-infinity" 

return "NaN" 

  

cdef double_str(double x): 

""" 

Convert a double to an informal string. 

""" 

if gsl_finite(x): 

return str(x) 

return double_repr(x)