File: singular.py

package info (click to toggle)
sagemath 7.4-9
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 108,312 kB
  • ctags: 72,147
  • sloc: python: 800,328; sh: 10,775; cpp: 7,154; ansic: 2,301; objc: 1,372; makefile: 889; lisp: 1
file content (2714 lines) | stat: -rw-r--r-- 89,224 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
2161
2162
2163
2164
2165
2166
2167
2168
2169
2170
2171
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
2192
2193
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206
2207
2208
2209
2210
2211
2212
2213
2214
2215
2216
2217
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
2231
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
2283
2284
2285
2286
2287
2288
2289
2290
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325
2326
2327
2328
2329
2330
2331
2332
2333
2334
2335
2336
2337
2338
2339
2340
2341
2342
2343
2344
2345
2346
2347
2348
2349
2350
2351
2352
2353
2354
2355
2356
2357
2358
2359
2360
2361
2362
2363
2364
2365
2366
2367
2368
2369
2370
2371
2372
2373
2374
2375
2376
2377
2378
2379
2380
2381
2382
2383
2384
2385
2386
2387
2388
2389
2390
2391
2392
2393
2394
2395
2396
2397
2398
2399
2400
2401
2402
2403
2404
2405
2406
2407
2408
2409
2410
2411
2412
2413
2414
2415
2416
2417
2418
2419
2420
2421
2422
2423
2424
2425
2426
2427
2428
2429
2430
2431
2432
2433
2434
2435
2436
2437
2438
2439
2440
2441
2442
2443
2444
2445
2446
2447
2448
2449
2450
2451
2452
2453
2454
2455
2456
2457
2458
2459
2460
2461
2462
2463
2464
2465
2466
2467
2468
2469
2470
2471
2472
2473
2474
2475
2476
2477
2478
2479
2480
2481
2482
2483
2484
2485
2486
2487
2488
2489
2490
2491
2492
2493
2494
2495
2496
2497
2498
2499
2500
2501
2502
2503
2504
2505
2506
2507
2508
2509
2510
2511
2512
2513
2514
2515
2516
2517
2518
2519
2520
2521
2522
2523
2524
2525
2526
2527
2528
2529
2530
2531
2532
2533
2534
2535
2536
2537
2538
2539
2540
2541
2542
2543
2544
2545
2546
2547
2548
2549
2550
2551
2552
2553
2554
2555
2556
2557
2558
2559
2560
2561
2562
2563
2564
2565
2566
2567
2568
2569
2570
2571
2572
2573
2574
2575
2576
2577
2578
2579
2580
2581
2582
2583
2584
2585
2586
2587
2588
2589
2590
2591
2592
2593
2594
2595
2596
2597
2598
2599
2600
2601
2602
2603
2604
2605
2606
2607
2608
2609
2610
2611
2612
2613
2614
2615
2616
2617
2618
2619
2620
2621
2622
2623
2624
2625
2626
2627
2628
2629
2630
2631
2632
2633
2634
2635
2636
2637
2638
2639
2640
2641
2642
2643
2644
2645
2646
2647
2648
2649
2650
2651
2652
2653
2654
2655
2656
2657
2658
2659
2660
2661
2662
2663
2664
2665
2666
2667
2668
2669
2670
2671
2672
2673
2674
2675
2676
2677
2678
2679
2680
2681
2682
2683
2684
2685
2686
2687
2688
2689
2690
2691
2692
2693
2694
2695
2696
2697
2698
2699
2700
2701
2702
2703
2704
2705
2706
2707
2708
2709
2710
2711
2712
2713
2714
r"""
Interface to Singular

AUTHORS:

- David Joyner and William Stein (2005): first version

- Martin Albrecht (2006-03-05): code so singular.[tab] and x =
  singular(...), x.[tab] includes all singular commands.

- Martin Albrecht (2006-03-06): This patch adds the equality symbol to
  singular. Also fix a problem in which " " as prompt means comparison
  will break all further communication with Singular.

- Martin Albrecht (2006-03-13): added current_ring() and
  current_ring_name()

- William Stein (2006-04-10): Fixed problems with ideal constructor

- Martin Albrecht (2006-05-18): added sage_poly.

- Simon King (2010-11-23): Reduce the overhead caused by waiting for
  the Singular prompt by doing garbage collection differently.

- Simon King (2011-06-06): Make conversion from Singular to Sage more flexible.

- Simon King (2015): Extend pickling capabilities.

Introduction
------------

This interface is extremely flexible, since it's exactly like
typing into the Singular interpreter, and anything that works there
should work here.

The Singular interface will only work if Singular is installed on
your computer; this should be the case, since Singular is included
with Sage. The interface offers three pieces of functionality:


#. ``singular_console()`` - A function that dumps you
   into an interactive command-line Singular session.

#. ``singular(expr, type='def')`` - Creation of a
   Singular object. This provides a Pythonic interface to Singular.
   For example, if ``f=singular(10)``, then
   ``f.factorize()`` returns the factorization of
   `10` computed using Singular.

#. ``singular.eval(expr)`` - Evaluation of arbitrary
   Singular expressions, with the result returned as a string.

Of course, there are polynomial rings and ideals in Sage as well
(often based on a C-library interface to Singular). One can convert
an object in the Singular interpreter interface to Sage by the
method ``sage()``.


Tutorial
--------

EXAMPLES: First we illustrate multivariate polynomial
factorization::

    sage: R1 = singular.ring(0, '(x,y)', 'dp')
    sage: R1
    polynomial ring, over a field, global ordering
    //   characteristic : 0
    //   number of vars : 2
    //        block   1 : ordering dp
    //                  : names    x y
    //        block   2 : ordering C
    sage: f = singular('9x16 - 18x13y2 - 9x12y3 + 9x10y4 - 18x11y2 + 36x8y4 + 18x7y5 - 18x5y6 + 9x6y4 - 18x3y6 - 9x2y7 + 9y8')
    sage: f
    9*x^16-18*x^13*y^2-9*x^12*y^3+9*x^10*y^4-18*x^11*y^2+36*x^8*y^4+18*x^7*y^5-18*x^5*y^6+9*x^6*y^4-18*x^3*y^6-9*x^2*y^7+9*y^8
    sage: f.parent()
    Singular

::

    sage: F = f.factorize(); F
    [1]:
       _[1]=9
       _[2]=x^6-2*x^3*y^2-x^2*y^3+y^4
       _[3]=-x^5+y^2
    [2]:
       1,1,2

::

    sage: F[1]
    9,
    x^6-2*x^3*y^2-x^2*y^3+y^4,
    -x^5+y^2
    sage: F[1][2]
    x^6-2*x^3*y^2-x^2*y^3+y^4

We can convert `f` and each exponent back to Sage objects
as well.

::

    sage: g = f.sage(); g
    9*x^16 - 18*x^13*y^2 - 9*x^12*y^3 + 9*x^10*y^4 - 18*x^11*y^2 + 36*x^8*y^4 + 18*x^7*y^5 - 18*x^5*y^6 + 9*x^6*y^4 - 18*x^3*y^6 - 9*x^2*y^7 + 9*y^8
    sage: F[1][2].sage()
    x^6 - 2*x^3*y^2 - x^2*y^3 + y^4
    sage: g.parent()
    Multivariate Polynomial Ring in x, y over Rational Field

This example illustrates polynomial GCD's::

    sage: R2 = singular.ring(0, '(x,y,z)', 'lp')
    sage: a = singular.new('3x2*(x+y)')
    sage: b = singular.new('9x*(y2-x2)')
    sage: g = a.gcd(b)
    sage: g
    x^2+x*y

This example illustrates computation of a Groebner basis::

    sage: R3 = singular.ring(0, '(a,b,c,d)', 'lp')
    sage: I = singular.ideal(['a + b + c + d', 'a*b + a*d + b*c + c*d', 'a*b*c + a*b*d + a*c*d + b*c*d', 'a*b*c*d - 1'])
    sage: I2 = I.groebner()
    sage: I2
    c^2*d^6-c^2*d^2-d^4+1,
    c^3*d^2+c^2*d^3-c-d,
    b*d^4-b+d^5-d,
    b*c-b*d^5+c^2*d^4+c*d-d^6-d^2,
    b^2+2*b*d+d^2,
    a+b+c+d

The following example is the same as the one in the Singular - Gap
interface documentation::

    sage: R  = singular.ring(0, '(x0,x1,x2)', 'lp')
    sage: I1 = singular.ideal(['x0*x1*x2 -x0^2*x2', 'x0^2*x1*x2-x0*x1^2*x2-x0*x1*x2^2', 'x0*x1-x0*x2-x1*x2'])
    sage: I2 = I1.groebner()
    sage: I2
    x1^2*x2^2,
    x0*x2^3-x1^2*x2^2+x1*x2^3,
    x0*x1-x0*x2-x1*x2,
    x0^2*x2-x0*x2^2-x1*x2^2
    sage: I2.sage()
    Ideal (x1^2*x2^2, x0*x2^3 - x1^2*x2^2 + x1*x2^3, x0*x1 - x0*x2 - x1*x2, x0^2*x2 - x0*x2^2 - x1*x2^2) of Multivariate Polynomial Ring in x0, x1, x2 over Rational Field


This example illustrates moving a polynomial from one ring to
another. It also illustrates calling a method of an object with an
argument.

::

    sage: R = singular.ring(0, '(x,y,z)', 'dp')
    sage: f = singular('x3+y3+(x-y)*x2y2+z2')
    sage: f
    x^3*y^2-x^2*y^3+x^3+y^3+z^2
    sage: R1 = singular.ring(0, '(x,y,z)', 'ds')
    sage: f = R.fetch(f)
    sage: f
    z^2+x^3+y^3+x^3*y^2-x^2*y^3

We can calculate the Milnor number of `f`::

    sage: _=singular.LIB('sing.lib')     # assign to _ to suppress printing
    sage: f.milnor()
    4

The Jacobian applied twice yields the Hessian matrix of
`f`, with which we can compute.

::

    sage: H = f.jacob().jacob()
    sage: H
    6*x+6*x*y^2-2*y^3,6*x^2*y-6*x*y^2,  0,
    6*x^2*y-6*x*y^2,  6*y+2*x^3-6*x^2*y,0,
    0,                0,                2
    sage: H.sage()
    [6*x + 6*x*y^2 - 2*y^3     6*x^2*y - 6*x*y^2                     0]
    [    6*x^2*y - 6*x*y^2 6*y + 2*x^3 - 6*x^2*y                     0]
    [                    0                     0                     2]
    sage: H.det()   # This is a polynomial in Singular
    72*x*y+24*x^4-72*x^3*y+72*x*y^3-24*y^4-48*x^4*y^2+64*x^3*y^3-48*x^2*y^4
    sage: H.det().sage()   # This is the corresponding polynomial in Sage
    72*x*y + 24*x^4 - 72*x^3*y + 72*x*y^3 - 24*y^4 - 48*x^4*y^2 + 64*x^3*y^3 - 48*x^2*y^4

The 1x1 and 2x2 minors::

    sage: H.minor(1)
    2,
    6*y+2*x^3-6*x^2*y,
    6*x^2*y-6*x*y^2,
    6*x^2*y-6*x*y^2,
    6*x+6*x*y^2-2*y^3
    sage: H.minor(2)
    12*y+4*x^3-12*x^2*y,
    12*x^2*y-12*x*y^2,
    12*x^2*y-12*x*y^2,
    12*x+12*x*y^2-4*y^3,
    -36*x*y-12*x^4+36*x^3*y-36*x*y^3+12*y^4+24*x^4*y^2-32*x^3*y^3+24*x^2*y^4

::

    sage: _=singular.eval('option(redSB)')
    sage: H.minor(1).groebner()
    1

Computing the Genus
-------------------

We compute the projective genus of ideals that define curves over
`\QQ`. It is *very important* to load the
``normal.lib`` library before calling the
``genus`` command, or you'll get an error message.

EXAMPLE::

    sage: singular.lib('normal.lib')
    sage: R = singular.ring(0,'(x,y)','dp')
    sage: i2 = singular.ideal('y9 - x2*(x-1)^9 + x')
    sage: i2.genus()
    40

Note that the genus can be much smaller than the degree::

    sage: i = singular.ideal('y9 - x2*(x-1)^9')
    sage: i.genus()
    0

An Important Concept
--------------------

AUTHORS:

- Neal Harris

The following illustrates an important concept: how Sage interacts
with the data being used and returned by Singular. Let's compute a
Groebner basis for some ideal, using Singular through Sage.

::

    sage: singular.lib('poly.lib')
    sage: singular.ring(32003, '(a,b,c,d,e,f)', 'lp')
            polynomial ring, over a field, global ordering
            //   characteristic : 32003
            //   number of vars : 6
            //        block   1 : ordering lp
            //                        : names    a b c d e f
            //        block   2 : ordering C
    sage: I = singular.ideal('cyclic(6)')
    sage: g = singular('groebner(I)')
    Traceback (most recent call last):
    ...
    TypeError: Singular error:
    ...

We restart everything and try again, but correctly.

::

    sage: singular.quit()
    sage: singular.lib('poly.lib'); R = singular.ring(32003, '(a,b,c,d,e,f)', 'lp')
    sage: I = singular.ideal('cyclic(6)')
    sage: I.groebner()
    f^48-2554*f^42-15674*f^36+12326*f^30-12326*f^18+15674*f^12+2554*f^6-1,
    ...

It's important to understand why the first attempt at computing a
basis failed. The line where we gave singular the input
'groebner(I)' was useless because Singular has no idea what 'I' is!
Although 'I' is an object that we computed with calls to Singular
functions, it actually lives in Sage. As a consequence, the name
'I' means nothing to Singular. When we called
``I.groebner()``, Sage was able to call the groebner
function on'I' in Singular, since 'I' actually means something to
Sage.

Long Input
----------

The Singular interface reads in even very long input (using files)
in a robust manner, as long as you are creating a new object.

::

    sage: t = '"%s"'%10^15000   # 15 thousand character string (note that normal Singular input must be at most 10000)
    sage: a = singular.eval(t)
    sage: a = singular(t)

TESTS:

We test an automatic coercion::

    sage: a = 3*singular('2'); a
    6
    sage: type(a)
    <class 'sage.interfaces.singular.SingularElement'>
    sage: a = singular('2')*3; a
    6
    sage: type(a)
    <class 'sage.interfaces.singular.SingularElement'>

Create a ring over GF(9) to check that ``gftables`` has been installed,
see :trac:`11645`::

    sage: singular.eval("ring testgf9 = (9,x),(a,b,c,d,e,f),(M((1,2,3,0)),wp(2,3),lp);")
    ''
"""

#*****************************************************************************
#       Copyright (C) 2005 David Joyner and William Stein
#
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 2 of the License, or
# (at your option) any later version.
#                  http://www.gnu.org/licenses/
#*****************************************************************************
from __future__ import print_function
from __future__ import absolute_import

import os
import re
import sys
import pexpect
from time import sleep

from .expect import Expect, ExpectElement, FunctionElement, ExpectFunction

from sage.interfaces.tab_completion import ExtraTabCompletion
from sage.structure.sequence import Sequence
from sage.structure.element import RingElement

import sage.rings.integer

from sage.misc.misc import get_verbose
from sage.misc.superseded import deprecation

from six import reraise as raise_

class SingularError(RuntimeError):
    """
    Raised if Singular printed an error message
    """
    pass


class Singular(ExtraTabCompletion, Expect):
    r"""
    Interface to the Singular interpreter.

    EXAMPLES: A Groebner basis example.

    ::

        sage: R = singular.ring(0, '(x0,x1,x2)', 'lp')
        sage: I = singular.ideal([ 'x0*x1*x2 -x0^2*x2', 'x0^2*x1*x2-x0*x1^2*x2-x0*x1*x2^2', 'x0*x1-x0*x2-x1*x2'])
        sage: I.groebner()
        x1^2*x2^2,
        x0*x2^3-x1^2*x2^2+x1*x2^3,
        x0*x1-x0*x2-x1*x2,
        x0^2*x2-x0*x2^2-x1*x2^2

    AUTHORS:

    - David Joyner and William Stein
    """
    def __init__(self, maxread=None, script_subdirectory=None,
                 logfile=None, server=None,server_tmpdir=None,
                 seed=None):
        """
        EXAMPLES::

            sage: singular == loads(dumps(singular))
            True
        """
        prompt = '> '
        Expect.__init__(self,
                        terminal_echo=False,
                        name = 'singular',
                        prompt = prompt,
                        # no tty, fine grained cputime()
                        # and do not display CTRL-C prompt
                        command = "Singular -t --ticks-per-sec 1000 --cntrlc=a",
                        server = server,
                        server_tmpdir = server_tmpdir,
                        script_subdirectory = script_subdirectory,
                        restart_on_ctrlc = True,
                        verbose_start = False,
                        logfile = logfile,
                        eval_using_file_cutoff=100 if os.uname()[0]=="SunOS" else 1000)
        self.__libs  = []
        self._prompt_wait = prompt
        self.__to_clear = []   # list of variable names that need to be cleared.
        self._seed = seed

    def set_seed(self,seed=None):
        """
        Sets the seed for singular interpeter.
        The seed should be an integer at least 1
        and not more than 30 bits.
        See
        http://www.singular.uni-kl.de/Manual/html/sing_19.htm#SEC26
        and
        http://www.singular.uni-kl.de/Manual/html/sing_283.htm#SEC323

        EXAMPLES::

            sage: s = Singular()
            sage: s.set_seed(1)
            1
            sage: [s.random(1,10) for i in range(5)]
            [8, 10, 4, 9, 1]
        """
        if seed is None:
            seed = self.rand_seed()
        self.eval('system("--random",%d)' % seed)
        self._seed = seed
        return seed

    def _start(self, alt_message=None):
        """
        EXAMPLES::

            sage: s = Singular()
            sage: s.is_running()
            False
            sage: s._start()
            sage: s.is_running()
            True
            sage: s.quit()
        """
        self.__libs = []
        Expect._start(self, alt_message)
        # Load some standard libraries.
        self.lib('general')   # assumed loaded by misc/constants.py

        # these options are required by the new coefficient rings
        # supported by Singular 3-1-0.
        self.option("redTail")
        self.option("redThrough")
        self.option("intStrategy")
        self._saved_options = self.option('get')
        # set random seed
        self.set_seed(self._seed)

    def __reduce__(self):
        """
        EXAMPLES::

            sage: singular.__reduce__()
            (<function reduce_load_Singular at 0x...>, ())
        """
        return reduce_load_Singular, ()

    def _equality_symbol(self):
        """
        EXAMPLES::

            sage: singular._equality_symbol()
            '=='
        """
        return '=='

    def _true_symbol(self):
        """
        EXAMPLES::

            sage: singular._true_symbol()
            '1'
        """
        return '1'

    def _false_symbol(self):
        """
        EXAMPLES::

            sage: singular._false_symbol()
            '0'
        """
        return '0'

    def _quit_string(self):
        """
        EXAMPLES::

            sage: singular._quit_string()
            'quit'
        """
        return 'quit'

    def _send_interrupt(self):
        """
        Send an interrupt to Singular. If needed, additional
        semi-colons are sent until we get back at the prompt.

        TESTS:

        The following works without restarting Singular::

            sage: a = singular(1)
            sage: _ = singular._expect.sendline('1+')  # unfinished input
            sage: try:
            ....:     alarm(0.5)
            ....:     singular._expect_expr('>')  # interrupt this
            ....: except KeyboardInterrupt:
            ....:     pass
            Control-C pressed.  Interrupting Singular. Please wait a few seconds...

        We can still access a::

            sage: 2*a
            2
        """
        # Work around for Singular bug
        # http://www.singular.uni-kl.de:8002/trac/ticket/727
        sleep(0.1)

        E = self._expect
        E.sendline(chr(3))
        for i in range(5):
            try:
                E.expect_upto(self._prompt, timeout=1.0)
                return
            except Exception:
                pass
            E.sendline(";")

    def _read_in_file_command(self, filename):
        r"""
        EXAMPLES::

            sage: singular._read_in_file_command('test')
            '< "...";'

            sage: filename = tmp_filename()
            sage: f = open(filename, 'w')
            sage: f.write('int x = 2;\n')
            sage: f.close()
            sage: singular.read(filename)
            sage: singular.get('x')
            '2'
        """
        return '< "%s";'%filename

    def eval(self, x, allow_semicolon=True, strip=True, **kwds):
        r"""
        Send the code x to the Singular interpreter and return the output
        as a string.

        INPUT:


        -  ``x`` - string (of code)

        -  ``allow_semicolon`` - default: False; if False then
           raise a TypeError if the input line contains a semicolon.

        -  ``strip`` - ignored


        EXAMPLES::

            sage: singular.eval('2 > 1')
            '1'
            sage: singular.eval('2 + 2')
            '4'

        if the verbosity level is `> 1` comments are also printed
        and not only returned.

        ::

            sage: r = singular.ring(0,'(x,y,z)','dp')
            sage: i = singular.ideal(['x^2','y^2','z^2'])
            sage: s = i.std()
            sage: singular.eval('hilb(%s)'%(s.name()))
            '// 1 t^0\n// -3 t^2\n// 3 t^4\n// -1 t^6\n\n// 1 t^0\n//
            3 t^1\n// 3 t^2\n// 1 t^3\n// dimension (affine) = 0\n//
            degree (affine) = 8'

        ::

            sage: set_verbose(1)
            sage: o = singular.eval('hilb(%s)'%(s.name()))
            //         1 t^0
            //        -3 t^2
            //         3 t^4
            //        -1 t^6
            //         1 t^0
            //         3 t^1
            //         3 t^2
            //         1 t^3
            // dimension (affine) = 0
            // degree (affine)  = 8

        This is mainly useful if this method is called implicitly. Because
        then intermediate results, debugging outputs and printed statements
        are printed

        ::

            sage: o = s.hilb()
            //         1 t^0
            //        -3 t^2
            //         3 t^4
            //        -1 t^6
            //         1 t^0
            //         3 t^1
            //         3 t^2
            //         1 t^3
            // dimension (affine) = 0
            // degree (affine)  = 8
            // ** right side is not a datum, assignment ignored
            ...

        rather than ignored

        ::

            sage: set_verbose(0)
            sage: o = s.hilb()
        """
        # Simon King:
        # In previous versions, the interface was first synchronised and then
        # unused variables were killed. This created a considerable overhead.
        # By trac ticket #10296, killing unused variables is now done inside
        # singular.set(). Moreover, it is not done by calling a separate _eval_line.
        # In that way, the time spent by waiting for the singular prompt is reduced.

        # Before #10296, it was possible that garbage collection occured inside
        # of _eval_line. But collection of the garbage would launch another call
        # to _eval_line. The result would have been a dead lock, that could only
        # be avoided by synchronisation. Since garbage collection is now done
        # without an additional call to _eval_line, synchronisation is not
        # needed anymore, saving even more waiting time for the prompt.

        # Uncomment the print statements below for low-level debugging of
        # code that involves the singular interfaces.  Everything goes
        # through here.

        x = str(x).rstrip().rstrip(';')
        x = x.replace("> ",">\t") #don't send a prompt  (added by Martin Albrecht)
        if not allow_semicolon and x.find(";") != -1:
            raise TypeError("singular input must not contain any semicolons:\n%s"%x)
        if len(x) == 0 or x[len(x) - 1] != ';':
            x += ';'

        s = Expect.eval(self, x, **kwds)

        if s.find("error") != -1 or s.find("Segment fault") != -1:
            raise SingularError('Singular error:\n%s'%s)

        if get_verbose() > 0:
            for line in s.splitlines():
                if line.startswith("//"):
                    print(line)
            return s
        else:
            return s

    def set(self, type, name, value):
        """
        Set the variable with given name to the given value.

        REMARK:

        If a variable in the Singular interface was previously marked for
        deletion, the actual deletion is done here, before the new variable
        is created in Singular.

        EXAMPLES::

            sage: singular.set('int', 'x', '2')
            sage: singular.get('x')
            '2'

        We test that an unused variable is only actually deleted if this method
        is called::

            sage: a = singular(3)
            sage: n = a.name()
            sage: del a
            sage: singular.eval(n)
            '3'
            sage: singular.set('int', 'y', '5')
            sage: singular.eval('defined(%s)'%n)
            '0'

        """
        cmd = ''.join('if(defined(%s)){kill %s;};'%(v,v) for v in self.__to_clear)
        cmd += '%s %s=%s;'%(type, name, value)
        self.__to_clear = []
        self.eval(cmd)

    def get(self, var):
        """
        Get string representation of variable named var.

        EXAMPLES::

            sage: singular.set('int', 'x', '2')
            sage: singular.get('x')
            '2'
        """
        return self.eval('print(%s);'%var)

    def clear(self, var):
        """
        Clear the variable named ``var``.

        EXAMPLES::

            sage: singular.set('int', 'x', '2')
            sage: singular.get('x')
            '2'
            sage: singular.clear('x')

        "Clearing the variable" means to allow to free the memory
        that it uses in the Singular sub-process. However, the
        actual deletion of the variable is only committed when
        the next element in the Singular interface is created::

            sage: singular.get('x')
            '2'
            sage: a = singular(3)
            sage: singular.get('x')
            '`x`'

        """
        # We add the variable to the list of vars to clear when we do an eval.
        # We queue up all the clears and do them at once to avoid synchronizing
        # the interface at the same time we do garbage collection, which can
        # lead to subtle problems.    This was Willem Jan's ideas, implemented
        # by William Stein.
        self.__to_clear.append(var)

    def _create(self, value, type='def'):
        """
        Creates a new variable in the Singular session and returns the name
        of that variable.

        EXAMPLES::

            sage: singular._create('2', type='int')
            'sage...'
            sage: singular.get(_)
            '2'
        """
        name = self._next_var_name()
        self.set(type, name, value)
        return name

    def __call__(self, x, type='def'):
        """
        Create a singular object X with given type determined by the string
        x. This returns var, where var is built using the Singular
        statement type var = ... x ... Note that the actual name of var
        could be anything, and can be recovered using X.name().

        The object X returned can be used like any Sage object, and wraps
        an object in self. The standard arithmetic operators work. Moreover
        if foo is a function then X.foo(y,z,...) calls foo(X, y, z, ...)
        and returns the corresponding object.

        EXAMPLES::

            sage: R = singular.ring(0, '(x0,x1,x2)', 'lp')
            sage: I = singular.ideal([ 'x0*x1*x2 -x0^2*x2', 'x0^2*x1*x2-x0*x1^2*x2-x0*x1*x2^2', 'x0*x1-x0*x2-x1*x2'])
            sage: I
             -x0^2*x2+x0*x1*x2,
            x0^2*x1*x2-x0*x1^2*x2-x0*x1*x2^2,
            x0*x1-x0*x2-x1*x2
            sage: type(I)
            <class 'sage.interfaces.singular.SingularElement'>
            sage: I.parent()
            Singular
        """
        if isinstance(x, SingularElement) and x.parent() is self:
            return x
        elif isinstance(x, ExpectElement):
            return self(x.sage())
        elif not isinstance(x, ExpectElement) and hasattr(x, '_singular_'):
            return x._singular_(self)

        # some convenient conversions
        if type in ("module","list") and isinstance(x,(list,tuple,Sequence)):
            x = str(x)[1:-1]

        return SingularElement(self, type, x, False)

    def _coerce_map_from_(self, S):
        """
        Return ``True`` if ``S`` admits a coercion map into the
        Singular interface.

        EXAMPLES::

            sage: singular._coerce_map_from_(ZZ)
            True
            sage: singular.coerce_map_from(ZZ)
            Call morphism:
              From: Integer Ring
              To:   Singular
            sage: singular.coerce_map_from(float)
        """
        # we want to implement this without coercing, since singular has state.
        if hasattr(S, 'an_element'):
            if hasattr(S.an_element(), '_singular_'):
                return True
            try:
                self._coerce_(S.an_element())
                return True
            except TypeError:
                pass
        elif S is int or S is long:
            return True
        return None

    def cputime(self, t=None):
        r"""
        Returns the amount of CPU time that the Singular session has used.
        If ``t`` is not None, then it returns the difference
        between the current CPU time and ``t``.

        EXAMPLES::

            sage: t = singular.cputime()
            sage: R = singular.ring(0, '(x0,x1,x2)', 'lp')
            sage: I = singular.ideal([ 'x0*x1*x2 -x0^2*x2', 'x0^2*x1*x2-x0*x1^2*x2-x0*x1*x2^2', 'x0*x1-x0*x2-x1*x2'])
            sage: gb = I.groebner()
            sage: singular.cputime(t) #random
            0.02
        """
        if t:
            return float(self.eval('timer-(%d)'%(int(1000*t))))/1000.0
        else:
            return float(self.eval('timer'))/1000.0

    ###################################################################
    # Singular libraries
    ###################################################################
    def lib(self, lib, reload=False):
        """
        Load the Singular library named lib.

        Note that if the library was already loaded during this session it
        is not reloaded unless the optional reload argument is True (the
        default is False).

        EXAMPLES::

            sage: singular.lib('sing.lib')
            sage: singular.lib('sing.lib', reload=True)
        """
        if lib[-4:] != ".lib":
            lib += ".lib"
        if not reload and lib in self.__libs:
            return
        self.eval('LIB "%s"'%lib)
        self.__libs.append(lib)

    LIB = lib
    load = lib

    ###################################################################
    # constructors
    ###################################################################
    def ideal(self, *gens):
        """
        Return the ideal generated by gens.

        INPUT:


        -  ``gens`` - list or tuple of Singular objects (or
           objects that can be made into Singular objects via evaluation)


        OUTPUT: the Singular ideal generated by the given list of gens

        EXAMPLES: A Groebner basis example done in a different way.

        ::

            sage: _ = singular.eval("ring R=0,(x0,x1,x2),lp")
            sage: i1 = singular.ideal([ 'x0*x1*x2 -x0^2*x2', 'x0^2*x1*x2-x0*x1^2*x2-x0*x1*x2^2', 'x0*x1-x0*x2-x1*x2'])
            sage: i1
            -x0^2*x2+x0*x1*x2,
            x0^2*x1*x2-x0*x1^2*x2-x0*x1*x2^2,
            x0*x1-x0*x2-x1*x2

        ::

            sage: i2 = singular.ideal('groebner(%s);'%i1.name())
            sage: i2
            x1^2*x2^2,
            x0*x2^3-x1^2*x2^2+x1*x2^3,
            x0*x1-x0*x2-x1*x2,
            x0^2*x2-x0*x2^2-x1*x2^2
        """
        if isinstance(gens, str):
            gens = self(gens)

        if isinstance(gens, SingularElement):
            return self(gens.name(), 'ideal')

        if not isinstance(gens, (list, tuple)):
            raise TypeError("gens (=%s) must be a list, tuple, string, or Singular element"%gens)

        if len(gens) == 1 and isinstance(gens[0], (list, tuple)):
            gens = gens[0]
        gens2 = []
        for g in gens:
            if not isinstance(g, SingularElement):
                gens2.append(self.new(g))
            else:
                gens2.append(g)
        return self(",".join([g.name() for g in gens2]), 'ideal')

    def list(self, x):
        r"""
        Creates a list in Singular from a Sage list ``x``.

        EXAMPLES::

            sage: singular.list([1,2])
            [1]:
               1
            [2]:
               2
        """
        return self(x, 'list')

    def matrix(self, nrows, ncols, entries=None):
        """
        EXAMPLES::

            sage: singular.lib("matrix")
            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(3,2,'1,2,3,4,5,6')
            sage: A
            1,2,
            3,4,
            5,6
            sage: A.gauss_col()
            2,-1,
            1,0,
            0,1

        AUTHORS:

        - Martin Albrecht (2006-01-14)
        """
        name = self._next_var_name()
        if entries is None:
            self.eval('matrix %s[%s][%s]'%(name, nrows, ncols))
        else:
            self.eval('matrix %s[%s][%s] = %s'%(name, nrows, ncols, entries))
        return SingularElement(self, None, name, True)

    def ring(self, char=0, vars='(x)', order='lp', check=True):
        r"""
        Create a Singular ring and makes it the current ring.

        INPUT:


        -  ``char`` - characteristic of the base ring (see
           examples below), which must be either 0, prime (!), or one of
           several special codes (see examples below).

        -  ``vars`` - a tuple or string that defines the
           variable names

        -  ``order`` - string - the monomial order (default:
           'lp')

        -  ``check`` - if True, check primality of the
           characteristic if it is an integer.


        OUTPUT: a Singular ring

        .. note::

           This function is *not* identical to calling the Singular
           ``ring`` function. In particular, it also attempts to
           "kill" the variable names, so they can actually be used
           without getting errors, and it sets printing of elements
           for this range to short (i.e., with \*'s and carets).

        EXAMPLES: We first declare `\QQ[x,y,z]` with degree reverse
        lexicographic ordering.

        ::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: R
            polynomial ring, over a field, global ordering
            //   characteristic : 0
            //   number of vars : 3
            //        block   1 : ordering dp
            //                  : names    x y z
            //        block   2 : ordering C

        ::

            sage: R1 = singular.ring(32003, '(x,y,z)', 'dp')
            sage: R2 = singular.ring(32003, '(a,b,c,d)', 'lp')

        This is a ring in variables named x(1) through x(10) over the
        finite field of order `7`::

            sage: R3 = singular.ring(7, '(x(1..10))', 'ds')

        This is a polynomial ring over the transcendental extension
        `\QQ(a)` of `\QQ`::

            sage: R4 = singular.ring('(0,a)', '(mu,nu)', 'lp')

        This is a ring over the field of single-precision floats::

            sage: R5 = singular.ring('real', '(a,b)', 'lp')

        This is over 50-digit floats::

            sage: R6 = singular.ring('(real,50)', '(a,b)', 'lp')
            sage: R7 = singular.ring('(complex,50,i)', '(a,b)', 'lp')

        To use a ring that you've defined, use the set_ring() method on
        the ring. This sets the ring to be the "current ring". For
        example,

        ::

            sage: R = singular.ring(7, '(a,b)', 'ds')
            sage: S = singular.ring('real', '(a,b)', 'lp')
            sage: singular.new('10*a')
            (1.000e+01)*a
            sage: R.set_ring()
            sage: singular.new('10*a')
            3*a
        """
        if len(vars) > 2:
            s = '; '.join(['if(defined(%s)>0){kill %s;};'%(x,x)
                           for x in vars[1:-1].split(',')])
            self.eval(s)

        if check and isinstance(char, (int, long, sage.rings.integer.Integer)):
            if char != 0:
                n = sage.rings.integer.Integer(char)
                if not n.is_prime():
                    raise ValueError("the characteristic must be 0 or prime")
        R = self('%s,%s,%s'%(char, vars, order), 'ring')
        self.eval('short=0')  # make output include *'s for multiplication for *THIS* ring.
        return R

    def string(self, x):
        """
        Creates a Singular string from a Sage string. Note that the Sage
        string has to be "double-quoted".

        EXAMPLES::

            sage: singular.string('"Sage"')
            Sage
        """
        return self(x, 'string')

    def set_ring(self, R):
        """
        Sets the current Singular ring to R.

        EXAMPLES::

            sage: R = singular.ring(7, '(a,b)', 'ds')
            sage: S = singular.ring('real', '(a,b)', 'lp')
            sage: singular.current_ring()
            polynomial ring, over a field, global ordering
            //   characteristic : 0 (real)
            //   number of vars : 2
            //        block   1 : ordering lp
            //                  : names    a b
            //        block   2 : ordering C
            sage: singular.set_ring(R)
            sage: singular.current_ring()
            polynomial ring, over a field, local/mixed ordering
            //   characteristic : 7
            //   number of vars : 2
            //        block   1 : ordering ds
            //                  : names    a b
            //        block   2 : ordering C
        """
        if not isinstance(R, SingularElement):
            raise TypeError("R must be a singular ring")
        self.eval("setring %s; short=0"%R.name(), allow_semicolon=True)

    setring = set_ring

    def current_ring_name(self):
        """
        Returns the Singular name of the currently active ring in
        Singular.

        OUTPUT: currently active ring's name

        EXAMPLES::

            sage: r = PolynomialRing(GF(127),3,'xyz')
            sage: r._singular_().name() == singular.current_ring_name()
            True
        """
        ringlist = self.eval("listvar(ring)").splitlines()
        p = re.compile("// ([a-zA-Z0-9_]*).*\[.*\].*\*.*") #do this in constructor?
        for line in ringlist:
            m = p.match(line)
            if m:
                return m.group(int(1))
        return None

    def current_ring(self):
        """
        Returns the current ring of the running Singular session.

        EXAMPLES::

            sage: r = PolynomialRing(GF(127),3,'xyz', order='invlex')
            sage: r._singular_()
            polynomial ring, over a field, global ordering
            //   characteristic : 127
            //   number of vars : 3
            //        block   1 : ordering rp
            //                  : names    x y z
            //        block   2 : ordering C
            sage: singular.current_ring()
            polynomial ring, over a field, global ordering
            //   characteristic : 127
            //   number of vars : 3
            //        block   1 : ordering rp
            //                  : names    x y z
            //        block   2 : ordering C
        """
        name = self.current_ring_name()
        if name:
            return self(name)
        else:
            return None

    def _tab_completion(self):
        """
         Return a list of all Singular commands.

         EXAMPLES::

             sage: singular._tab_completion()
             ['exteriorPower',
              ...
              'stdfglm']
         """
        p = re.compile("// *([a-z0-9A-Z_]*).*") #compiles regular expression
        proclist = self.eval("listvar(proc)").splitlines()
        return [p.match(line).group(int(1)) for line in proclist]

    def console(self):
        """
        EXAMPLES::

            sage: singular_console() #not tested
                                 SINGULAR                             /  Development
             A Computer Algebra System for Polynomial Computations   /   version 3-0-4
                                                                   0<
                 by: G.-M. Greuel, G. Pfister, H. Schoenemann        \   Nov 2007
            FB Mathematik der Universitaet, D-67653 Kaiserslautern    \
        """
        singular_console()

    def version(self):
        """
        EXAMPLES:
        """
        return singular_version()

    def _function_class(self):
        """
        EXAMPLES::

            sage: singular._function_class()
            <class 'sage.interfaces.singular.SingularFunction'>
        """
        return SingularFunction

    def _function_element_class(self):
        """
        EXAMPLES::

            sage: singular._function_element_class()
            <class 'sage.interfaces.singular.SingularFunctionElement'>
        """
        return SingularFunctionElement

    def option(self, cmd=None, val=None):
        """
        Access to Singular's options as follows:

        Syntax: option() Returns a string of all defined options.

        Syntax: option( 'option_name' ) Sets an option. Note to disable an
        option, use the prefix no.

        Syntax: option( 'get' ) Returns an intvec of the state of all
        options.

        Syntax: option( 'set', intvec_expression ) Restores the state of
        all options from an intvec (produced by option('get')).

        EXAMPLES::

            sage: singular.option()
            //options: redefine loadLib usage prompt
            sage: singular.option('get')
            0,
            10321
            sage: old_options = _
            sage: singular.option('noredefine')
            sage: singular.option()
            //options: loadLib usage prompt
            sage: singular.option('set', old_options)
            sage: singular.option('get')
            0,
            10321
        """
        if cmd is None:
            return SingularFunction(self,"option")()
        elif cmd == "get":
            #return SingularFunction(self,"option")("\"get\"")
            return self(self.eval("option(get)"),"intvec")
        elif cmd == "set":
            if not isinstance(val,SingularElement):
                raise TypeError("singular.option('set') needs SingularElement as second parameter")
            #SingularFunction(self,"option")("\"set\"",val)
            self.eval("option(set,%s)"%val.name())
        else:
            SingularFunction(self,"option")("\""+str(cmd)+"\"")

    def _keyboard_interrupt(self):
        print("Interrupting %s..." % self)
        try:
            self._expect.sendline(chr(4))
        except pexpect.ExceptionPexpect as msg:
            raise pexpect.ExceptionPexpect("THIS IS A BUG -- PLEASE REPORT. This should never happen.\n" + msg)
        self._start()
        raise KeyboardInterrupt("Restarting %s (WARNING: all variables defined in previous session are now invalid)" % self)

    
class SingularElement(ExtraTabCompletion, ExpectElement):
    
    def __init__(self, parent, type, value, is_name=False):
        """
        EXAMPLES::

            sage: a = singular(2)
            sage: loads(dumps(a))
            2
        """
        RingElement.__init__(self, parent)
        if parent is None: return
        if not is_name:
            try:
                self._name = parent._create( value, type)
            # Convert SingularError to TypeError for
            # coercion to work properly.
            except SingularError as x:
                self._session_number = -1
                raise_(TypeError, x, sys.exc_info()[2])
            except BaseException:
                self._session_number = -1
                raise
        else:
            self._name = value
        self._session_number = parent._session_number

    def __repr__(self):
        r"""
        Return string representation of ``self``.

        EXAMPLE::

            sage: r = singular.ring(0,'(x,y)','dp')
            sage: singular(0)
            0
            sage: singular('x') # indirect doctest
            x
            sage: singular.matrix(2,2)
            0,0,
            0,0
            sage: singular.matrix(2,2,"(25/47*x^2*y^4 + 63/127*x + 27)^3,y,0,1")
            15625/103823*x^6*y.., y,
            0,                    1

        Note that the output is truncated

        ::

            sage: M= singular.matrix(2,2,"(25/47*x^2*y^4 + 63/127*x + 27)^3,y,0,1")
            sage: M.rename('T')
            sage: M
            T[1,1],y,
            0,         1

        if ``self`` has a custom name, it is used to print the
        matrix, rather than abbreviating its contents
        """
        try:
            self._check_valid()
        except ValueError:
            return '(invalid object -- defined in terms of closed session)'
        try:
            if self._get_using_file:
                s = self.parent().get_using_file(self._name)
        except AttributeError:
            s = self.parent().get(self._name)
        if self._name in s:
            if hasattr(self, '__custom_name'):
                s =  s.replace(self._name, self.__dict__['__custom_name'])
            elif self.type() == 'matrix':
                s = self.parent().eval('pmat(%s,20)'%(self.name()))
        return s

    def __copy__(self):
        r"""
        Returns a copy of ``self``.

        EXAMPLES::

            sage: R=singular.ring(0,'(x,y)','dp')
            sage: M=singular.matrix(3,3,'0,0,-x, 0,y,0, x*y,0,0')
            sage: N=copy(M)
            sage: N[1,1]=singular('x+y')
            sage: N
            x+y,0,-x,
            0,  y,0,
            x*y,0,0
            sage: M
            0,  0,-x,
            0,  y,0,
            x*y,0,0
            sage: L=R.ringlist()
            sage: L[4]=singular.ideal('x**2-5')
            sage: Q=L.ring()
            sage: otherR=singular.ring(5,'(x)','dp')
            sage: cpQ=copy(Q)
            sage: cpQ.set_ring()
            sage: cpQ
            polynomial ring, over a field, global ordering
            //   characteristic : 0
            //   number of vars : 2
            //        block   1 : ordering dp
            //                  : names    x y
            //        block   2 : ordering C
            // quotient ring from ideal
            _[1]=x^2-5
            sage: R.fetch(M)
            0,  0,-x,
            0,  y,0,
            x*y,0,0
        """
        if (self.type()=='ring') or (self.type()=='qring'):
            # Problem: singular has no clean method to produce
            # a copy of a ring/qring. We use ringlist, but this
            # is only possible if we make self the active ring,
            # use ringlist, and switch back to the previous
            # base ring.
            br=self.parent().current_ring()
            self.set_ring()
            OUT = (self.ringlist()).ring()
            br.set_ring()
            return OUT
        else:
            return self.parent()(self.name())

    def __len__(self):
        """
        Returns the size of this Singular element.

        EXAMPLES::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(2,2)
            sage: len(A)
            4
        """
        return int(self.size())

    def __setitem__(self, n, value):
        """
        Set the n-th element of self to x.

        INPUT:


        -  ``n`` - an integer *or* a 2-tuple (for setting
           matrix elements)

        -  ``value`` - anything (is coerced to a Singular
           object if it is not one already)


        OUTPUT: Changes elements of self.

        EXAMPLES::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(2,2)
            sage: A
            0,0,
            0,0
            sage: A[1,1] = 5
            sage: A
            5,0,
            0,0
            sage: A[1,2] = '5*x + y + z3'
            sage: A
            5,z^3+5*x+y,
            0,0
        """
        P = self.parent()
        if not isinstance(value, SingularElement):
            value = P(value)
        if isinstance(n, tuple):
            if len(n) != 2:
                raise ValueError("If n (=%s) is a tuple, it must be a 2-tuple"%n)
            x, y = n
            P.eval('%s[%s,%s] = %s'%(self.name(), x, y, value.name()))
        else:
            P.eval('%s[%s] = %s'%(self.name(), n, value.name()))

    def __nonzero__(self):
        """
        Returns True if this Singular element is not zero.

        EXAMPLES::

            sage: singular(0).__nonzero__()
            False
            sage: singular(1).__nonzero__()
            True
        """
        P = self.parent()
        return P.eval('%s == 0'%self.name()) == '0'

    def sage_polystring(self):
        r"""
        If this Singular element is a polynomial, return a string
        representation of this polynomial that is suitable for evaluation
        in Python. Thus \* is used for multiplication and \*\* for
        exponentiation. This function is primarily used internally.

        The short=0 option *must* be set for the parent ring or this
        function will not work as expected. This option is set by default
        for rings created using ``singular.ring`` or set using
        ``ring_name.set_ring()``.

        EXAMPLES::

            sage: R = singular.ring(0,'(x,y)')
            sage: f = singular('x^3 + 3*y^11 + 5')
            sage: f
            x^3+3*y^11+5
            sage: f.sage_polystring()
            'x**3+3*y**11+5'
        """
        return str(self).replace('^','**')

    def sage_global_ring(self):
        """
        Return the current basering in Singular as a polynomial ring or quotient ring.

        EXAMPLE::

            sage: singular.eval('ring r1 = (9,x),(a,b,c,d,e,f),(M((1,2,3,0)),wp(2,3),lp)')
            ''
            sage: R = singular('r1').sage_global_ring()
            sage: R
            Multivariate Polynomial Ring in a, b, c, d, e, f over Finite Field in x of size 3^2
            sage: R.term_order()
            Block term order with blocks:
            (Matrix term order with matrix
            [1 2]
            [3 0],
             Weighted degree reverse lexicographic term order with weights (2, 3),
             Lexicographic term order of length 2)

        ::

            sage: singular.eval('ring r2 = (0,x),(a,b,c),dp')
            ''
            sage: singular('r2').sage_global_ring()
            Multivariate Polynomial Ring in a, b, c over Fraction Field of Univariate Polynomial Ring in x over Rational Field

        ::

            sage: singular.eval('ring r3 = (3,z),(a,b,c),dp')
            ''
            sage: singular.eval('minpoly = 1+z+z2+z3+z4')
            ''
            sage: singular('r3').sage_global_ring()
            Multivariate Polynomial Ring in a, b, c over Finite Field in z of size 3^4

        Real and complex fields in both Singular and Sage are defined with a precision.
        The precision in Singular is given in terms of digits, but in Sage it is given
        in terms of bits. So, the digit precision is internally converted to a reasonable
        bit precision::

            sage: singular.eval('ring r4 = (real,20),(a,b,c),dp')
            ''
            sage: singular('r4').sage_global_ring()
            Multivariate Polynomial Ring in a, b, c over Real Field with 70 bits of precision

        The case of complex coefficients is not fully supported, yet, since
        the generator of a complex field in Sage is always called "I"::

            sage: singular.eval('ring r5 = (complex,15,j),(a,b,c),dp')
            ''
            sage: R = singular('r5').sage_global_ring(); R
            Multivariate Polynomial Ring in a, b, c over Complex Field with 54 bits of precision
            sage: R.base_ring()('j')
            Traceback (most recent call last):
            ...
            NameError: name 'j' is not defined
            sage: R.base_ring()('I')
            1.00000000000000*I

        An example where the base ring is a polynomial ring over an extension of the rational field::

            sage: singular.eval('ring r7 = (0,a), (x,y), dp')
            ''
            sage: singular.eval('minpoly = a2 + 1')
            ''
            sage: singular('r7').sage_global_ring()
            Multivariate Polynomial Ring in x, y over Number Field in a with defining polynomial a^2 + 1

        In our last example, the base ring is a quotient ring::

            sage: singular.eval('ring r6 = (9,a), (x,y,z),lp')
            ''
            sage: Q = singular('std(ideal(x^2,x+y^2+z^3))', type='qring')
            sage: Q.sage_global_ring()
            Quotient of Multivariate Polynomial Ring in x, y, z over Finite Field in a of size 3^2 by the ideal (y^4 - y^2*z^3 + z^6, x + y^2 + z^3)

        AUTHOR:

        - Simon King (2011-06-06)

        """
        # extract the ring of coefficients
        singular = self.parent()
        charstr = singular.eval('charstr(basering)').split(',',1)
        from sage.all import ZZ
        is_extension = len(charstr)==2
        if charstr[0]=='integer':
            br = ZZ
            is_extension = False
        elif charstr[0]=='0':
            from sage.all import QQ
            br = QQ
        elif charstr[0]=='real':
            from sage.all import RealField, ceil, log
            prec = singular.eval('ringlist(basering)[1][2][1]')
            br = RealField(ceil((ZZ(prec)+1)/log(2,10)))
            is_extension = False
        elif charstr[0]=='complex':
            from sage.all import ComplexField, ceil, log
            prec = singular.eval('ringlist(basering)[1][2][1]')
            br = ComplexField(ceil((ZZ(prec)+1)/log(2,10)))
            is_extension = False
        else:
            # it ought to be a finite field
            q = ZZ(charstr[0])
            from sage.all import GF
            if q.is_prime():
                br = GF(q)
            else:
                br = GF(q,charstr[1])
                # Singular has no extension of a non-prime field
                is_extension = False

        # We have the base ring of the base ring. But is it
        # an extension?
        if is_extension:
            minpoly = singular.eval('minpoly')
            if minpoly == '0':
                from sage.all import Frac
                BR = Frac(br[charstr[1]])
            else:
                is_short = singular.eval('short')
                if is_short != '0':
                    singular.eval('short=0')
                    minpoly = ZZ[charstr[1]](singular.eval('minpoly'))
                    singular.eval('short=%s'%is_short)
                else:
                    minpoly = ZZ[charstr[1]](minpoly)
                BR = br.extension(minpoly,names=charstr[1])
        else:
            BR = br

        # Now, we form the polynomial ring over BR with the given variables,
        # using Singular's term order
        from sage.rings.polynomial.term_order import termorder_from_singular
        from sage.all import PolynomialRing
        # Meanwhile Singulars quotient rings are also of 'ring' type, not 'qring' as it was in the past.
        # To find out if a singular ring is a quotient ring or not checking for ring type does not help
        # and instead of that we we check if the quotient ring is zero or not:
        if (singular.eval('ideal(basering)==0')=='1'):
            return PolynomialRing(BR, names=singular.eval('varstr(basering)'), order=termorder_from_singular(singular))
        P = PolynomialRing(BR, names=singular.eval('varstr(basering)'), order=termorder_from_singular(singular))
        return P.quotient(singular('ringlist(basering)[4]')._sage_(P), names=singular.eval('varstr(basering)'))

    def sage_poly(self, R=None, kcache=None):
        """
        Returns a Sage polynomial in the ring r matching the provided poly
        which is a singular polynomial.

        INPUT:


        -  ``R`` - (default: None); an optional polynomial ring.
           If it is provided, then you have to make sure that it
           matches the current singular ring as, e.g., returned by
           singular.current_ring(). By default, the output of
           :meth:`sage_global_ring` is used.

        -  ``kcache`` - (default: None); an optional dictionary
           for faster finite field lookups, this is mainly useful for finite
           extension fields


        OUTPUT: MPolynomial

        EXAMPLES::

            sage: R = PolynomialRing(GF(2^8,'a'),2,'xy')
            sage: f=R('a^20*x^2*y+a^10+x')
            sage: f._singular_().sage_poly(R)==f
            True
            sage: R = PolynomialRing(GF(2^8,'a'),1,'x')
            sage: f=R('a^20*x^3+x^2+a^10')
            sage: f._singular_().sage_poly(R)==f
            True

        ::

            sage: P.<x,y> = PolynomialRing(QQ, 2)
            sage: f = x*y**3 - 1/9 * x + 1; f
            x*y^3 - 1/9*x + 1
            sage: singular(f)
            x*y^3-1/9*x+1
            sage: P(singular(f))
            x*y^3 - 1/9*x + 1

        TESTS::

            sage: singular.eval('ring r = (3,z),(a,b,c),dp')
            ''
            sage: singular.eval('minpoly = 1+z+z2+z3+z4')
            ''
            sage: p = singular('z^4*a^3+z^2*a*b*c')
            sage: p.sage_poly()
            (-z^3 - z^2 - z - 1)*a^3 + (z^2)*a*b*c
            sage: singular('z^4')
            (-z3-z2-z-1)

        AUTHORS:

        - Martin Albrecht (2006-05-18)
        - Simon King (2011-06-06): Deal with Singular's short polynomial representation,
          automatic construction of a polynomial ring, if it is not explicitly given.

        .. note::

           For very simple polynomials
           ``eval(SingularElement.sage_polystring())`` is faster than
           SingularElement.sage_poly(R), maybe we should detect the
           crossover point (in dependence of the string length) and
           choose an appropriate conversion strategy
        """
        # TODO: Refactor imports to move this to the top
        from sage.rings.polynomial.multi_polynomial_ring import MPolynomialRing_polydict
        from sage.rings.polynomial.multi_polynomial_libsingular import MPolynomialRing_libsingular
        from sage.rings.polynomial.multi_polynomial_element import MPolynomial_polydict
        from sage.rings.polynomial.polynomial_ring import is_PolynomialRing
        from sage.rings.polynomial.polydict import PolyDict,ETuple
        from sage.rings.polynomial.polynomial_singular_interface import can_convert_to_singular
        from sage.rings.quotient_ring import QuotientRing_generic
        from sage.rings.quotient_ring_element import QuotientRingElement

        ring_is_fine = False
        if R is None:
            ring_is_fine = True
            R = self.sage_global_ring()

        sage_repr = {}
        k = R.base_ring()

        variable_str = "*".join(R.variable_names())

        # This returns a string which looks like a list where the first
        # half of the list is filled with monomials occurring in the
        # Singular polynomial and the second half filled with the matching
        # coefficients.
        #
        # Our strategy is to split the monomials at "*" to get the powers
        # in the single variables and then to split the result to get
        # actual exponent.
        #
        # So e.g. ['x^3*y^3','a'] get's split to
        # [[['x','3'],['y','3']],'a']. We may do this quickly,
        # as we know what to expect.

        is_short = self.parent().eval('short')
        if is_short!='0':
            self.parent().eval('short=0')
            if isinstance(R, MPolynomialRing_libsingular):
                out = R(self)
                self.parent().eval('short=%s'%is_short)
                return out
            singular_poly_list = self.parent().eval("string(coef(%s,%s))"%(\
                    self.name(),variable_str)).split(",")
            self.parent().eval('short=%s'%is_short)
        else:
            if isinstance(R, MPolynomialRing_libsingular):
                return R(self)
            singular_poly_list = self.parent().eval("string(coef(%s,%s))"%(\
                    self.name(),variable_str)).split(",")

        # Directly treat constants
        if singular_poly_list[0] in ['1', '(1.000e+00)']:
            return R(singular_poly_list[1])

        coeff_start = len(singular_poly_list) // 2

        # Singular 4 puts parentheses around floats and sign outside them
        charstr = self.parent().eval('charstr(basering)').split(',',1)
        if charstr[0] in ['real', 'complex']:
              for i in range(coeff_start, 2*coeff_start):
                  singular_poly_list[i] = singular_poly_list[i].replace('(','').replace(')','')

        if isinstance(R,(MPolynomialRing_polydict,QuotientRing_generic)) and (ring_is_fine or can_convert_to_singular(R)):
            # we need to lookup the index of a given variable represented
            # through a string
            var_dict = dict(zip(R.variable_names(),range(R.ngens())))

            ngens = R.ngens()

            for i in range(coeff_start):
                exp = dict()
                monomial = singular_poly_list[i]

                if monomial!="1":
                    variables = [var.split("^") for var in monomial.split("*") ]
                    for e in variables:
                        var = e[0]
                        if len(e)==int(2):
                            power = int(e[1])
                        else:
                            power=1
                        exp[var_dict[var]]=power

                if kcache is None:
                    sage_repr[ETuple(exp,ngens)]=k(singular_poly_list[coeff_start+i])
                else:
                    elem = singular_poly_list[coeff_start+i]
                    if elem not in kcache:
                        kcache[elem] = k( elem )
                    sage_repr[ETuple(exp,ngens)]= kcache[elem]

            p = MPolynomial_polydict(R,PolyDict(sage_repr,force_int_exponents=False,force_etuples=False))
            if isinstance(R, MPolynomialRing_polydict):
                return p
            else:
                return QuotientRingElement(R,p,reduce=False)

        elif is_PolynomialRing(R) and (ring_is_fine or can_convert_to_singular(R)):

            sage_repr = [0]*int(self.deg()+1)

            for i in range(coeff_start):
                monomial = singular_poly_list[i]
                exp = int(0)

                if monomial!="1":
                    term =  monomial.split("^")
                    if len(term)==int(2):
                        exp = int(term[1])
                    else:
                        exp = int(1)

                if kcache is None:
                    sage_repr[exp] = k(singular_poly_list[coeff_start+i])
                else:
                    elem = singular_poly_list[coeff_start+i]
                    if elem not in kcache:
                        kcache[elem] = k( elem )
                    sage_repr[ exp ]= kcache[elem]

            return R(sage_repr)

        else:
            raise TypeError("Cannot coerce %s into %s"%(self,R))

    def sage_matrix(self, R, sparse=True):
        """
        Returns Sage matrix for self

        INPUT:

        -  ``R`` - (default: None); an optional ring, over which
           the resulting matrix is going to be defined.
           By default, the output of :meth:`sage_global_ring` is used.

        - ``sparse`` - (default: True); determines whether the
          resulting matrix is sparse or not.

        EXAMPLES::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(2,2)
            sage: A.sage_matrix(ZZ)
            [0 0]
            [0 0]
            sage: A.sage_matrix(RDF)
            [0.0 0.0]
            [0.0 0.0]
        """
        from sage.matrix.constructor import Matrix
        nrows, ncols = int(self.nrows()),int(self.ncols())

        if R is None:
            R = self.sage_global_ring()
            A = Matrix(R, nrows, ncols, sparse=sparse)
            #this is slow
            for x in range(nrows):
                for y in range(ncols):
                    A[x,y]=self[x+1,y+1].sage_poly(R)
            return A

        A = Matrix(R, nrows, ncols, sparse=sparse)
        #this is slow
        for x in range(nrows):
            for y in range(ncols):
                A[x,y]=R(self[x+1,y+1])

        return A

    def _sage_(self, R=None):
        r"""
        Convert self to Sage.

        EXAMPLES::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(2,2)
            sage: A.sage(ZZ)   # indirect doctest
            [0 0]
            [0 0]
            sage: A = random_matrix(ZZ,3,3); A
            [ -8   2   0]
            [  0   1  -1]
            [  2   1 -95]
            sage: As = singular(A); As
            -8     2     0
            0     1    -1
            2     1   -95
            sage: As.sage()
            [ -8   2   0]
            [  0   1  -1]
            [  2   1 -95]

        ::

            sage: singular.eval('ring R = integer, (x,y,z),lp')
            '// ** redefining R (ring R = integer, (x,y,z),lp;)'
            sage: I = singular.ideal(['x^2','y*z','z+x'])
            sage: I.sage()
            Ideal (x^2, y*z, x + z) of Multivariate Polynomial Ring in x, y, z over Integer Ring

        ::

            sage: singular('ringlist(basering)').sage()
            [['integer'], ['x', 'y', 'z'], [['lp', (1, 1, 1)], ['C', (0)]], Ideal (0) of Multivariate Polynomial Ring in x, y, z over Integer Ring]

        ::

            sage: singular.eval('ring r10 = (9,a), (x,y,z),lp')
            ''
            sage: singular.eval('setring R')
            ''
            sage: singular('r10').sage()
            Multivariate Polynomial Ring in x, y, z over Finite Field in a of size 3^2

        Note that the current base ring has not been changed by asking for another ring::

            sage: singular('basering')
            polynomial ring, over a domain, global ordering
            //   coeff. ring is : integer
            //   number of vars : 3
            //        block   1 : ordering lp
            //                  : names    x y z
            //        block   2 : ordering C

        ::

            sage: singular.eval('setring r10')
            ''
            sage: Q = singular('std(ideal(x^2,x+y^2+z^3))', type='qring')
            sage: Q.sage()
            Quotient of Multivariate Polynomial Ring in x, y, z over Finite Field in a of size 3^2 by the ideal (y^4 - y^2*z^3 + z^6, x + y^2 + z^3)
            sage: singular('x^2+y').sage()
            x^2 + y
            sage: singular('x^2+y').sage().parent()
            Quotient of Multivariate Polynomial Ring in x, y, z over Finite Field in a of size 3^2 by the ideal (y^4 - y^2*z^3 + z^6, x + y^2 + z^3)

        Test that :trac:`18848` is fixed::

            sage: singular(5).sage()
            5
            sage: type(singular(int(5)).sage())
            <type 'sage.rings.integer.Integer'>

        """
        typ = self.type()
        if typ=='poly':
            return self.sage_poly(R)
        elif typ=='int':
            return sage.rings.integer.Integer(repr(self))
        elif typ == 'module':
            return self.sage_matrix(R,sparse=True)
        elif typ == 'matrix':
            return self.sage_matrix(R,sparse=False)
        elif typ == 'list':
            return [ f._sage_(R) for f in self ]
        elif typ == 'intvec':
            from sage.modules.free_module_element import vector
            return vector([sage.rings.integer.Integer(str(e)) for e in self])
        elif typ == 'intmat':
            from sage.matrix.constructor import matrix
            from sage.rings.integer_ring import ZZ
            A =  matrix(ZZ, int(self.nrows()), int(self.ncols()))
            for i in xrange(A.nrows()):
                for j in xrange(A.ncols()):
                    A[i,j] = sage.rings.integer.Integer(str(self[i+1,j+1]))
            return A
        elif typ == 'string':
            return repr(self)
        elif typ == 'ideal':
            R = R or self.sage_global_ring()
            return R.ideal([p.sage_poly(R) for p in self])
        elif typ in ['ring', 'qring']:
            br = singular('basering')
            self.set_ring()
            R = self.sage_global_ring()
            br.set_ring()
            return R
        raise NotImplementedError("Coercion of this datatype not implemented yet")

    def is_string(self):
        """
        Tell whether this element is a string.

        EXAMPLES::

            sage: singular('"abc"').is_string()
            True
            sage: singular('1').is_string()
            False

        """
        return self.type() == 'string'

    def set_ring(self):
        """
        Sets the current ring in Singular to be self.

        EXAMPLES::

            sage: R = singular.ring(7, '(a,b)', 'ds')
            sage: S = singular.ring('real', '(a,b)', 'lp')
            sage: singular.current_ring()
            polynomial ring, over a field, global ordering
            //   characteristic : 0 (real)
            //   number of vars : 2
            //        block   1 : ordering lp
            //                  : names    a b
            //        block   2 : ordering C
            sage: R.set_ring()
            sage: singular.current_ring()
            polynomial ring, over a field, local/mixed ordering
            //   characteristic : 7
            //   number of vars : 2
            //        block   1 : ordering ds
            //                  : names    a b
            //        block   2 : ordering C
        """
        self.parent().set_ring(self)


    def sage_flattened_str_list(self):
        """
        EXAMPLES::

            sage: R=singular.ring(0,'(x,y)','dp')
            sage: RL = R.ringlist()
            sage: RL.sage_flattened_str_list()
            ['0', 'x', 'y', 'dp', '1,1', 'C', '0', '_[1]=0']
        """
        s = str(self)
        c = '\[[0-9]*\]:'
        r = re.compile(c)
        s = r.sub('',s).strip()
        return s.split()

    def sage_structured_str_list(self):
        r"""
        If self is a Singular list of lists of Singular elements, returns
        corresponding Sage list of lists of strings.

        EXAMPLES::

            sage: R=singular.ring(0,'(x,y)','dp')
            sage: RL=R.ringlist()
            sage: RL
            [1]:
               0
            [2]:
               [1]:
                  x
               [2]:
                  y
            [3]:
               [1]:
                  [1]:
                     dp
                  [2]:
                     1,1
               [2]:
                  [1]:
                     C
                  [2]:
                     0
            [4]:
               _[1]=0
            sage: RL.sage_structured_str_list()
            ['0', ['x', 'y'], [['dp', '1,\n1 '], ['C', '0 ']], '0']
        """
        if not (self.type()=='list'):
            return str(self)
        return [X.sage_structured_str_list() for X in self]

    def _tab_completion(self):
        """
        Returns the possible tab-completions for self. In this case, we
        just return all the tab completions for the Singular object.

        EXAMPLES::

            sage: R = singular.ring(0,'(x,y)','dp')
            sage: R._tab_completion()
            ['exteriorPower',
             ...
             'stdfglm']
        """
        return self.parent()._tab_completion()

    def type(self):
        """
        Returns the internal type of this element.

        EXAMPLES::

            sage: R = PolynomialRing(GF(2^8,'a'),2,'x')
            sage: R._singular_().type()
            'ring'
            sage: fs = singular('x0^2','poly')
            sage: fs.type()
            'poly'
        """
        # singular reports // $varname $type $stuff
        p = re.compile("// [\w]+ (\w+) [\w]*")
        m = p.match(self.parent().eval("type(%s)"%self.name()))
        return m.group(1)

    def __iter__(self):
        """
        EXAMPLES::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(2,2)
            sage: list(iter(A))
            [[0], [0]]
            sage: A[1,1] = 1; A[1,2] = 2
            sage: A[2,1] = 3; A[2,2] = 4
            sage: list(iter(A))
            [[1,3], [2,4]]
        """
        if self.type()=='matrix':
            l = self.ncols()
        else:
            l = len(self)
        for i in range(1, l+1):
            yield self[i]

    def _singular_(self):
        """
        EXAMPLES::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(2,2)
            sage: A._singular_() is A
            True
        """
        return self

    def attrib(self, name, value=None):
        """
        Get and set attributes for self.

        INPUT:


        -  ``name`` - string to choose the attribute

        -  ``value`` - boolean value or None for reading,
           (default:None)


        VALUES: isSB - the standard basis property is set by all commands
        computing a standard basis like groebner, std, stdhilb etc.; used
        by lift, dim, degree, mult, hilb, vdim, kbase isHomog - the weight
        vector for homogeneous or quasihomogeneous ideals/modules isCI -
        complete intersection property isCM - Cohen-Macaulay property rank
        - set the rank of a module (see nrows) withSB - value of type
        ideal, resp. module, is std withHilb - value of type intvec is
        hilb(_,1) (see hilb) withRes - value of type list is a free
        resolution withDim - value of type int is the dimension (see dim)
        withMult - value of type int is the multiplicity (see mult)

        EXAMPLE::

            sage: P.<x,y,z> = PolynomialRing(QQ)
            sage: I = Ideal([z^2, y*z, y^2, x*z, x*y, x^2])
            sage: Ibar = I._singular_()
            sage: Ibar.attrib('isSB')
            0
            sage: singular.eval('vdim(%s)'%Ibar.name()) # sage7 name is random
            // ** sage7 is no standard basis
            4
            sage: Ibar.attrib('isSB',1)
            sage: singular.eval('vdim(%s)'%Ibar.name())
            '4'
        """
        if value is None:
            return int(self.parent().eval('attrib(%s,"%s")'%(self.name(),name)))
        else:
            self.parent().eval('attrib(%s,"%s",%d)'%(self.name(),name,value))

class SingularFunction(ExpectFunction):
    def _sage_doc_(self):
        """
        EXAMPLES::

            sage: 'groebner' in singular.groebner._sage_doc_()
            True
        """
        if not nodes:
            generate_docstring_dictionary()

        prefix = \
"""
This function is an automatically generated pexpect wrapper around the Singular
function '%s'.

EXAMPLE::

    sage: groebner = singular.groebner
    sage: P.<x, y> = PolynomialRing(QQ)
    sage: I = P.ideal(x^2-y, y+x)
    sage: groebner(singular(I))
    x+y,
    y^2-y
"""%(self._name,)
        prefix2 = \
"""

The Singular documentation for '%s' is given below.
"""%(self._name,)

        try:
            return prefix + prefix2 + nodes[node_names[self._name]]
        except KeyError:
            return prefix

class SingularFunctionElement(FunctionElement):
    def _sage_doc_(self):
        r"""
        EXAMPLES::

            sage: R = singular.ring(0, '(x,y,z)', 'dp')
            sage: A = singular.matrix(2,2)
            sage: 'matrix_expression' in A.nrows._sage_doc_()
            True
        """
        if not nodes:
            generate_docstring_dictionary()
        try:
            return nodes[node_names[self._name]]
        except KeyError:
            return ""

def is_SingularElement(x):
    r"""
    Returns True is x is of type ``SingularElement``.

    EXAMPLES::

        sage: from sage.interfaces.singular import is_SingularElement
        sage: is_SingularElement(singular(2))
        True
        sage: is_SingularElement(2)
        False
    """
    return isinstance(x, SingularElement)

# This is only for backwards compatibility, in order to be able
# to unpickle the invalid objects that are in the pickle jar.
def reduce_load():
    """
    This is for backwards compatibility only.

    To be precise, it only serves at unpickling the invalid
    singular elements that are stored in the pickle jar.

    EXAMPLES::

        sage: from sage.interfaces.singular import reduce_load
        sage: reduce_load()
        doctest:...: DeprecationWarning: This function is only used to unpickle invalid objects
        See http://trac.sagemath.org/18848 for details.
        (invalid object -- defined in terms of closed session)

    By :trac:`18848`, pickling actually often works::

        sage: loads(dumps(singular.ring()))
        polynomial ring, over a field, global ordering
        //   characteristic : 0
        //   number of vars : 1
        //        block   1 : ordering lp
        //                  : names    x
        //        block   2 : ordering C

    """
    deprecation(18848, "This function is only used to unpickle invalid objects")
    return SingularElement(None, None, None)

nodes = {}
node_names = {}

def generate_docstring_dictionary():
    """
    Generate global dictionaries which hold the docstrings for
    Singular functions.

    EXAMPLE::

        sage: from sage.interfaces.singular import generate_docstring_dictionary
        sage: generate_docstring_dictionary()
    """
    from sage.env import SAGE_LOCAL

    global nodes
    global node_names

    nodes.clear()
    node_names.clear()

    singular_docdir = "/usr/share/doc/singular-doc/"

    new_node = re.compile("File: singular\.hlp,  Node: ([^,]*),.*")
    new_lookup = re.compile("\* ([^:]*):*([^.]*)\..*")

    L, in_node, curr_node = [], False, None

    for line in open(singular_docdir + "singular.hlp"):
        m = re.match(new_node,line)
        if m:
            # a new node starts
            in_node = True
            nodes[curr_node] = "".join(L)
            L = []
            curr_node, = m.groups()
        elif in_node: # we are in a node
           L.append(line)
        else:
           m = re.match(new_lookup, line)
           if m:
               a,b = m.groups()
               node_names[a] = b.strip()

        if line == "6 Index\n":
            in_node = False

    nodes[curr_node] = "".join(L) # last node

def get_docstring(name):
    """
    Return the docstring for the function ``name``.

    INPUT:

    - ``name`` - a Singular function name

    EXAMPLE::

        sage: from sage.interfaces.singular import get_docstring
        sage: 'groebner' in get_docstring('groebner')
        True
        sage: 'standard.lib' in get_docstring('groebner')
        True

    """
    if not nodes:
        generate_docstring_dictionary()
    try:
        return nodes[node_names[name]]
    except KeyError:
        return ""

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

singular = Singular()

def reduce_load_Singular():
    """
    EXAMPLES::

        sage: from sage.interfaces.singular import reduce_load_Singular
        sage: reduce_load_Singular()
        Singular
    """
    return singular


def singular_console():
    """
    Spawn a new Singular command-line session.

    EXAMPLES::

        sage: singular_console() #not tested
                             SINGULAR                             /  Development
         A Computer Algebra System for Polynomial Computations   /   version 3-0-4
                                                               0<
             by: G.-M. Greuel, G. Pfister, H. Schoenemann        \   Nov 2007
        FB Mathematik der Universitaet, D-67653 Kaiserslautern    \
    """
    from sage.repl.rich_output.display_manager import get_display_manager
    if not get_display_manager().is_in_terminal():
        raise RuntimeError('Can use the console only in the terminal. Try %%singular magics instead.')
    os.system('Singular')


def singular_version():
    """
    Returns the version of Singular being used.

    EXAMPLES:
    """
    return singular.eval('system("--version");')



class SingularGBLogPrettyPrinter:
    """
    A device which prints Singular Groebner basis computation logs
    more verbatim.
    """
    rng_chng = re.compile("\[\d+:\d+\]")# [m:n] internal ring change to
                                        # poly representation with
                                        # exponent bound m and n words in
                                        # exponent vector
    new_elem = re.compile("s")          # found a new element of the standard basis
    red_zero = re.compile("-")          # reduced a pair/S-polynomial to 0
    red_post = re.compile("\.")         # postponed a reduction of a pair/S-polynomial
    cri_hilb = re.compile("h")          # used Hilbert series criterion
    hig_corn = re.compile("H\(\d+\)")   # found a 'highest corner' of degree d, no need to consider higher degrees
    num_crit = re.compile("\(\d+\)")    # n critical pairs are still to be reduced
    red_num =  re.compile("\(S:\d+\)")  # doing complete reduction of n elements
    deg_lead = re.compile("\d+")        # the degree of the leading terms is currently d

    # SlimGB
    red_para = re.compile("M\[(\d+),(\d+)\]") # parallel reduction of n elements with m non-zero output elements
    red_betr = re.compile("b")                # exchange of a reductor by a 'better' one
    non_mini = re.compile("e")                # a new reductor with non-minimal leading term

    crt_lne1 = re.compile("product criterion:(\d+) chain criterion:(\d+)")
    crt_lne2 = re.compile("NF:(\d+) product criterion:(\d+), ext_product criterion:(\d+)")

    pat_sync = re.compile("1\+(\d+);")

    global_pattern = re.compile("(\[\d+:\d+\]|s|-|\.|h|H\(\d+\)|\(\d+\)|\(S:\d+\)|\d+|M\[\d+,[b,e]*\d+\]|b|e).*")

    def __init__(self, verbosity=1):
        """
        Construct a new Singular Groebner Basis log pretty printer.

        INPUT:

        - ``verbosity`` - how much information should be printed
          (between 0 and 3)

        EXAMPLE::

            sage: from sage.interfaces.singular import SingularGBLogPrettyPrinter
            sage: s0 = SingularGBLogPrettyPrinter(verbosity=0)
            sage: s1 = SingularGBLogPrettyPrinter(verbosity=1)
            sage: s0.write("[1:2]12")

            sage: s1.write("[1:2]12")
            Leading term degree: 12.
        """
        self.verbosity = verbosity

        self.curr_deg = 0 # current degree
        self.max_deg = 0  # maximal degree in total

        self.nf = 0 # number of normal forms computed (SlimGB only)
        self.prod = 0 # number of S-polynomials discarded using product criterion
        self.ext_prod = 0 # number of S-polynomials discarded using extended product criterion
        self.chain = 0 # number of S-polynomials discarded using chain criterion

        self.storage = "" # stores incomplete strings
        self.sync = None # should we expect a sync integer?

    def write(self, s):
        """
        EXAMPLE::

            sage: from sage.interfaces.singular import SingularGBLogPrettyPrinter
            sage: s3 = SingularGBLogPrettyPrinter(verbosity=3)
            sage: s3.write("(S:1337)")
            Performing complete reduction of 1337 elements.
            sage: s3.write("M[389,12]")
            Parallel reduction of 389 elements with 12 non-zero output elements.
        """
        verbosity = self.verbosity

        if self.storage:
            s = self.storage + s
            self.storage = ""

        for line in s.splitlines():
            # deal with the Sage <-> Singular syncing code
            match = re.match(SingularGBLogPrettyPrinter.pat_sync,line)
            if match:
                self.sync = int(match.groups()[0])
                continue

            if self.sync and line == "%d"%(self.sync+1):
                self.sync = None
                continue

            if line.endswith(";"):
                continue
            if line.startswith(">"):
                continue

            if line.startswith("std") or line.startswith("slimgb"):
                continue

            # collect stats returned about avoided reductions to zero
            match = re.match(SingularGBLogPrettyPrinter.crt_lne1,line)
            if match:
                self.prod,self.chain = map(int,re.match(SingularGBLogPrettyPrinter.crt_lne1,line).groups())
                self.storage = ""
                continue
            match = re.match(SingularGBLogPrettyPrinter.crt_lne2,line)
            if match:
                self.nf,self.prod,self.ext_prod = map(int,re.match(SingularGBLogPrettyPrinter.crt_lne2,line).groups())
                self.storage = ""
                continue

            while line:
                match = re.match(SingularGBLogPrettyPrinter.global_pattern, line)
                if not match:
                    self.storage = line
                    line = None
                    continue

                token, = match.groups()
                line = line[len(token):]

                if re.match(SingularGBLogPrettyPrinter.rng_chng,token):
                    continue

                elif re.match(SingularGBLogPrettyPrinter.new_elem,token) and verbosity >= 3:
                    print("New element found.")

                elif re.match(SingularGBLogPrettyPrinter.red_zero,token) and verbosity >= 2:
                    print("Reduction to zero.")

                elif re.match(SingularGBLogPrettyPrinter.red_post, token) and verbosity >= 2:
                    print("Reduction postponed.")

                elif re.match(SingularGBLogPrettyPrinter.cri_hilb, token) and verbosity >= 2:
                    print("Hilber series criterion applied.")

                elif re.match(SingularGBLogPrettyPrinter.hig_corn, token) and verbosity >= 1:
                    print("Maximal degree found: %s" % token)

                elif re.match(SingularGBLogPrettyPrinter.num_crit, token) and verbosity >= 1:
                    print("Leading term degree: %2d. Critical pairs: %s."%(self.curr_deg,token[1:-1]))

                elif re.match(SingularGBLogPrettyPrinter.red_num, token) and verbosity >= 3:
                    print("Performing complete reduction of %s elements."%token[3:-1])

                elif re.match(SingularGBLogPrettyPrinter.deg_lead, token):
                    if verbosity >= 1:
                        print("Leading term degree: %2d." % int(token))
                    self.curr_deg = int(token)
                    if self.max_deg < self.curr_deg:
                        self.max_deg = self.curr_deg

                elif re.match(SingularGBLogPrettyPrinter.red_para, token) and verbosity >= 3:
                    m,n = re.match(SingularGBLogPrettyPrinter.red_para,token).groups()
                    print("Parallel reduction of %s elements with %s non-zero output elements." % (m, n))

                elif re.match(SingularGBLogPrettyPrinter.red_betr, token) and verbosity >= 3:
                    print("Replaced reductor by 'better' one.")

                elif re.match(SingularGBLogPrettyPrinter.non_mini, token) and verbosity >= 2:
                    print("New reductor with non-minimal leading term found.")

    def flush(self):
        """
        EXAMPLE::

            sage: from sage.interfaces.singular import SingularGBLogPrettyPrinter
            sage: s3 = SingularGBLogPrettyPrinter(verbosity=3)
            sage: s3.flush()
        """
        sys.stdout.flush()

