summaryrefslogtreecommitdiff
path: root/mpfr.texi
blob: bf8e402940c846b9ce3f8ccbd614fbddeeca8c2f (plain)
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
\input texinfo    @c -*-texinfo-*-
@c %**start of header
@setfilename mpfr.info
@documentencoding ISO-8859-1
@set VERSION 2.3.0-dev
@set UPDATED-MONTH February 2007
@settitle MPFR @value{VERSION}
@synindex tp fn
@iftex
@afourpaper
@end iftex
@comment %**end of header

@copying
This manual documents how to install and use the Multiple Precision
Floating-Point Reliable Library, version @value{VERSION}.

Copyright 1991, 1993, 1994, 1995, 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007 Free Software Foundation, Inc.

Permission is granted to copy, distribute and/or modify this document under
the terms of the GNU Free Documentation License, Version 1.1 or any later
version published by the Free Software Foundation; with no Invariant Sections,
with the Front-Cover Texts being ``A GNU Manual'', and with the Back-Cover
Texts being ``You have freedom to copy and modify this GNU Manual, like GNU
software''.  A copy of the license is included in @ref{GNU Free Documentation
License}.
@end copying


@c  Texinfo version 4.2 or up will be needed to process this file.
@c
@c  A suitable texinfo.tex is supplied, a newer one should work
@c  equally well.
@c
@c  The edition number is in the VERSION variable above and should be
@c  updated where appropriate.  Also, update the month and year in
@c  UPDATED-MONTH.


@dircategory GNU libraries
@direntry
* mpfr: (mpfr).                 Multiple Precision Floating-Point Reliable Library.
@end direntry

@c  html <meta name=description content="...">
@documentdescription
How to install and use MPFR, a library for reliable multiple precision floating-point arithmetic, version @value{VERSION}.
@end documentdescription

@c smallbook
@finalout
@setchapternewpage on

@ifnottex
@node Top, Copying, (dir), (dir)
@top MPFR
@end ifnottex

@iftex
@titlepage
@title MPFR
@subtitle The Multiple Precision Floating-Point Reliable Library
@subtitle Edition @value{VERSION}
@subtitle @value{UPDATED-MONTH}

@author The MPFR team
@email{mpfr@@loria.fr}

@c Include the Distribution inside the titlepage so
@c that headings are turned off.

@tex
\global\parindent=0pt
\global\parskip=8pt
\global\baselineskip=13pt
@end tex

@page
@vskip 0pt plus 1filll
@end iftex

@insertcopying
@ifnottex
@sp 1
@end ifnottex

@iftex
@end titlepage
@headings double
@end iftex

@c  Don't bother with contents for html, the menus seem adequate.
@ifnothtml
@contents
@end ifnothtml

@menu
* Copying::                     MPFR Copying Conditions (LGPL).
* Introduction to MPFR::        Brief introduction to MPFR.
* Installing MPFR::             How to configure and compile the MPFR library.
* Reporting Bugs::              How to usefully report bugs.
* MPFR Basics::                 What every MPFR user should now.
* MPFR Interface::              MPFR functions and macros.
* Contributors::
* References::
* GNU Free Documentation License::
* Concept Index::
* Function Index::
@end menu


@c  @m{T,N} is $T$ in tex or @math{N} otherwise.  This is an easy way to give
@c  different forms for math in tex and info.  Commas in N or T don't work,
@c  but @C{} can be used instead.  \, works in info but not in tex.
@iftex
@macro m {T,N}
@tex$\T\$@end tex
@end macro
@end iftex
@ifnottex
@macro m {T,N}
@math{\N\}
@end macro
@end ifnottex

@c  Usage: @GMPabs{x}
@c  Give either |x| in tex, or abs(x) in info or html.
@tex
\gdef\GMPabs#1{|#1|}
@end tex
@ifnottex
@macro GMPabs {X}
@abs{}(\X\)
@end macro
@end ifnottex

@c  Usage: @GMPtimes{}
@c  Give either \times or the word "times".
@tex
\gdef\GMPtimes{\times}
@end tex
@ifnottex
@macro GMPtimes
times
@end macro
@end ifnottex

@c  New math operators.
@c  @abs{} can be used in both tex and info, or just \abs in tex.
@tex
\gdef\abs{\mathop{\rm abs}}
@end tex
@ifnottex
@macro abs
abs
@end macro
@end ifnottex

@c  @times{} made available as a "*" in info and html (already works in tex).
@ifnottex
@macro times
*
@end macro
@end ifnottex

@c  Math operators already available in tex, made available in info too.
@c  For example @log{} can be used in both tex and info.
@ifnottex
@macro le
<=
@end macro
@macro ge
>=
@end macro
@macro ne
<>
@end macro
@macro log
log
@end macro
@end ifnottex

@c  @pom{} definition
@tex
\gdef\pom{\ifmmode\pm\else$\pm$\fi}
@end tex
@ifnottex
@macro pom
±
@end macro
@end ifnottex

@node Copying, Introduction to MPFR, Top, Top
@comment  node-name, next, previous,  up
@unnumbered MPFR Copying Conditions
@cindex Copying conditions
@cindex Conditions for copying MPFR

This library is @dfn{free}; this means that everyone is free to use it and
free to redistribute it on a free basis.  The library is not in the public
domain; it is copyrighted and there are restrictions on its distribution, but
these restrictions are designed to permit everything that a good cooperating
citizen would want to do.  What is not allowed is to try to prevent others
from further sharing any version of this library that they might get from
you.@refill

Specifically, we want to make sure that you have the right to give away copies
of the library, that you receive source code or else can get it if you want
it, that you can change this library or use pieces of it in new free programs,
and that you know you can do these things.@refill

To make sure that everyone has such rights, we have to forbid you to deprive
anyone else of these rights.  For example, if you distribute copies of the
MPFR library, you must give the recipients all the rights that you have.  You
must make sure that they, too, receive or can get the source code.  And you
must tell them their rights.@refill

Also, for our own protection, we must make certain that everyone finds out
that there is no warranty for the MPFR library.  If it is modified by
someone else and passed on, we want their recipients to know that what they
have is not what we distributed, so that any problems introduced by others
will not reflect on our reputation.@refill

The precise conditions of the license for the MPFR library are found in the
Lesser General Public License that accompanies the source code.
See the file COPYING.LIB.@refill

@node Introduction to MPFR, Installing MPFR, Copying, Top
@comment  node-name,  next,  previous,  up
@chapter Introduction to MPFR


MPFR is a portable library written in C for arbitrary precision arithmetic
on floating-point numbers. It is based on the GNU MP library.
It aims to extend the class of floating-point numbers provided by the
GNU MP library by a precise semantics. The main differences
with the @code{mpf} class from GNU MP are:

@itemize @bullet
@item the @code{mpfr} code is portable, i.e.@: the result of any operation
does not depend (or should not) on the machine word size
@code{mp_bits_per_limb} (32 or 64 on most machines);
@item the precision in bits can be set exactly to any valid value
for each variable (including very small precision);
@item @code{mpfr} provides the four rounding modes from the IEEE 754-1985
standard.
@end itemize

In particular, with a precision of 53 bits, @code{mpfr} should be able
to exactly reproduce all computations with double-precision machine
floating-point
numbers (@code{double} type in C), except the default exponent range
is much wider and subnormal numbers are not implemented but can be emulated.

This version of MPFR is released under the GNU Lesser General Public
License.
It is permitted to link MPFR to non-free programs, as long as when
distributing them the MPFR source code and a means to re-link with a
modified MPFR library is provided.

@section How to use this Manual

Everyone should read @ref{MPFR Basics}.  If you need to install the library
yourself, you need to read @ref{Installing MPFR}, too.

The rest of the manual can be used for later reference, although it is
probably a good idea to glance through it.

@node Installing MPFR, Reporting Bugs, Introduction to MPFR, Top
@comment  node-name,  next,  previous,  up
@chapter Installing MPFR
@cindex Installation

@section How to install

Here are the steps needed to install the library on Unix systems
(more details are provided in the @file{INSTALL} file):

@enumerate
@item
To build MPFR, you first have to install GNU MP
(version 4.1 or higher) on your computer.
You need a C compiler, preferably GCC, but any reasonable compiler should
work.  And you need a standard Unix @samp{make} program, plus some other
standard Unix utility programs.

@item
In the MPFR build directory, type
@samp{./configure}

This will prepare the build and setup the options according to your system.
If you get error messages, you might check that you use the same compiler
and compile options as for GNU MP (see the @file{INSTALL} file).

@item
@samp{make}

This will compile MPFR, and create a library archive file @file{libmpfr.a}.
A dynamic library may be produced too (see configure).

@item
@samp{make check}

This will make sure MPFR was built correctly.
If you get error messages, please
report this to @samp{mpfr@@loria.fr}.  (@xref{Reporting Bugs}, for
information on what to include in useful bug reports.)

@item
@samp{make install}

This will copy the files @file{mpfr.h} and @file{mpf2mpfr.h} to the directory
@file{/usr/local/include}, the file @file{libmpfr.a} to the directory
@file{/usr/local/lib}, and the file @file{mpfr.info} to the directory
@file{/usr/local/share/info} (or if you passed the @samp{--prefix} option to
 @file{configure}, using the prefix directory given as argument to
@samp{--prefix} instead of @file{/usr/local}).
@end enumerate

@section Other make targets

There are some other useful make targets:

@itemize @bullet
@item
@samp{mpfr.info} or @samp{info}

Create an info version of the manual, in @file{mpfr.info}.

@item
@samp{mpfr.dvi} or @samp{dvi}

Create a DVI version of the manual, in @file{mpfr.dvi}.

@item
@samp{mpfr.ps}

Create a Postscript version of the manual, in @file{mpfr.ps}.

@c @item
@c @samp{html}
@c Create a HTML version of the manual, in @file{mpfr.html}.

@item
@samp{clean}

Delete all object files and archive files, but not the configuration files.

@item
@samp{distclean}

Delete all files not included in the distribution.

@item
@samp{uninstall}

Delete all files copied by @samp{make install}.
@end itemize


@section Known Build Problems

MPFR suffers from all bugs from the GNU MP library, plus many more.

