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

2915

2916

2917

2918

2919

2920

2921

2922

2923

2924

2925

2926

2927

2928

2929

2930

2931

2932

2933

2934

2935

2936

2937

2938

2939

2940

2941

2942

2943

2944

2945

2946

2947

2948

2949

2950

2951

2952

2953

2954

2955

2956

2957

2958

2959

2960

2961

2962

2963

2964

2965

2966

2967

2968

2969

2970

2971

2972

2973

2974

2975

2976

2977

2978

2979

2980

2981

2982

2983

2984

2985

2986

2987

2988

2989

2990

2991

2992

2993

2994

2995

2996

2997

2998

2999

3000

3001

3002

3003

3004

3005

3006

3007

3008

3009

3010

3011

3012

3013

3014

3015

3016

3017

3018

3019

3020

3021

3022

3023

3024

3025

3026

3027

3028

3029

3030

3031

3032

3033

3034

3035

3036

3037

3038

3039

3040

3041

3042

3043

3044

3045

3046

3047

3048

3049

3050

3051

3052

3053

3054

3055

3056

3057

3058

3059

3060

3061

3062

3063

3064

3065

3066

3067

3068

3069

3070

3071

3072

3073

3074

3075

3076

3077

3078

3079

3080

3081

3082

3083

3084

3085

3086

3087

3088

3089

3090

3091

3092

3093

3094

3095

3096

3097

3098

3099

3100

3101

3102

3103

3104

3105

3106

3107

3108

3109

3110

3111

3112

3113

3114

3115

3116

3117

3118

3119

3120

3121

3122

3123

3124

3125

3126

3127

3128

3129

3130

3131

3132

3133

3134

3135

3136

3137

3138

3139

3140

3141

3142

3143

3144

3145

3146

3147

3148

3149

3150

3151

3152

3153

3154

3155

3156

3157

3158

3159

3160

3161

3162

3163

3164

3165

3166

3167

3168

3169

3170

3171

3172

3173

3174

3175

3176

3177

3178

3179

3180

3181

3182

3183

3184

3185

3186

3187

3188

3189

3190

3191

3192

3193

3194

3195

3196

3197

3198

3199

3200

3201

3202

3203

3204

3205

3206

3207

3208

3209

3210

3211

3212

3213

3214

3215

3216

3217

3218

3219

3220

3221

3222

3223

3224

3225

3226

3227

3228

3229

3230

3231

3232

3233

3234

3235

3236

3237

3238

3239

3240

3241

3242

3243

3244

3245

3246

3247

3248

3249

3250

3251

3252

3253

3254

3255

3256

3257

3258

3259

3260

3261

3262

3263

3264

3265

3266

3267

3268

3269

3270

3271

3272

3273

3274

3275

3276

3277

3278

3279

3280

3281

3282

3283

3284

3285

3286

3287

3288

3289

3290

3291

3292

3293

3294

3295

3296

3297

3298

3299

3300

3301

3302

3303

3304

3305

3306

3307

3308

3309

3310

3311

3312

3313

3314

3315

3316

3317

3318

3319

3320

3321

3322

3323

3324

3325

3326

3327

3328

3329

3330

3331

3332

3333

3334

3335

3336

3337

3338

3339

3340

3341

3342

3343

3344

3345

3346

3347

3348

3349

3350

3351

3352

3353

3354

3355

3356

3357

3358

3359

3360

3361

3362

3363

3364

3365

3366

3367

3368

3369

3370

3371

3372

3373

3374

3375

3376

3377

3378

3379

3380

3381

3382

3383

3384

3385

3386

3387

3388

3389

3390

3391

3392

3393

3394

3395

3396

3397

3398

3399

3400

3401

3402

3403

3404

3405

3406

3407

3408

3409

3410

3411

3412

3413

3414

3415

3416

3417

3418

3419

3420

3421

3422

3423

3424

3425

3426

3427

3428

3429

3430

3431

3432

3433

3434

3435

3436

3437

3438

3439

3440

3441

3442

3443

3444

3445

3446

3447

3448

3449

3450

3451

3452

3453

3454

3455

3456

3457

3458

3459

3460

3461

3462

3463

3464

3465

3466

3467

3468

3469

3470

3471

3472

3473

3474

3475

3476

3477

3478

3479

3480

3481

3482

3483

3484

3485

3486

3487

3488

3489

3490

3491

3492

3493

3494

3495

3496

3497

3498

3499

3500

3501

3502

3503

3504

3505

3506

3507

3508

3509

3510

3511

3512

3513

3514

3515

3516

3517

3518

3519

3520

3521

3522

3523

3524

3525

3526

3527

3528

3529

3530

3531

3532

3533

3534

3535

3536

3537

3538

3539

3540

3541

3542

3543

3544

3545

3546

3547

3548

3549

3550

3551

3552

3553

3554

3555

3556

3557

3558

3559

3560

3561

3562

3563

3564

3565

3566

3567

3568

3569

3570

3571

3572

3573

3574

3575

3576

# -*- coding: utf-8 

r""" 

Arbitrary precision real balls using Arb 

  

This is a binding to the `Arb library <http://fredrikj.net/arb/>`_ for ball 

arithmetic. It may be useful to refer to its documentation for more details. 

  

Parts of the documentation for this module are copied or adapted from 

Arb's own documentation, licenced under the GNU General Public License 

version 2, or later. 

  

.. SEEALSO:: 

  

- :mod:`Complex balls using Arb <sage.rings.complex_arb>` 

- :mod:`Real intervals using MPFI <sage.rings.real_mpfi>` 

  

Data Structure 

============== 

  

Ball arithmetic, also known as mid-rad interval arithmetic, is an extension of 

floating-point arithmetic in which an error bound is attached to each variable. 

This allows doing rigorous computations over the real numbers, while avoiding 

the overhead of traditional (inf-sup) interval arithmetic at high precision, 

and eliminating much of the need for time-consuming and bug-prone manual error 

analysis associated with standard floating-point arithmetic. 

  

Sage :class:`RealBall` objects wrap Arb objects of type ``arb_t``. A real 

ball represents a ball over the real numbers, that is, an interval `[m-r,m+r]` 

where the midpoint `m` and the radius `r` are (extended) real numbers:: 

  

sage: RBF(pi) 

[3.141592653589793 +/- 5.61e-16] 

sage: RBF(pi).mid(), RBF(pi).rad() 

(3.14159265358979, 4.4408921e-16) 

  

The midpoint is represented as an arbitrary-precision floating-point number 

with arbitrary-precision exponent. The radius is a floating-point number with 

fixed-precision mantissa and arbitrary-precision exponent. :: 

  

sage: RBF(2)^(2^100) 

[2.285367694229514e+381600854690147056244358827360 +/- 2.98e+381600854690147056244358827344] 

  

:class:`RealBallField` objects (the parents of real balls) model the field of 

real numbers represented by balls on which computations are carried out with a 

certain precision:: 

  

sage: RBF 

Real ball field with 53 bits of precision 

  

It is possible to construct a ball whose parent is the real ball field with 

precision `p` but whose midpoint does not fit on `p` bits. However, the results 

of operations involving such a ball will (usually) be rounded to its parent's 

precision:: 

  

sage: RBF(factorial(50)).mid(), RBF(factorial(50)).rad() 

(3.0414093201713378043612608166064768844377641568961e64, 0.00000000) 

sage: (RBF(factorial(50)) + 0).mid() 

3.04140932017134e64 

  

Comparison 

========== 

  

.. WARNING:: 

  

In accordance with the semantics of Arb, identical :class:`RealBall` 

objects are understood to give permission for algebraic simplification. 

This assumption is made to improve performance. For example, setting ``z = 

x*x`` may set `z` to a ball enclosing the set `\{t^2 : t \in x\}` and not 

the (generally larger) set `\{tu : t \in x, u \in x\}`. 

  

Two elements are equal if and only if they are exact and equal (in spite of the 

above warning, inexact balls are not considered equal to themselves):: 

  

sage: a = RBF(1) 

sage: b = RBF(1) 

sage: a is b 

False 

sage: a == a 

True 

sage: a == b 

True 

  

:: 

  

sage: a = RBF(1/3) 

sage: b = RBF(1/3) 

sage: a.is_exact() 

False 

sage: b.is_exact() 

False 

sage: a is b 

False 

sage: a == a 

False 

sage: a == b 

False 

  

A ball is non-zero in the sense of comparison if and only if it does not 

contain zero. :: 

  

sage: a = RBF(RIF(-0.5, 0.5)) 

sage: a != 0 

False 

sage: b = RBF(1/3) 

sage: b != 0 

True 

  

However, ``bool(b)`` returns ``False`` for a ball ``b`` only if ``b`` is exactly 

zero:: 

  

sage: bool(a) 

True 

sage: bool(b) 

True 

sage: bool(RBF.zero()) 

False 

  

A ball ``left`` is less than a ball ``right`` if all elements of 

``left`` are less than all elements of ``right``. :: 

  

sage: a = RBF(RIF(1, 2)) 

sage: b = RBF(RIF(3, 4)) 

sage: a < b 

True 

sage: a <= b 

True 

sage: a > b 

False 

sage: a >= b 

False 

sage: a = RBF(RIF(1, 3)) 

sage: b = RBF(RIF(2, 4)) 

sage: a < b 

False 

sage: a <= b 

False 

sage: a > b 

False 

sage: a >= b 

False 

  

Comparisons with Sage symbolic infinities work with some limitations:: 

  

sage: -infinity < RBF(1) < +infinity 

True 

sage: -infinity < RBF(infinity) 

True 

sage: RBF(infinity) < infinity 

False 

sage: RBF(NaN) < infinity 

Traceback (most recent call last): 

... 

ValueError: infinite but not with +/- phase 

sage: 1/RBF(0) <= infinity 

Traceback (most recent call last): 

... 

ValueError: infinite but not with +/- phase 

  

Comparisons between elements of real ball fields, however, support special 

values and should be preferred:: 

  

sage: RBF(NaN) < RBF(infinity) 

False 

sage: 1/RBF(0) <= RBF(infinity) 

True 

  

TESTS:: 

  

sage: (RBF(pi) * identity_matrix(QQ, 3)).parent() 

Full MatrixSpace of 3 by 3 dense matrices over Real ball field 

with 53 bits of precision 

  

sage: polygen(RBF, 'x')^3 

x^3 

  

:: 

  

sage: SR.coerce(RBF(0.42)) 

[0.4200000000000000 +/- 1.56e-17] 

sage: RBF(0.42) + SR(1) 

[1.420000000000000 +/- 2.94e-16] 

sage: _.parent() 

Symbolic Ring 

  

Classes and Methods 

=================== 

""" 

  

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

# Copyright (C) 2014 Clemens Heuberger <clemens.heuberger@aau.at> 

# 2017 Vincent Delecroix <20100.delecroix@gmail.com> 

# 

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

# it under the terms of the GNU General Public License as published by 

# the Free Software Foundation, either version 2 of the License, or 

# (at your option) any later version. 

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

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

from __future__ import absolute_import 

  

from cysignals.signals cimport sig_on, sig_str, sig_off 

  

from cpython.float cimport PyFloat_AS_DOUBLE 

from cpython.int cimport PyInt_AS_LONG 

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

from libc.stdlib cimport abort 

  

from sage.libs.arb.arb cimport * 

from sage.libs.arb.arf cimport * 

from sage.libs.arb.arf cimport * 

from sage.libs.arb.mag cimport * 

from sage.libs.flint.flint cimport flint_free 

from sage.libs.flint.fmpz cimport * 

from sage.libs.flint.fmpq cimport * 

from sage.libs.gmp.mpz cimport * 

from sage.libs.mpfi cimport * 

from sage.libs.mpfr cimport * 

from sage.libs.mpfr cimport MPFR_RNDN, MPFR_RNDU, MPFR_RNDD, MPFR_RNDZ 

  

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

from sage.rings.ring cimport Field 

from sage.rings.integer cimport Integer 

from sage.rings.rational cimport Rational 

from sage.rings.real_double cimport RealDoubleElement 

from sage.rings.real_mpfr cimport RealField_class, RealField, RealNumber 

  

import operator 

  

import sage.categories.fields 

  

from sage.rings.integer_ring import ZZ 

from sage.rings.rational_field import QQ 

from sage.rings.real_mpfi import RealIntervalField, RealIntervalField_class 

from sage.structure.unique_representation import UniqueRepresentation 

from sage.cpython.string cimport char_to_str 

  

cdef void mpfi_to_arb(arb_t target, const mpfi_t source, const long precision): 

""" 

Convert an MPFI interval to an Arb ball. 

  

INPUT: 

  

- ``target`` -- an ``arb_t``. 

  

- ``source`` -- an ``mpfi_t``. 

  

- ``precision`` -- an integer `\ge 2`. 

  

TESTS:: 

  

sage: RBF(RIF(infinity)).endpoints() 

(+infinity, +infinity) 

sage: RBF(RIF(-infinity)).endpoints() 

(-infinity, -infinity) 

sage: RIF(RBF(infinity)).endpoints() 

(+infinity, +infinity) 

sage: RIF(RBF(-infinity)).endpoints() 

(-infinity, -infinity) 

""" 

cdef mpfr_t left 

cdef mpfr_t right 

  

mpfr_init2(left, precision) 

mpfr_init2(right, precision) 

  

if _do_sig(precision): sig_on() 

mpfi_get_left(left, source) 

mpfi_get_right(right, source) 

arb_set_interval_mpfr(target, left, right, precision) 

# Work around weakness of arb_set_interval_mpfr(tgt, inf, inf) 

if mpfr_equal_p(left, right): 

mag_zero(arb_radref(target)) 

if _do_sig(precision): sig_off() 

  

mpfr_clear(left) 

mpfr_clear(right) 

  

cdef int arb_to_mpfi(mpfi_t target, arb_t source, const long precision) except -1: 

""" 

Convert an Arb ball to an MPFI interval. 

  

INPUT: 

  

- ``target`` -- an ``mpfi_t``. 

  

- ``source`` -- an ``arb_t``. 

  

- ``precision`` -- an integer `\ge 2`. 

  

EXAMPLES:: 

  

sage: RIF(RBF(2)**(2**100)) # indirect doctest 

Traceback (most recent call last): 

... 

ArithmeticError: Error converting arb to mpfi. Overflow? 

""" 

cdef mpfr_t left 

cdef mpfr_t right 

  

mpfr_init2(left, precision) 

mpfr_init2(right, precision) 

  

try: 

sig_on() 

arb_get_interval_mpfr(left, right, source) 

mpfi_interv_fr(target, left, right) 

sig_off() 

except RuntimeError: 

raise ArithmeticError("Error converting arb to mpfi. Overflow?") 

finally: 

mpfr_clear(left) 

mpfr_clear(right) 

  

class RealBallField(UniqueRepresentation, Field): 

r""" 

An approximation of the field of real numbers using mid-rad intervals, also 

known as balls. 

  

INPUT: 

  

- ``precision`` -- an integer `\ge 2`. 

  

EXAMPLES:: 

  

sage: RBF = RealBallField() # indirect doctest 

sage: RBF(1) 

1.000000000000000 

  

:: 

  

sage: (1/2*RBF(1)) + AA(sqrt(2)) - 1 + polygen(QQ, 'x') 

x + [0.914213562373095 +/- 4.10e-16] 

  

TESTS:: 

  

sage: RBF.bracket(RBF(1/2), RBF(1/3)) 

[+/- 5.56e-17] 

sage: RBF.cardinality() 

+Infinity 

sage: RBF.cartesian_product(QQ).an_element()**2 

([1.440000000000000 +/- 4.98e-16], 1/4) 

sage: RBF.coerce_embedding() is None 

True 

sage: RBF['x'].gens_dict_recursive() 

{'x': x} 

sage: RBF.is_finite() 

False 

sage: RBF.is_zero() 

False 

sage: RBF.one() 

1.000000000000000 

sage: RBF.zero() 

0 

  

sage: NF.<sqrt2> = QuadraticField(2, embedding=AA(2).sqrt()) 

sage: a = (sqrt2 - 1)^1000 

sage: RBF(a) 

[1.676156872756536e-383 +/- 4.39e-399] 

""" 

Element = RealBall 

  

@staticmethod 

def __classcall__(cls, long precision=53): 

r""" 

Normalize the arguments for caching. 

  

TESTS:: 

  

sage: RealBallField(53) is RealBallField() is RBF 

True 

""" 

return super(RealBallField, cls).__classcall__(cls, precision) 

  

def __init__(self, long precision=53): 

r""" 

Initialize the real ball field. 

  

INPUT: 

  

- ``precision`` -- an integer `\ge 2`. 

  

EXAMPLES:: 

  

sage: RBF = RealBallField() 

sage: RBF(1) 

1.000000000000000 

sage: RealBallField(0) 

Traceback (most recent call last): 

... 

ValueError: precision must be at least 2 

sage: RealBallField(1) 

Traceback (most recent call last): 

... 

ValueError: precision must be at least 2 

  

TESTS:: 

  

sage: RBF.base() 

Real ball field with 53 bits of precision 

sage: RBF.base_ring() 

Real ball field with 53 bits of precision 

  

""" 

if precision < 2: 

raise ValueError("precision must be at least 2") 

Field.__init__(self, 

base_ring=self, 

category=sage.categories.fields.Fields().Infinite()) 

self._prec = precision 

from sage.rings.real_lazy import RLF 

self._populate_coercion_lists_([ZZ, QQ], convert_method_name='_arb_') 

  

def _repr_(self): 

r""" 

String representation of ``self``. 

  

EXAMPLES:: 

  

sage: RealBallField() 

Real ball field with 53 bits of precision 

sage: RealBallField(106) 

Real ball field with 106 bits of precision 

""" 

return "Real ball field with {} bits of precision".format(self._prec) 

  

def _coerce_map_from_(self, other): 

r""" 

Parents that canonically coerce into real ball fields include: 

  

- some exact or lazy parents representing subsets of the reals, such as 

``ZZ``, ``QQ``, ``AA``, and ``RLF``; 

  

- real ball fields with a larger precision. 

  

TESTS:: 

  

sage: RealBallField().has_coerce_map_from(RealBallField(54)) 

True 

sage: RealBallField().has_coerce_map_from(RealBallField(52)) 

False 

sage: RealBallField().has_coerce_map_from(RIF) 

False 

sage: RealBallField().has_coerce_map_from(SR) 

False 

sage: RealBallField().has_coerce_map_from(RR) 

False 

sage: K = QuadraticField(2, embedding=AA(2).sqrt()) 

sage: RBF.has_coerce_map_from(K) 

True 

sage: RBF.has_coerce_map_from(QuadraticField(2, embedding=None)) 

False 

sage: RBF.has_coerce_map_from(QuadraticField(-2)) 

False 

  

Check that the map goes through the ``_arb_`` method:: 

  

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

Conversion via _arb_ method map: 

... 

sage: RBF.convert_map_from(QuadraticField(2)) 

Conversion via _arb_ method map: 

... 

""" 

if isinstance(other, RealBallField): 

return other._prec >= self._prec 

  

from sage.rings.qqbar import AA 

from sage.rings.real_lazy import RLF 

if other in [AA, RLF]: 

return True 

  

from sage.rings.number_field.number_field_base import is_NumberField 

if is_NumberField(other): 

emb = other.coerce_embedding() 

return emb is not None and self.has_coerce_map_from(emb.codomain()) 

  

def _element_constructor_(self, mid=None, rad=None): 

""" 

Convert ``mid`` to an element of this real ball field, perhaps 

non-canonically. 

  

In addition to the inputs supported by :meth:`RealBall.__init__`, 

anything that is convertible to a real interval can also be used to 

construct a real ball:: 

  

sage: RBF(RIF(0, 1)) # indirect doctest 

[+/- 1.01] 

sage: RBF(1) 

1.000000000000000 

sage: RBF(x) 

Traceback (most recent call last): 

... 

TypeError: unable to convert x to a RealBall 

  

Various symbolic constants can be converted without going through real 

intervals. (This is faster and yields tighter error bounds.) :: 

  

sage: RBF(e) 

[2.718281828459045 +/- 5.35e-16] 

sage: RBF(pi) 

[3.141592653589793 +/- 5.61e-16] 

  

Symbolic expressions are parsed :: 

  

sage: RBF(4*zeta(3)) 

[4.808227612638377 +/- 9.42e-16] 

""" 

try: 

return self.element_class(self, mid, rad) 

except TypeError: 

pass 

try: 

return self.element_class(self, mid.pyobject(), rad) 

except (AttributeError, TypeError): 

pass 

try: 

return mid.operator()(*[self(operand) for operand in mid.operands()]) 

except (AttributeError, TypeError): 

pass 

try: 

mid = RealIntervalField(self._prec)(mid) 

return self.element_class(self, mid, rad) 

except TypeError: 

pass 

raise TypeError("unable to convert {!r} to a RealBall".format(mid)) 

  

def _repr_option(self, key): 

""" 

Declare that real balls print atomically. 

  

TESTS:: 

  

sage: RBF._repr_option('element_is_atomic') 

True 

sage: RBF['x']([-2,-2,-2/3]) 

[-0.666666666666667 +/- 4.82e-16]*x^2 - 2.000000000000000*x 

- 2.000000000000000 

sage: RBF._repr_option('element_is_atomic_typo') 

Traceback (most recent call last): 

... 

KeyError: 'element_is_atomic_typo' 

""" 

if key == 'element_is_atomic': 

return True 

  

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

  

def gens(self): 

r""" 

EXAMPLES:: 

  

sage: RBF.gens() 

(1.000000000000000,) 

sage: RBF.gens_dict() 

{'1.000000000000000': 1.000000000000000} 

""" 

return (self.one(),) 

  

def construction(self): 

""" 

Return the construction of a real ball field as a completion of the 

rationals. 

  

EXAMPLES:: 

  

sage: RBF = RealBallField(42) 

sage: functor, base = RBF.construction() 

sage: functor, base 

(Completion[+Infinity, prec=42], Rational Field) 

sage: functor(base) is RBF 

True 

""" 

from sage.categories.pushout import CompletionFunctor 

functor = CompletionFunctor(sage.rings.infinity.Infinity, 

self._prec, 

{'type': 'Ball'}) 

return functor, QQ 

  

def complex_field(self): 

""" 

Return the complex ball field with the same precision. 

  

EXAMPLES:: 

  

sage: from sage.rings.complex_arb import ComplexBallField 

sage: RBF.complex_field() 

Complex ball field with 53 bits of precision 

sage: RealBallField(3).algebraic_closure() 

Complex ball field with 3 bits of precision 

""" 

from sage.rings.complex_arb import ComplexBallField 

return ComplexBallField(self._prec) 

  

algebraic_closure = complex_field 

  

def precision(self): 

""" 

Return the bit precision used for operations on elements of this field. 

  

EXAMPLES:: 

  

sage: RealBallField().precision() 

53 

""" 

return self._prec 

  

def is_exact(self): 

""" 

Real ball fields are not exact. 

  

EXAMPLES:: 

  

sage: RealBallField().is_exact() 

False 

""" 

return False 

  

def is_finite(self): 

""" 

Real ball fields are infinite. 

  

They already specify it via their category, but we currently need to 

re-implement this method due to the legacy implementation in 

:class:`sage.rings.ring.Ring`. 

  

EXAMPLES:: 

  

sage: RealBallField().is_finite() 

False 

""" 

return False 

  

def characteristic(self): 

""" 

Real ball fields have characteristic zero. 

  

EXAMPLES:: 

  

sage: RealBallField().characteristic() 

0 

""" 

return 0 

  

def some_elements(self): 

""" 

Real ball fields contain exact balls, inexact balls, infinities, and 

more. 

  

EXAMPLES:: 

  

sage: RBF.some_elements() 

[1.000000000000000, 

[0.3333333333333333 +/- 7.04e-17], 

[-4.733045976388941e+363922934236666733021124 +/- 3.46e+363922934236666733021108], 

[+/- inf], 

[+/- inf], 

nan] 

""" 

import sage.symbolic.constants 

return [self(1), self(1)/3, 

-self(2)**(Integer(2)**80), 

self(sage.rings.infinity.Infinity), ~self(0), 

self.element_class(self, sage.symbolic.constants.NotANumber())] 

  

def _sum_of_products(self, terms): 

r""" 

Compute a sum of product of real balls without creating temporary 

Python objects 

  

The input objects should be real balls, but need not belong to this 

parent. The computation is performed at the precision of this parent. 

  

EXAMPLES:: 

  

sage: Pol.<x> = RealBallField(1000)[] 

sage: pol = (x + 1/3)^100 

sage: RBF._sum_of_products((c, c) for c in pol) 

[6.3308767660842e+23 +/- 4.59e+9] 

  

TESTS:: 

  

sage: RBF._sum_of_products([]) 

0 

sage: RBF._sum_of_products([[]]) 

1.000000000000000 

sage: RBF._sum_of_products([["a"]]) 

Traceback (most recent call last): 

... 

TypeError: Cannot convert str to sage.rings.real_arb.RealBall 

""" 

cdef RealBall res = RealBall.__new__(RealBall) 

cdef RealBall factor 

cdef arb_t tmp 

res._parent = self 

arb_zero(res.value) 

arb_init(tmp) 

try: 

for term in terms: 

arb_one(tmp) 

for factor in term: 

arb_mul(tmp, tmp, factor.value, self._prec) 

arb_add(res.value, res.value, tmp, self._prec) 

finally: 

arb_clear(tmp) 

return res 

# Constants 

  

def pi(self): 

r""" 

Return a ball enclosing `\pi`. 

  

EXAMPLES:: 

  

sage: RBF.pi() 

[3.141592653589793 +/- 5.61e-16] 

sage: RealBallField(128).pi() 

[3.1415926535897932384626433832795028842 +/- 1.65e-38] 

""" 

cdef RealBall res = RealBall.__new__(RealBall) 

res._parent = self 

if _do_sig(self._prec): sig_on() 

arb_const_pi(res.value, self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

  

def log2(self): 

r""" 

Return a ball enclosing `\log(2)`. 

  

EXAMPLES:: 

  

sage: RBF.log2() 

[0.693147180559945 +/- 3.98e-16] 

sage: RealBallField(128).log2() 

[0.69314718055994530941723212145817656807 +/- 7.70e-39] 

""" 

cdef RealBall res = RealBall.__new__(RealBall) 

res._parent = self 

if _do_sig(self._prec): sig_on() 

arb_const_log2(res.value, self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

  

def euler_constant(self): 

r""" 

Return a ball enclosing the Euler constant. 

  

EXAMPLES:: 

  

sage: RBF.euler_constant() 

[0.577215664901533 +/- 3.57e-16] 

sage: RealBallField(128).euler_constant() 

[0.57721566490153286060651209008240243104 +/- 3.25e-39] 

""" 

cdef RealBall res = RealBall.__new__(RealBall) 

res._parent = self 

if _do_sig(self._prec): sig_on() 

arb_const_euler(res.value, self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

  

def catalan_constant(self): 

r""" 

Return a ball enclosing the Catalan constant. 

  

EXAMPLES:: 

  

sage: RBF.catalan_constant() 

[0.915965594177219 +/- 1.23e-16] 

sage: RealBallField(128).catalan_constant() 

[0.91596559417721901505460351493238411077 +/- 6.37e-39] 

""" 

cdef RealBall res = RealBall.__new__(RealBall) 

res._parent = self 

if _do_sig(self._prec): sig_on() 

arb_const_catalan(res.value, self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

  

# Ball functions of non-ball arguments 

  

def sinpi(self, x): 

""" 

Return a ball enclosing `\sin(\pi x)`. 

  

This works even if ``x`` itself is not a ball, and may be faster or 

more accurate where ``x`` is a rational number. 

  

EXAMPLES:: 

  

sage: RBF.sinpi(1) 

0 

sage: RBF.sinpi(1/3) 

[0.866025403784439 +/- 5.15e-16] 

sage: RBF.sinpi(1 + 2^(-100)) 

[-2.478279624546525e-30 +/- 5.90e-46] 

  

.. SEEALSO:: :meth:`~sage.rings.real_arb.RealBall.sin` 

  

TESTS:: 

  

sage: RBF.sinpi(RLF(sqrt(2))) 

[-0.96390253284988 +/- 4.11e-15] 

""" 

cdef RealBall res, x_as_ball 

cdef Rational x_as_Rational 

cdef fmpq_t tmpq 

res = self.element_class(self) 

try: 

x_as_Rational = QQ.coerce(x) 

try: 

if _do_sig(self._prec): sig_on() 

fmpq_init(tmpq) 

fmpq_set_mpq(tmpq, x_as_Rational.value) 

arb_sin_pi_fmpq(res.value, tmpq, self._prec) 

if _do_sig(self._prec): sig_off() 

finally: 

fmpq_clear(tmpq) 

return res 

except TypeError: 

pass 

x_as_ball = self.coerce(x) 

if _do_sig(self._prec): sig_on() 

arb_sin_pi(res.value, x_as_ball.value, self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

  

def cospi(self, x): 

""" 

Return a ball enclosing `\cos(\pi x)`. 

  

This works even if ``x`` itself is not a ball, and may be faster or 

more accurate where ``x`` is a rational number. 

  

EXAMPLES:: 

  

sage: RBF.cospi(1) 

-1.000000000000000 

sage: RBF.cospi(1/3) 

0.5000000000000000 

  

.. SEEALSO:: :meth:`~sage.rings.real_arb.RealBall.cos` 

  

TESTS:: 

  

sage: RBF.cospi(RLF(sqrt(2))) 

[-0.26625534204142 +/- 5.38e-15] 

""" 

cdef RealBall res, x_as_ball 

cdef Rational x_as_Rational 

cdef fmpq_t tmpq 

res = self.element_class(self) 

try: 

x_as_Rational = QQ.coerce(x) 

try: 

if _do_sig(self._prec): sig_on() 

fmpq_init(tmpq) 

fmpq_set_mpq(tmpq, x_as_Rational.value) 

arb_cos_pi_fmpq(res.value, tmpq, self._prec) 

if _do_sig(self._prec): sig_off() 

finally: 

fmpq_clear(tmpq) 

return res 

except TypeError: 

pass 

x_as_ball = self.coerce(x) 

if _do_sig(self._prec): sig_on() 

arb_cos_pi(res.value, x_as_ball.value, self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

  

def gamma(self, x): 

""" 

Return a ball enclosing the gamma function of ``x``. 

  

This works even if ``x`` itself is not a ball, and may be more 

efficient in the case where ``x`` is an integer or a rational number. 

  

EXAMPLES:: 

  

sage: RBF.gamma(5) 

24.00000000000000 

sage: RBF.gamma(10**20) 

[+/- 5.92e+1956570551809674821757] 

sage: RBF.gamma(1/3) 

[2.678938534707747 +/- 8.99e-16] 

sage: RBF.gamma(-5) 

nan 

  

.. SEEALSO:: :meth:`~sage.rings.real_arb.RealBall.gamma` 

  

TESTS:: 

  

sage: RBF.gamma(RLF(pi)) 

[2.2880377953400 +/- 4.29e-14] 

""" 

cdef RealBall res 

cdef Integer x_as_Integer 

cdef Rational x_as_Rational 

cdef fmpz_t tmpz 

cdef fmpq_t tmpq 

res = self.element_class(self) 

try: 

x_as_Integer = ZZ.coerce(x) 

try: 

if _do_sig(self._prec): sig_on() 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, x_as_Integer.value) 

arb_gamma_fmpz(res.value, tmpz, self._prec) 

if _do_sig(self._prec): sig_off() 

finally: 

fmpz_clear(tmpz) 

return res 

except TypeError: 

pass 

try: 

x_as_Rational = QQ.coerce(x) 

try: 

if _do_sig(self._prec): sig_on() 

fmpq_init(tmpq) 

fmpq_set_mpq(tmpq, x_as_Rational.value) 

arb_gamma_fmpq(res.value, tmpq, self._prec) 

if _do_sig(self._prec): sig_off() 

finally: 

fmpq_clear(tmpq) 

return res 

except TypeError: 

pass 

return self.coerce(x).gamma() 

  

def zeta(self, s): 

""" 

Return a ball enclosing the Riemann zeta function of ``s``. 

  

This works even if ``s`` itself is not a ball, and may be more 

efficient in the case where ``s`` is an integer. 

  

EXAMPLES:: 

  

sage: RBF.zeta(3) # abs tol 5e-16 

[1.202056903159594 +/- 2.87e-16] 

sage: RBF.zeta(1) 

nan 

sage: RBF.zeta(1/2) 

[-1.460354508809587 +/- 1.94e-16] 

  

.. SEEALSO:: :meth:`~sage.rings.real_arb.RealBall.zeta` 

""" 

cdef RealBall res 

cdef Integer s_as_Integer 

try: 

s_as_Integer = ZZ.coerce(s) 

if mpz_fits_ulong_p(s_as_Integer.value): 

res = self.element_class(self) 

if _do_sig(self._prec): sig_on() 

arb_zeta_ui(res.value, mpz_get_ui(s_as_Integer.value), self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

except TypeError: 

pass 

return self.coerce(s).zeta() 

  

def bernoulli(self, n): 

""" 

Return a ball enclosing the ``n``-th Bernoulli number. 

  

EXAMPLES:: 

  

sage: [RBF.bernoulli(n) for n in range(4)] 

[1.000000000000000, -0.5000000000000000, [0.1666666666666667 +/- 7.04e-17], 0] 

sage: RBF.bernoulli(2**20) 

[-1.823002872104961e+5020717 +/- 7.16e+5020701] 

sage: RBF.bernoulli(2**1000) 

Traceback (most recent call last): 

... 

ValueError: argument too large 

  

TESTS:: 

  

sage: RBF.bernoulli(2r) 

[0.1666666666666667 +/- 7.04e-17] 

sage: RBF.bernoulli(2/3) 

Traceback (most recent call last): 

... 

TypeError: no canonical coercion from Rational Field to Integer Ring 

sage: RBF.bernoulli(-1) 

Traceback (most recent call last): 

... 

ValueError: expected a nonnegative index 

""" 

cdef RealBall res 

cdef Integer n_as_Integer = ZZ.coerce(n) 

if mpz_fits_ulong_p(n_as_Integer.value): 

res = self.element_class(self) 

if _do_sig(self._prec): sig_on() 

arb_bernoulli_ui(res.value, mpz_get_ui(n_as_Integer.value), self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

elif n_as_Integer < 0: 

raise ValueError("expected a nonnegative index") 

else: 

# TODO: Fall back to a Sage implementation in this case? 

raise ValueError("argument too large") 

  

def fibonacci(self, n): 

""" 

Return a ball enclosing the ``n``-th Fibonacci number. 

  

EXAMPLES:: 

  

sage: [RBF.fibonacci(n) for n in range(7)] 

[0, 

1.000000000000000, 

1.000000000000000, 

2.000000000000000, 

3.000000000000000, 

5.000000000000000, 

8.000000000000000] 

sage: RBF.fibonacci(-2) 

-1.000000000000000 

sage: RBF.fibonacci(10**20) 

[3.78202087472056e+20898764024997873376 +/- 4.01e+20898764024997873361] 

""" 

cdef fmpz_t tmpz 

cdef RealBall res = self.element_class(self) 

cdef Integer n_as_Integer = ZZ.coerce(n) 

try: 

if _do_sig(self._prec): sig_on() 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, n_as_Integer.value) 

arb_fib_fmpz(res.value, tmpz, self._prec) 

if _do_sig(self._prec): sig_off() 

finally: 

fmpz_clear(tmpz) 

return res 

  

def bell_number(self, n): 

""" 

Return a ball enclosing the ``n``-th Bell number. 

  

EXAMPLES:: 

  

sage: [RBF.bell_number(n) for n in range(7)] 

[1.000000000000000, 

1.000000000000000, 

2.000000000000000, 

5.000000000000000, 

15.00000000000000, 

52.00000000000000, 

203.0000000000000] 

sage: RBF.bell_number(-1) 

Traceback (most recent call last): 

... 

ValueError: expected a nonnegative index 

sage: RBF.bell_number(10**20) 

[5.38270113176282e+1794956117137290721328 +/- 5.44e+1794956117137290721313] 

""" 

cdef fmpz_t tmpz 

cdef RealBall res = self.element_class(self) 

cdef Integer n_as_Integer = ZZ.coerce(n) 

if n_as_Integer < 0: 

raise ValueError("expected a nonnegative index") 

try: 

if _do_sig(self._prec): sig_on() 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, n_as_Integer.value) 

arb_bell_fmpz(res.value, tmpz, self._prec) 

if _do_sig(self._prec): sig_off() 

finally: 

fmpz_clear(tmpz) 

return res 

  

def double_factorial(self, n): 

""" 

Return a ball enclosing the ``n``-th double factorial. 

  

EXAMPLES:: 

  

sage: [RBF.double_factorial(n) for n in range(7)] 

[1.000000000000000, 

1.000000000000000, 

2.000000000000000, 

3.000000000000000, 

8.000000000000000, 

15.00000000000000, 

48.00000000000000] 

sage: RBF.double_factorial(2**20) 

[1.4483729903e+2928836 +/- 8.96e+2928825] 

sage: RBF.double_factorial(2**1000) 

Traceback (most recent call last): 

... 

ValueError: argument too large 

sage: RBF.double_factorial(-1) 

Traceback (most recent call last): 

... 

ValueError: expected a nonnegative index 

  

""" 

cdef RealBall res 

cdef Integer n_as_Integer = ZZ.coerce(n) 

if mpz_fits_ulong_p(n_as_Integer.value): 

res = self.element_class(self) 

if _do_sig(self._prec): sig_on() 

arb_doublefac_ui(res.value, mpz_get_ui(n_as_Integer.value), self._prec) 

if _do_sig(self._prec): sig_off() 

return res 

elif n_as_Integer < 0: 

raise ValueError("expected a nonnegative index") 

else: 

# TODO: Fall back to a Sage implementation in this case? 

raise ValueError("argument too large") 

  

def maximal_accuracy(self): 

r""" 

Return the relative accuracy of exact elements measured in bits. 

  

OUTPUT: 

  

An integer. 

  

EXAMPLES:: 

  

sage: RBF.maximal_accuracy() 

9223372036854775807 # 64-bit 

2147483647 # 32-bit 

  

.. SEEALSO:: :meth:`RealBall.accuracy` 

""" 

return ARF_PREC_EXACT 

  

cdef inline bint _do_sig(long prec): 

""" 

Whether signal handlers should be installed for calls to arb. 

  

TESTS:: 

  

sage: _ = RealBallField()(1).psi() # indirect doctest 

sage: _ = RealBallField(1500)(1).psi() 

""" 

return (prec > 1000) 

  

cdef inline long prec(RealBall ball): 

return ball._parent._prec 

  

cdef class RealBall(RingElement): 

""" 

Hold one ``arb_t`` of the `Arb library 

<http://fredrikj.net/arb/>`_ 

  

EXAMPLES:: 

  

sage: a = RealBallField()(RIF(1)) # indirect doctest 

sage: b = a.psi() 

sage: b 

[-0.577215664901533 +/- 3.85e-16] 

sage: RIF(b) 

-0.577215664901533? 

""" 

  

def __cinit__(self): 

""" 

Initialize the parent and allocate memory. 

  

EXAMPLES:: 

  

sage: RealBallField()(RIF(1)) # indirect doctest 

1.000000000000000 

""" 

arb_init(self.value) 

  

def __dealloc__(self): 

""" 

Deallocate memory of the encapsulated value. 

  

EXAMPLES:: 

  

sage: a = RealBallField()(RIF(1)) # indirect doctest 

sage: del a 

""" 

arb_clear(self.value) 

  

def __init__(self, parent, mid=None, rad=None): 

""" 

Initialize the :class:`RealBall`. 

  

INPUT: 

  

- ``parent`` -- a :class:`RealBallField`. 

  

- ``mid`` (optional) -- ball midpoint, see examples below. If omitted, 

initialize the ball to zero, ignoring the ``rad`` argument. 

  

- ``rad`` (optional) -- a :class:`RealNumber` or a Python float, ball 

radius. If the midpoint is not exactly representable in 

floating-point, the radius is adjusted to account for the roundoff 

error. 

  

EXAMPLES:: 

  

sage: RBF = RealBallField() 

sage: RBF() 

0 

  

One can create exact real balls using elements of various exact parents, 

or using floating-point numbers:: 

  

sage: RBF(3) 

3.000000000000000 

sage: RBF(3r) 

3.000000000000000 

sage: RBF(1/3) 

[0.3333333333333333 +/- 7.04e-17] 

sage: RBF(3.14) 

[3.140000000000000 +/- 1.25e-16] 

  

:: 

  

sage: RBF(3, 0.125) 

[3e+0 +/- 0.126] 

sage: RBF(pi, 0.125r) 

[3e+0 +/- 0.267] 

  

:: 

  

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

sage: RBF(1/5 + sqrt2/2) 

[0.907106781186547 +/- 7.06e-16] 

  

Note that integers and floating-point numbers are ''not'' rounded to 

the parent's precision:: 

  

sage: b = RBF(11111111111111111111111111111111111111111111111); b 

[1.111111111111111e+46 +/- 1.12e+30] 

sage: b.mid().exact_rational() 

11111111111111111111111111111111111111111111111 

  

Similarly, converting a real ball from one real ball field to another 

(with a different precision) only changes the way it is displayed and 

the precision of operations involving it, not the actual representation 

of its center:: 

  

sage: RBF100 = RealBallField(100) 

sage: b100 = RBF100(1/3); b100 

[0.333333333333333333333333333333 +/- 4.65e-31] 

sage: b53 = RBF(b100); b53 

[0.3333333333333333 +/- 3.34e-17] 

sage: RBF100(b53) 

[0.333333333333333333333333333333 +/- 4.65e-31] 

  

Special values are supported:: 

  

sage: RBF(oo).mid(), RBF(-oo).mid(), RBF(unsigned_infinity).mid() 

(+infinity, -infinity, 0.000000000000000) 

sage: RBF(NaN) 

nan 

  

Strings can be given as input. Strings must contain decimal 

floating-point literals. A valid string must consist of a midpoint, 

a midpoint and a radius separated by "+/-", or just a 

radius prefixed by "+/-". Optionally, the whole string can be enclosed 

in square brackets. In general, the string representation of a 

real ball as returned by ``str()`` can be parsed back (the result 

will be larger than the original ball if rounding occurs). 

A few examples:: 

  

sage: RBF("1.1") 

[1.100000000000000 +/- 3.56e-16] 

sage: RBF(str(RBF("1.1"))) 

[1.100000000000000 +/- 7.12e-16] 

sage: RBF("3.25") 

3.250000000000000 

sage: RBF("-3.1 +/- 1e-10") 

[-3.100000000 +/- 1.01e-10] 

sage: RBF("[+/-1]") 

[+/- 1.01] 

sage: RBF("inf +/- inf") 

[+/- inf] 

  

.. SEEALSO:: :meth:`RealBallField._element_constructor_` 

  

TESTS:: 

  

sage: from sage.rings.real_arb import RealBall 

sage: RealBall(RBF, sage.symbolic.constants.Pi()) # abs tol 1e-16 

[3.141592653589793 +/- 5.62e-16] 

sage: RealBall(RBF, sage.symbolic.constants.Log2()) # abs tol 1e-16 

[0.693147180559945 +/- 4.06e-16] 

sage: RealBall(RBF, sage.symbolic.constants.Catalan()) 

[0.915965594177219 +/- 1.23e-16] 

sage: RealBall(RBF, sage.symbolic.constants.Khinchin()) 

[2.685452001065306 +/- 6.82e-16] 

sage: RealBall(RBF, sage.symbolic.constants.Glaisher()) 

[1.282427129100623 +/- 6.02e-16] 

sage: RealBall(RBF, sage.symbolic.constants.e) 

[2.718281828459045 +/- 5.35e-16] 

sage: RealBall(RBF, sage.symbolic.constants.EulerGamma()) 

[0.577215664901533 +/- 3.57e-16] 

sage: RBF("1 +/- 0.001") 

[1.00 +/- 1.01e-3] 

sage: RBF("2.3e10000000000000000000000 +/- 0.00005e10000000000000000000000") 

[2.3000e+10000000000000000000000 +/- 5.01e+9999999999999999999995] 

sage: RBF("0.3 +/- 0.2 +/- 0.1") 

Traceback (most recent call last): 

... 

ValueError: unsupported string format 

  

sage: NF.<a> = QuadraticField(2, embedding=AA(2).sqrt()) 

sage: RBF.coerce(a) 

[1.414213562373095 +/- 2.99e-16] 

sage: NF.<a> = QuadraticField(2, embedding=-AA(2).sqrt()) 

sage: RBF.coerce(a) 

[-1.414213562373095 +/- 2.99e-16] 

sage: NF.<a> = QuadraticField(2, embedding=None) 

sage: RBF.coerce(a) 

Traceback (most recent call last): 

... 

TypeError: no canonical coercion ... 

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

sage: RBF(QQi(3)) 

3.000000000000000 

sage: RBF(i) 

Traceback (most recent call last): 

... 

ValueError: nonzero imaginary part 

sage: RBF.coerce(QQi(3)) 

Traceback (most recent call last): 

... 

TypeError: no canonical coercion... 

""" 

cdef fmpz_t tmpz 

cdef fmpq_t tmpq 

cdef arf_t tmpr 

cdef mag_t tmpm 

  

Element.__init__(self, parent) 

  

if mid is None: 

return 

  

elif isinstance(mid, RealBall): 

arb_set(self.value, (<RealBall> mid).value) # no rounding! 

elif isinstance(mid, int): 

arb_set_si(self.value, PyInt_AS_LONG(mid)) # no rounding! 

elif isinstance(mid, Integer): 

if _do_sig(prec(self)): sig_on() 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, (<Integer> mid).value) 

arb_set_fmpz(self.value, tmpz) # no rounding! 

fmpz_clear(tmpz) 

if _do_sig(prec(self)): sig_off() 

elif isinstance(mid, Rational): 

if _do_sig(prec(self)): sig_on() 

fmpq_init(tmpq) 

fmpq_set_mpq(tmpq, (<Rational> mid).value) 

arb_set_fmpq(self.value, tmpq, prec(self)) 

fmpq_clear(tmpq) 

if _do_sig(prec(self)): sig_off() 

elif isinstance(mid, RealNumber): 

if _do_sig(prec(self)): sig_on() 

arf_init(tmpr) 

arf_set_mpfr(tmpr, (<RealNumber> mid).value) 

arb_set_arf(self.value, tmpr) # no rounding! 

arf_clear(tmpr) 

if _do_sig(prec(self)): sig_off() 

elif isinstance(mid, RealIntervalFieldElement): 

mpfi_to_arb(self.value, 

(<RealIntervalFieldElement> mid).value, 

prec(self)) 

elif isinstance(mid, str): 

if arb_set_str(self.value, mid, prec(self)) != 0: 

raise ValueError("unsupported string format") 

else: 

# the initializers that trigger imports 

import sage.symbolic.constants 

if isinstance(mid, sage.rings.infinity.AnInfinity): 

if isinstance(mid, sage.rings.infinity.PlusInfinity): 

arb_pos_inf(self.value) 

elif isinstance(mid, sage.rings.infinity.MinusInfinity): 

arb_neg_inf(self.value) 

else: 

arb_zero_pm_inf(self.value) 

elif isinstance(mid, sage.symbolic.constants.Constant): 

if _do_sig(prec(self)): sig_on() 

try: 

if isinstance(mid, sage.symbolic.constants.NotANumber): 

arb_indeterminate(self.value) 

elif isinstance(mid, sage.symbolic.constants.Pi): 

arb_const_pi(self.value, prec(self)) 

elif isinstance(mid, sage.symbolic.constants.Log2): 

arb_const_log2(self.value, prec(self)) 

elif isinstance(mid, sage.symbolic.constants.Catalan): 

arb_const_catalan(self.value, prec(self)) 

elif isinstance(mid, sage.symbolic.constants.Khinchin): 

arb_const_khinchin(self.value, prec(self)) 

elif isinstance(mid, sage.symbolic.constants.Glaisher): 

arb_const_glaisher(self.value, prec(self)) 

elif isinstance(mid, sage.symbolic.constants.EulerGamma): 

arb_const_euler(self.value, prec(self)) 

else: 

raise TypeError("unsupported constant") 

finally: 

if _do_sig(prec(self)): sig_off() 

elif isinstance(mid, sage.symbolic.constants_c.E): 

if _do_sig(prec(self)): sig_on() 

arb_const_e(self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

else: 

raise TypeError("unsupported midpoint type") 

  

if rad is not None: 

mag_init(tmpm) 

if isinstance(rad, RealNumber): 

arf_init(tmpr) 

arf_set_mpfr(tmpr, (<RealNumber> rad).value) 

arf_get_mag(tmpm, tmpr) 

arf_clear(tmpr) 

elif isinstance(rad, float): 

mag_set_d(tmpm, PyFloat_AS_DOUBLE(rad)) 

else: 

raise TypeError("rad should be a RealNumber or a Python float") 

mag_add(arb_radref(self.value), arb_radref(self.value), tmpm) 

mag_clear(tmpm) 

  

def __hash__(self): 

""" 

TESTS:: 

  

sage: hash(RealBallField(10)(1)) == hash(RealBallField(20)(1)) 

True 

sage: hash(RBF(1/3)) == hash(RBF(1/3, rad=.1)) 

False 

sage: vals = [0, 1, 3/4, 5/8, 7/8, infinity, 'nan', 2^1000 - 1] 

sage: len({hash(RBF(v)) for v in vals}) == len(vals) 

True 

""" 

cdef arf_t mid = arb_midref(self.value) 

cdef fmpz_t mant, expo 

fmpz_init(mant) 

fmpz_init(expo) 

arf_get_fmpz_2exp(mant, expo, mid) 

cdef long h = ( 

fmpz_fdiv_ui(mant, 1073741789) 

^ fmpz_fdiv_ui(expo, 2**30) 

^ (arf_abs_bound_lt_2exp_si(mid) << 10) 

^ arb_rel_error_bits(self.value) << 20) 

fmpz_clear(expo) 

fmpz_clear(mant) 

return h 

  

def _repr_(self): 

""" 

Return a string representation of ``self``. 

  

OUTPUT: 

  

A string. 

  

EXAMPLES:: 

  

sage: RealBallField()(RIF(1.9, 2)) 

[2e+0 +/- 0.101] 

""" 

cdef char* c_result 

  

c_result = arb_get_str(self.value, (prec(self) * 31) // 100, 0) 

try: 

py_string = char_to_str(c_result) 

finally: 

flint_free(c_result) 

  

return py_string 

  

# Conversions 

  

cpdef RealIntervalFieldElement _real_mpfi_(self, RealIntervalField_class parent): 

""" 

Return a :mod:`real interval <sage.rings.real_mpfi>` containing this ball. 

  

OUTPUT: 

  

A :class:`~sage.rings.real_mpfi.RealIntervalFieldElement`. 

  

EXAMPLES:: 

  

sage: a = RealBallField()(RIF(2)) 

sage: RIF(a) # indirect doctest 

2 

""" 

cdef RealIntervalFieldElement result 

result = parent(0) 

arb_to_mpfi(result.value, self.value, prec(self)) 

return result 

  

def _integer_(self, _): 

""" 

Check that this ball contains a single integer and return that integer. 

  

EXAMPLES:: 

  

sage: ZZ(RBF(1, rad=0.1r)) 

1 

sage: ZZ(RBF(1, rad=1.0r)) 

Traceback (most recent call last): 

... 

ValueError: [+/- 2.01] does not contain a unique integer 

sage: ZZ(RBF(pi)) 

Traceback (most recent call last): 

... 

ValueError: [3.141592653589793 +/- 5.61e-16] does not contain a unique integer 

  

""" 

cdef Integer res 

cdef fmpz_t tmp 

fmpz_init(tmp) 

try: 

if arb_get_unique_fmpz(tmp, self.value): 

res = Integer.__new__(Integer) 

fmpz_get_mpz(res.value, tmp) 

else: 

raise ValueError("{} does not contain a unique integer".format(self)) 

finally: 

fmpz_clear(tmp) 

return res 

  

def _rational_(self): 

""" 

Check that this ball contains a single rational number and return that 

number. 

  

EXAMPLES:: 

  

sage: QQ(RBF(123456/2^12)) 

1929/64 

sage: QQ(RBF(1/3)) 

Traceback (most recent call last): 

... 

ValueError: [0.3333333333333333 +/- 7.04e-17] does not contain a unique rational number 

""" 

if arb_is_exact(self.value): 

return self.mid().exact_rational() 

else: 

raise ValueError("{} does not contain a unique rational number".format(self)) 

  

def _mpfr_(self, RealField_class field): 

""" 

Convert this real ball to a real number. 

  

This attempts to do something sensible for all rounding modes, as 

illustrated below. 

  

EXAMPLES:: 

  

sage: mypi = RBF(pi) 

sage: RR(mypi) 

3.14159265358979 

sage: Reals(rnd='RNDU')(mypi) 

3.14159265358980 

sage: Reals(rnd='RNDD')(mypi) 

3.14159265358979 

sage: Reals(rnd='RNDZ')(mypi) 

3.14159265358979 

sage: Reals(rnd='RNDZ')(-mypi) 

-3.14159265358979 

sage: Reals(rnd='RNDU')(-mypi) 

-3.14159265358979 

  

:: 

  

sage: b = RBF(RIF(-1/2, 1)) 

sage: RR(b) 

0.250000000000000 

sage: Reals(rnd='RNDU')(b) 

1.00000000093133 

sage: Reals(rnd='RNDD')(b) 

-0.500000000931323 

sage: Reals(rnd='RNDZ')(b) 

0.250000000000000 

""" 

cdef RealNumber left, mid, right 

cdef long prec = field.precision() 

cdef int sl, sr 

if (field.rnd == MPFR_RNDN or 

field.rnd == MPFR_RNDZ and arb_contains_zero(self.value)): 

mid = RealNumber(field, None) 

sig_str("unable to convert to MPFR (exponent out of range?)") 

arf_get_mpfr(mid.value, arb_midref(self.value), field.rnd) 

sig_off() 

return mid 

else: 

left = RealNumber(field, None) 

right = RealNumber(field, None) 

sig_str("unable to convert to MPFR (exponent out of range?)") 

arb_get_interval_mpfr(left.value, right.value, self.value) 

sig_off() 

if field.rnd == MPFR_RNDD: 

return left 

elif field.rnd == MPFR_RNDU: 

return right 

elif field.rnd == MPFR_RNDZ: 

sl, sr = mpfr_sgn(left.value), mpfr_sgn(left.value) 

if sr > 0 and sl > 0: 

return left 

elif sr < 0 and sl < 0: 

return right 

else: 

return field(0) 

raise ValueError("unknown rounding mode") 

  

def __float__(self): 

""" 

Convert ``self`` to a ``float``. 

  

EXAMPLES:: 

  

sage: float(RBF(1)) 

1.0 

""" 

return float(self.n(prec(self))) 

  

def __complex__(self): 

""" 

Convert ``self`` to a ``complex``. 

  

EXAMPLES:: 

  

sage: complex(RBF(1)) 

(1+0j) 

""" 

return complex(self.n(prec(self))) 

  

# Center and radius, absolute value, endpoints 

  

def mid(self): 

""" 

Return the center of this ball. 

  

EXAMPLES:: 

  

sage: RealBallField(16)(1/3).mid() 

0.3333 

sage: RealBallField(16)(1/3).mid().parent() 

Real Field with 16 bits of precision 

sage: RealBallField(16)(RBF(1/3)).mid().parent() 

Real Field with 53 bits of precision 

sage: RBF('inf').mid() 

+infinity 

  

:: 

  

sage: b = RBF(2)^(2^1000) 

sage: b.mid() 

Traceback (most recent call last): 

... 

RuntimeError: unable to convert to MPFR (exponent out of range?) 

  

.. SEEALSO:: :meth:`rad`, :meth:`squash` 

""" 

cdef long mid_prec = max(arb_bits(self.value), prec(self)) 

if mid_prec < MPFR_PREC_MIN: 

mid_prec = MPFR_PREC_MIN 

cdef RealField_class mid_field = RealField(mid_prec) 

return self._mpfr_(mid_field) 

  

center = mid 

  

def rad(self): 

""" 

Return the radius of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1/3).rad() 

5.5511151e-17 

sage: RBF(1/3).rad().parent() 

Real Field with 30 bits of precision 

  

.. SEEALSO:: :meth:`mid`, :meth:`rad_as_ball`, :meth:`diameter` 

  

TESTS:: 

  

sage: (RBF(1, rad=.1) << (2^64)).rad() 

Traceback (most recent call last): 

... 

RuntimeError: unable to convert the radius to MPFR (exponent out of range?) 

""" 

# Should we return a real number with rounding towards +∞ (or away from 

# zero if/when implemented)? 

cdef RealField_class rad_field = RealField(MAG_BITS) 

cdef RealNumber rad = RealNumber(rad_field, None) 

cdef arf_t tmp 

arf_init(tmp) 

sig_str("unable to convert the radius to MPFR (exponent out of range?)") 

arf_set_mag(tmp, arb_radref(self.value)) 

if arf_get_mpfr(rad.value, tmp, MPFR_RNDN): 

abort() 

sig_off() 

arf_clear(tmp) 

return rad 

  

def diameter(self): 

r""" 

Return the diameter of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1/3).diameter() 

1.1102230e-16 

sage: RBF(1/3).diameter().parent() 

Real Field with 30 bits of precision 

sage: RBF(RIF(1.02, 1.04)).diameter() 

0.020000000 

  

.. SEEALSO:: :meth:`rad`, :meth:`rad_as_ball`, :meth:`mid` 

""" 

return 2 * self.rad() 

  

def squash(self): 

""" 

Return an exact ball with the same center as this ball. 

  

EXAMPLES:: 

  

sage: mid = RealBallField(16)(1/3).squash() 

sage: mid 

[0.3333 +/- 2.83e-5] 

sage: mid.is_exact() 

True 

sage: mid.parent() 

Real ball field with 16 bits of precision 

  

.. SEEALSO:: :meth:`mid`, :meth:`rad_as_ball` 

""" 

cdef RealBall res = self._new() 

arf_set(arb_midref(res.value), arb_midref(self.value)) 

mag_zero(arb_radref(res.value)) 

return res 

  

def rad_as_ball(self): 

""" 

Return an exact ball with center equal to the radius of this ball. 

  

EXAMPLES:: 

  

sage: rad = RBF(1/3).rad_as_ball() 

sage: rad 

[5.55111512e-17 +/- 3.13e-26] 

sage: rad.is_exact() 

True 

sage: rad.parent() 

Real ball field with 30 bits of precision 

  

.. SEEALSO:: :meth:`squash`, :meth:`rad` 

""" 

cdef RealBall res = self._parent.element_class(RealBallField(MAG_BITS)) 

arf_set_mag(arb_midref(res.value), arb_radref(self.value)) 

mag_zero(arb_radref(res.value)) 

return res 

  

def __abs__(self): 

""" 

Return the absolute value of this ball. 

  

EXAMPLES:: 

  

sage: RBF(-1/3).abs() # indirect doctest 

[0.3333333333333333 +/- 7.04e-17] 

sage: abs(RBF(-1)) 

1.000000000000000 

""" 

cdef RealBall r = self._new() 

arb_abs(r.value, self.value) 

return r 

  

def below_abs(self, test_zero=False): 

""" 

Return a lower bound for the absolute value of this ball. 

  

INPUT: 

  

- ``test_zero`` (boolean, default ``False``) -- if ``True``, 

make sure that the returned lower bound is positive, raising 

an error if the ball contains zero. 

  

OUTPUT: 

  

A ball with zero radius 

  

EXAMPLES:: 

  

sage: RealBallField(8)(1/3).below_abs() 

[0.33 +/- 7.82e-5] 

sage: b = RealBallField(8)(1/3).below_abs() 

sage: b 

[0.33 +/- 7.82e-5] 

sage: b.is_exact() 

True 

sage: QQ(b) 

169/512 

  

sage: RBF(0).below_abs() 

0 

sage: RBF(0).below_abs(test_zero=True) 

Traceback (most recent call last): 

... 

ValueError: ball contains zero 

  

.. SEEALSO:: :meth:`above_abs` 

""" 

cdef RealBall res = self._new() 

arb_get_abs_lbound_arf(arb_midref(res.value), self.value, prec(self)) 

if test_zero and arb_contains_zero(res.value): 

assert arb_contains_zero(self.value) 

raise ValueError("ball contains zero") 

return res 

  

def above_abs(self): 

""" 

Return an upper bound for the absolute value of this ball. 

  

OUTPUT: 

  

A ball with zero radius 

  

EXAMPLES:: 

  

sage: b = RealBallField(8)(1/3).above_abs() 

sage: b 

[0.33 +/- 3.99e-3] 

sage: b.is_exact() 

True 

sage: QQ(b) 

171/512 

  

.. SEEALSO:: :meth:`below_abs` 

""" 

cdef RealBall res = self._new() 

arb_get_abs_ubound_arf(arb_midref(res.value), self.value, prec(self)) 

return res 

  

def upper(self, rnd=None): 

""" 

Return the right endpoint of this ball, rounded upwards. 

  

INPUT: 

  

- ``rnd`` (string) -- rounding mode for the parent of the result (does 

not affect its value!), see 

:meth:`sage.rings.real_mpfi.RealIntervalFieldElement.upper` 

  

OUTPUT: 

  

A real number. 

  

EXAMPLES:: 

  

sage: RBF(-1/3).upper() 

-0.333333333333333 

sage: RBF(-1/3).upper().parent() 

Real Field with 53 bits of precision and rounding RNDU 

  

.. SEEALSO:: 

  

:meth:`lower`, :meth:`endpoints` 

""" 

# naive and slow 

return self._real_mpfi_(RealIntervalField(prec(self))).upper(rnd) 

  

def lower(self, rnd=None): 

""" 

Return the right endpoint of this ball, rounded downwards. 

  

INPUT: 

  

- ``rnd`` (string) -- rounding mode for the parent of the result (does 

not affect its value!), see 

:meth:`sage.rings.real_mpfi.RealIntervalFieldElement.lower` 

  

OUTPUT: 

  

A real number. 

  

EXAMPLES:: 

  

sage: RBF(-1/3).lower() 

-0.333333333333334 

sage: RBF(-1/3).lower().parent() 

Real Field with 53 bits of precision and rounding RNDD 

  

.. SEEALSO:: :meth:`upper`, :meth:`endpoints` 

""" 

# naive and slow 

return self._real_mpfi_(RealIntervalField(prec(self))).lower(rnd) 

  

def endpoints(self, rnd=None): 

""" 

Return the endpoints of this ball, rounded outwards. 

  

INPUT: 

  

- ``rnd`` (string) -- rounding mode for the parent of the resulting 

floating-point numbers (does not affect their values!), see 

:meth:`sage.rings.real_mpfi.RealIntervalFieldElement.upper` 

  

OUTPUT: 

  

A pair of real numbers. 

  

EXAMPLES:: 

  

sage: RBF(-1/3).endpoints() 

(-0.333333333333334, -0.333333333333333) 

  

.. SEEALSO:: :meth:`lower`, :meth:`upper` 

""" 

# naive and slow 

return self._real_mpfi_(RealIntervalField(prec(self))).endpoints(rnd) 

  

def union(self, other): 

r""" 

Return a ball containing the convex hull of ``self`` and ``other``. 

  

EXAMPLES:: 

  

sage: RBF(0).union(1).endpoints() 

(-9.31322574615479e-10, 1.00000000093133) 

""" 

cdef RealBall my_other = self._parent.coerce(other) 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_union(res.value, self.value, my_other.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

# Precision and accuracy 

  

def round(self): 

""" 

Return a copy of this ball with center rounded to the precision of the 

parent. 

  

EXAMPLES: 

  

It is possible to create balls whose midpoint is more precise that 

their parent's nominal precision (see :mod:`~sage.rings.real_arb` for 

more information):: 

  

sage: b = RBF(pi.n(100)) 

sage: b.mid() 

3.141592653589793238462643383 

  

The ``round()`` method rounds such a ball to its parent's precision:: 

  

sage: b.round().mid() 

3.14159265358979 

  

.. SEEALSO:: :meth:`trim` 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_set_round(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def accuracy(self): 

""" 

Return the effective relative accuracy of this ball measured in bits. 

  

The accuracy is defined as the difference between the position of the 

top bit in the midpoint and the top bit in the radius, minus one. 

The result is clamped between plus/minus 

:meth:`~RealBallField.maximal_accuracy`. 

  

EXAMPLES:: 

  

sage: RBF(pi).accuracy() 

51 

sage: RBF(1).accuracy() == RBF.maximal_accuracy() 

True 

sage: RBF(NaN).accuracy() == -RBF.maximal_accuracy() 

True 

  

.. SEEALSO:: :meth:`~RealBallField.maximal_accuracy` 

""" 

return arb_rel_accuracy_bits(self.value) 

  

def trim(self): 

""" 

Return a trimmed copy of this ball. 

  

Round ``self`` to a number of bits equal to the :meth:`accuracy` of 

``self`` (as indicated by its radius), plus a few guard bits. The 

resulting ball is guaranteed to contain ``self``, but is more economical 

if ``self`` has less than full accuracy. 

  

EXAMPLES:: 

  

sage: b = RBF(0.11111111111111, rad=.001) 

sage: b.mid() 

0.111111111111110 

sage: b.trim().mid() 

0.111111104488373 

  

.. SEEALSO:: :meth:`round` 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_trim(res.value, self.value) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def add_error(self, ampl): 

""" 

Increase the radius of this ball by (an upper bound on) ``ampl``. 

  

If ``ampl`` is negative, the radius is unchanged. 

  

INPUT: 

  

- ``ampl`` -- A real ball (or an object that can be coerced to a real 

ball). 

  

OUTPUT: 

  

A new real ball. 

  

EXAMPLES:: 

  

sage: err = RBF(10^-16) 

sage: RBF(1).add_error(err) 

[1.000000000000000 +/- 1.01e-16] 

  

TESTS:: 

  

sage: RBF(1).add_error(-1) 

1.000000000000000 

sage: RBF(0).add_error(RBF(1, rad=2.)).endpoints() 

(-3.00000000745059, 3.00000000745059) 

""" 

cdef RealBall res = self._new() 

cdef RealBall my_ampl = self._parent(ampl) 

if my_ampl < 0: 

my_ampl = self._parent.zero() 

arb_set(res.value, self.value) 

arb_add_error(res.value, my_ampl.value) 

return res 

  

# Comparisons and predicates 

  

def is_zero(self): 

""" 

Return ``True`` iff the midpoint and radius of this ball are both zero. 

  

EXAMPLES:: 

  

sage: RBF = RealBallField() 

sage: RBF(0).is_zero() 

True 

sage: RBF(RIF(-0.5, 0.5)).is_zero() 

False 

  

.. SEEALSO:: :meth:`is_nonzero` 

""" 

return arb_is_zero(self.value) 

  

def is_nonzero(self): 

""" 

Return ``True`` iff zero is not contained in the interval represented 

by this ball. 

  

.. NOTE:: 

  

This method is not the negation of :meth:`is_zero`: it only 

returns ``True`` if zero is known not to be contained in the ball. 

  

Use ``bool(b)`` (or, equivalently, ``not b.is_zero()``) to check if 

a ball ``b`` **may** represent a nonzero number (for instance, to 

determine the “degree” of a polynomial with ball coefficients). 

  

EXAMPLES:: 

  

sage: RBF = RealBallField() 

sage: RBF(pi).is_nonzero() 

True 

sage: RBF(RIF(-0.5, 0.5)).is_nonzero() 

False 

  

.. SEEALSO:: :meth:`is_zero` 

""" 

return arb_is_nonzero(self.value) 

  

def __nonzero__(self): 

""" 

Return ``True`` iff this ball is not the zero ball, i.e. if it its 

midpoint and radius are not both zero. 

  

This is the preferred way, for instance, to determine the “degree” of a 

polynomial with ball coefficients. 

  

.. WARNING:: 

  

A “nonzero” ball in the sense of this method may represent the 

value zero. Use :meth:`is_nonzero` to check that a real number 

represented by a ``RealBall`` object is known to be nonzero. 

  

EXAMPLES:: 

  

sage: bool(RBF(0)) # indirect doctest 

False 

sage: bool(RBF(1/3)) 

True 

sage: bool(RBF(RIF(-0.5, 0.5))) 

True 

""" 

return not arb_is_zero(self.value) 

  

def is_exact(self): 

""" 

Return ``True`` iff the radius of this ball is zero. 

  

EXAMPLES:: 

  

sage: RBF = RealBallField() 

sage: RBF(1).is_exact() 

True 

sage: RBF(RIF(0.1, 0.2)).is_exact() 

False 

""" 

return arb_is_exact(self.value) 

  

cpdef _richcmp_(left, right, int op): 

""" 

Compare ``left`` and ``right``. 

  

For more information, see :mod:`sage.rings.real_arb`. 

  

EXAMPLES:: 

  

sage: RBF = RealBallField() 

sage: a = RBF(1) 

sage: b = RBF(1) 

sage: a is b 

False 

sage: a == b 

True 

sage: a = RBF(1/3) 

sage: a.is_exact() 

False 

sage: b = RBF(1/3) 

sage: b.is_exact() 

False 

sage: a == b 

False 

  

TESTS: 

  

Balls whose intersection consists of one point:: 

  

sage: a = RBF(RIF(1, 2)) 

sage: b = RBF(RIF(2, 4)) 

sage: a < b 

False 

sage: a > b 

False 

sage: a <= b 

False 

sage: a >= b 

False 

sage: a == b 

False 

sage: a != b 

False 

  

Balls with non-trivial intersection:: 

  

sage: a = RBF(RIF(1, 4)) 

sage: a = RBF(RIF(2, 5)) 

sage: a < b 

False 

sage: a <= b 

False 

sage: a > b 

False 

sage: a >= b 

False 

sage: a == b 

False 

sage: a != b 

False 

  

One ball contained in another:: 

  

sage: a = RBF(RIF(1, 4)) 

sage: b = RBF(RIF(2, 3)) 

sage: a < b 

False 

sage: a <= b 

False 

sage: a > b 

False 

sage: a >= b 

False 

sage: a == b 

False 

sage: a != b 

False 

  

Disjoint balls:: 

  

sage: a = RBF(1/3) 

sage: b = RBF(1/2) 

sage: a < b 

True 

sage: a <= b 

True 

sage: a > b 

False 

sage: a >= b 

False 

sage: a == b 

False 

sage: a != b 

True 

  

Exact elements:: 

  

sage: a = RBF(2) 

sage: b = RBF(2) 

sage: a.is_exact() 

True 

sage: b.is_exact() 

True 

sage: a < b 

False 

sage: a <= b 

True 

sage: a > b 

False 

sage: a >= b 

True 

sage: a == b 

True 

sage: a != b 

False 

  

Special values:: 

  

sage: inf = RBF(+infinity) 

sage: other_inf = RBF(+infinity, 42.r) 

sage: neg_inf = RBF(-infinity) 

sage: extended_line = 1/RBF(0) 

sage: exact_nan = inf - inf 

sage: exact_nan.mid(), exact_nan.rad() 

(NaN, 0.00000000) 

sage: other_exact_nan = inf - inf 

  

:: 

  

sage: exact_nan == exact_nan, exact_nan <= exact_nan, exact_nan >= exact_nan 

(False, False, False) 

sage: exact_nan != exact_nan, exact_nan < exact_nan, exact_nan > exact_nan 

(False, False, False) 

sage: from operator import eq, ne, le, lt, ge, gt 

sage: ops = [eq, ne, le, lt, ge, gt] 

sage: any(op(exact_nan, other_exact_nan) for op in ops) 

False 

sage: any(op(exact_nan, b) for op in ops for b in [RBF(1), extended_line, inf, neg_inf]) 

False 

  

:: 

  

sage: neg_inf < a < inf and inf > a > neg_inf 

True 

sage: neg_inf <= b <= inf and inf >= b >= neg_inf 

True 

sage: neg_inf <= extended_line <= inf and inf >= extended_line >= neg_inf 

True 

sage: neg_inf < extended_line or extended_line < inf 

False 

sage: inf > extended_line or extended_line > neg_inf 

False 

  

:: 

  

sage: all(b <= b == b >= b and not (b < b or b != b or b > b) 

....: for b in [inf, neg_inf, other_inf]) 

True 

sage: any(b1 == b2 for b1 in [inf, neg_inf, a, extended_line] 

....: for b2 in [inf, neg_inf, a, extended_line] 

....: if not b1 is b2) 

False 

sage: all(b1 != b2 and not b1 == b2 

....: for b1 in [inf, neg_inf, a] 

....: for b2 in [inf, neg_inf, a] 

....: if not b1 is b2) 

True 

sage: neg_inf <= -other_inf == neg_inf == -other_inf < other_inf == inf <= other_inf 

True 

sage: any(inf < b or b > inf 

....: for b in [inf, other_inf, a, extended_line]) 

False 

sage: any(inf <= b or b >= inf for b in [a, extended_line]) 

False 

""" 

cdef RealBall lt, rt 

cdef arb_t difference 

  

lt = left 

rt = right 

  

if op == Py_EQ: 

return arb_eq(lt.value, rt.value) 

elif op == Py_NE: 

return arb_ne(lt.value, rt.value) 

elif op == Py_GT: 

return arb_gt(lt.value, rt.value) 

elif op == Py_LT: 

return arb_lt(lt.value, rt.value) 

elif op == Py_GE: 

return arb_ge(lt.value, rt.value) 

elif op == Py_LE: 

return arb_le(lt.value, rt.value) 

  

def min(self, *others): 

""" 

Return a ball containing the minimum of this ball and the 

remaining arguments. 

  

EXAMPLES:: 

  

sage: RBF(1, rad=.5).min(0) 

0 

  

sage: RBF(0, rad=2.).min(RBF(0, rad=1.)).endpoints() 

(-2.00000000651926, 1.00000000465662) 

  

sage: RBF(infinity).min(3, 1/3) 

[0.3333333333333333 +/- 7.04e-17] 

  

Note that calls involving NaNs try to return a number when possible. 

This is consistent with IEEE-754-2008 but may be surprising. :: 

  

sage: RBF('nan').min(0) 

0 

sage: RBF('nan').min(RBF('nan')) 

nan 

  

.. SEEALSO:: :meth:`max` 

  

TESTS:: 

  

sage: RBF(0).min() 

0 

sage: RBF(infinity).min().rad() 

0.00000000 

""" 

iv = self._real_mpfi_(RealIntervalField(prec(self))) 

my_others = [self._parent.coerce(x) for x in others] 

return self._parent(iv.min(*my_others)) 

  

def max(self, *others): 

""" 

Return a ball containing the maximum of this ball and the 

remaining arguments. 

  

EXAMPLES:: 

  

sage: RBF(-1, rad=.5).max(0) 

0 

  

sage: RBF(0, rad=2.).max(RBF(0, rad=1.)).endpoints() 

(-1.00000000465662, 2.00000000651926) 

  

sage: RBF(-infinity).max(-3, 1/3) 

[0.3333333333333333 +/- 7.04e-17] 

  

Note that calls involving NaNs try to return a number when possible. 

This is consistent with IEEE-754-2008 but may be surprising. :: 

  

sage: RBF('nan').max(0) 

0 

sage: RBF('nan').max(RBF('nan')) 

nan 

  

.. SEEALSO:: :meth:`min` 

  

TESTS:: 

  

sage: RBF(0).max() 

0 

""" 

iv = self._real_mpfi_(RealIntervalField(prec(self))) 

my_others = [self._parent.coerce(x) for x in others] 

return self._parent(iv.max(*my_others)) 

  

def is_finite(self): 

""" 

Return True iff the midpoint and radius of this ball are both 

finite floating-point numbers, i.e. not infinities or NaN. 

  

EXAMPLES:: 

  

sage: (RBF(2)^(2^1000)).is_finite() 

True 

sage: RBF(oo).is_finite() 

False 

""" 

return arb_is_finite(self.value) 

  

def identical(self, RealBall other): 

""" 

Return True iff ``self`` and ``other`` are equal as balls, i.e. 

have both the same midpoint and radius. 

  

Note that this is not the same thing as testing whether both ``self`` 

and ``other`` certainly represent the same real number, unless either 

``self`` or ``other`` is exact (and neither contains NaN). To test 

whether both operands might represent the same mathematical quantity, 

use :meth:`overlaps` or :meth:`contains`, depending on the 

circumstance. 

  

EXAMPLES:: 

  

sage: RBF(1).identical(RBF(3)-RBF(2)) 

True 

sage: RBF(1, rad=0.25r).identical(RBF(1, rad=0.25r)) 

True 

sage: RBF(1).identical(RBF(1, rad=0.25r)) 

False 

""" 

return arb_equal(self.value, other.value) 

  

def overlaps(self, RealBall other): 

""" 

Return True iff ``self`` and ``other`` have some point in common. 

  

If either ``self`` or ``other`` contains NaN, this method always 

returns nonzero (as a NaN could be anything, it could in particular 

contain any number that is included in the other operand). 

  

EXAMPLES:: 

  

sage: RBF(pi).overlaps(RBF(pi) + 2**(-100)) 

True 

sage: RBF(pi).overlaps(RBF(3)) 

False 

""" 

return arb_overlaps(self.value, other.value) 

  

def contains_exact(self, other): 

""" 

Return ``True`` *iff* the given number (or ball) ``other`` is contained 

in the interval represented by ``self``. 

  

If ``self`` contains NaN, this function always returns ``True`` (as 

it could represent anything, and in particular could represent all the 

points included in ``other``). If ``other`` contains NaN and ``self`` 

does not, it always returns ``False``. 

  

Use ``other in self`` for a test that works for a wider range of inputs 

but may return false negatives. 

  

EXAMPLES:: 

  

sage: b = RBF(1) 

sage: b.contains_exact(1) 

True 

sage: b.contains_exact(QQ(1)) 

True 

sage: b.contains_exact(1.) 

True 

sage: b.contains_exact(b) 

True 

  

:: 

  

sage: RBF(1/3).contains_exact(1/3) 

True 

sage: RBF(sqrt(2)).contains_exact(sqrt(2)) 

Traceback (most recent call last): 

... 

TypeError: unsupported type: <type 'sage.symbolic.expression.Expression'> 

  

TESTS:: 

  

sage: b.contains_exact(1r) 

True 

  

""" 

cdef fmpz_t tmpz 

cdef fmpq_t tmpq 

if _do_sig(prec(self)): sig_on() 

try: 

if isinstance(other, RealBall): 

res = arb_contains(self.value, (<RealBall> other).value) 

elif isinstance(other, int): 

res = arb_contains_si(self.value, PyInt_AS_LONG(other)) 

elif isinstance(other, Integer): 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, (<Integer> other).value) 

res = arb_contains_fmpz(self.value, tmpz) 

fmpz_clear(tmpz) 

elif isinstance(other, Rational): 

fmpq_init(tmpq) 

fmpq_set_mpq(tmpq, (<Rational> other).value) 

res = arb_contains_fmpq(self.value, tmpq) 

fmpq_clear(tmpq) 

elif isinstance(other, RealNumber): 

res = arb_contains_mpfr(self.value, (<RealNumber> other).value) 

else: 

raise TypeError("unsupported type: " + str(type(other))) 

finally: 

if _do_sig(prec(self)): sig_off() 

return res 

  

def __contains__(self, other): 

""" 

Return True if ``other`` can be verified to be contained in ``self``. 

  

The test is done using interval arithmetic with a precision determined 

by the parent of ``self`` and may return false negatives. 

  

EXAMPLES:: 

  

sage: sqrt(2) in RBF(sqrt(2)) 

True 

  

A false negative:: 

  

sage: sqrt(2) in RBF(RealBallField(100)(sqrt(2))) 

False 

  

.. SEEALSO:: :meth:`contains_exact` 

""" 

return self.contains_exact(self._parent(other)) 

  

def contains_zero(self): 

""" 

Return ``True`` iff this ball contains zero. 

  

EXAMPLES:: 

  

sage: RBF(0).contains_zero() 

True 

sage: RBF(RIF(-1, 1)).contains_zero() 

True 

sage: RBF(1/3).contains_zero() 

False 

""" 

return arb_contains_zero(self.value) 

  

def contains_integer(self): 

""" 

Return ``True`` iff this ball contains any integer. 

  

EXAMPLES:: 

  

sage: RBF(3.1, 0.1).contains_integer() 

True 

sage: RBF(3.1, 0.05).contains_integer() 

False 

""" 

return arb_contains_int(self.value) 

  

def is_negative_infinity(self): 

""" 

Return ``True`` if this ball is the point -∞. 

  

EXAMPLES:: 

  

sage: RBF(-infinity).is_negative_infinity() 

True 

""" 

return (arf_is_neg_inf(arb_midref(self.value)) 

and mag_is_finite(arb_radref(self.value))) 

  

def is_positive_infinity(self): 

""" 

Return ``True`` if this ball is the point +∞. 

  

EXAMPLES:: 

  

sage: RBF(infinity).is_positive_infinity() 

True 

""" 

return (arf_is_pos_inf(arb_midref(self.value)) 

and mag_is_finite(arb_radref(self.value))) 

  

def is_infinity(self): 

""" 

Return ``True`` if this ball contains or may represent a point at 

infinity. 

  

This is the exact negation of :meth:`is_finite`, used in comparisons 

with Sage symbolic infinities. 

  

.. WARNING:: 

  

Contrary to the usual convention, a return value of True does 

not imply that all points of the ball satisfy the predicate. 

This is due to the way comparisons with symbolic infinities work in 

sage. 

  

EXAMPLES:: 

  

sage: RBF(infinity).is_infinity() 

True 

sage: RBF(-infinity).is_infinity() 

True 

sage: RBF(NaN).is_infinity() 

True 

sage: (~RBF(0)).is_infinity() 

True 

sage: RBF(42, rad=1.r).is_infinity() 

False 

""" 

return not self.is_finite() 

  

def is_NaN(self): 

""" 

Return ``True`` if this ball is not-a-number. 

  

EXAMPLES:: 

  

sage: RBF(NaN).is_NaN() 

True 

sage: RBF(-5).gamma().is_NaN() 

True 

sage: RBF(infinity).is_NaN() 

False 

sage: RBF(42, rad=1.r).is_NaN() 

False 

""" 

return arf_is_nan(arb_midref(self.value)) 

  

# Arithmetic 

  

def __neg__(self): 

""" 

Return the opposite of this ball. 

  

EXAMPLES:: 

  

sage: -RBF(1/3) 

[-0.3333333333333333 +/- 7.04e-17] 

""" 

cdef RealBall res = self._new() 

arb_neg(res.value, self.value) 

return res 

  

def __invert__(self): 

""" 

Return the inverse of this ball. 

  

The result is guaranteed to contain the inverse of any point of the 

input ball. 

  

EXAMPLES:: 

  

sage: ~RBF(5) 

[0.2000000000000000 +/- 4.45e-17] 

sage: ~RBF(0) 

[+/- inf] 

sage: RBF(RIF(-0.1,0.1)) 

[+/- 0.101] 

  

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_inv(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

cpdef _add_(self, other): 

""" 

Return the sum of two balls, rounded to the ambient field's precision. 

  

The resulting ball is guaranteed to contain the sums of any two points 

of the respective input balls. 

  

EXAMPLES:: 

  

sage: RBF(1) + RBF(1/3) 

[1.333333333333333 +/- 5.37e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_add(res.value, self.value, (<RealBall> other).value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

cpdef _sub_(self, other): 

""" 

Return the difference of two balls, rounded to the ambient field's 

precision. 

  

The resulting ball is guaranteed to contain the differences of any two 

points of the respective input balls. 

  

EXAMPLES:: 

  

sage: RBF(1) - RBF(1/3) 

[0.666666666666667 +/- 5.37e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_sub(res.value, self.value, (<RealBall> other).value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

cpdef _mul_(self, other): 

""" 

Return the product of two balls, rounded to the ambient field's 

precision. 

  

The resulting ball is guaranteed to contain the products of any two 

points of the respective input balls. 

  

EXAMPLES:: 

  

sage: RBF(-2) * RBF(1/3) 

[-0.666666666666667 +/- 4.82e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_mul(res.value, self.value, (<RealBall> other).value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

cpdef _div_(self, other): 

""" 

Return the quotient of two balls, rounded to the ambient field's 

precision. 

  

The resulting ball is guaranteed to contain the quotients of any two 

points of the respective input balls. 

  

EXAMPLES:: 

  

sage: RBF(pi)/RBF(e) 

[1.155727349790922 +/- 8.43e-16] 

sage: RBF(2)/RBF(0) 

[+/- inf] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_div(res.value, self.value, (<RealBall> other).value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def __pow__(base, expo, _): 

""" 

EXAMPLES:: 

  

sage: RBF(e)^17 

[24154952.7535753 +/- 9.30e-8] 

sage: RBF(e)^(-1) 

[0.367879441171442 +/- 4.50e-16] 

sage: RBF(e)^(1/2) 

[1.648721270700128 +/- 4.96e-16] 

sage: RBF(e)^RBF(pi) 

[23.1406926327793 +/- 9.16e-14] 

  

:: 

  

sage: RBF(-1)^(1/3) 

nan 

sage: RBF(0)^(-1) 

[+/- inf] 

sage: RBF(-e)**RBF(pi) 

nan 

  

TESTS:: 

  

sage: RBF(e)**(2r) 

[7.38905609893065 +/- 4.68e-15] 

sage: RBF(e)**(-1r) 

[0.367879441171442 +/- 4.50e-16] 

""" 

cdef fmpz_t tmpz 

if not isinstance(base, RealBall): 

return sage.structure.element.bin_op(base, expo, operator.pow) 

cdef RealBall self = base 

cdef RealBall res = self._new() 

if isinstance(expo, int) and expo > 0: 

if _do_sig(prec(self)): sig_on() 

arb_pow_ui(res.value, self.value, PyInt_AS_LONG(expo), prec(self)) 

if _do_sig(prec(self)): sig_off() 

elif isinstance(expo, Integer): 

if _do_sig(prec(self)): sig_on() 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, (<Integer> expo).value) 

arb_pow_fmpz(res.value, self.value, tmpz, prec(self)) 

fmpz_clear(tmpz) 

if _do_sig(prec(self)): sig_off() 

elif isinstance(expo, RealBall): 

if _do_sig(prec(self)): sig_on() 

arb_pow(res.value, self.value, (<RealBall> expo).value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

else: 

return sage.structure.element.bin_op(base, expo, operator.pow) 

return res 

  

def sqrt(self): 

""" 

Return the square root of this ball. 

  

EXAMPLES:: 

  

sage: RBF(2).sqrt() 

[1.414213562373095 +/- 2.99e-16] 

sage: RBF(-1/3).sqrt() 

nan 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_sqrt(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def sqrtpos(self): 

""" 

Return the square root of this ball, assuming that it represents a 

nonnegative number. 

  

Any negative numbers in the input interval are discarded. 

  

EXAMPLES:: 

  

sage: RBF(2).sqrtpos() 

[1.414213562373095 +/- 2.99e-16] 

sage: RBF(-1/3).sqrtpos() 

0 

sage: RBF(0, rad=2.r).sqrtpos() 

[+/- 1.42] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_sqrtpos(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def rsqrt(self): 

""" 

Return the reciprocal square root of ``self``. 

  

At high precision, this is faster than computing a square root. 

  

EXAMPLES:: 

  

sage: RBF(2).rsqrt() 

[0.707106781186547 +/- 5.73e-16] 

sage: RBF(0).rsqrt() 

nan 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_rsqrt(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def sqrt1pm1(self): 

""" 

Return `\sqrt{1+\mathrm{self}}-1`, computed accurately when ``self`` is 

close to zero. 

  

EXAMPLES:: 

  

sage: eps = RBF(10^(-20)) 

sage: (1 + eps).sqrt() - 1 

[+/- 1.12e-16] 

sage: eps.sqrt1pm1() 

[5.00000000000000e-21 +/- 2.54e-36] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_sqrt1pm1(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

# Floor, ceil, etc. 

  

def floor(self): 

""" 

Return the floor of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1000+1/3, rad=1.r).floor() 

[1.00e+3 +/- 1.01] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_floor(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def ceil(self): 

""" 

Return the ceil of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1000+1/3, rad=1.r).ceil() 

[1.00e+3 +/- 2.01] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_ceil(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def __lshift__(val, shift): 

r""" 

If ``val`` is a ``RealBall`` and ``shift`` is an integer, return the 

ball obtained by shifting the center and radius of ``val`` to the left 

by ``shift`` bits. 

  

INPUT: 

  

- ``shift`` -- integer, may be negative. 

  

EXAMPLES:: 

  

sage: RBF(1/3) << 2 # indirect doctest 

[1.333333333333333 +/- 4.82e-16] 

sage: RBF(1) << -1 

0.5000000000000000 

  

TESTS:: 

  

sage: RBF(1) << (2^100) 

[2.285367694229514e+381600854690147056244358827360 +/- 2.98e+381600854690147056244358827344] 

sage: RBF(1) << (-2^100) 

[4.375663498372584e-381600854690147056244358827361 +/- 9.57e-381600854690147056244358827378] 

  

sage: "a" << RBF(1/3) 

Traceback (most recent call last): 

... 

TypeError: unsupported operand type(s) for <<: 'str' and 'RealBall' 

sage: RBF(1) << RBF(1/3) 

Traceback (most recent call last): 

... 

TypeError: shift should be an integer 

""" 

cdef fmpz_t tmpz 

# the RealBall might be shift, not val 

if not isinstance(val, RealBall): 

raise TypeError("unsupported operand type(s) for <<: '{}' and '{}'" 

.format(type(val).__name__, type(shift).__name__)) 

cdef RealBall self = val 

cdef RealBall res = self._new() 

if isinstance(shift, int): 

arb_mul_2exp_si(res.value, self.value, PyInt_AS_LONG(shift)) 

elif isinstance(shift, Integer): 

sig_on() 

fmpz_init(tmpz) 

fmpz_set_mpz(tmpz, (<Integer> shift).value) 

arb_mul_2exp_fmpz(res.value, self.value, tmpz) 

fmpz_clear(tmpz) 

sig_off() 

else: 

raise TypeError("shift should be an integer") 

return res 

  

def __rshift__(val, shift): 

r""" 

If ``val`` is a ``RealBall`` and ``shift`` is an integer, return the 

ball obtained by shifting the center and radius of ``val`` to the right 

by ``shift`` bits. 

  

INPUT: 

  

- ``shift`` -- integer, may be negative. 

  

EXAMPLES:: 

  

sage: RBF(4) >> 2 

1.000000000000000 

sage: RBF(1/3) >> -2 

[1.333333333333333 +/- 4.82e-16] 

  

TESTS:: 

  

sage: "a" >> RBF(1/3) 

Traceback (most recent call last): 

... 

TypeError: unsupported operand type(s) for >>: 'str' and 'RealBall' 

""" 

# the RealBall might be shift, not val 

if isinstance(val, RealBall): 

return val << (-shift) 

else: 

raise TypeError("unsupported operand type(s) for >>: '{}' and '{}'" 

.format(type(val).__name__, type(shift).__name__)) 

  

# Elementary functions 

  

def log(self, base=None): 

""" 

Return the logarithm of this ball. 

  

INPUT: 

  

- ``base`` (optional, positive real ball or number) -- if ``None``, 

return the natural logarithm ``ln(self)``, otherwise, return the 

general logarithm ``ln(self)/ln(base)`` 

  

EXAMPLES:: 

  

sage: RBF(3).log() 

[1.098612288668110 +/- 6.63e-16] 

sage: RBF(3).log(2) 

[1.584962500721156 +/- 7.53e-16] 

sage: log(RBF(5), 2) 

[2.32192809488736 +/- 3.04e-15] 

  

sage: RBF(-1/3).log() 

nan 

sage: RBF(3).log(-1) 

nan 

sage: RBF(2).log(0) 

nan 

""" 

cdef RealBall cst 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_log(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

if base is not None: 

cst = self._parent.coerce(base).log() 

if _do_sig(prec(self)): sig_on() 

arb_div(res.value, res.value, cst.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def log1p(self): 

""" 

Return ``log(1 + self)``, computed accurately when ``self`` is close to 

zero. 

  

EXAMPLES:: 

  

sage: eps = RBF(1e-30) 

sage: (1 + eps).log() 

[+/- 2.23e-16] 

sage: eps.log1p() 

[1.00000000000000e-30 +/- 2.68e-46] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_log1p(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def exp(self): 

""" 

Return the exponential of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).exp() 

[2.718281828459045 +/- 5.41e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_exp(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def expm1(self): 

""" 

Return ``exp(self) - 1``, computed accurately when ``self`` is close to 

zero. 

  

EXAMPLES:: 

  

sage: eps = RBF(1e-30) 

sage: exp(eps) - 1 

[+/- 3.16e-30] 

sage: eps.expm1() 

[1.000000000000000e-30 +/- 8.34e-47] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_expm1(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def sin(self): 

""" 

Return the sine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(pi).sin() # abs tol 1e-16 

[+/- 5.69e-16] 

  

.. SEEALSO:: :meth:`~sage.rings.real_arb.RealBallField.sinpi` 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_sin(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def cos(self): 

""" 

Return the cosine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(pi).cos() # abs tol 1e-16 

[-1.00000000000000 +/- 6.69e-16] 

  

.. SEEALSO:: :meth:`~sage.rings.real_arb.RealBallField.cospi` 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_cos(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def tan(self): 

""" 

Return the tangent of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).tan() 

[1.557407724654902 +/- 3.26e-16] 

sage: RBF(pi/2).tan() 

[+/- inf] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_tan(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def cot(self): 

""" 

Return the cotangent of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).cot() 

[0.642092615934331 +/- 4.79e-16] 

sage: RBF(pi).cot() 

[+/- inf] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_cot(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def arcsin(self): 

""" 

Return the arcsine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).arcsin() 

[1.570796326794897 +/- 6.65e-16] 

sage: RBF(1, rad=.125r).arcsin() 

nan 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_asin(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def arccos(self): 

""" 

Return the arccosine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).arccos() 

0 

sage: RBF(1, rad=.125r).arccos() 

nan 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_acos(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def arctan(self): 

""" 

Return the arctangent of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).arctan() 

[0.785398163397448 +/- 3.91e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_atan(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def sinh(self): 

""" 

Return the hyperbolic sine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).sinh() 

[1.175201193643801 +/- 6.18e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_sinh(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def cosh(self): 

""" 

Return the hyperbolic cosine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).cosh() 

[1.543080634815244 +/- 5.28e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_cosh(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def tanh(self): 

""" 

Return the hyperbolic tangent of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).tanh() 

[0.761594155955765 +/- 2.81e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_tanh(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def coth(self): 

""" 

Return the hyperbolic cotangent of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).coth() 

[1.313035285499331 +/- 4.97e-16] 

sage: RBF(0).coth() 

[+/- inf] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_coth(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def arcsinh(self): 

""" 

Return the inverse hyperbolic sine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(1).arcsinh() 

[0.881373587019543 +/- 1.87e-16] 

sage: RBF(0).arcsinh() 

0 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_asinh(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def arccosh(self): 

""" 

Return the inverse hyperbolic cosine of this ball. 

  

EXAMPLES:: 

  

sage: RBF(2).arccosh() 

[1.316957896924817 +/- 6.61e-16] 

sage: RBF(1).arccosh() 

0 

sage: RBF(0).arccosh() 

nan 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_acosh(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def arctanh(self): 

""" 

Return the inverse hyperbolic tangent of this ball. 

  

EXAMPLES:: 

  

sage: RBF(0).arctanh() 

0 

sage: RBF(1/2).arctanh() 

[0.549306144334055 +/- 3.32e-16] 

sage: RBF(1).arctanh() 

nan 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_atanh(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

# Special functions 

  

def gamma(self): 

""" 

Return the image of this ball by the Euler Gamma function. 

  

For integer and rational arguments, 

:meth:`~sage.rings.real_arb.RealBallField.gamma` may be faster. 

  

EXAMPLES:: 

  

sage: RBF(1/2).gamma() 

[1.772453850905516 +/- 3.41e-16] 

  

.. SEEALSO:: 

:meth:`~sage.rings.real_arb.RealBallField.gamma` 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_gamma(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def log_gamma(self): 

""" 

Return the image of this ball by the logarithmic Gamma function. 

  

The complex branch structure is assumed, so if ``self`` <= 0, the result 

is an indeterminate interval. 

  

EXAMPLES:: 

  

sage: RBF(1/2).log_gamma() 

[0.572364942924700 +/- 2.67e-16] 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_lgamma(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def rgamma(self): 

""" 

Return the image of this ball by the function 1/Γ, avoiding division by 

zero at the poles of the gamma function. 

  

EXAMPLES:: 

  

sage: RBF(-1).rgamma() 

0 

sage: RBF(3).rgamma() 

0.5000000000000000 

""" 

cdef RealBall res = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_rgamma(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def rising_factorial(self, n): 

""" 

Return the ``n``-th rising factorial of this ball. 

  

The `n`-th rising factorial of `x` is equal to `x (x+1) \cdots (x+n-1)`. 

  

For real `n`, it is a quotient of gamma functions. 

  

EXAMPLES:: 

  

sage: RBF(1).rising_factorial(5) 

120.0000000000000 

sage: RBF(1/2).rising_factorial(1/3) 

[0.63684988431797 +/- 5.71e-15] 

""" 

cdef RealBall result = self._new() 

cdef RealBall my_n = self._parent.coerce(n) 

if _do_sig(prec(self)): sig_on() 

arb_rising(result.value, self.value, my_n.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return result 

  

cpdef RealBall psi(self): 

""" 

Compute the digamma function with argument self. 

  

EXAMPLES:: 

  

sage: RBF(1).psi() 

[-0.577215664901533 +/- 3.85e-16] 

""" 

  

cdef RealBall result = self._new() 

if _do_sig(prec(self)): sig_on() 

arb_digamma(result.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return result 

  

def zeta(self, a=None): 

""" 

Return the image of this ball by the Hurwitz zeta function. 

  

For ``a = 1`` (or ``a = None``), this computes the Riemann zeta function. 

  

Use :meth:`RealBallField.zeta` to compute the Riemann zeta function of 

a small integer without first converting it to a real ball. 

  

EXAMPLES:: 

  

sage: RBF(-1).zeta() 

[-0.0833333333333333 +/- 4.26e-17] 

sage: RBF(-1).zeta(1) 

[-0.0833333333333333 +/- 4.26e-17] 

sage: RBF(-1).zeta(2) 

[-1.083333333333333 +/- 4.96e-16] 

""" 

cdef RealBall a_ball 

cdef RealBall res = self._new() 

if a is None: 

if _do_sig(prec(self)): sig_on() 

arb_zeta(res.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

else: 

a_ball = self._parent.coerce(a) 

if _do_sig(prec(self)): sig_on() 

arb_hurwitz_zeta(res.value, self.value, a_ball.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def polylog(self, s): 

""" 

Return the polylogarithm `\operatorname{Li}_s(\mathrm{self})`. 

  

EXAMPLES:: 

  

sage: polylog(0, -1) 

-1/2 

sage: RBF(-1).polylog(0) 

[-0.50000000000000 +/- 1.29e-15] 

sage: polylog(1, 1/2) 

-log(1/2) 

sage: RBF(1/2).polylog(1) 

[0.6931471805599 +/- 5.02e-14] 

sage: RBF(1/3).polylog(1/2) 

[0.44210883528067 +/- 6.7...e-15] 

sage: RBF(1/3).polylog(RLF(pi)) 

[0.34728895057225 +/- 5.51e-15] 

  

TESTS:: 

  

sage: RBF(1/3).polylog(2r) 

[0.366213229977063 +/- 5.85e-16] 

""" 

cdef RealBall s_as_ball 

cdef Integer s_as_Integer 

cdef RealBall res = self._new() 

try: 

s_as_Integer = ZZ.coerce(s) 

if mpz_fits_slong_p(s_as_Integer.value): 

if _do_sig(prec(self)): sig_on() 

arb_polylog_si(res.value, mpz_get_si(s_as_Integer.value), self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

except TypeError: 

pass 

s_as_ball = self._parent.coerce(s) 

if _do_sig(prec(self)): sig_on() 

arb_polylog(res.value, s_as_ball.value, self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

def chebyshev_T(self, n): 

""" 

Evaluate the Chebyshev polynomial of the first kind ``T_n`` at this 

ball. 

  

EXAMPLES:: 

  

sage: RBF(pi).chebyshev_T(0) 

1.000000000000000 

sage: RBF(pi).chebyshev_T(1) # abs tol 1e-16 

[3.141592653589793 +/- 5.62e-16] 

sage: RBF(pi).chebyshev_T(10**20) 

Traceback (most recent call last): 

... 

ValueError: index too large 

sage: RBF(pi).chebyshev_T(-1) 

Traceback (most recent call last): 

... 

ValueError: expected a nonnegative index 

""" 

cdef RealBall res = self._new() 

cdef Integer n_as_Integer = ZZ.coerce(n) 

if mpz_fits_ulong_p(n_as_Integer.value): 

if _do_sig(prec(self)): sig_on() 

arb_chebyshev_t_ui(res.value, mpz_get_ui(n_as_Integer.value), self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

elif n_as_Integer < 0: 

raise ValueError("expected a nonnegative index") 

else: 

raise ValueError("index too large") 

  

def chebyshev_U(self, n): 

""" 

Evaluate the Chebyshev polynomial of the second kind ``U_n`` at this 

ball. 

  

EXAMPLES:: 

  

sage: RBF(pi).chebyshev_U(0) 

1.000000000000000 

sage: RBF(pi).chebyshev_U(1) 

[6.28318530717959 +/- 4.66e-15] 

sage: RBF(pi).chebyshev_U(10**20) 

Traceback (most recent call last): 

... 

ValueError: index too large 

sage: RBF(pi).chebyshev_U(-1) 

Traceback (most recent call last): 

... 

ValueError: expected a nonnegative index 

""" 

cdef RealBall res = self._new() 

cdef Integer n_as_Integer = ZZ.coerce(n) 

if mpz_fits_ulong_p(n_as_Integer.value): 

if _do_sig(prec(self)): sig_on() 

arb_chebyshev_u_ui(res.value, mpz_get_ui(n_as_Integer.value), self.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

elif n_as_Integer < 0: 

raise ValueError("expected a nonnegative index") 

else: 

raise ValueError("index too large") 

  

def agm(self, other): 

""" 

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

  

EXAMPLES:: 

  

sage: RBF(1).agm(1) 

1.000000000000000 

sage: RBF(sqrt(2)).agm(1)^(-1) 

[0.83462684167407 +/- 3.9...e-15] 

""" 

cdef RealBall other_as_ball 

cdef RealBall res = self._new() 

other_as_ball = self._parent.coerce(other) 

if _do_sig(prec(self)): sig_on() 

arb_agm(res.value, self.value, other_as_ball.value, prec(self)) 

if _do_sig(prec(self)): sig_off() 

return res 

  

RBF = RealBallField()