class SingularGBDefaultContext:
    """
    Within this context all Singular Groebner basis calculations are
    reduced automatically.

    AUTHORS:

    - Martin Albrecht
    - Simon King
    """
    def __init__(self, singular=None):
        """
        Within this context all Singular Groebner basis calculations
        are reduced automatically.

        INPUT:

        -  ``singular`` - Singular instance (default: default instance)

        EXAMPLE::

            sage: from sage.interfaces.singular import SingularGBDefaultContext
            sage: P.<a,b,c> = PolynomialRing(QQ,3, order='lex')
            sage: I = sage.rings.ideal.Katsura(P,3)
            sage: singular.option('noredTail')
            sage: singular.option('noredThrough')
            sage: Is = I._singular_()
            sage: gb = Is.groebner()
            sage: gb
            84*c^4-40*c^3+c^2+c,
            7*b+210*c^3-79*c^2+3*c,
            a+2*b+2*c-1

        ::

            sage: with SingularGBDefaultContext(): rgb = Is.groebner()
            sage: rgb
            84*c^4-40*c^3+c^2+c,
            7*b+210*c^3-79*c^2+3*c,
            7*a-420*c^3+158*c^2+8*c-7

        Note that both bases are Groebner bases because they have
        pairwise prime leading monomials but that the monic version of
        the last element in ``rgb`` is smaller than the last element
        of ``gb`` with respect to the lexicographical term ordering. ::

            sage: (7*a-420*c^3+158*c^2+8*c-7)/7 < (a+2*b+2*c-1)
            True

        .. note::

           This context is used automatically internally whenever a
           Groebner basis is computed so the user does not need to use
           it manually.
        """
        if singular is None:
            from sage.interfaces.all import singular as singular_default
            singular = singular_default
        self.singular = singular

    def __enter__(self):
        """
        EXAMPLE::

            sage: from sage.interfaces.singular import SingularGBDefaultContext
            sage: P.<a,b,c> = PolynomialRing(QQ,3, order='lex')
            sage: I = sage.rings.ideal.Katsura(P,3)
            sage: singular.option('noredTail')
            sage: singular.option('noredThrough')
            sage: Is = I._singular_()
            sage: with SingularGBDefaultContext(): rgb = Is.groebner()
            sage: rgb
            84*c^4-40*c^3+c^2+c,
            7*b+210*c^3-79*c^2+3*c,
            7*a-420*c^3+158*c^2+8*c-7
        """
        from sage.interfaces.singular import SingularError
        try:
            self.bck_degBound = int(self.singular.eval('degBound'))
        except SingularError:
            self.bck_degBound = int(0)
        try:
            self.bck_multBound = int(self.singular.eval('multBound'))
        except SingularError:
            self.bck_multBound = int(0)
        self.o = self.singular.option("get")
        self.singular.option('set',self.singular._saved_options)
        self.singular.option("redSB")
        self.singular.option("redTail")
        try:
            self.singular.eval('degBound=0')
        except SingularError:
            pass
        try:
            self.singular.eval('multBound=0')
        except SingularError:
            pass

    def __exit__(self, typ, value, tb):
        """
        EXAMPLE::

            sage: from sage.interfaces.singular import SingularGBDefaultContext
            sage: P.<a,b,c> = PolynomialRing(QQ,3, order='lex')
            sage: I = sage.rings.ideal.Katsura(P,3)
            sage: singular.option('noredTail')
            sage: singular.option('noredThrough')
            sage: Is = I._singular_()
            sage: with SingularGBDefaultContext(): rgb = Is.groebner()
            sage: rgb
            84*c^4-40*c^3+c^2+c,
            7*b+210*c^3-79*c^2+3*c,
            7*a-420*c^3+158*c^2+8*c-7
        """
        from sage.interfaces.singular import SingularError
        self.singular.option("set",self.o)
        try:
            self.singular.eval('degBound=%d'%self.bck_degBound)
        except SingularError:
            pass
        try:
            self.singular.eval('multBound=%d'%self.bck_multBound)
        except SingularError:
            pass

def singular_gb_standard_options(func):
    r"""
    Decorator to force a reduced Singular groebner basis.

    TESTS::

        sage: P.<a,b,c,d,e> = PolynomialRing(GF(127))
        sage: J = sage.rings.ideal.Cyclic(P).homogenize()
        sage: from sage.misc.sageinspect import sage_getsource
        sage: "basis" in sage_getsource(J.interreduced_basis) #indirect doctest
        True

    The following tests against a bug that was fixed in :trac:`11298`::

        sage: from sage.misc.sageinspect import sage_getsourcelines, sage_getargspec
        sage: P.<x,y> = QQ[]
        sage: I = P*[x,y]
        sage: sage_getargspec(I.interreduced_basis)
        ArgSpec(args=['self'], varargs=None, keywords=None, defaults=None)
        sage: sage_getsourcelines(I.interreduced_basis)
        (['    @singular_gb_standard_options\n',
          '    @libsingular_gb_standard_options\n',
          '    def interreduced_basis(self):\n', '
          ...
          '        return self.basis.reduced()\n'], ...)

    .. note::

       This decorator is used automatically internally so the user
       does not need to use it manually.
    """
    from sage.misc.decorators import sage_wraps
    @sage_wraps(func)
    def wrapper(*args, **kwds):
        with SingularGBDefaultContext():
            return func(*args, **kwds)
    return wrapper