Please report other problems to @samp{mpfr@@loria.fr}.
@xref{Reporting Bugs}.
Some bug fixes are available on the MPFR web page
@url{http://www.mpfr.org/}.

@section Getting the Latest Version of MPFR

The latest version of MPFR is available from @url{http://www.mpfr.org/}.

@node Reporting Bugs, MPFR Basics, Installing MPFR, Top
@comment  node-name,  next,  previous,  up
@chapter Reporting Bugs
@cindex Reporting bugs

If you think you have found a bug in the MPFR library, first have a look on the
MPFR web page @url{http://www.mpfr.org/}: perhaps this bug is already known,
in which case you may find there a workaround for it.
Otherwise, please investigate
and report it. We have made this library available to you, and it is not to ask
too much from you, to ask you to report the bugs that you find.

There are a few things you should think about when you put your bug report
together.

You have to send us a test case that makes it possible for us to reproduce the
bug.  Include instructions on how to run the test case.

You also have to explain what is wrong; if you get a crash, or if the results
printed are incorrect and in that case, in what way.

Please include compiler version information in your bug report. This can
be extracted using @samp{cc -V} on some machines, or, if you're using gcc,
@samp{gcc -v}. Also, include the output from @samp{uname -a} and the MPFR
version (the GMP version may be useful too).

If your bug report is good, we will do our best to help you to get a corrected
version of the library; if the bug report is poor, we won't do anything about
it (aside of chiding you to send better bug reports).

Send your bug report to: @samp{mpfr@@loria.fr}.

If you think something in this manual is unclear, or downright incorrect, or if
the language needs to be improved, please send a note to the same address.

@node MPFR Basics, MPFR Interface, Reporting Bugs, Top
@comment  node-name,  next,  previous,  up
@chapter MPFR Basics


@cindex @file{mpfr.h}
All declarations needed to use MPFR are collected in the include file
@file{mpfr.h}.  It is designed to work with both C and C++ compilers.
You should include that file in any program using the MPFR library:
@verbatim
#include <mpfr.h>
@end verbatim

@section Nomenclature and Types

@cindex Floating-point number
@tindex @code{mpfr_t}
@noindent
A @dfn{floating-point number} or @dfn{float} for short, is an arbitrary
precision mantissa with a limited precision exponent. The C data type
for such objects is @code{mpfr_t}. A floating-point number can have
three special values: Not-a-Number (NaN) or plus or minus Infinity. NaN
represents an uninitialized object, the result of an invalid operation
(like 0 divided by 0), or a value that cannot be determined (like
+Infinity minus +Infinity). Moreover, like in the IEEE 754-1985 standard,
zero is signed, i.e.@: there are both +0 and @minus{}0; the behavior
is the same as in the IEEE 754-1985 standard and it is generalized to
the other functions supported by MPFR.


@cindex Precision
@tindex @code{mp_prec_t}
@noindent
The @dfn{precision} is the number of bits used to represent the mantissa
of a floating-point number;
the corresponding C data type is @code{mp_prec_t}.
The precision can be any integer between @code{MPFR_PREC_MIN} and
@code{MPFR_PREC_MAX}. In the current implementation, @code{MPFR_PREC_MIN}
is equal to 2.

@cindex Rounding Modes
@tindex @code{mp_rnd_t}
@noindent
The @dfn{rounding mode} specifies the way to round the result of a
floating-point operation, in case the exact result can not be represented
exactly in the destination mantissa;
the corresponding C data type is @code{mp_rnd_t}.

@cindex Limb
@c @tindex @code{mp_limb_t}
@noindent
A @dfn{limb} means the part of a multi-precision number that fits in a single
word.  (We chose this word because a limb of the human body is analogous to a
digit, only larger, and containing several digits.)  Normally a limb contains
32 or 64 bits.  The C data type for a limb is @code{mp_limb_t}.

@section Function Classes

There is only one class of functions in the MPFR library:

@enumerate
@item
Functions for floating-point arithmetic, with names beginning with
@code{mpfr_}.  The associated type is @code{mpfr_t}.
@end enumerate


@section MPFR Variable Conventions

As a general rule, all MPFR functions expect output arguments before input
arguments.  This notation is based on an analogy with the assignment operator.

MPFR allows you to use the same variable for both input and output in the same
expression.  For example, the main function for floating-point multiplication,
@code{mpfr_mul}, can be used like this: @code{mpfr_mul (x, x, x, rnd_mode)}.
This
computes the square of @var{x} with rounding mode @code{rnd_mode}
and puts the result back in @var{x}.

Before you can assign to an MPFR variable, you need to initialize it by calling
one of the special initialization functions.  When you're done with a
variable, you need to clear it out, using one of the functions for that
purpose.

A variable should only be initialized once, or at least cleared out between
each initialization.  After a variable has been initialized, it may be
assigned to any number of times.

For efficiency reasons, avoid to initialize and clear out a variable in loops.
Instead, initialize it before entering the loop, and clear it out after the
loop has exited.

You don't need to be concerned about allocating additional space for MPFR
variables, since any variable has a mantissa of fixed size.
Hence unless you change its precision, or clear and reinitialize it,
a floating-point variable will have the same allocated space during all its
life.

@section Rounding modes

The following four rounding modes are supported:

@itemize @bullet
@item @code{GMP_RNDN}: round to nearest
@item @code{GMP_RNDZ}: round towards zero
@item @code{GMP_RNDU}: round towards plus infinity
@item @code{GMP_RNDD}: round towards minus infinity
@end itemize

The @samp{round to nearest} mode works as in the IEEE 754-1985 standard: in
case the number to be rounded lies exactly in the middle of two representable
numbers, it is rounded to the one with the least significant bit set to zero.
For example, the number 5/2, which is represented by (10.1) in binary, is
rounded to (10.0)=2 with a precision of two bits, and not to (11.0)=3.
This rule avoids the @dfn{drift} phenomenon mentioned by Knuth in volume 2
of The Art of Computer Programming (Section 4.2.2).

Most MPFR functions take as first argument the destination variable, as
second and following arguments the input variables, as last argument a
rounding mode, and have a return value of type @code{int}, called the
@dfn{ternary value}. The value stored in the destination variable is
correctly rounded, i.e.@: MPFR behaves as if it computed the result with
an infinite precision, then rounded it to the precision of this variable.
The input variables are regarded as exact (in particular, their precision
does not affect the result).

As a consequence, in case of a non-zero real rounded result, the error
on the result is less or equal to 1/2 ulp (unit in the last place) of
the target in the rounding to nearest mode, and less than 1 ulp of the
target in the directed rounding modes (a ulp is the weight of the least
significant represented bit of the target after rounding).
@c Since subnormals are not supported, we must take into account the ulp of
@c the rounded result, not the one of the exact result, for full generality.

Unless documented otherwise, functions returning an @code{int} return
a ternary value.
If the ternary value is zero, it means that the value stored in the
destination variable is the exact result of the corresponding mathematical
function. If the ternary value is positive (resp.@: negative), it means
the value stored in the destination variable is greater (resp.@: lower)
than the exact result. For example with the @code{GMP_RNDU} rounding mode,
the ternary value is usually positive, except when the result is exact, in
which case it is zero. In the case of an infinite result, it is considered
as inexact when it was obtained by overflow, and exact otherwise. A NaN
result (Not-a-Number) always corresponds to an exact return value.
The opposite of a returned ternary value is guaranteed to be representable
in an @code{int}.

Unless documented otherwise, functions returning a @code{1}
(or any other value specified in this manual)
for special cases (like @code{acos(0)}) should return an overflow or
an underflow if @code{1} is not representable in the current exponent range.

@section Exceptions

MPFR supports 5 exception types:

@itemize @bullet

@item Underflow:
An underflow occurs when the exact result of a function is a non-zero
real number and the result obtained after the rounding, assuming an
unbounded exponent range (for the rounding), has an exponent smaller
than the minimum exponent of the current range. In the round-to-nearest
mode, the halfway case is rounded toward zero.

Note: This is not the single definition of the underflow. MPFR chooses
to consider the underflow after rounding. The underflow before rounding
can also be defined. For instance, consider a function that has the
exact result @m{7 \times 2^{e-4}, 7 multiplied by two to the power
@var{e}@minus{}4}, where @var{e} is the smallest exponent (for a
mantissa between 1/2 and 1) in the current
range, with a 2-bit target precision and rounding towards plus infinity.
The exact result has the exponent @var{e}@minus{}1. With the underflow
before rounding, such a function call would yield an underflow, as
@var{e}@minus{}1 is outside the current exponent range. However, MPFR
first considers the rounded result assuming an unbounded exponent range.
The exact result cannot be represented exactly in precision 2, and here,
it is rounded to @m{0.5 @times 2^e, 0.5 times 2 to @var{e}}, which is
representable in the current exponent range. As a consequence, this will
not yield an underflow in MPFR.

@item Overflow:
An overflow occurs when the exact result of a function is a non-zero
real number and the result obtained after the rounding, assuming an
unbounded exponent range (for the rounding), has an exponent larger
than the maximum exponent of the current range. In the round-to-nearest
mode, the result is infinite.

@item Invalid (NaN):
An invalid (or NaN) exception occurs when the result of a function is
a NaN.
@c NaN is defined above. So, we don't say anything more.

@item Inexact:
An inexact exception occurs when the result of a function cannot be
represented exactly and must be rounded.

@item Range error:
A range exception occurs when a function that does not return a MPFR
number (such as comparisons and conversions to an integer) has an
invalid result.

@end itemize

MPFR has a global flag for each exception, which can be cleared, set
or tested by functions described in @ref{Exceptions}.

@node MPFR Interface, Contributors, MPFR Basics, Top
@comment  node-name,  next,  previous,  up
@chapter MPFR Interface
@cindex Floating-point functions
@cindex Float functions

The floating-point functions expect arguments of type @code{mpfr_t}.

The MPFR floating-point functions have an interface that is similar to the
GNU MP
integer functions.  The function prefix for floating-point operations is
@code{mpfr_}.

There is one significant characteristic of floating-point numbers that has
motivated a difference between this function class and other GNU MP function
classes: the inherent inexactness of floating-point arithmetic.  The user has
to specify the precision for each variable.  A computation that assigns a
variable will take place with the precision of the assigned variable; the
cost of that computation should not depend from the
precision of variables used as input (on average).

@cindex Precision
The semantics of a calculation in MPFR is specified as follows: Compute the
requested operation exactly (with ``infinite accuracy''), and round the result
to the precision of the destination variable, with the given rounding mode.
The MPFR floating-point functions are intended to be a smooth extension
of the IEEE 754-1985 arithmetic. The results obtained on one computer should
not differ from the results obtained on a computer with a different word size.

@cindex Accuracy
MPFR does not keep track of the accuracy of a computation. This is left
to the user or to a higher layer.
As a consequence, if two variables are used to store
only a few significant bits, and their product is stored in a variable with large
precision, then MPFR will still compute the result with full precision.

The value of errno may be set to non-zero by any MPFR function or macro,
whether or not there is an error.

@menu
* Initialization Functions::
* Assignment Functions::
* Combined Initialization and Assignment Functions::
* Conversion Functions::
* Basic Arithmetic Functions::
* Comparison Functions::
* Special Functions::
* Input and Output Functions::
* Integer Related Functions::
* Miscellaneous Functions::
* Rounding Modes::
* Exceptions::
* Advanced Functions::
* Compatibility with MPF::
* Custom Interface::
* Internals::
@end menu

@node Initialization Functions, Assignment Functions, MPFR Interface, MPFR Interface
@comment  node-name,  next,  previous,  up
@cindex Initialization functions
@section Initialization Functions

An @code{mpfr_t} object must be initialized before storing the first value in
it.  The functions @code{mpfr_init} and @code{mpfr_init2} are used for that
purpose.

@deftypefun void mpfr_init2 (mpfr_t @var{x}, mp_prec_t @var{prec})
Initialize @var{x}, set its precision to be @strong{exactly}
@var{prec} bits and its value to NaN. (Warning: the corresponding
@code{mpf} functions initialize to zero instead.)

Normally, a variable should be initialized once only or at
least be cleared, using @code{mpfr_clear}, between initializations.
To change the precision of a variable which has already been initialized,
use @code{mpfr_set_prec}.
The precision @var{prec} must be an integer between @code{MPFR_PREC_MIN} and
@code{MPFR_PREC_MAX} (otherwise the behavior is undefined).
@end deftypefun

@deftypefun void mpfr_clear (mpfr_t @var{x})
Free the space occupied by @var{x}.  Make sure to call this function for all
@code{mpfr_t} variables when you are done with them.
@end deftypefun

@deftypefun void mpfr_init (mpfr_t @var{x})
Initialize @var{x} and set its value to NaN.

Normally, a variable should be initialized once only
or at least be cleared, using @code{mpfr_clear}, between initializations.  The
precision of @var{x} is the default precision, which can be changed
by a call to @code{mpfr_set_default_prec}.
@end deftypefun

@deftypefun void mpfr_set_default_prec (mp_prec_t @var{prec})
Set the default precision to be @strong{exactly} @var{prec} bits.  The
precision of a variable means the number of bits used to store its mantissa.
All
subsequent calls to @code{mpfr_init} will use this precision, but previously
initialized variables are unaffected.
This default precision is set to 53 bits initially.
The precision can be any integer between @code{MPFR_PREC_MIN} and
@code{MPFR_PREC_MAX}.
@end deftypefun

@deftypefun mp_prec_t mpfr_get_default_prec (void)
Return the default MPFR precision in bits.
@end deftypefun

@need 2000
Here is an example on how to initialize floating-point variables:

@example
@{
  mpfr_t x, y;
  mpfr_init (x);			/* use default precision */
  mpfr_init2 (y, 256);		/* precision @emph{exactly} 256 bits */
  @dots{}
  /* When the program is about to exit, do ... */
  mpfr_clear (x);
  mpfr_clear (y);
@}
@end example

The following functions are useful for changing the precision during a
calculation.  A typical use would be for adjusting the precision gradually in
iterative algorithms like Newton-Raphson, making the computation precision
closely match the actual accurate part of the numbers.

@deftypefun void mpfr_set_prec (mpfr_t @var{x}, mp_prec_t @var{prec})
Reset the precision of @var{x} to be @strong{exactly} @var{prec} bits,
and set its value to NaN.
The previous value stored in @var{x} is lost. It is equivalent to
a call to @code{mpfr_clear(x)} followed by a call to
@code{mpfr_init2(x, prec)}, but more efficient as no allocation is done in
case the current allocated space for the mantissa of @var{x} is enough.
The precision @var{prec} can be any integer between @code{MPFR_PREC_MIN} and
@code{MPFR_PREC_MAX}.

In case you want to keep the previous value stored in @var{x},
use @code{mpfr_prec_round} instead.
@end deftypefun

@deftypefun mp_prec_t mpfr_get_prec (mpfr_t @var{x})
Return the precision actually used for assignments of @var{x}, i.e.@: the
number of bits used to store its mantissa.
@end deftypefun

@node Assignment Functions, Combined Initialization and Assignment Functions, Initialization Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@cindex Assignment functions
@section Assignment Functions

These functions assign new values to already initialized floats
(@pxref{Initialization Functions}). When using any functions using
@code{intmax_t}, you must include @code{<stdint.h>} or @code{<inttypes.h>}
before @file{mpfr.h}, to allow @file{mpfr.h} to define prototypes for
these functions.

@deftypefun int mpfr_set (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_ui (mpfr_t @var{rop}, unsigned long int @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_si (mpfr_t @var{rop}, long int @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_uj (mpfr_t @var{rop}, uintmax_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_sj (mpfr_t @var{rop}, intmax_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_d (mpfr_t @var{rop}, double @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_ld (mpfr_t @var{rop}, long double @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_decimal64 (mpfr_t @var{rop}, _Decimal64 @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_z (mpfr_t @var{rop}, mpz_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_q (mpfr_t @var{rop}, mpq_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_f (mpfr_t @var{rop}, mpf_t @var{op}, mp_rnd_t @var{rnd})
Set the value of @var{rop} from @var{op}, rounded
towards the given direction @var{rnd}.
Note that the input 0 is converted to +0 by @code{mpfr_set_ui},
@code{mpfr_set_si}, @code{mpfr_set_sj}, @code{mpfr_set_uj},
@code{mpfr_set_z}, @code{mpfr_set_q} and
@code{mpfr_set_f}, regardless of the rounding mode.
If the system doesn't support the IEEE-754 standard, @code{mpfr_set_d},
@code{mpfr_set_ld} and 
@code{mpfr_set_decimal64} might not preserve the signed zeros.
The @code{mpfr_set_decimal64} function is built only with the configure
option @samp{--enable-decimal-float}, and when the compiler or
system provides the @code{_Decimal64} data type.
@code{mpfr_set_q} might not be able to work if the numerator (or the
denominator) can not be representable as a @code{mpfr_t}.

Note: If you want to store a floating-point constant to a @code{mpfr_t},
you should use @code{mpfr_set_str} (or one of the MPFR constant functions,
such as @code{mpfr_const_pi} for @m{\pi,Pi}) instead of @code{mpfr_set_d},
@code{mpfr_set_ld} or @code{mpfr_set_decimal64}.
Otherwise the floating-point constant will be first
converted into a reduced-precision (e.g., 53-bit) binary number before
MPFR can work with it.
@end deftypefun

@deftypefun int mpfr_set_ui_2exp (mpfr_t @var{rop}, unsigned long int @var{op}, mp_exp_t @var{e}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_si_2exp (mpfr_t @var{rop}, long int @var{op}, mp_exp_t @var{e}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_uj_2exp (mpfr_t @var{rop}, uintmax_t @var{op}, intmax_t @var{e}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_set_sj_2exp (mpfr_t @var{rop}, intmax_t @var{op}, intmax_t @var{e}, mp_rnd_t @var{rnd})
Set the value of @var{rop} from @m{@var{op} \times 2^e, @var{op} multiplied by
two to the power @var{e}}, rounded towards the given direction @var{rnd}.
Note that the input 0 is converted to +0.
@end deftypefun

@deftypefun int mpfr_set_str (mpfr_t @var{rop}, const char *@var{s}, int @var{base}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the whole string @var{s} in base @var{base},
rounded in the direction @var{rnd}.
See the documentation of @code{mpfr_strtofr} for a detailed description
of the valid string formats.
@c Additionally, special values
@c @code{@@NaN@@}, @code{@@Inf@@}, @code{+@@Inf@@} and @code{-@@Inf@@},
@c all case insensitive, without leading whitespace and possibly followed by
@c other characters, are accepted too (it may change).
This function returns 0 if the entire string up to the final null character
is a valid number in base @var{base}; otherwise it returns @minus{}1, and
@var{rop} may have changed.
@end deftypefun

@deftypefun int mpfr_strtofr (mpfr_t @var{rop}, const char *@var{nptr}, char **@var{endptr}, int @var{base}, mp_rnd_t @var{rnd})

Read a floating-point number from a string @var{nptr} in base @var{base},
rounded in the direction @var{rnd}. If successful, the
result is stored in @var{rop} and @code{*@var{endptr}} points to the
character just after those parsed. If @var{str} doesn't start with a
valid number then @var{rop} is set to zero and the value of @var{nptr}
is stored in the location referenced by @var{endptr}.

Parsing follows the standard C @code{strtod} function. This means optional
leading whitespace, an optional @code{+} or @code{-}, mantissa digits with
an optional decimal point, and an
optional exponent consisting of an @code{e} or @code{E} (if
@math{@var{base} @le{} 10}) or @code{@@}, an optional sign, and digits.
The decimal point can be either the one defined by the current locale or
the period (the first one is accepted for consistency with the C standard
and the practice, the second one is accepted to allow the programmer to
provide MPFR numbers from strings in a way that does not depend on the
current locale).
A hexadecimal mantissa can be given with a leading @code{0x} or @code{0X}, in
which case @code{p} or @code{P} may introduce an optional binary exponent,
indicating the power of 2 by which the mantissa is to be scaled. A binary
mantissa can be given with a leading @code{0b} or @code{0B}, in which case
@code{e}, @code{E}, @code{p}, @code{P} or @code{@@} may introduce the
binary exponent. The exponent is always written in base 10.

In addition, @code{infinity}, @code{inf} (if @math{@var{base} @le{} 10})
or @code{@@inf@@} with an optional sign, or @code{nan},
@code{nan(n-char-sequence)} (if @math{@var{base} @le{} 10}), @code{@@nan@@}
or @code{@@nan@@(n-char-sequence)} all case insensitive (as Latin letters),
can be given. A @code{n-char-sequence} is a non-empty string containing
only digits, Latin letters and the underscore (0, 1, 2, ..., 9, a, b, ...,
z, A, B, ..., Z, _).

There must be at least one digit in the mantissa for the number to
be valid. If an exponent has no digits it's ignored and parsing
stops after the mantissa. If an @code{0x}, @code{0X}, @code{0b} or
@code{0B} is not followed by hexadecimal/binary digits, parsing stops
after the first @code{0}:
the subject sequence is defined as the longest initial
subsequence of the input string, starting with the first
non-white-space character, that is of the expected form.
The subject sequence contains no characters if the input
string is not of the expected form.

Note that in the hex format the exponent @code{P} represents a power of 2,
whereas @code{@@} represents a power of the base (i.e. 16).

If the argument @var{base} is different from 0, it must be in the range 2
to 36.
@c For bases up to 36,
Case is ignored; uppercase and lowercase letters have the same
value.
@c ; for bases 37 to 62, uppercase letters represent the usual
@c  10..35, while lowercase letters represent 36..61.

If @code{base} is 0, then it tries to identify the used base: if the
mantissa begins with the @code{0x} prefix, it assumes that @var{base} is 16.
If it begins with @code{0b}, it assumes that @var{base} is 2. Otherwise, it
assumes it is 10.

It returns a usual ternary value.
If @var{endptr} is not a null pointer, a pointer to the character after
the last character used in the conversion is stored in the location
referenced by @var{endptr}.

@end deftypefun

@deftypefun void mpfr_set_inf (mpfr_t @var{x}, int @var{sign})
@deftypefunx void mpfr_set_nan (mpfr_t @var{x})
Set the variable @var{x} to infinity or NaN (Not-a-Number) respectively.
In @code{mpfr_set_inf}, @var{x} is set to plus infinity iff @var{sign} is
nonnegative.
@end deftypefun

@deftypefun void mpfr_swap (mpfr_t @var{x}, mpfr_t @var{y})
Swap the values @var{x} and @var{y} efficiently. Warning: the
precisions are exchanged too; in case the precisions are different,
@code{mpfr_swap} is thus not equivalent to three @code{mpfr_set} calls
using a third auxiliary variable.
@end deftypefun

@node Combined Initialization and Assignment Functions, Conversion Functions, Assignment Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@cindex Combined initialization and assignment functions
@section Combined Initialization and Assignment Functions

@deftypefn Macro int mpfr_init_set (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefnx Macro int mpfr_init_set_ui (mpfr_t @var{rop}, unsigned long int @var{op}, mp_rnd_t @var{rnd})
@deftypefnx Macro int mpfr_init_set_si (mpfr_t @var{rop}, signed long int @var{op}, mp_rnd_t @var{rnd})
@deftypefnx Macro int mpfr_init_set_d (mpfr_t @var{rop}, double @var{op}, mp_rnd_t @var{rnd})
@deftypefnx Macro int mpfr_init_set_ld (mpfr_t @var{rop}, long double @var{op}, mp_rnd_t @var{rnd})
@deftypefnx Macro int mpfr_init_set_z (mpfr_t @var{rop}, mpz_t @var{op}, mp_rnd_t @var{rnd})
@deftypefnx Macro int mpfr_init_set_q (mpfr_t @var{rop}, mpq_t @var{op}, mp_rnd_t @var{rnd})
@deftypefnx Macro int mpfr_init_set_f (mpfr_t @var{rop}, mpf_t @var{op}, mp_rnd_t @var{rnd})
Initialize @var{rop} and set its value from @var{op}, rounded in the direction
@var{rnd}.
The precision of @var{rop} will be taken from the active default precision,
as set by @code{mpfr_set_default_prec}.
@end deftypefn

@deftypefun int mpfr_init_set_str (mpfr_t @var{x}, const char *@var{s}, int @var{base}, mp_rnd_t @var{rnd})
Initialize @var{x} and set its value from
the string @var{s} in base @var{base},
rounded in the direction @var{rnd}.
See @code{mpfr_set_str}.
@end deftypefun

@node Conversion Functions, Basic Arithmetic Functions, Combined Initialization and Assignment Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@cindex Conversion functions
@section Conversion Functions

@deftypefun double mpfr_get_d (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx {long double} mpfr_get_ld (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx _Decimal64 mpfr_get_decimal64 (mpfr_t @var{op}, mp_rnd_t @var{rnd})
Convert @var{op} to a @code{double} (respectively @code{_Decimal64} or
@code{long double}), using the rounding mode @var{rnd}.
If the system doesn't support IEEE 754 standard, these functions
might not preserve the signed zeros.
The @code{mpfr_get_decimal64} function is built only with the configure
option @samp{--enable-decimal-float}, and when the compiler or
system provides the @code{_Decimal64} data type.
@end deftypefun

@deftypefun double mpfr_get_d_2exp (long *@var{exp}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx {long double} mpfr_get_ld_2exp (long *@var{exp}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Return @var{d} and set @var{exp} such that @math{0.5@le{}@GMPabs{@var{d}}<1}
and @m{@var{d}\times 2^{exp}, @var{d} times 2 raised to @var{exp}} equals
@var{op} rounded to double (resp. long double)
precision, using the given rounding mode.
@comment See ISO C standard, frexp function.
If @var{op} is zero, then a zero of the same sign (or an unsigned zero,
if the implementation does not have signed zeros) is returned, and
@var{exp} is set to 0.
If @var{op} is NaN or an infinity, then the corresponding double precision
(resp. long-double precision)
value is returned, and @var{exp} is undefined.
@end deftypefun

@deftypefun long mpfr_get_si (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx {unsigned long} mpfr_get_ui (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx intmax_t mpfr_get_sj (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx uintmax_t mpfr_get_uj (mpfr_t @var{op}, mp_rnd_t @var{rnd})
Convert @var{op} to a @code{long}, an @code{unsigned long},
an @code{intmax_t} or an @code{uintmax_t} (respectively) after rounding
it with respect to @var{rnd}.
If @var{op} is NaN, the result is undefined.
If @var{op} is too big for the return type, it returns the maximum
or the minimum of the corresponding C type, depending on the direction
of the overflow. The flag erange is set too.
See also @code{mpfr_fits_slong_p}, @code{mpfr_fits_ulong_p},
@code{mpfr_fits_intmax_p} and @code{mpfr_fits_uintmax_p}.
@end deftypefun

@deftypefun mp_exp_t mpfr_get_z_exp (mpz_t @var{rop}, mpfr_t @var{op})
Put the scaled mantissa of @var{op} (regarded as an integer, with the
precision of @var{op}) into @var{rop}, and return the exponent @var{exp}
(which may be outside the current exponent range) such that @var{op}
exactly equals
@ifnottex
@var{rop} multiplied by two exponent @var{exp}.
@end ifnottex
@tex
$rop \times 2^{\rm exp}$.
@end tex
If the exponent is not representable in the @code{mp_exp_t} type, the
behavior is undefined.
@end deftypefun

@deftypefun void mpfr_get_z (mpz_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Convert @var{op} to a @code{mpz_t}, after rounding it with respect to
@var{rnd}. If @var{op} is NaN or Inf, the result is undefined.
@end deftypefun

@deftypefun int mpfr_get_f (mpf_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Convert @var{op} to a @code{mpf_t}, after rounding it with respect to
@var{rnd}. Return zero iff no error occurred,
in particular a non-zero value is returned if
@var{op} is NaN or Inf, which do not exist in @code{mpf}.
@end deftypefun

@deftypefun {char *} mpfr_get_str (char *@var{str}, mp_exp_t *@var{expptr}, int @var{b}, size_t @var{n}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Convert @var{op} to a string of digits in base @var{b}, with rounding in
the direction @var{rnd}, where @var{n} is either zero (see below) or the
number of significant digits; in the latter case, @var{n} must be greater
or equal to 2. The base may vary from 2 to 36.

The generated string is a fraction, with an implicit radix point immediately
to the left of the first digit.  For example, the number @minus{}3.1416 would
be returned as "@minus{}31416" in the string and 1 written at @var{expptr}.
If @var{rnd} is to nearest, and @var{op} is exactly in the middle of two
possible outputs, the one with an even last digit is chosen
(for an odd base, this may not correspond to an even mantissa).

If @var{n} is zero, the number of digits of the mantissa is chosen
large enough so that re-reading the printed value with the same precision,
assuming both output and input use rounding to nearest, will recover
the original value of @var{op}.
More precisely, in most cases, the chosen precision of @var{str} is
the minimal precision depending on @var{n} and @var{b} only that
satisfies the above property, i.e.,
@ifnottex
m = 1 + ceil(@var{n}*log(2)/log(@var{b})),
@end ifnottex
@tex
$m = 1 + \lceil n {\log 2 \over \log b} \rceil$,
@end tex
but in some very rare cases, it might be @math{m+1}.

If @var{str} is a null pointer, space for the mantissa is allocated using
the current allocation function, and a pointer to the string is returned.
To free the returned string, you must use @code{mpfr_free_str}.

If @var{str} is not a null pointer, it should point to a block of storage
large enough for the mantissa, i.e., at least @code{max(@var{n} + 2, 7)}.
The extra two bytes are for a possible minus sign, and for the terminating null
character.

If the input number is an ordinary number, the exponent is written through
the pointer @var{expptr} (the current minimal exponent for 0).

A pointer to the string is returned, unless there is an error, in which
case a null pointer is returned.
@end deftypefun

@deftypefun void mpfr_free_str (char *@var{str})
Free a string allocated by @code{mpfr_get_str} using the current unallocation
function (preliminary interface).
The block is assumed to be @code{strlen(@var{str})+1} bytes.
For more information about how it is done:
@ifnothtml
@pxref{Custom Allocation,,, gmp,GNU MP}.
@end ifnothtml
@ifhtml
see Custom Allocation (GNU MP).
@end ifhtml
@end deftypefun

@deftypefun int mpfr_fits_ulong_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_fits_slong_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_fits_uint_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_fits_sint_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_fits_ushort_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_fits_sshort_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_fits_intmax_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_fits_uintmax_p (mpfr_t @var{op}, mp_rnd_t @var{rnd})
Return non-zero if @var{op} would fit in the respective C data type, when
rounded to an integer in the direction @var{rnd}.
@end deftypefun

@node Basic Arithmetic Functions, Comparison Functions, Conversion Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@section Basic Arithmetic Functions
@cindex Basic arithmetic functions
@cindex Float arithmetic functions
@cindex Arithmetic functions

@deftypefun int mpfr_add (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_add_ui (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_add_si (mpfr_t @var{rop}, mpfr_t @var{op1}, long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_add_z (mpfr_t @var{rop}, mpfr_t @var{op1}, mpz_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_add_q (mpfr_t @var{rop}, mpfr_t @var{op1}, mpq_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to @math{@var{op1} + @var{op2}} rounded in the direction
@var{rnd}. For types having no signed zero, it is considered unsigned
(i.e. (+0) + 0 = (+0) and (@minus{}0) + 0 = (@minus{}0)).
@end deftypefun

@deftypefun int mpfr_sub (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_ui_sub (mpfr_t @var{rop}, unsigned long int @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_sub_ui (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_si_sub (mpfr_t @var{rop}, long int @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_sub_si (mpfr_t @var{rop}, mpfr_t @var{op1}, long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_sub_z (mpfr_t @var{rop}, mpfr_t @var{op1}, mpz_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_sub_q (mpfr_t @var{rop}, mpfr_t @var{op1}, mpq_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to @math{@var{op1} - @var{op2}} rounded in the direction
@var{rnd}. For types having no signed zero, it is considered unsigned
(i.e. (+0) @minus{} 0 = (+0), (@minus{}0) @minus{} 0 = (@minus{}0),
0 @minus{} (+0) = (@minus{}0) and 0 @minus{} (@minus{}0) = (+0)).
@end deftypefun

@deftypefun int mpfr_mul (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_mul_ui (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_mul_si (mpfr_t @var{rop}, mpfr_t @var{op1}, long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_mul_z (mpfr_t @var{rop}, mpfr_t @var{op1}, mpz_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_mul_q (mpfr_t @var{rop}, mpfr_t @var{op1}, mpq_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to @math{@var{op1} @GMPtimes{} @var{op2}} rounded in the
direction @var{rnd}.
When a result is zero, its sign is the product of the signs of the operands
(for types having no signed zero, it is considered positive).
@end deftypefun

@deftypefun int mpfr_sqr (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to @m{@var{op}^{2}, the square of @var{op}}
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_div (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_ui_div (mpfr_t @var{rop}, unsigned long int @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_div_ui (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_si_div (mpfr_t @var{rop}, long int @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_div_si (mpfr_t @var{rop}, mpfr_t @var{op1}, long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_div_z (mpfr_t @var{rop}, mpfr_t @var{op1}, mpz_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_div_q (mpfr_t @var{rop}, mpfr_t @var{op1}, mpq_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to @math{@var{op1}/@var{op2}} rounded in the direction @var{rnd}.
When a result is zero, its sign is the product of the signs of the operands
(for types having no signed zero, it is considered positive).
@end deftypefun

@deftypefun int mpfr_sqrt (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_sqrt_ui (mpfr_t @var{rop}, unsigned long int @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to @m{\sqrt{@var{op}}, the square root of @var{op}}
rounded in the direction @var{rnd}. Return @minus{}0 if @var{op} is
@minus{}0 (to be consistent with the IEEE 754-1985 standard).
Set @var{rop} to NaN if @var{op} is negative.
@end deftypefun

@deftypefun int mpfr_cbrt (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_root (mpfr_t @var{rop}, mpfr_t @var{op}, unsigned long int @var{k}, mp_rnd_t @var{rnd})
Set @var{rop} to the cubic root (resp. the @var{k}th root)
of @var{op} rounded in the direction @var{rnd}.
An odd (resp. even) root of a negative number (including @minus{}Inf)
returns a negative number (resp. NaN).
The @var{k}th root of @minus{}0 is defined to be @minus{}0,
whatever the parity of @var{k}.
@end deftypefun

@deftypefun int mpfr_pow (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_pow_ui (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_pow_si (mpfr_t @var{rop}, mpfr_t @var{op1}, long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_pow_z (mpfr_t @var{rop}, mpfr_t @var{op1}, mpz_t @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_ui_pow_ui (mpfr_t @var{rop}, unsigned long int @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_ui_pow (mpfr_t @var{rop}, unsigned long int @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to @m{@var{op1}^{op2}, @var{op1} raised to @var{op2}},
rounded in the direction @var{rnd}.
Special values are currently handled as described in the ISO C99 standard
for the @code{pow} function (note this may change in future versions):
@itemize @bullet
@item @code{pow(@pom{}0, @var{y})} returns plus or minus infinity for @var{y} a negative odd integer.
@item @code{pow(@pom{}0, @var{y})} returns plus infinity for @var{y} negative and not an odd integer.
@item @code{pow(@pom{}0, @var{y})} returns plus or minus zero for @var{y} a positive odd integer.
@item @code{pow(@pom{}0, @var{y})} returns plus zero for @var{y} positive and not an odd integer.
@item @code{pow(-1, @pom{}inf)} returns 1.
@item @code{pow(+1, @var{y})} returns 1 for any @var{y}, even a NaN.
@item @code{pow(@var{x}, @var{y})} returns NaN for finite negative @var{x} and finite non-integer @var{y}.
@item @code{pow(@var{x}, -inf)} returns plus infinity for @math{0 < @GMPabs{x} < 1}, and plus zero for @math{@GMPabs{x} > 1}.
@item @code{pow(@var{x}, +inf)} returns plus zero for @math{0 < @GMPabs{x} < 1}, and plus infinity for @math{@GMPabs{x} > 1}.
@item @code{pow(-inf, @var{y})} returns minus zero for @var{y} a negative odd integer.
@item @code{pow(-inf, @var{y})} returns plus zero for @var{y} negative and not an odd integer.
@item @code{pow(-inf, @var{y})} returns minus infinity for @var{y} a positive odd integer.
@item @code{pow(-inf, @var{y})} returns plus infinity for @var{y} positive and not an odd integer.
@item @code{pow(+inf, @var{y})} returns plus zero for @var{y} negative, and plus infinity for @var{y} positive.
@end itemize
@end deftypefun

@deftypefun int mpfr_neg (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to @math{-@var{op}} rounded in the direction @var{rnd}.
Just changes the sign if @var{rop} and @var{op} are the same variable.
@end deftypefun

@deftypefun int mpfr_abs (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the absolute value of @var{op},
rounded in the direction @var{rnd}.
Just changes the sign if @var{rop} and @var{op} are the same variable.
@end deftypefun

@deftypefun int mpfr_dim (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to the positive difference of @var{op1} and @var{op2}, i.e.,
@math{@var{op1} - @var{op2}} rounded in the direction @var{rnd}
if @math{@var{op1} > @var{op2}}, and +0 otherwise.
Returns NaN when @var{op1} or @var{op2} is NaN.
@end deftypefun

@deftypefun int mpfr_mul_2ui (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_mul_2si (mpfr_t @var{rop}, mpfr_t @var{op1}, long int @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to @m{@var{op1} \times 2^{op2}, @var{op1} times 2 raised
to @var{op2}}
rounded in the direction @var{rnd}. Just increases the exponent by @var{op2}
when @var{rop} and @var{op1} are identical.
@end deftypefun

@deftypefun int mpfr_div_2ui (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_div_2si (mpfr_t @var{rop}, mpfr_t @var{op1}, long int @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to @m{@var{op1}/2^{op2}, @var{op1} divided by 2 raised
to @var{op2}}
rounded in the direction @var{rnd}. Just decreases the exponent by @var{op2}
when @var{rop} and @var{op1} are identical.
@end deftypefun

@node Comparison Functions, Special Functions, Basic Arithmetic Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@section Comparison Functions
@cindex Float comparisons functions
@cindex Comparison functions

@deftypefun int mpfr_cmp (mpfr_t @var{op1}, mpfr_t @var{op2})
@deftypefunx int mpfr_cmp_ui (mpfr_t @var{op1}, unsigned long int @var{op2})
@deftypefunx int mpfr_cmp_si (mpfr_t @var{op1}, signed long int @var{op2})
@deftypefunx int mpfr_cmp_d (mpfr_t @var{op1}, double @var{op2})
@deftypefunx int mpfr_cmp_ld (mpfr_t @var{op1}, long double @var{op2})
@deftypefunx int mpfr_cmp_z (mpfr_t @var{op1}, mpz_t @var{op2})
@deftypefunx int mpfr_cmp_q (mpfr_t @var{op1}, mpq_t @var{op2})
@deftypefunx int mpfr_cmp_f (mpfr_t @var{op1}, mpf_t @var{op2})
Compare @var{op1} and @var{op2}.  Return a positive value if @math{@var{op1} >
@var{op2}}, zero if @math{@var{op1} = @var{op2}}, and a negative value if
@math{@var{op1} < @var{op2}}.
Both @var{op1} and @var{op2} are considered to their full own precision,
which may differ.
If one of the operands is NaN (Not-a-Number), return zero and set
the erange flag.

Note: These functions may be useful to distinguish the three possible cases.
If you need to distinguish two cases only, it is recommended to use the
predicate functions (e.g., @code{mpfr_equal_p} for the equality) described
below; they behave like the IEEE-754 comparisons, in particular when one
or both arguments are NaN. But only floating-point numbers can be compared
(you may need to do a conversion first).
@end deftypefun

@deftypefun int mpfr_cmp_ui_2exp (mpfr_t @var{op1}, unsigned long int @var{op2}, mp_exp_t @var{e})
@deftypefunx int mpfr_cmp_si_2exp (mpfr_t @var{op1}, long int @var{op2}, mp_exp_t @var{e})
Compare @var{op1} and @m{@var{op2} \times 2^e, @var{op2} multiplied by two to
the power @var{e}}. Similar as above.
@end deftypefun

@deftypefun int mpfr_cmpabs (mpfr_t @var{op1}, mpfr_t @var{op2})
Compare @math{|@var{op1}|} and @math{|@var{op2}|}.  Return a positive value if
@math{|@var{op1}| > |@var{op2}|}, zero if @math{|@var{op1}| = |@var{op2}|}, and
a negative value if @math{|@var{op1}| < |@var{op2}|}.
If one of the operands is NaN (Not-a-Number), return zero and set
the erange flag.
@end deftypefun

@deftypefun int mpfr_nan_p (mpfr_t @var{op})
@deftypefunx int mpfr_inf_p (mpfr_t @var{op})
@deftypefunx int mpfr_number_p (mpfr_t @var{op})
@deftypefunx int mpfr_zero_p (mpfr_t @var{op})
Return non-zero if @var{op} is respectively Not-a-Number (NaN),
an infinity, an ordinary number (i.e.@: neither Not-a-Number nor
an infinity) or zero. Return zero otherwise.
@end deftypefun

@deftypefn Macro int mpfr_sgn (mpfr_t @var{op})
Return a positive value if @math{@var{op} > 0}, zero if @math{@var{op} = 0},
and a negative value if @math{@var{op} < 0}.
Its result is undefined when @var{op} is NaN (Not-a-Number).
@end deftypefn

@deftypefun int mpfr_greater_p (mpfr_t @var{op1}, mpfr_t @var{op2})
Return non-zero if @math{@var{op1} > @var{op2}}, zero otherwise.
@end deftypefun

@deftypefun int mpfr_greaterequal_p (mpfr_t @var{op1}, mpfr_t @var{op2})
Return non-zero if @math{@var{op1} @ge{} @var{op2}}, zero otherwise.
@end deftypefun

@deftypefun int mpfr_less_p (mpfr_t @var{op1}, mpfr_t @var{op2})
Return non-zero if @math{@var{op1} < @var{op2}}, zero otherwise.
@end deftypefun

@deftypefun int mpfr_lessequal_p (mpfr_t @var{op1}, mpfr_t @var{op2})
Return non-zero if @math{@var{op1} @le{} @var{op2}}, zero otherwise.
@end deftypefun

@deftypefun int mpfr_lessgreater_p (mpfr_t @var{op1}, mpfr_t @var{op2})
Return non-zero if @math{@var{op1} < @var{op2}} or
@math{@var{op1} > @var{op2}} (i.e.@: neither @var{op1}, nor @var{op2} is
NaN, and @math{@var{op1} @ne{} @var{op2}}), zero otherwise (i.e.@: @var{op1}
and/or @var{op2} are NaN, or @math{@var{op1} = @var{op2}}).
@end deftypefun

@deftypefun int mpfr_equal_p (mpfr_t @var{op1}, mpfr_t @var{op2})
Return non-zero if @math{@var{op1} = @var{op2}}, zero otherwise
(i.e.@: @var{op1} and/or @var{op2} are NaN, or
@math{@var{op1} @ne{} @var{op2}}).
@end deftypefun

@deftypefun int mpfr_unordered_p (mpfr_t @var{op1}, mpfr_t @var{op2})
Return non-zero if @var{op1} or @var{op2} is a NaN (i.e.@: they cannot be
compared), zero otherwise.
@end deftypefun

@node Special Functions, Input and Output Functions, Comparison Functions, MPFR Interface
@section Special Functions
@cindex Special functions

All those functions, except explicitly stated, return zero for an
exact return value, a positive value for a return value larger than the
exact result, and a negative value otherwise.

@deftypefun int mpfr_log (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_log2 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_log10 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the natural logarithm of @var{op},
@m{\log_2 @var{op}, log2(@var{op})} or
@m{\log_{10} @var{op}, log10(@var{op})}, respectively,
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_exp (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_exp2 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_exp10 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the exponential of @var{op},
 to @m{2^{op}, 2 power of @var{op}}
or to @m{10^{op}, 10 power of @var{op}}, respectively,
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_cos (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_sin (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_tan (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the cosine of @var{op}, sine of @var{op},
tangent of @var{op}, rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_sec (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_csc (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_cot (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the secant of @var{op}, cosecant of @var{op},
cotangent of @var{op}, rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_sin_cos (mpfr_t @var{sop}, mpfr_t @var{cop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set simultaneously @var{sop} to the sine of @var{op} and
                   @var{cop} to the cosine of @var{op},
rounded in the direction @var{rnd} with the corresponding precisions of
@var{sop} and @var{cop}.
Return 0 iff both results are exact.
@end deftypefun

@deftypefun int mpfr_acos (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_asin (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_atan (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the arc-cosine, arc-sine or arc-tangent of @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_atan2 (mpfr_t @var{rop}, mpfr_t @var{y}, mpfr_t @var{x}, mp_rnd_t @var{rnd})
Set @var{rop} to the arc-tangent2 of @var{y} and @var{x},
rounded in the direction @var{rnd}:
if @code{x > 0}, @code{atan2(y, x) = atan (y/x)};
if @code{x < 0}, @code{atan2(y, x) = sign(y)*(Pi - atan (@GMPabs{y/x}))}.

@code{atan2(y, 0)} does not raise any floating-point exception.
Special values are currently handled as described in the ISO C99 standard
for the @code{atan2} function (note this may change in future versions):
@itemize @bullet
@item @code{atan2(+0, -0)} returns @m{+\pi,+Pi}.
@item @code{atan2(-0, -0)} returns @m{-\pi,-Pi}.
@item @code{atan2(+0, +0)} returns +0.
@item @code{atan2(-0, +0)} returns @minus{}0.
@item @code{atan2(+0, x)} returns @m{+\pi,+Pi} for @math{x < 0}.
@item @code{atan2(-0, x)} returns @m{-\pi,-Pi} for @math{x < 0}.
@item @code{atan2(+0, x)} returns +0 for @math{x > 0}.
@item @code{atan2(-0, x)} returns @minus{}0 for @math{x > 0}.
@item @code{atan2(y, 0)} returns @m{-\pi/2,-Pi/2} for @math{y < 0}.
@item @code{atan2(y, 0)} returns @m{+\pi/2,+Pi/2} for @math{y > 0}.
@item @code{atan2(+Inf, -Inf)} returns @m{+3*\pi/4,+3*Pi/4}.
@item @code{atan2(-Inf, -Inf)} returns @m{-3*\pi/4,-3*Pi/4}.
@item @code{atan2(+Inf, +Inf)} returns @m{+\pi/4,+Pi/4}.
@item @code{atan2(-Inf, +Inf)} returns @m{-\pi/4,-Pi/4}.
@item @code{atan2(+Inf, x)} returns @m{+\pi/2,+Pi/2} for finite @math{x}.
@item @code{atan2(-Inf, x)} returns @m{-\pi/2,-Pi/2} for finite @math{x}.
@item @code{atan2(y, -Inf)} returns @m{+\pi,+Pi} for finite @math{y > 0}.
@item @code{atan2(y, -Inf)} returns @m{-\pi,-Pi} for finite @math{y < 0}.
@item @code{atan2(y, +Inf)} returns +0 for finite @math{y > 0}.
@item @code{atan2(y, +Inf)} returns @minus{}0 for finite @math{y < 0}.
@end itemize
@end deftypefun

@deftypefun int mpfr_cosh (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_sinh (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_tanh (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the hyperbolic cosine, sine or tangent of @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_sech (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_csch (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_coth (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the hyperbolic secant of @var{op}, cosecant of @var{op},
cotangent of @var{op}, rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_acosh (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_asinh (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_atanh (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the inverse hyperbolic cosine, sine or tangent of @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_fac_ui (mpfr_t @var{rop}, unsigned long int  @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the factorial of the @code{unsigned long int} @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_log1p (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the logarithm of one plus @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_expm1 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the exponential of @var{op} minus one,
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_eint (mpfr_t @var{y}, mpfr_t @var{x}, mp_rnd_t @var{rnd})
Set @var{y} to the exponential integral of @var{x},
rounded in the direction @var{rnd}.
For positive @var{x},
the exponential integral is the sum of Euler's constant, of the logarithm
of @var{x}, and of the sum for k from 1 to infinity of
@ifnottex
@var{x} to the power k, divided by k and factorial(k).
@end ifnottex
@tex
$x^k/k/k!$.
@end tex
For negative @var{x}, the returned value is NaN.
@end deftypefun

@deftypefun int mpfr_gamma (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the Gamma function on @var{op}, rounded in the
direction @var{rnd}. When @var{op} is a negative integer, NaN is returned.
@end deftypefun

@deftypefun int mpfr_lngamma (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the logarithm of the Gamma function on @var{op},
rounded in the direction @var{rnd}.
When @math{@minus{}2@var{k}@minus{}1 @le{} @var{x} @le{} @minus{}2@var{k}},
@var{k} being a non-negative integer, NaN is returned.
Note that on these negative values, this function is different from the
corresponding function defined in some languages, such as @code{lgamma}
in ISO C99; a new function @code{mpfr_lgamma} will be added in a future
MPFR version to match this alternative definition.
@end deftypefun

@deftypefun int mpfr_zeta (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_zeta_ui (mpfr_t @var{rop}, unsigned long @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the Riemann Zeta function on @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_erf (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the error function on @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_erfc (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the complementary error function on @var{op},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_j0 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_j1 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_jn (mpfr_t @var{rop}, long @var{n}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the first order Bessel function of order 0, 1
and @var{n} on @var{op}, rounded in the direction @var{rnd}. When @var{op} is
NaN, @var{rop} is always set to NaN. When @var{op} is plus or minus Infinity,
@var{rop} is set to +0. When @var{op} is zero, and @var{n} is not zero,
@var{rop} is +0 or @minus{}0 depending on the parity and sign of @var{n},
and the sign of @var{op}.
@end deftypefun

@deftypefun int mpfr_y0 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_y1 (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_yn (mpfr_t @var{rop}, long @var{n}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the value of the second order Bessel function of order 0, 1
and @var{n} on @var{op}, rounded in the direction @var{rnd}. When @var{op} is
NaN or negative,
@var{rop} is always set to NaN. When @var{op} is +Inf,
@var{rop} is +0. When @var{op} is zero,
@var{rop} is +Inf or @minus{}Inf depending on the parity and sign of @var{n}.
@end deftypefun

@deftypefun int mpfr_fma (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mpfr_t @var{op3}, mp_rnd_t @var{rnd})
Set @var{rop} to @math{@var{op1} @GMPtimes{} @var{op2} + @var{op3}},
rounded in the direction @var{rnd}.
@end deftypefun

@deftypefun int mpfr_agm (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to the arithmetic-geometric mean of @var{op1} and @var{op2},
rounded in the direction @var{rnd}.
The arithmetic-geometric mean is the common limit of the sequences
u[n] and v[n], where u[0]=@var{op1}, v[0]=@var{op2}, u[n+1] is the
arithmetic mean of u[n] and v[n], and v[n+1] is the geometric mean of
u[n] and v[n].
If any operand is negative, the return value is NaN.
@end deftypefun

@deftypefun int mpfr_hypot (mpfr_t @var{rop}, mpfr_t @var{x}, mpfr_t @var{y}, mp_rnd_t @var{rnd})
Set @var{rop} to the Euclidean norm of @var{x} and @var{y},
i.e.
@ifnottex
the square root of the sum of the squares of @var{x} and @var{y},
@end ifnottex
@tex
$\sqrt{x^2+y^2}$,
@end tex
rounded in the direction @var{rnd}.
Special values are currently handled as described in Section F.9.4.3 of
the ISO C99 standard, for the @code{hypot} function (note this may change
in future versions): If @var{x} or @var{y} is an infinity, then plus
infinity is returned in @var{rop}, even if the other number is NaN.
@end deftypefun

@deftypefun int mpfr_const_log2 (mpfr_t @var{rop}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_const_pi (mpfr_t @var{rop}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_const_euler (mpfr_t @var{rop}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_const_catalan (mpfr_t @var{rop}, mp_rnd_t @var{rnd})
Set @var{rop} to the logarithm of 2, the value of @m{\pi,Pi},
of Euler's constant 0.577@dots{}, of Catalan's constant 0.915@dots{},
respectively, rounded in the direction
@var{rnd}. These functions cache the computed values to avoid other
calculations if a lower or equal precision is requested. To free these caches,
use @code{mpfr_free_cache}.
@end deftypefun

@deftypefun void mpfr_free_cache (void)
Free the cache used by the functions computing constants if needed
(currently
@code{mpfr_const_log2}, @code{mpfr_const_pi},
@code{mpfr_const_euler} and @code{mpfr_const_catalan}).
@end deftypefun

@deftypefun int mpfr_sum (mpfr_t @var{rop}, mpfr_ptr const @var{tab}[], unsigned long @var{n}, mp_rnd_t @var{rnd})
Set @var{ret} to the sum of all elements of @var{tab} whose size is @var{n},
rounded in the direction @var{rnd}. Warning, @var{tab} is a table of pointers
to mpfr_t, not a table of mpfr_t (preliminary interface). The returned
@code{int} value is zero when the computed value is the exact value,
and non-zero when this cannot be guaranteed, without giving the
direction of the error as the other functions do.
@end deftypefun

@node Input and Output Functions, Integer Related Functions, Special Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@section Input and Output Functions
@cindex Float input and output functions
@cindex Input functions
@cindex Output functions
@cindex I/O functions

This section describes functions that perform input from an input/output
stream, and functions that output to an input/output stream.
Passing a null pointer for a @var{stream} argument to any of
these functions will make them read from @code{stdin} and write to
@code{stdout}, respectively.

When using any of these functions, you must include the @code{<stdio.h>}
standard header before @file{mpfr.h}, to allow @file{mpfr.h} to define
prototypes for these functions.

@deftypefun size_t mpfr_out_str (FILE *@var{stream}, int @var{base}, size_t @var{n}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Output @var{op} on stream @var{stream}, as a string of digits in
base @var{base}, rounded in the direction @var{rnd}.
The base may vary from 2 to 36.  Print @var{n} significant digits exactly,
or if @var{n} is 0, enough digits so that @var{op} can be read back
exactly (see @code{mpfr_get_str}).

In addition to the significant digits, a decimal point (defined by the
current locale) at the right of the
first digit and a trailing exponent in base 10, in the form @samp{eNNN},
are printed. If @var{base} is greater than 10, @samp{@@} will be used
instead of @samp{e} as exponent delimiter.

Return the number of bytes written, or if an error occurred, return 0.
@end deftypefun

@deftypefun size_t mpfr_inp_str (mpfr_t @var{rop}, FILE *@var{stream}, int @var{base}, mp_rnd_t @var{rnd})
Input a string in base @var{base} from stream @var{stream},
rounded in the direction @var{rnd}, and put the
read float in @var{rop}.
@c The argument @var{base} must be in the range 2 to 36.

@c The string is of the form @samp{M@@N} or, if the
@c base is 10 or less, alternatively @samp{MeN} or @samp{MEN}, or, if the base
@c is 16, alternatively @samp{MpB} or @samp{MPB}.
@c @samp{M} is the mantissa in the specified base, @samp{N} is the exponent
@c written in decimal for the specified base, and in base 16, @samp{B} is the
@c binary exponent written in decimal (i.e.@: it indicates the power of 2 by
@c which the mantissa is to be scaled).
This function reads a word (defined as a sequence of characters between
whitespace) and parses it using @code{mpfr_set_str} (it may change).
See the documentation of @code{mpfr_strtofr} for a detailed description
of the valid string formats.
@c Special values can be read as follows (the case does not matter):
@c @code{@@NaN@@}, @code{@@Inf@@}, @code{+@@Inf@@} and @code{-@@Inf@@},
@c possibly followed by other characters; if the base is smaller or equal
@c to 16, the following strings are accepted too: @code{NaN}, @code{Inf},
@c @code{+Inf} and @code{-Inf}.

Return the number of bytes read, or if an error occurred, return 0.
@end deftypefun

@c @deftypefun void mpfr_inp_raw (mpfr_t @var{float}, FILE *@var{stream})
@c Input from stdio stream @var{stream} in the format written by
@c @code{mpfr_out_raw}, and put the result in @var{float}.
@c @end deftypefun

@node Integer Related Functions, Miscellaneous Functions, Input and Output Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@section Integer Related Functions
@cindex Integer Related Functions

@deftypefun int mpfr_rint (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_ceil (mpfr_t @var{rop}, mpfr_t @var{op})
@deftypefunx int mpfr_floor (mpfr_t @var{rop}, mpfr_t @var{op})
@deftypefunx int mpfr_round (mpfr_t @var{rop}, mpfr_t @var{op})
@deftypefunx int mpfr_trunc (mpfr_t @var{rop}, mpfr_t @var{op})
Set @var{rop} to @var{op} rounded to an integer.
@code{mpfr_rint} rounds to the nearest representable integer in the
given rounding mode, @code{mpfr_ceil} rounds
to the next higher or equal representable integer, @code{mpfr_floor} to
the next lower or equal representable integer, @code{mpfr_round} to the
nearest representable integer, rounding halfway cases away from zero,
and @code{mpfr_trunc} to the next representable integer towards zero.

The returned value is zero when the result is exact, positive when it is
greater than the original value of @var{op}, and negative when it is smaller.
More precisely, the returned value is 0 when @var{op} is an integer
representable in @var{rop}, 1 or @minus{}1 when @var{op} is an integer
that is not representable in @var{rop}, 2 or @minus{}2 when @var{op} is
not an integer.

Note that @code{mpfr_round} is different from @code{mpfr_rint} called with
the rounding to the nearest mode (where halfway cases are rounded to an even
integer or mantissa). Note also that no double rounding is performed; for
instance, 4.5 (100.1 in binary) is rounded by @code{mpfr_round} to 4 (100
in binary) in 2-bit precision, though @code{round(4.5)} is equal to 5 and
5 (101 in binary) is rounded to 6 (110 in binary) in 2-bit precision.
@end deftypefun

@deftypefun int mpfr_rint_ceil (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_rint_floor (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_rint_round (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_rint_trunc (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to @var{op} rounded to an integer.
@code{mpfr_rint_ceil} rounds to the next higher or equal integer,
@code{mpfr_rint_floor} to the next lower or equal integer,
@code{mpfr_rint_round} to the nearest integer, rounding halfway cases away
from zero, and @code{mpfr_rint_trunc} to the next integer towards zero.
If the result is not representable, it is rounded in the direction @var{rnd}.
The returned value is the ternary value associated with the considered
round-to-integer function (regarded in the same way as any other
mathematical function).
@end deftypefun

@deftypefun int mpfr_frac (mpfr_t @var{rop}, mpfr_t @var{op}, mp_rnd_t @var{rnd})
Set @var{rop} to the fractional part of @var{op}, having the same sign as
@var{op}, rounded in the direction @var{rnd} (unlike in @code{mpfr_rint},
@var{rnd} affects only how the exact fractional part is rounded, not how
the fractional part is generated).
@end deftypefun

@deftypefun int mpfr_integer_p (mpfr_t @var{op})
Return non-zero iff @var{op} is an integer.
@end deftypefun

@node Miscellaneous Functions, Rounding Modes, Integer Related Functions, MPFR Interface
@comment  node-name,  next,  previous,  up
@section Miscellaneous Functions
@cindex Miscellaneous float functions

@deftypefun void mpfr_nexttoward (mpfr_t @var{x}, mpfr_t @var{y})
If @var{x} or @var{y} is NaN, set @var{x} to NaN. Otherwise, if @var{x}
is different from @var{y}, replace @var{x} by the next floating-point
number (with the precision of @var{x} and the current exponent range)
in the direction of @var{y}, if there is one
(the infinite values are seen as the smallest and largest floating-point
numbers). If the result is zero, it keeps the same sign. No underflow or
overflow is generated.
@end deftypefun

@deftypefun void mpfr_nextabove (mpfr_t @var{x})
Equivalent to @code{mpfr_nexttoward} where @var{y} is plus infinity.
@end deftypefun

@deftypefun void mpfr_nextbelow (mpfr_t @var{x})
Equivalent to @code{mpfr_nexttoward} where @var{y} is minus infinity.
@end deftypefun

@deftypefun int mpfr_min (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to the minimum of @var{op1} and @var{op2}. If @var{op1}
and @var{op2} are both NaN, then @var{rop} is set to NaN. If @var{op1}
or @var{op2} is NaN, then @var{rop} is set to the numeric value. If
@var{op1} and @var{op2} are zeros of different signs, then @var{rop}
is set to @minus{}0.
@end deftypefun

@deftypefun int mpfr_max (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
Set @var{rop} to the maximum of @var{op1} and @var{op2}. If @var{op1}
and @var{op2} are both NaN, then @var{rop} is set to NaN. If @var{op1}
or @var{op2} is NaN, then @var{rop} is set to the numeric value. If
@var{op1} and @var{op2} are zeros of different signs, then @var{rop}
is set to +0.
@end deftypefun

@deftypefun int mpfr_urandomb (mpfr_t @var{rop}, gmp_randstate_t @var{state})
Generate a uniformly distributed random float in the interval
@math{0 @le{} @var{rop} < 1}.
Return 0, unless the exponent is not in the current exponent range, in
which case @var{rop} is set to NaN and a non-zero value is returned. The
second argument is a @code{gmp_randstate_t} structure which should be
created using the GMP @code{gmp_randinit} function, see the GMP manual.
@end deftypefun

@deftypefun void mpfr_random (mpfr_t @var{rop})
Generate a uniformly distributed random float in the interval
@math{0 @le{} @var{rop} < 1}.
This function is deprecated; @code{mpfr_urandomb} should be used instead.
@end deftypefun

@deftypefun void mpfr_random2 (mpfr_t @var{rop}, mp_size_t @var{size}, mp_exp_t @var{exp})
Generate a random float of at most @var{size} limbs, with long strings of
zeros and ones in the binary representation. The exponent of the number is in
the interval @minus{}@var{exp} to @var{exp}.
This function is useful for
testing functions and algorithms, since this kind of random numbers have
proven to be more likely to trigger corner-case bugs.
Negative random numbers are generated when @var{size} is negative.
Put +0 in @var{rop} when size if zero. The internal state of the default
pseudorandom number generator is modified by a call to this function (the
same one as GMP if MPFR was built using @samp{--with-gmp-build}).
@end deftypefun

@deftypefun mp_exp_t mpfr_get_exp (mpfr_t @var{x})
Get the exponent of @var{x}, assuming that @var{x} is a non-zero ordinary
number. The behavior for NaN, Infinity or Zero is undefined.
@end deftypefun

@deftypefun int mpfr_set_exp (mpfr_t @var{x}, mp_exp_t @var{e})
Set the exponent of @var{x} if @var{e} is in the current exponent range,
and return 0 (even if @var{x} is not a non-zero ordinary number);
otherwise, return a non-zero value.
@end deftypefun


@deftypefun {const char *} mpfr_get_version (void)
Return the MPFR version, as a null-terminated string.
@end deftypefun

@defmac MPFR_VERSION
@defmacx MPFR_VERSION_MAJOR
@defmacx MPFR_VERSION_MINOR
@defmacx MPFR_VERSION_PATCHLEVEL
@defmacx MPFR_VERSION_STRING
@code{MPFR_VERSION} is the version of MPFR as a preprocessing constant.
@code{MPFR_VERSION_MAJOR}, @code{MPFR_VERSION_MINOR} and
@code{MPFR_VERSION_PATCHLEVEL} are respectively the major, minor and patch
level of MPFR version, as preprocessing constants.
@code{MPFR_VERSION_STRING} is the version as a string constant, which can
be compared to the result of @code{mpfr_get_version} to check at run time
the header file and library used match:
@example
if (strcmp (mpfr_get_version (), MPFR_VERSION_STRING))
   fprintf (stderr, "Error, header and library files do not match\n");
@end example
@end defmac

@deftypefn Macro long MPFR_VERSION_NUM (@var{major}, @var{minor}, @var{patchlevel})
Create an integer in the same format as used by @code{MPFR_VERSION} from the
given @var{major}, @var{minor} and @var{patchlevel}.
Here is an example of how to check the MPFR version at compile time:
@example
#if (!defined(MPFR_VERSION) || (MPFR_VERSION<MPFR_VERSION_NUM(2,1,0)))
# error "Wrong MPFR version."
#endif
@end example
@end deftypefn

@node Rounding Modes, Exceptions, Miscellaneous Functions, MPFR Interface
@section Rounding Modes
@cindex Rounding Modes

@deftypefun void mpfr_set_default_rounding_mode (mp_rnd_t @var{rnd})
Set the default rounding mode to @var{rnd}.
The default rounding mode is to nearest initially.
@end deftypefun

@deftypefun mp_rnd_t mpfr_get_default_rounding_mode (void)
Get the default rounding mode.
@end deftypefun

@deftypefun int mpfr_prec_round (mpfr_t @var{x}, mp_prec_t @var{prec}, mp_rnd_t @var{rnd})
Round @var{x} according to @var{rnd} with precision @var{prec}, which
must be an integer between @code{MPFR_PREC_MIN} and @code{MPFR_PREC_MAX}
(otherwise the behavior is undefined).
If @var{prec} is greater or equal to the precision of @var{x}, then new
space is allocated for the mantissa, and it is filled with zeros.
Otherwise, the mantissa is rounded to precision @var{prec} with the given
direction. In both cases, the precision of @var{x} is changed to @var{prec}.
@end deftypefun

@deftypefun int mpfr_round_prec (mpfr_t @var{x}, mp_rnd_t @var{rnd}, mp_prec_t @var{prec})
[This function is obsolete. Please use @code{mpfr_prec_round} instead.]
@end deftypefun

@deftypefun {const char *} mpfr_print_rnd_mode (mp_rnd_t @var{rnd})
Return the input string (GMP_RNDD, GMP_RNDU, GMP_RNDN, GMP_RNDZ)
corresponding to the rounding mode @var{rnd} or a null pointer if
@var{rnd} is an invalid rounding mode.
@end deftypefun


@node Exceptions, Advanced Functions, Rounding Modes, MPFR Interface
@comment  node-name,  next,  previous,  up
@section Exceptions
@cindex Exceptions

@deftypefun mp_exp_t mpfr_get_emin (void)
@deftypefunx mp_exp_t mpfr_get_emax (void)
Return the (current) smallest and largest exponents allowed for a
floating-point variable. The smallest positive value of a floating-point
variable is @m{1/2 \times 2^{\rm emin}, one half times 2 raised to the
smallest exponent} and the largest value has the form @m{(1 - \varepsilon)
\times 2^{\rm emax}, (1 - epsilon) times 2 raised to the largest exponent}.
@end deftypefun

@deftypefun int mpfr_set_emin (mp_exp_t @var{exp})
@deftypefunx int mpfr_set_emax (mp_exp_t @var{exp})
Set the smallest and largest exponents allowed for a floating-point variable.
Return a non-zero value when @var{exp} is not in the range accepted by the
implementation (in that case the smallest or largest exponent is not changed),
and zero otherwise.
If the user changes the exponent range, it is her/his responsibility to check
that all current floating-point variables are in the new allowed range
(for example using @code{mpfr_check_range}), otherwise the subsequent
behavior will be undefined, in the sense of the ISO C standard.
@c It is also her/his responsibility to check that @m {emin <= emax}.
@end deftypefun

@deftypefun mp_exp_t mpfr_get_emin_min (void)
@deftypefunx mp_exp_t mpfr_get_emin_max (void)
@deftypefunx mp_exp_t mpfr_get_emax_min (void)
@deftypefunx mp_exp_t mpfr_get_emax_max (void)
Return the minimum and maximum of the smallest and largest exponents
allowed for @code{mpfr_set_emin} and @code{mpfr_set_emax}. These values
are implementation dependent; it is possible to create a non
portable program by writing @code{mpfr_set_emax(mpfr_get_emax_max())}
and @code{mpfr_set_emin(mpfr_get_emin_min())} since the values
of the smallest and largest exponents become implementation dependent.
@end deftypefun

@deftypefun int mpfr_check_range (mpfr_t @var{x}, int @var{t}, mp_rnd_t @var{rnd})
This function forces @var{x} to be in the current range of acceptable
values, @var{t} being the current ternary value: negative if @var{x}
is smaller than the exact value, positive if @var{x} is larger than
the exact value and zero if @var{x} is exact (before the call). It
generates an underflow or an overflow if the exponent of @var{x} is
outside the current allowed range; the value of @var{t} may be used
to avoid a double rounding. This function returns zero if the rounded
result is equal to the exact one, a positive value if the rounded
result is larger than the exact one, a negative value if the rounded
result is smaller than the exact one. Note that unlike most functions,
the result is compared to the exact one, not the input value @var{x},
i.e.@: the ternary value is propagated.
@end deftypefun

@deftypefun int mpfr_subnormalize (mpfr_t @var{x}, int @var{t}, mp_rnd_t @var{rnd})
This function rounds @var{x} emulating subnormal number arithmetic:
if @var{x} is outside the subnormal exponent range, it just propagates the
ternary value @var{t}; otherwise, it rounds @var{x} to precision
@code{EXP(x)-emin+1} according to rounding mode @var{rnd} and previous
ternary value @var{t}, avoiding double rounding problems.
@code{PREC(x)} is not modified by this function.
@var{rnd} and @var{t} must be the used rounding mode for computing @var{x}
and the returned ternary value when computing @var{x}.
The subnormal exponent range is from @code{emin} to @code{emin+PREC(x)-1}.
This functions assumes that @code{emax-emin >= PREC(x)}.
Note that unlike most functions, the result is compared to the exact one,
not the input value @var{x}, i.e.@: the ternary value is propagated.
This is a preliminary interface.
@end deftypefun

This is an example of how to emulate double IEEE-754 arithmetic
using MPFR:

@example
@{
  mpfr_t xa, xb;
  int i;
  volatile double a, b;

  mpfr_set_default_prec (53);
  mpfr_set_emin (-1073);
  mpfr_set_emax (1024);

  mpfr_init (xa); mpfr_init (xb);

  b = 34.3; mpfr_set_d (xb, b, GMP_RNDN);
  a = 0x1.1235P-1021; mpfr_set_d (xa, a, GMP_RNDN);

  a /= b;
  i = mpfr_div (xa, xa, xb, GMP_RNDN);
  i = mpfr_subnormalize (xa, i, GMP_RNDN);

  mpfr_clear (xa); mpfr_clear (xb);
@}
@end example

Warning: this emulates a double IEEE-754 arithmetic with correct rounding
in the subnormal range, which may not be the case for your hardware.

@deftypefun void mpfr_clear_underflow (void)
@deftypefunx void mpfr_clear_overflow (void)
@deftypefunx void mpfr_clear_nanflag (void)
@deftypefunx void mpfr_clear_inexflag (void)
@deftypefunx void mpfr_clear_erangeflag (void)
Clear the underflow, overflow, invalid, inexact and erange flags.
@end deftypefun

@deftypefun void mpfr_set_underflow (void)
@deftypefunx void mpfr_set_overflow (void)
@deftypefunx void mpfr_set_nanflag (void)
@deftypefunx void mpfr_set_inexflag (void)
@deftypefunx void mpfr_set_erangeflag (void)
Set the underflow, overflow, invalid, inexact and erange flags.
@end deftypefun

@deftypefun void mpfr_clear_flags (void)
Clear all global flags (underflow, overflow, inexact, invalid, erange).
@end deftypefun

@deftypefun int mpfr_underflow_p (void)
@deftypefunx int mpfr_overflow_p (void)
@deftypefunx int mpfr_nanflag_p (void)
@deftypefunx int mpfr_inexflag_p (void)
@deftypefunx int mpfr_erangeflag_p (void)
Return the corresponding (underflow, overflow, invalid, inexact, erange)
flag, which is non-zero iff the flag is set.
@end deftypefun

@node Advanced Functions, Compatibility with MPF, Exceptions, MPFR Interface
@comment  node-name,  next,  previous,  up
@section Advanced Functions
@cindex Advanced Functions

All the given interfaces are preliminary. They might change incompatibly in
future revisions.

@defmac MPFR_DECL_INIT (@var{name}, @var{prec})
This macro declares @var{name} as an automatic variable of type @code{mpfr_t},
initializes it and sets its precision to be @strong{exactly} @var{prec} bits
and its value to NaN. @var{name} must be a valid identifier.
You must use this macro in the declaration section.
This macro is much faster than using @code{mpfr_init2} but has some
drawbacks:

@itemize @bullet
@item You @strong{must not} call @code{mpfr_clear} with variables
created with this macro (The storage is allocated at the point of declaration
and deallocated when the brace-level is exited.).
@item You @strong{can not} change their precision.
@item You @strong{should not} create variables with huge precision with this macro.
@item Your compiler must support @samp{Non-Constant Initializers} (standard
in C++ and ISO C99) and @samp{Token Pasting}
(standard in ISO C89). If @var{prec} is not a compiler constant, your compiler
must support @samp{Variable-length automatic arrays} (standard in ISO C99).
@samp{GCC 2.95.3} supports all these features. If you compile your program
with gcc in c89 mode and with @samp{-pedantic}, you may want to define the
@code{MPFR_USE_EXTENSION} macro to avoid warnings due to the
@code{MPFR_DECL_INIT} implementation.
@end itemize
@end defmac

@deftypefun void mpfr_inits (mpfr_t @var{x}, ...)
Initialize all the @code{mpfr_t} variables of the given @code{va_list},
set their precision to be the default precision and their value to NaN.
See @code{mpfr_init} for more details.
The @code{va_list} is assumed to be composed only of type @code{mpfr_t}.
It begins from @var{x}. It ends when it encounters a null pointer.
@end deftypefun


@deftypefun void mpfr_inits2 (mp_prec_t @var{prec}, mpfr_t @var{x}, ...)
Initialize all the @code{mpfr_t} variables of the given @code{va_list},
set their precision to be @strong{exactly}
@var{prec} bits and their value to NaN.
See @code{mpfr_init2} for more details.
The @code{va_list} is assumed to be composed only of type @code{mpfr_t}.
It begins from @var{x}. It ends when it encounters a null pointer.
@end deftypefun

@deftypefun void mpfr_clears (mpfr_t @var{x}, ...)
Free the space occupied by all the @code{mpfr_t} variables of the given
@code{va_list}. See @code{mpfr_clear} for more details.
The @code{va_list} is assumed to be composed only of type @code{mpfr_t}.
It begins from @var{x}. It ends when it encounters a null pointer.
@end deftypefun

Here is an example of how to use multiple initialization functions:

@example
@{
  mpfr_t x, y, z, t;
  mpfr_inits2 (256, x, y, z, t, (void *) 0);
  @dots{}
  mpfr_clears (x, y, z, t, (void *) 0);
@}
@end example


@node Compatibility with MPF, Custom Interface, Advanced Functions, MPFR Interface
@section Compatibility with MPF
@cindex Compatibility with MPF

A header file @file{mpf2mpfr.h} is included in the distribution of MPFR for
compatibility with the GNU MP class MPF.
After inserting the following two lines after the @code{#include <gmp.h>}
line,
@verbatim
#include <mpfr.h>
#include <mpf2mpfr.h>
@end verbatim
@noindent
any program written for
MPF can be compiled directly with MPFR without any changes.
All operations are then performed with the default MPFR rounding mode,
which can be reset with @code{mpfr_set_default_rounding_mode}.

Warning: the @code{mpf_init} and @code{mpf_init2} functions initialize
to zero, whereas the corresponding @code{mpfr} functions initialize to NaN:
this is useful to detect uninitialized values, but is slightly incompatible
with @code{mpf}.

@deftypefun void mpfr_set_prec_raw (mpfr_t @var{x}, mp_prec_t @var{prec})
Reset the precision of @var{x} to be @strong{exactly} @var{prec} bits.
The only difference with @code{mpfr_set_prec} is that @var{prec} is assumed to
be small enough so that the mantissa fits into the current allocated memory
space for @var{x}. Otherwise the behavior is undefined.
@end deftypefun

@deftypefun int mpfr_eq (mpfr_t @var{op1}, mpfr_t @var{op2}, unsigned long int @var{op3})
Return non-zero if @var{op1} and @var{op2} are both non-zero ordinary
numbers with the same exponent and the same first @var{op3} bits, both
zero, or both infinities of the same sign. Return zero otherwise. This
function is defined for compatibility with @code{mpf}. Do not use it if
you want to know whether two numbers are close to each other; for instance,
1.011111 and 1.100000 are regarded as different for any value of @var{op3}
larger than 1.
@end deftypefun

@deftypefun void mpfr_reldiff (mpfr_t @var{rop}, mpfr_t @var{op1}, mpfr_t @var{op2}, mp_rnd_t @var{rnd})
Compute the relative difference between @var{op1} and @var{op2}
and store the result in @var{rop}.
This function does not guarantee the correct rounding on the relative
difference; it just computes @math{|@var{op1}-@var{op2}|/@var{op1}}, using the
rounding mode @var{rnd} for all operations and the precision of @var{rop}.
@end deftypefun

@deftypefun int mpfr_mul_2exp (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
@deftypefunx int mpfr_div_2exp (mpfr_t @var{rop}, mpfr_t @var{op1}, unsigned long int @var{op2}, mp_rnd_t @var{rnd})
See @code{mpfr_mul_2ui} and @code{mpfr_div_2ui}. These functions are only kept
for compatibility with MPF.
@end deftypefun


@node Custom Interface, Internals, Compatibility with MPF, MPFR Interface
@section Custom Interface
@cindex Custom Interface

Some applications use a stack to handle the memory and their objects.
However, the MPFR memory design is not well suited for such a thing. So that
such applications are able to use MPFR, an auxiliary memory interface has
been created: the Custom Interface.

The following interface allows them to use MPFR in two ways:
@itemize
@item Either they directly store the MPFR FP number as a @code{mpfr_t}
on the stack.
@item Either they store their own representation of a FP number on the
stack and construct a new temporary @code{mpfr_t} each time it is needed.
@end itemize
Nothing has to be done to destroy the FP numbers except garbaging the used
memory: all the memory stuff (allocating, destroying, garbaging) is kept to
the application.

Each function is this interface is also implemented as a macro for
efficiency reasons.

Note 1: MPFR functions may still initialize temporary FP numbers
using standard mpfr_init. See Custom Allocation (GNU MP).

Note 2: MPFR functions may use the cached functions (mpfr_const_pi for
example), even if they are not explicitly called. You have to call
@code{mpfr_free_cache} each time you garbage the memory iff mpfr_init, through
GMP Custom Allocation, allocates its memory on the application stack.

Note 3: This interface is preliminary.

@deftypefun size_t mpfr_custom_get_size (mp_prec_t @var{prec})
Return the needed size in bytes to store the mantissa of a FP number
of precision @var{prec}.
@end deftypefun

@deftypefun void mpfr_custom_init (void *@var{mantissa}, mp_prec_t @var{prec})
Initialize a mantissa of precision @var{prec}.
@var{mantissa} must be an area of @code{mpfr_custom_get_size (prec)} bytes
at least and be suitably aligned for an array of @code{mp_limb_t}.
@end deftypefun

@deftypefun void mpfr_custom_init_set (mpfr_t @var{x}, int @var{kind}, mp_exp_t @var{exp}, mp_prec_t @var{prec}, void *@var{mantissa})
Perform a dummy initialization of a @code{mpfr_t} and set it to:
@itemize
@item if @code{ABS(kind) == MPFR_NAN_LIND}, @var{x} is set to NaN;
@item if @code{ABS(kind) == MPFR_INF_KIND}, @var{x} is set to the infinity
of sign @code{sign(kind)};
@item if @code{ABS(kind) == MPFR_ZERO_KIND}, @var{x} is set to the zero of
sign @code{sign(kind)};
@item if @code{ABS(kind) == MPFR_REGULAR_KIND}, @var{x} is set to a regular
number: @code{x = sign(kind)*mantissa*2^exp}
@end itemize
In all cases, it uses @var{mantissa} directly for further computing
involving @var{x}. It will not allocate anything.
A FP number initialized with this function cannot be resized using
@code{mpfr_set_prec}, or cleared using @code{mpfr_clear}!
@var{mantissa} must have been initialized with @code{mpfr_custom_init}
using the same precision @var{prec}.
@end deftypefun

@deftypefun int mpfr_custom_get_kind (mpfr_t @var{x})
Return the current kind of a @code{mpfr_t} as used by
@code{mpfr_custom_init_set}.
The behavior of this function for any @code{mpfr_t} not initialized
with @code{mpfr_custom_init_set} is undefined.
@end deftypefun

@deftypefun void *mpfr_custom_get_mantissa (mpfr_t @var{x})
Return a pointer to the mantissa used by a @code{mpfr_t} initialized with
@code{mpfr_custom_init_set}.
The behavior of this function for any @code{mpfr_t} not initialized
with @code{mpfr_custom_init_set} is undefined.
@end deftypefun

@deftypefun mp_exp_t mpfr_custom_get_exp (mpfr_t @var{x})
Return the exponent of @var{x}, assuming that @var{x} is a non-zero ordinary
number. The return value for NaN, Infinity or Zero is unspecified but doesn't
produced any trap.
The behavior of this function for any @code{mpfr_t} not initialized
with @code{mpfr_custom_init_set} is undefined.
@end deftypefun

@deftypefun void mpfr_custom_move (mpfr_t @var{x}, void *@var{new_position})
Inform MPFR that the mantissa has moved due to a garbage collect
and update its new position to @code{new_position}.
However the application has to move the mantissa and the @code{mpfr_t} itself.
The behavior of this function for any @code{mpfr_t} not initialized
with @code{mpfr_custom_init_set} is undefined.
@end deftypefun

See the test suite for examples.

@node Internals,  , Custom Interface, MPFR Interface
@section Internals
@cindex Internals

The following types and
functions were mainly designed for the implementation of @code{mpfr},
but may be useful for users too.
However no upward compatibility is guaranteed.
You may need to include @file{mpfr-impl.h} to use them.

The @code{mpfr_t} type consists of four fields.

@itemize @bullet

@item The @code{_mpfr_prec} field is used to store the precision of
the variable (in bits); this is not less than @code{MPFR_PREC_MIN}.

@item The @code{_mpfr_sign} field is used to store the sign of the variable.

@item The @code{_mpfr_exp} field stores the exponent.
An exponent of 0 means a radix point just above the most significant
limb.  Non-zero values @math{n} are a multiplier @math{2^n} relative to that
point.
A NaN, an Infinity and a Zero are indicated by a special value of the exponent.

@item Finally, the @code{_mpfr_d} is a pointer to the limbs, least
significant limbs stored first.
The number of limbs in use is controlled by @code{_mpfr_prec}, namely
ceil(@code{_mpfr_prec}/@code{mp_bits_per_limb}).
Non-singular values always have the most significant bit of the most
significant limb set to 1.  When the precision does not correspond to a
whole number of limbs, the excess bits at the low end of the data are zero.

@end itemize

@c @deftypefun int mpfr_add_one_ulp (mpfr_t @var{x}, mp_rnd_t @var{rnd})
@c Add one unit in last place (ulp) to @var{x} if @var{x} is finite
@c and positive, subtract one ulp if @var{x} is finite and negative;
@c otherwise, @var{x} is not changed.
@c The return value is zero unless an overflow occurs, in which case the
@c @code{mpfr_add_one_ulp} function behaves like a conventional addition.
@c @end deftypefun

@c @deftypefun int mpfr_sub_one_ulp (mpfr_t @var{x}, mp_rnd_t @var{rnd})
@c Subtract one ulp to @var{x} if @var{x} is finite and positive, add one
@c ulp if @var{x} is finite and negative; otherwise, @var{x} is not changed.
@c The return value is zero unless an underflow occurs, in which case the
@c @code{mpfr_sub_one_ulp} function behaves like a conventional subtraction.
@c @end deftypefun

@deftypefun int mpfr_can_round (mpfr_t @var{b}, mp_exp_t @var{err}, mp_rnd_t @var{rnd1}, mp_rnd_t @var{rnd2}, mp_prec_t @var{prec})
Assuming @var{b} is an approximation of an unknown number
@var{x} in the direction @var{rnd1} with error at most two to the power
E(b)-@var{err} where E(b) is the exponent of @var{b}, return a non-zero
value if one is able to round correctly @var{x} to precision
@var{prec} with the direction @var{rnd2},
and 0 otherwise (including for NaN and Inf).
This function @strong{does not modify} its arguments.
@end deftypefun

@deftypefun double mpfr_get_d1 (mpfr_t @var{op})
Convert @var{op} to a @code{double}, using the default MPFR rounding mode
(see function @code{mpfr_set_default_rounding_mode}). This function is
obsolete.
@end deftypefun


@c @deftypefun void mpfr_set_str_binary (mpfr_t @var{x}, const char *@var{s})
@c Set @var{x} to the value of the binary number in string @var{s}, which has to
@c be of the
@c form +/-xxxx.xxxxxxEyy. The exponent is read in decimal, but is interpreted
@c as the power of two to be multiplied by the mantissa.
@c The mantissa length of @var{s} has to be less or equal to the precision of
@c @var{x}, otherwise an error occurs.
@c If @var{s} starts with @code{N}, it is interpreted as NaN (Not-a-Number);
@c if it starts with @code{I} after the sign, it is interpreted as infinity,
@c with the corresponding sign.
@c @end deftypefun

@c @deftypefun void mpfr_print_binary (mpfr_t @var{float})
@c Output @var{float} on stdout
@c in raw binary format (the exponent is written in decimal, yet).
@c @end deftypefun

@node Contributors, References, MPFR Interface, Top
@comment  node-name,  next,  previous,  up
@unnumbered Contributors

The main developers consist of Guillaume Hanrot, Vincent Lef@`evre,
Patrick P@'elissier and Paul Zimmermann.

We would like to thank Jean-Michel Muller and Joris van der Hoeven for very
fruitful discussions at the beginning of that project, Torbj@"orn Granlund
and Kevin Ryde for their help about design issues,
and Nathalie Revol for her careful reading of a previous version of
this documentation.
Kevin Ryde did a tremendous job for the portability of MPFR,
and integrating it into GMP 4.x;
alas the GMP developers decided in January 2004 not to include MPFR any more.

Sylvie Boldo from ENS-Lyon, France,
contributed the functions @code{mpfr_agm} and @code{mpfr_log}.
Emmanuel Jeandel, from ENS-Lyon too,
contributed the generic hypergeometric code in
@code{generic.c}, as well as the @code{mpfr_exp3},
a first implementation of the sine and cosine,
and improved versions of
@code{mpfr_const_log2} and @code{mpfr_const_pi}.
Mathieu Dutour contributed the functions @code{mpfr_atan} and @code{mpfr_asin},
and a previous version of @code{mpfr_gamma};
David Daney contributed the hyperbolic and inverse hyperbolic functions,
the base-2 exponential, and the factorial function. Fabrice Rouillier
contributed the original version of @file{mul_ui.c}, the @file{gmp_op.c}
file, and helped to the Windows porting.
Jean-Luc R@'emy contributed the @code{mpfr_zeta} code.
Ludovic Meunier helped in the design of the @code{mpfr_erf} code.
Damien Stehl@'e contributed the @code{mpfr_get_ld_2exp} function.

The development of the MPFR library would not have been possible without the
continuous support of INRIA, and of the LORIA and LIP laboratories.
In particular the main authors were or are members of the
PolKA, Spaces, Cacao project-teams at LORIA (Nancy, France)
and of the Arenaire project-team at LIP (Lyon, France).
The development of MPFR was also supported by a grant
(202F0659 00 MPN 121) from the Conseil R@'egional de Lorraine in 2002.

@node References, GNU Free Documentation License, Contributors, Top
@comment  node-name,  next,  previous,  up
@unnumbered References

@itemize @bullet

@item
Torbj@"orn Granlund, "GNU MP: The GNU Multiple Precision Arithmetic Library",
  version 4.1.2, 2002.

@item
IEEE standard for binary floating-point arithmetic, Technical Report
ANSI-IEEE Standard 754-1985, New York, 1985.
Approved March 21, 1985: IEEE Standards Board; approved July 26,
  1985: American National Standards Institute, 18 pages.

@item
Donald E. Knuth, "The Art of Computer Programming", vol 2,
"Seminumerical Algorithms", 2nd edition, Addison-Wesley, 1981.

@item
Jean-Michel Muller, "Elementary Functions, Algorithms and Implementation",
Birkhauser, Boston, 1997.

@end itemize


@node GNU Free Documentation License, Concept Index, References, Top
@appendix GNU Free Documentation License
@cindex GNU Free Documentation License
@include fdl.texi


@node Concept Index, Function Index, GNU Free Documentation License, Top
@comment  node-name,  next,  previous,  up
@unnumbered Concept Index
@printindex cp

@node Function Index,  , Concept Index, Top
@comment  node-name,  next,  previous,  up
@unnumbered Function and Type Index
@printindex fn

@bye

@c Local variables:
@c fill-column: 78
@c End: