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


import contextlib
import warnings

from llvmlite import ir

import numpy as np
import operator

from numba.core.imputils import (lower_builtin, impl_ret_borrowed,
                                    impl_ret_new_ref, impl_ret_untracked)
from numba.core.typing import signature
from numba.core.extending import intrinsic, overload, register_jitable
from numba.core import types, cgutils
from numba.core.errors import TypingError, NumbaTypeError, \
    NumbaPerformanceWarning
from .arrayobj import make_array, _empty_nd_impl, array_copy
from numba.np import numpy_support as np_support

ll_char = ir.IntType(8)
ll_char_p = ll_char.as_pointer()
ll_void_p = ll_char_p
ll_intc = ir.IntType(32)
ll_intc_p = ll_intc.as_pointer()
intp_t = cgutils.intp_t
ll_intp_p = intp_t.as_pointer()


# fortran int type, this needs to match the F_INT C declaration in
# _lapack.c and is present to accommodate potential future 64bit int
# based LAPACK use.
F_INT_nptype = np.int32
F_INT_nbtype = types.int32

# BLAS kinds as letters
_blas_kinds = {
    types.float32: 's',
    types.float64: 'd',
    types.complex64: 'c',
    types.complex128: 'z',
}


def get_blas_kind(dtype, func_name="<BLAS function>"):
    kind = _blas_kinds.get(dtype)
    if kind is None:
        raise TypeError("unsupported dtype for %s()" % (func_name,))
    return kind


def ensure_blas():
    try:
        import scipy.linalg.cython_blas
    except ImportError:
        raise ImportError("scipy 0.16+ is required for linear algebra")


def ensure_lapack():
    try:
        import scipy.linalg.cython_lapack
    except ImportError:
        raise ImportError("scipy 0.16+ is required for linear algebra")


def make_constant_slot(context, builder, ty, val):
    const = context.get_constant_generic(builder, ty, val)
    return cgutils.alloca_once_value(builder, const)


class _BLAS:
    """
    Functions to return type signatures for wrapped
    BLAS functions.
    """

    def __init__(self):
        ensure_blas()

    @classmethod
    def numba_xxnrm2(cls, dtype):
        rtype = getattr(dtype, "underlying_float", dtype)
        sig = types.intc(types.char,             # kind
                         types.intp,             # n
                         types.CPointer(dtype),  # x
                         types.intp,             # incx
                         types.CPointer(rtype))  # returned

        return types.ExternalFunction("numba_xxnrm2", sig)

    @classmethod
    def numba_xxgemm(cls, dtype):
        sig = types.intc(
            types.char,             # kind
            types.char,             # transa
            types.char,             # transb
            types.intp,             # m
            types.intp,             # n
            types.intp,             # k
            types.CPointer(dtype),  # alpha
            types.CPointer(dtype),  # a
            types.intp,             # lda
            types.CPointer(dtype),  # b
            types.intp,             # ldb
            types.CPointer(dtype),  # beta
            types.CPointer(dtype),  # c
            types.intp              # ldc
        )
        return types.ExternalFunction("numba_xxgemm", sig)


class _LAPACK:
    """
    Functions to return type signatures for wrapped
    LAPACK functions.
    """

    def __init__(self):
        ensure_lapack()

    @classmethod
    def numba_xxgetrf(cls, dtype):
        sig = types.intc(types.char,                   # kind
                         types.intp,                   # m
                         types.intp,                   # n
                         types.CPointer(dtype),        # a
                         types.intp,                   # lda
                         types.CPointer(F_INT_nbtype)  # ipiv
                         )
        return types.ExternalFunction("numba_xxgetrf", sig)

    @classmethod
    def numba_ez_xxgetri(cls, dtype):
        sig = types.intc(types.char,                   # kind
                         types.intp,                   # n
                         types.CPointer(dtype),        # a
                         types.intp,                   # lda
                         types.CPointer(F_INT_nbtype)  # ipiv
                         )
        return types.ExternalFunction("numba_ez_xxgetri", sig)

    @classmethod
    def numba_ez_rgeev(cls, dtype):
        sig = types.intc(types.char,             # kind
                         types.char,             # jobvl
                         types.char,             # jobvr
                         types.intp,             # n
                         types.CPointer(dtype),  # a
                         types.intp,             # lda
                         types.CPointer(dtype),  # wr
                         types.CPointer(dtype),  # wi
                         types.CPointer(dtype),  # vl
                         types.intp,             # ldvl
                         types.CPointer(dtype),  # vr
                         types.intp              # ldvr
                         )
        return types.ExternalFunction("numba_ez_rgeev", sig)

    @classmethod
    def numba_ez_cgeev(cls, dtype):
        sig = types.intc(types.char,             # kind
                         types.char,             # jobvl
                         types.char,             # jobvr
                         types.intp,             # n
                         types.CPointer(dtype),  # a
                         types.intp,             # lda
                         types.CPointer(dtype),  # w
                         types.CPointer(dtype),  # vl
                         types.intp,             # ldvl
                         types.CPointer(dtype),  # vr
                         types.intp              # ldvr
                         )
        return types.ExternalFunction("numba_ez_cgeev", sig)

    @classmethod
    def numba_ez_xxxevd(cls, dtype):
        wtype = getattr(dtype, "underlying_float", dtype)
        sig = types.intc(types.char,             # kind
                         types.char,             # jobz
                         types.char,             # uplo
                         types.intp,             # n
                         types.CPointer(dtype),  # a
                         types.intp,             # lda
                         types.CPointer(wtype),  # w
                         )
        return types.ExternalFunction("numba_ez_xxxevd", sig)

    @classmethod
    def numba_xxpotrf(cls, dtype):
        sig = types.intc(types.char,             # kind
                         types.char,             # uplo
                         types.intp,             # n
                         types.CPointer(dtype),  # a
                         types.intp              # lda
                         )
        return types.ExternalFunction("numba_xxpotrf", sig)

    @classmethod
    def numba_ez_gesdd(cls, dtype):
        stype = getattr(dtype, "underlying_float", dtype)
        sig = types.intc(
            types.char,             # kind
            types.char,             # jobz
            types.intp,             # m
            types.intp,             # n
            types.CPointer(dtype),  # a
            types.intp,             # lda
            types.CPointer(stype),  # s
            types.CPointer(dtype),  # u
            types.intp,             # ldu
            types.CPointer(dtype),  # vt
            types.intp              # ldvt
        )

        return types.ExternalFunction("numba_ez_gesdd", sig)

    @classmethod
    def numba_ez_geqrf(cls, dtype):
        sig = types.intc(
            types.char,             # kind
            types.intp,             # m
            types.intp,             # n
            types.CPointer(dtype),  # a
            types.intp,             # lda
            types.CPointer(dtype),  # tau
        )
        return types.ExternalFunction("numba_ez_geqrf", sig)

    @classmethod
    def numba_ez_xxgqr(cls, dtype):
        sig = types.intc(
            types.char,             # kind
            types.intp,             # m
            types.intp,             # n
            types.intp,             # k
            types.CPointer(dtype),  # a
            types.intp,             # lda
            types.CPointer(dtype),  # tau
        )
        return types.ExternalFunction("numba_ez_xxgqr", sig)

    @classmethod
    def numba_ez_gelsd(cls, dtype):
        rtype = getattr(dtype, "underlying_float", dtype)
        sig = types.intc(
            types.char,                 # kind
            types.intp,                 # m
            types.intp,                 # n
            types.intp,                 # nrhs
            types.CPointer(dtype),      # a
            types.intp,                 # lda
            types.CPointer(dtype),      # b
            types.intp,                 # ldb
            types.CPointer(rtype),      # S
            types.float64,              # rcond
            types.CPointer(types.intc)  # rank
        )
        return types.ExternalFunction("numba_ez_gelsd", sig)

    @classmethod
    def numba_xgesv(cls, dtype):
        sig = types.intc(
            types.char,                    # kind
            types.intp,                    # n
            types.intp,                    # nhrs
            types.CPointer(dtype),         # a
            types.intp,                    # lda
            types.CPointer(F_INT_nbtype),  # ipiv
            types.CPointer(dtype),         # b
            types.intp                     # ldb
        )
        return types.ExternalFunction("numba_xgesv", sig)


@contextlib.contextmanager
def make_contiguous(context, builder, sig, args):
    """
    Ensure that all array arguments are contiguous, if necessary by
    copying them.
    A new (sig, args) tuple is yielded.
    """
    newtys = []
    newargs = []
    copies = []
    for ty, val in zip(sig.args, args):
        if not isinstance(ty, types.Array) or ty.layout in 'CF':
            newty, newval = ty, val
        else:
            newty = ty.copy(layout='C')
            copysig = signature(newty, ty)
            newval = array_copy(context, builder, copysig, (val,))
            copies.append((newty, newval))
        newtys.append(newty)
        newargs.append(newval)
    yield signature(sig.return_type, *newtys), tuple(newargs)
    for ty, val in copies:
        context.nrt.decref(builder, ty, val)


def check_c_int(context, builder, n):
    """
    Check whether *n* fits in a C `int`.
    """
    _maxint = 2**31 - 1

    def impl(n):
        if n > _maxint:
            raise OverflowError("array size too large to fit in C int")

    context.compile_internal(builder, impl,
                             signature(types.none, types.intp), (n,))


def check_blas_return(context, builder, res):
    """
    Check the integer error return from one of the BLAS wrappers in
    _helperlib.c.
    """
    with builder.if_then(cgutils.is_not_null(builder, res), likely=False):
        # Those errors shouldn't happen, it's easier to just abort the process
        pyapi = context.get_python_api(builder)
        pyapi.gil_ensure()
        pyapi.fatal_error("BLAS wrapper returned with an error")


def check_lapack_return(context, builder, res):
    """
    Check the integer error return from one of the LAPACK wrappers in
    _helperlib.c.
    """
    with builder.if_then(cgutils.is_not_null(builder, res), likely=False):
        # Those errors shouldn't happen, it's easier to just abort the process
        pyapi = context.get_python_api(builder)
        pyapi.gil_ensure()
        pyapi.fatal_error("LAPACK wrapper returned with an error")


def call_xxdot(context, builder, conjugate, dtype,
               n, a_data, b_data, out_data):
    """
    Call the BLAS vector * vector product function for the given arguments.
    """
    fnty = ir.FunctionType(ir.IntType(32),
                           [ll_char, ll_char, intp_t,    # kind, conjugate, n
                            ll_void_p, ll_void_p, ll_void_p,  # a, b, out
                            ])
    fn = cgutils.get_or_insert_function(builder.module, fnty, "numba_xxdot")

    kind = get_blas_kind(dtype)
    kind_val = ir.Constant(ll_char, ord(kind))
    conjugate = ir.Constant(ll_char, int(conjugate))

    res = builder.call(fn, (kind_val, conjugate, n,
                            builder.bitcast(a_data, ll_void_p),
                            builder.bitcast(b_data, ll_void_p),
                            builder.bitcast(out_data, ll_void_p)))
    check_blas_return(context, builder, res)


def call_xxgemv(context, builder, do_trans,
                m_type, m_shapes, m_data, v_data, out_data):
    """
    Call the BLAS matrix * vector product function for the given arguments.
    """
    fnty = ir.FunctionType(ir.IntType(32),
                           [ll_char, ll_char,                 # kind, trans
                            intp_t, intp_t,                   # m, n
                            ll_void_p, ll_void_p, intp_t,     # alpha, a, lda
                            ll_void_p, ll_void_p, ll_void_p,  # x, beta, y
                            ])
    fn = cgutils.get_or_insert_function(builder.module, fnty, "numba_xxgemv")

    dtype = m_type.dtype
    alpha = make_constant_slot(context, builder, dtype, 1.0)
    beta = make_constant_slot(context, builder, dtype, 0.0)

    if m_type.layout == 'F':
        m, n = m_shapes
        lda = m_shapes[0]
    else:
        n, m = m_shapes
        lda = m_shapes[1]

    kind = get_blas_kind(dtype)
    kind_val = ir.Constant(ll_char, ord(kind))
    trans = ir.Constant(ll_char, ord('t') if do_trans else ord('n'))

    res = builder.call(fn, (kind_val, trans, m, n,
                            builder.bitcast(alpha, ll_void_p),
                            builder.bitcast(m_data, ll_void_p), lda,
                            builder.bitcast(v_data, ll_void_p),
                            builder.bitcast(beta, ll_void_p),
                            builder.bitcast(out_data, ll_void_p)))
    check_blas_return(context, builder, res)


def call_xxgemm(context, builder,
                x_type, x_shapes, x_data,
                y_type, y_shapes, y_data,
                out_type, out_shapes, out_data):
    """
    Call the BLAS matrix * matrix product function for the given arguments.
    """
    fnty = ir.FunctionType(ir.IntType(32),
                           [ll_char,                       # kind
                            ll_char, ll_char,              # transa, transb
                            intp_t, intp_t, intp_t,        # m, n, k
                            ll_void_p, ll_void_p, intp_t,  # alpha, a, lda
                            ll_void_p, intp_t, ll_void_p,  # b, ldb, beta
                            ll_void_p, intp_t,             # c, ldc
                            ])
    fn = cgutils.get_or_insert_function(builder.module, fnty, "numba_xxgemm")

    m, k = x_shapes
    _k, n = y_shapes
    dtype = x_type.dtype
    alpha = make_constant_slot(context, builder, dtype, 1.0)
    beta = make_constant_slot(context, builder, dtype, 0.0)

    trans = ir.Constant(ll_char, ord('t'))
    notrans = ir.Constant(ll_char, ord('n'))

    def get_array_param(ty, shapes, data):
        return (
            # Transpose if layout different from result's
            notrans if ty.layout == out_type.layout else trans,
            # Size of the inner dimension in physical array order
            shapes[1] if ty.layout == 'C' else shapes[0],
            # The data pointer, unit-less
            builder.bitcast(data, ll_void_p),
        )

    transa, lda, data_a = get_array_param(y_type, y_shapes, y_data)
    transb, ldb, data_b = get_array_param(x_type, x_shapes, x_data)
    _, ldc, data_c = get_array_param(out_type, out_shapes, out_data)

    kind = get_blas_kind(dtype)
    kind_val = ir.Constant(ll_char, ord(kind))

    res = builder.call(fn, (kind_val, transa, transb, n, m, k,
                            builder.bitcast(alpha, ll_void_p), data_a, lda,
                            data_b, ldb, builder.bitcast(beta, ll_void_p),
                            data_c, ldc))
    check_blas_return(context, builder, res)


def dot_2_mm(context, builder, sig, args):
    """
    np.dot(matrix, matrix)
    """
    def dot_impl(a, b):
        m, k = a.shape
        _k, n = b.shape
        if k == 0:
            return np.zeros((m, n), a.dtype)
        out = np.empty((m, n), a.dtype)
        return np.dot(a, b, out)

    res = context.compile_internal(builder, dot_impl, sig, args)
    return impl_ret_new_ref(context, builder, sig.return_type, res)


def dot_2_vm(context, builder, sig, args):
    """
    np.dot(vector, matrix)
    """
    def dot_impl(a, b):
        m, = a.shape
        _m, n = b.shape
        if m == 0:
            return np.zeros((n, ), a.dtype)
        out = np.empty((n, ), a.dtype)
        return np.dot(a, b, out)

    res = context.compile_internal(builder, dot_impl, sig, args)
    return impl_ret_new_ref(context, builder, sig.return_type, res)


def dot_2_mv(context, builder, sig, args):
    """
    np.dot(matrix, vector)
    """
    def dot_impl(a, b):
        m, n = a.shape
        _n, = b.shape
        if n == 0:
            return np.zeros((m, ), a.dtype)
        out = np.empty((m, ), a.dtype)
        return np.dot(a, b, out)

    res = context.compile_internal(builder, dot_impl, sig, args)
    return impl_ret_new_ref(context, builder, sig.return_type, res)


def dot_2_vv(context, builder, sig, args, conjugate=False):
    """
    np.dot(vector, vector)
    np.vdot(vector, vector)
    """
    aty, bty = sig.args
    dtype = sig.return_type
    a = make_array(aty)(context, builder, args[0])
    b = make_array(bty)(context, builder, args[1])
    n, = cgutils.unpack_tuple(builder, a.shape)

    def check_args(a, b):
        m, = a.shape
        n, = b.shape
        if m != n:
            raise ValueError("incompatible array sizes for np.dot(a, b) "
                             "(vector * vector)")

    context.compile_internal(builder, check_args,
                             signature(types.none, *sig.args), args)
    check_c_int(context, builder, n)

    out = cgutils.alloca_once(builder, context.get_value_type(dtype))
    call_xxdot(context, builder, conjugate, dtype, n, a.data, b.data, out)
    return builder.load(out)


@overload(np.dot)
def dot_2(left, right):
    """
    np.dot(a, b)
    """
    return dot_2_impl('np.dot()', left, right)


@overload(operator.matmul)
def matmul_2(left, right):
    """
    a @ b
    """
    return dot_2_impl("'@'", left, right)


def dot_2_impl(name, left, right):
    if isinstance(left, types.Array) and isinstance(right, types.Array):
        @intrinsic
        def _impl(typingcontext, left, right):
            ndims = (left.ndim, right.ndim)

            def _dot2_codegen(context, builder, sig, args):
                ensure_blas()

                with make_contiguous(context, builder, sig, args) as (sig, args):
                    if ndims == (2, 2):
                        return dot_2_mm(context, builder, sig, args)
                    elif ndims == (2, 1):
                        return dot_2_mv(context, builder, sig, args)
                    elif ndims == (1, 2):
                        return dot_2_vm(context, builder, sig, args)
                    elif ndims == (1, 1):
                        return dot_2_vv(context, builder, sig, args)
                    else:
                        raise AssertionError('unreachable')

            if left.dtype != right.dtype:
                raise TypingError(
                    "%s arguments must all have the same dtype" % name)

            if ndims == (2, 2):
                return_type = types.Array(left.dtype, 2, 'C')
            elif ndims == (2, 1) or ndims == (1, 2):
                return_type = types.Array(left.dtype, 1, 'C')
            elif ndims == (1, 1):
                return_type = left.dtype
            else:
                raise TypingError(("%s: inputs must have compatible "
                                   "dimensions") % name)
            return signature(return_type, left, right), _dot2_codegen

        if left.layout not in 'CF' or right.layout not in 'CF':
            warnings.warn(
                "%s is faster on contiguous arrays, called on %s" % (
                    name, (left, right),), NumbaPerformanceWarning)

        return lambda left, right: _impl(left, right)


@overload(np.vdot)
def vdot(left, right):
    """
    np.vdot(a, b)
    """
    if isinstance(left, types.Array) and isinstance(right, types.Array):
        @intrinsic
        def _impl(typingcontext, left, right):
            def codegen(context, builder, sig, args):
                ensure_blas()

                with make_contiguous(context, builder, sig, args) as\
                        (sig, args):
                    return dot_2_vv(context, builder, sig, args, conjugate=True)

            if left.ndim != 1 or right.ndim != 1:
                raise TypingError("np.vdot() only supported on 1-D arrays")

            if left.dtype != right.dtype:
                raise TypingError(
                    "np.vdot() arguments must all have the same dtype")
            return signature(left.dtype, left, right), codegen

        if left.layout not in 'CF' or right.layout not in 'CF':
            warnings.warn(
                "np.vdot() is faster on contiguous arrays, called on %s"
                % ((left, right),), NumbaPerformanceWarning)

        return lambda left, right: _impl(left, right)


def dot_3_vm_check_args(a, b, out):
    m, = a.shape
    _m, n = b.shape
    if m != _m:
        raise ValueError("incompatible array sizes for "
                         "np.dot(a, b) (vector * matrix)")
    if out.shape != (n,):
        raise ValueError("incompatible output array size for "
                         "np.dot(a, b, out) (vector * matrix)")


def dot_3_mv_check_args(a, b, out):
    m, _n = a.shape
    n, = b.shape
    if n != _n:
        raise ValueError("incompatible array sizes for np.dot(a, b) "
                         "(matrix * vector)")
    if out.shape != (m,):
        raise ValueError("incompatible output array size for "
                         "np.dot(a, b, out) (matrix * vector)")


def dot_3_vm(context, builder, sig, args):
    """
    np.dot(vector, matrix, out)
    np.dot(matrix, vector, out)
    """
    xty, yty, outty = sig.args
    assert outty == sig.return_type
    dtype = xty.dtype

    x = make_array(xty)(context, builder, args[0])
    y = make_array(yty)(context, builder, args[1])
    out = make_array(outty)(context, builder, args[2])
    x_shapes = cgutils.unpack_tuple(builder, x.shape)
    y_shapes = cgutils.unpack_tuple(builder, y.shape)
    out_shapes = cgutils.unpack_tuple(builder, out.shape)
    if xty.ndim < yty.ndim:
        # Vector * matrix
        # Asked for x * y, we will compute y.T * x
        mty = yty
        m_shapes = y_shapes
        v_shape = x_shapes[0]
        lda = m_shapes[1]
        do_trans = yty.layout == 'F'
        m_data, v_data = y.data, x.data
        check_args = dot_3_vm_check_args
    else:
        # Matrix * vector
        # We will compute x * y
        mty = xty
        m_shapes = x_shapes
        v_shape = y_shapes[0]
        lda = m_shapes[0]
        do_trans = xty.layout == 'C'
        m_data, v_data = x.data, y.data
        check_args = dot_3_mv_check_args

    context.compile_internal(builder, check_args,
                             signature(types.none, *sig.args), args)
    for val in m_shapes:
        check_c_int(context, builder, val)

    zero = context.get_constant(types.intp, 0)
    both_empty = builder.icmp_signed('==', v_shape, zero)
    matrix_empty = builder.icmp_signed('==', lda, zero)
    is_empty = builder.or_(both_empty, matrix_empty)
    with builder.if_else(is_empty, likely=False) as (empty, nonempty):
        with empty:
            cgutils.memset(builder, out.data,
                           builder.mul(out.itemsize, out.nitems), 0)
        with nonempty:
            call_xxgemv(context, builder, do_trans, mty, m_shapes, m_data,
                        v_data, out.data)

    return impl_ret_borrowed(context, builder, sig.return_type,
                             out._getvalue())


def dot_3_mm(context, builder, sig, args):
    """
    np.dot(matrix, matrix, out)
    """
    xty, yty, outty = sig.args
    assert outty == sig.return_type
    dtype = xty.dtype

    x = make_array(xty)(context, builder, args[0])
    y = make_array(yty)(context, builder, args[1])
    out = make_array(outty)(context, builder, args[2])
    x_shapes = cgutils.unpack_tuple(builder, x.shape)
    y_shapes = cgutils.unpack_tuple(builder, y.shape)
    out_shapes = cgutils.unpack_tuple(builder, out.shape)
    m, k = x_shapes
    _k, n = y_shapes

    # The only case Numpy supports
    assert outty.layout == 'C'

    def check_args(a, b, out):
        m, k = a.shape
        _k, n = b.shape
        if k != _k:
            raise ValueError("incompatible array sizes for np.dot(a, b) "
                             "(matrix * matrix)")
        if out.shape != (m, n):
            raise ValueError("incompatible output array size for "
                             "np.dot(a, b, out) (matrix * matrix)")

    context.compile_internal(builder, check_args,
                             signature(types.none, *sig.args), args)

    check_c_int(context, builder, m)
    check_c_int(context, builder, k)
    check_c_int(context, builder, n)

    x_data = x.data
    y_data = y.data
    out_data = out.data

    # If eliminated dimension is zero, set all entries to zero and return
    zero = context.get_constant(types.intp, 0)
    both_empty = builder.icmp_signed('==', k, zero)
    x_empty = builder.icmp_signed('==', m, zero)
    y_empty = builder.icmp_signed('==', n, zero)
    is_empty = builder.or_(both_empty, builder.or_(x_empty, y_empty))
    with builder.if_else(is_empty, likely=False) as (empty, nonempty):
        with empty:
            cgutils.memset(builder, out.data,
                           builder.mul(out.itemsize, out.nitems), 0)
        with nonempty:
            # Check if any of the operands is really a 1-d vector represented
            # as a (1, k) or (k, 1) 2-d array.  In those cases, it is pessimal
            # to call the generic matrix * matrix product BLAS function.
            one = context.get_constant(types.intp, 1)
            is_left_vec = builder.icmp_signed('==', m, one)
            is_right_vec = builder.icmp_signed('==', n, one)

            with builder.if_else(is_right_vec) as (r_vec, r_mat):
                with r_vec:
                    with builder.if_else(is_left_vec) as (v_v, m_v):
                        with v_v:
                            # V * V
                            call_xxdot(context, builder, False, dtype,
                                       k, x_data, y_data, out_data)
                        with m_v:
                            # M * V
                            do_trans = xty.layout == outty.layout
                            call_xxgemv(context, builder, do_trans,
                                        xty, x_shapes, x_data, y_data, out_data)
                with r_mat:
                    with builder.if_else(is_left_vec) as (v_m, m_m):
                        with v_m:
                            # V * M
                            do_trans = yty.layout != outty.layout
                            call_xxgemv(context, builder, do_trans,
                                        yty, y_shapes, y_data, x_data, out_data)
                        with m_m:
                            # M * M
                            call_xxgemm(context, builder,
                                        xty, x_shapes, x_data,
                                        yty, y_shapes, y_data,
                                        outty, out_shapes, out_data)

    return impl_ret_borrowed(context, builder, sig.return_type,
                             out._getvalue())


@overload(np.dot)
def dot_3(left, right, out):
    """
    np.dot(a, b, out)
    """
    if (isinstance(left, types.Array) and isinstance(right, types.Array) and
            isinstance(out, types.Array)):
        @intrinsic
        def _impl(typingcontext, left, right, out):
            def codegen(context, builder, sig, args):
                ensure_blas()

                with make_contiguous(context, builder, sig, args) as (sig,
                                                                      args):
                    ndims = set(x.ndim for x in sig.args[:2])
                    if ndims == {2}:
                        return dot_3_mm(context, builder, sig, args)
                    elif ndims == {1, 2}:
                        return dot_3_vm(context, builder, sig, args)
                    else:
                        raise AssertionError('unreachable')
            if left.dtype != right.dtype or left.dtype != out.dtype:
                raise TypingError(
                    "np.dot() arguments must all have the same dtype")

            return signature(out, left, right, out), codegen

        if left.layout not in 'CF' or right.layout not in 'CF' or out.layout\
            not in 'CF':
            warnings.warn(
                "np.vdot() is faster on contiguous arrays, called on %s"
                % ((left, right),), NumbaPerformanceWarning)

        return lambda left, right, out: _impl(left, right, out)


fatal_error_func = types.ExternalFunction("numba_fatal_error", types.intc())


@register_jitable
def _check_finite_matrix(a):
    for v in np.nditer(a):
        if not np.isfinite(v.item()):
            raise np.linalg.LinAlgError(
                "Array must not contain infs or NaNs.")


def _check_linalg_matrix(a, func_name, la_prefix=True):
    # la_prefix is present as some functions, e.g. np.trace()
    # are documented under "linear algebra" but aren't in the
    # module
    prefix = "np.linalg" if la_prefix else "np"
    interp = (prefix, func_name)
    # Unpack optional type
    if isinstance(a, types.Optional):
        a = a.type
    if not isinstance(a, types.Array):
        msg = "%s.%s() only supported for array types" % interp
        raise TypingError(msg, highlighting=False)
    if not a.ndim == 2:
        msg = "%s.%s() only supported on 2-D arrays." % interp
        raise TypingError(msg, highlighting=False)
    if not isinstance(a.dtype, (types.Float, types.Complex)):
        msg = "%s.%s() only supported on "\
            "float and complex arrays." % interp
        raise TypingError(msg, highlighting=False)


def _check_homogeneous_types(func_name, *types):
    t0 = types[0].dtype
    for t in types[1:]:
        if t.dtype != t0:
            msg = "np.linalg.%s() only supports inputs that have homogeneous dtypes." % func_name
            raise TypingError(msg, highlighting=False)


def _copy_to_fortran_order():
    pass


@overload(_copy_to_fortran_order)
def ol_copy_to_fortran_order(a):
    # This function copies the array 'a' into a new array with fortran order.
    # This exists because the copy routines don't take order flags yet.
    F_layout = a.layout == 'F'
    A_layout = a.layout == 'A'
    def impl(a):
        if F_layout:
            # it's F ordered at compile time, just copy
            acpy = np.copy(a)
        elif A_layout:
            # decide based on runtime value
            flag_f = a.flags.f_contiguous
            if flag_f:
                # it's already F ordered, so copy but in a round about way to
                # ensure that the copy is also F ordered
                acpy = np.copy(a.T).T
            else:
                # it's something else ordered, so let asfortranarray deal with
                # copying and making it fortran ordered
                acpy = np.asfortranarray(a)
        else:
            # it's C ordered at compile time, asfortranarray it.
            acpy = np.asfortranarray(a)
        return acpy
    return impl


@register_jitable
def _inv_err_handler(r):
    if r != 0:
        if r < 0:
            fatal_error_func()
            assert 0   # unreachable
        if r > 0:
            raise np.linalg.LinAlgError(
                "Matrix is singular to machine precision.")

@register_jitable
def _dummy_liveness_func(a):
    """pass a list of variables to be preserved through dead code elimination"""
    return a[0]


@overload(np.linalg.inv)
def inv_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "inv")

    numba_xxgetrf = _LAPACK().numba_xxgetrf(a.dtype)

    numba_xxgetri = _LAPACK().numba_ez_xxgetri(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "inv"))

    def inv_impl(a):
        n = a.shape[-1]
        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        if n == 0:
            return acpy

        ipiv = np.empty(n, dtype=F_INT_nptype)

        r = numba_xxgetrf(kind, n, n, acpy.ctypes, n, ipiv.ctypes)
        _inv_err_handler(r)

        r = numba_xxgetri(kind, n, acpy.ctypes, n, ipiv.ctypes)
        _inv_err_handler(r)

        # help liveness analysis
        _dummy_liveness_func([acpy.size, ipiv.size])
        return acpy

    return inv_impl


@register_jitable
def _handle_err_maybe_convergence_problem(r):
    if r != 0:
        if r < 0:
            fatal_error_func()
            assert 0   # unreachable
        if r > 0:
            raise ValueError("Internal algorithm failed to converge.")


def _check_linalg_1_or_2d_matrix(a, func_name, la_prefix=True):
    # la_prefix is present as some functions, e.g. np.trace()
    # are documented under "linear algebra" but aren't in the
    # module
    prefix = "np.linalg" if la_prefix else "np"
    interp = (prefix, func_name)
    # checks that a matrix is 1 or 2D
    if not isinstance(a, types.Array):
        raise TypingError("%s.%s() only supported for array types "
                          % interp)
    if not a.ndim <= 2:
        raise TypingError("%s.%s() only supported on 1 and 2-D arrays "
                          % interp)
    if not isinstance(a.dtype, (types.Float, types.Complex)):
        raise TypingError("%s.%s() only supported on "
                          "float and complex arrays." % interp)


@overload(np.linalg.cholesky)
def cho_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "cholesky")

    numba_xxpotrf = _LAPACK().numba_xxpotrf(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "cholesky"))
    UP = ord('U')
    LO = ord('L')

    def cho_impl(a):
        n = a.shape[-1]
        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        # The output is allocated in C order
        out = a.copy()

        if n == 0:
            return out

        # Pass UP since xxpotrf() operates in F order
        # The semantics ensure this works fine
        # (out is really its Hermitian in F order, but UP instructs
        #  xxpotrf to compute the Hermitian of the upper triangle
        #  => they cancel each other)
        r = numba_xxpotrf(kind, UP, n, out.ctypes, n)
        if r != 0:
            if r < 0:
                fatal_error_func()
                assert 0   # unreachable
            if r > 0:
                raise np.linalg.LinAlgError(
                    "Matrix is not positive definite.")
        # Zero out upper triangle, in F order
        for col in range(n):
            out[:col, col] = 0
        return out

    return cho_impl

@overload(np.linalg.eig)
def eig_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "eig")

    numba_ez_rgeev = _LAPACK().numba_ez_rgeev(a.dtype)
    numba_ez_cgeev = _LAPACK().numba_ez_cgeev(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "eig"))

    JOBVL = ord('N')
    JOBVR = ord('V')

    def real_eig_impl(a):
        """
        eig() implementation for real arrays.
        """
        n = a.shape[-1]
        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        ldvl = 1
        ldvr = n
        wr = np.empty(n, dtype=a.dtype)
        wi = np.empty(n, dtype=a.dtype)
        vl = np.empty((n, ldvl), dtype=a.dtype)
        vr = np.empty((n, ldvr), dtype=a.dtype)

        if n == 0:
            return (wr, vr.T)

        r = numba_ez_rgeev(kind,
                            JOBVL,
                            JOBVR,
                            n,
                            acpy.ctypes,
                            n,
                            wr.ctypes,
                            wi.ctypes,
                            vl.ctypes,
                            ldvl,
                            vr.ctypes,
                            ldvr)
        _handle_err_maybe_convergence_problem(r)

        # By design numba does not support dynamic return types, however,
        # Numpy does. Numpy uses this ability in the case of returning
        # eigenvalues/vectors of a real matrix. The return type of
        # np.linalg.eig(), when operating on a matrix in real space
        # depends on the values present in the matrix itself (recalling
        # that eigenvalues are the roots of the characteristic polynomial
        # of the system matrix, which will by construction depend on the
        # values present in the system matrix). As numba cannot handle
        # the case of a runtime decision based domain change relative to
        # the input type, if it is required numba raises as below.
        if np.any(wi):
            raise ValueError(
                "eig() argument must not cause a domain change.")

        # put these in to help with liveness analysis,
        # `.ctypes` doesn't keep the vars alive
        _dummy_liveness_func([acpy.size, vl.size, vr.size, wr.size, wi.size])
        return (wr, vr.T)

    def cmplx_eig_impl(a):
        """
        eig() implementation for complex arrays.
        """
        n = a.shape[-1]
        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        ldvl = 1
        ldvr = n
        w = np.empty(n, dtype=a.dtype)
        vl = np.empty((n, ldvl), dtype=a.dtype)
        vr = np.empty((n, ldvr), dtype=a.dtype)

        if n == 0:
            return (w, vr.T)

        r = numba_ez_cgeev(kind,
                            JOBVL,
                            JOBVR,
                            n,
                            acpy.ctypes,
                            n,
                            w.ctypes,
                            vl.ctypes,
                            ldvl,
                            vr.ctypes,
                            ldvr)
        _handle_err_maybe_convergence_problem(r)

        # put these in to help with liveness analysis,
        # `.ctypes` doesn't keep the vars alive
        _dummy_liveness_func([acpy.size, vl.size, vr.size, w.size])
        return (w, vr.T)

    if isinstance(a.dtype, types.scalars.Complex):
        return cmplx_eig_impl
    else:
        return real_eig_impl

@overload(np.linalg.eigvals)
def eigvals_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "eigvals")

    numba_ez_rgeev = _LAPACK().numba_ez_rgeev(a.dtype)
    numba_ez_cgeev = _LAPACK().numba_ez_cgeev(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "eigvals"))

    JOBVL = ord('N')
    JOBVR = ord('N')

    def real_eigvals_impl(a):
        """
        eigvals() implementation for real arrays.
        """
        n = a.shape[-1]
        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        ldvl = 1
        ldvr = 1
        wr = np.empty(n, dtype=a.dtype)

        if n == 0:
            return wr

        wi = np.empty(n, dtype=a.dtype)

        # not referenced but need setting for MKL null check
        vl = np.empty((1), dtype=a.dtype)
        vr = np.empty((1), dtype=a.dtype)

        r = numba_ez_rgeev(kind,
                            JOBVL,
                            JOBVR,
                            n,
                            acpy.ctypes,
                            n,
                            wr.ctypes,
                            wi.ctypes,
                            vl.ctypes,
                            ldvl,
                            vr.ctypes,
                            ldvr)
        _handle_err_maybe_convergence_problem(r)

        # By design numba does not support dynamic return types, however,
        # Numpy does. Numpy uses this ability in the case of returning
        # eigenvalues/vectors of a real matrix. The return type of
        # np.linalg.eigvals(), when operating on a matrix in real space
        # depends on the values present in the matrix itself (recalling
        # that eigenvalues are the roots of the characteristic polynomial
        # of the system matrix, which will by construction depend on the
        # values present in the system matrix). As numba cannot handle
        # the case of a runtime decision based domain change relative to
        # the input type, if it is required numba raises as below.
        if np.any(wi):
            raise ValueError(
                "eigvals() argument must not cause a domain change.")

        # put these in to help with liveness analysis,
        # `.ctypes` doesn't keep the vars alive
        _dummy_liveness_func([acpy.size, vl.size, vr.size, wr.size, wi.size])
        return wr

    def cmplx_eigvals_impl(a):
        """
        eigvals() implementation for complex arrays.
        """
        n = a.shape[-1]
        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        ldvl = 1
        ldvr = 1
        w = np.empty(n, dtype=a.dtype)

        if n == 0:
            return w

        vl = np.empty((1), dtype=a.dtype)
        vr = np.empty((1), dtype=a.dtype)

        r = numba_ez_cgeev(kind,
                            JOBVL,
                            JOBVR,
                            n,
                            acpy.ctypes,
                            n,
                            w.ctypes,
                            vl.ctypes,
                            ldvl,
                            vr.ctypes,
                            ldvr)
        _handle_err_maybe_convergence_problem(r)

        # put these in to help with liveness analysis,
        # `.ctypes` doesn't keep the vars alive
        _dummy_liveness_func([acpy.size, vl.size, vr.size, w.size])
        return w

    if isinstance(a.dtype, types.scalars.Complex):
        return cmplx_eigvals_impl
    else:
        return real_eigvals_impl

@overload(np.linalg.eigh)
def eigh_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "eigh")

    # convert typing floats to numpy floats for use in the impl
    w_type = getattr(a.dtype, "underlying_float", a.dtype)
    w_dtype = np_support.as_dtype(w_type)

    numba_ez_xxxevd = _LAPACK().numba_ez_xxxevd(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "eigh"))

    JOBZ = ord('V')
    UPLO = ord('L')

    def eigh_impl(a):
        n = a.shape[-1]

        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        w = np.empty(n, dtype=w_dtype)

        if n == 0:
            return (w, acpy)

        r = numba_ez_xxxevd(kind,  # kind
                            JOBZ,  # jobz
                            UPLO,  # uplo
                            n,  # n
                            acpy.ctypes,  # a
                            n,  # lda
                            w.ctypes  # w
                            )
        _handle_err_maybe_convergence_problem(r)

        # help liveness analysis
        _dummy_liveness_func([acpy.size, w.size])
        return (w, acpy)

    return eigh_impl

@overload(np.linalg.eigvalsh)
def eigvalsh_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "eigvalsh")

    # convert typing floats to numpy floats for use in the impl
    w_type = getattr(a.dtype, "underlying_float", a.dtype)
    w_dtype = np_support.as_dtype(w_type)

    numba_ez_xxxevd = _LAPACK().numba_ez_xxxevd(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "eigvalsh"))

    JOBZ = ord('N')
    UPLO = ord('L')

    def eigvalsh_impl(a):
        n = a.shape[-1]

        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        w = np.empty(n, dtype=w_dtype)

        if n == 0:
            return w

        r = numba_ez_xxxevd(kind,  # kind
                            JOBZ,  # jobz
                            UPLO,  # uplo
                            n,  # n
                            acpy.ctypes,  # a
                            n,  # lda
                            w.ctypes  # w
                            )
        _handle_err_maybe_convergence_problem(r)

        # help liveness analysis
        _dummy_liveness_func([acpy.size, w.size])
        return w

    return eigvalsh_impl

@overload(np.linalg.svd)
def svd_impl(a, full_matrices=1):
    ensure_lapack()

    _check_linalg_matrix(a, "svd")

    # convert typing floats to numpy floats for use in the impl
    s_type = getattr(a.dtype, "underlying_float", a.dtype)
    s_dtype = np_support.as_dtype(s_type)

    numba_ez_gesdd = _LAPACK().numba_ez_gesdd(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "svd"))

    JOBZ_A = ord('A')
    JOBZ_S = ord('S')

    def svd_impl(a, full_matrices=1):
        n = a.shape[-1]
        m = a.shape[-2]

        if n == 0 or m == 0:
            raise np.linalg.LinAlgError("Arrays cannot be empty")

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        ldu = m
        minmn = min(m, n)

        if full_matrices:
            JOBZ = JOBZ_A
            ucol = m
            ldvt = n
        else:
            JOBZ = JOBZ_S
            ucol = minmn
            ldvt = minmn

        u = np.empty((ucol, ldu), dtype=a.dtype)
        s = np.empty(minmn, dtype=s_dtype)
        vt = np.empty((n, ldvt), dtype=a.dtype)

        r = numba_ez_gesdd(
            kind,  # kind
            JOBZ,  # jobz
            m,  # m
            n,  # n
            acpy.ctypes,  # a
            m,  # lda
            s.ctypes,  # s
            u.ctypes,  # u
            ldu,  # ldu
            vt.ctypes,  # vt
            ldvt          # ldvt
        )
        _handle_err_maybe_convergence_problem(r)

        # help liveness analysis
        _dummy_liveness_func([acpy.size, vt.size, u.size, s.size])
        return (u.T, s, vt.T)

    return svd_impl


@overload(np.linalg.qr)
def qr_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "qr")

    # Need two functions, the first computes R, storing it in the upper
    # triangle of A with the below diagonal part of A containing elementary
    # reflectors needed to construct Q. The second turns the below diagonal
    # entries of A into Q, storing Q in A (creates orthonormal columns from
    # the elementary reflectors).

    numba_ez_geqrf = _LAPACK().numba_ez_geqrf(a.dtype)
    numba_ez_xxgqr = _LAPACK().numba_ez_xxgqr(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "qr"))

    def qr_impl(a):
        n = a.shape[-1]
        m = a.shape[-2]

        if n == 0 or m == 0:
            raise np.linalg.LinAlgError("Arrays cannot be empty")

        _check_finite_matrix(a)

        # copy A as it will be destroyed
        q = _copy_to_fortran_order(a)

        lda = m

        minmn = min(m, n)
        tau = np.empty((minmn), dtype=a.dtype)

        ret = numba_ez_geqrf(
            kind,  # kind
            m,  # m
            n,  # n
            q.ctypes,  # a
            m,  # lda
            tau.ctypes  # tau
        )
        if ret < 0:
            fatal_error_func()
            assert 0   # unreachable

        # pull out R, this is transposed because of Fortran
        r = np.zeros((n, minmn), dtype=a.dtype).T

        # the triangle in R
        for i in range(minmn):
            for j in range(i + 1):
                r[j, i] = q[j, i]

        # and the possible square in R
        for i in range(minmn, n):
            for j in range(minmn):
                r[j, i] = q[j, i]

        ret = numba_ez_xxgqr(
            kind,  # kind
            m,  # m
            minmn,  # n
            minmn,  # k
            q.ctypes,  # a
            m,  # lda
            tau.ctypes  # tau
        )
        _handle_err_maybe_convergence_problem(ret)

        # help liveness analysis
        _dummy_liveness_func([tau.size, q.size])
        return (q[:, :minmn], r)

    return qr_impl


# helpers and jitted specialisations required for np.linalg.lstsq
# and np.linalg.solve. These functions have "system" in their name
# as a differentiator.

def _system_copy_in_b(bcpy, b, nrhs):
    """
    Correctly copy 'b' into the 'bcpy' scratch space.
    """
    raise NotImplementedError


@overload(_system_copy_in_b)
def _system_copy_in_b_impl(bcpy, b, nrhs):
    if b.ndim == 1:
        def oneD_impl(bcpy, b, nrhs):
            bcpy[:b.shape[-1], 0] = b
        return oneD_impl
    else:
        def twoD_impl(bcpy, b, nrhs):
            bcpy[:b.shape[-2], :nrhs] = b
        return twoD_impl


def _system_compute_nrhs(b):
    """
    Compute the number of right hand sides in the system of equations
    """
    raise NotImplementedError


@overload(_system_compute_nrhs)
def _system_compute_nrhs_impl(b):
    if b.ndim == 1:
        def oneD_impl(b):
            return 1
        return oneD_impl
    else:
        def twoD_impl(b):
            return b.shape[-1]
        return twoD_impl


def _system_check_dimensionally_valid(a, b):
    """
    Check that AX=B style system input is dimensionally valid.
    """
    raise NotImplementedError


@overload(_system_check_dimensionally_valid)
def _system_check_dimensionally_valid_impl(a, b):
    ndim = b.ndim
    if ndim == 1:
        def oneD_impl(a, b):
            am = a.shape[-2]
            bm = b.shape[-1]
            if am != bm:
                raise np.linalg.LinAlgError(
                    "Incompatible array sizes, system is not dimensionally valid.")
        return oneD_impl
    else:
        def twoD_impl(a, b):
            am = a.shape[-2]
            bm = b.shape[-2]
            if am != bm:
                raise np.linalg.LinAlgError(
                    "Incompatible array sizes, system is not dimensionally valid.")
        return twoD_impl


def _system_check_non_empty(a, b):
    """
    Check that AX=B style system input is not empty.
    """
    raise NotImplementedError


@overload(_system_check_non_empty)
def _system_check_non_empty_impl(a, b):
    ndim = b.ndim
    if ndim == 1:
        def oneD_impl(a, b):
            am = a.shape[-2]
            an = a.shape[-1]
            bm = b.shape[-1]
            if am == 0 or bm == 0 or an == 0:
                raise np.linalg.LinAlgError('Arrays cannot be empty')
        return oneD_impl
    else:
        def twoD_impl(a, b):
            am = a.shape[-2]
            an = a.shape[-1]
            bm = b.shape[-2]
            bn = b.shape[-1]
            if am == 0 or bm == 0 or an == 0 or bn == 0:
                raise np.linalg.LinAlgError('Arrays cannot be empty')
        return twoD_impl


def _lstsq_residual(b, n, nrhs):
    """
    Compute the residual from the 'b' scratch space.
    """
    raise NotImplementedError


@overload(_lstsq_residual)
def _lstsq_residual_impl(b, n, nrhs):
    ndim = b.ndim
    dtype = b.dtype
    real_dtype = np_support.as_dtype(getattr(dtype, "underlying_float", dtype))

    if ndim == 1:
        if isinstance(dtype, (types.Complex)):
            def cmplx_impl(b, n, nrhs):
                res = np.empty((1,), dtype=real_dtype)
                res[0] = np.sum(np.abs(b[n:, 0])**2)
                return res
            return cmplx_impl
        else:
            def real_impl(b, n, nrhs):
                res = np.empty((1,), dtype=real_dtype)
                res[0] = np.sum(b[n:, 0]**2)
                return res
            return real_impl
    else:
        assert ndim == 2
        if isinstance(dtype, (types.Complex)):
            def cmplx_impl(b, n, nrhs):
                res = np.empty((nrhs), dtype=real_dtype)
                for k in range(nrhs):
                    res[k] = np.sum(np.abs(b[n:, k])**2)
                return res
            return cmplx_impl
        else:
            def real_impl(b, n, nrhs):
                res = np.empty((nrhs), dtype=real_dtype)
                for k in range(nrhs):
                    res[k] = np.sum(b[n:, k]**2)
                return res
            return real_impl


def _lstsq_solution(b, bcpy, n):
    """
    Extract 'x' (the lstsq solution) from the 'bcpy' scratch space.
    Note 'b' is only used to check the system input dimension...
    """
    raise NotImplementedError


@overload(_lstsq_solution)
def _lstsq_solution_impl(b, bcpy, n):
    if b.ndim == 1:
        def oneD_impl(b, bcpy, n):
            return bcpy.T.ravel()[:n]
        return oneD_impl
    else:
        def twoD_impl(b, bcpy, n):
            return bcpy[:n, :].copy()
        return twoD_impl


@overload(np.linalg.lstsq)
def lstsq_impl(a, b, rcond=-1.0):
    ensure_lapack()

    _check_linalg_matrix(a, "lstsq")

    # B can be 1D or 2D.
    _check_linalg_1_or_2d_matrix(b, "lstsq")

    _check_homogeneous_types("lstsq", a, b)

    np_dt = np_support.as_dtype(a.dtype)
    nb_dt = a.dtype

    # convert typing floats to np floats for use in the impl
    r_type = getattr(nb_dt, "underlying_float", nb_dt)
    real_dtype = np_support.as_dtype(r_type)

    # lapack solver
    numba_ez_gelsd = _LAPACK().numba_ez_gelsd(a.dtype)

    kind = ord(get_blas_kind(nb_dt, "lstsq"))

    # The following functions select specialisations based on
    # information around 'b', a lot of this effort is required
    # as 'b' can be either 1D or 2D, and then there are
    # some optimisations available depending on real or complex
    # space.

    def lstsq_impl(a, b, rcond=-1.0):
        n = a.shape[-1]
        m = a.shape[-2]
        nrhs = _system_compute_nrhs(b)

        # check the systems have no inf or NaN
        _check_finite_matrix(a)
        _check_finite_matrix(b)

        # check the system is not empty
        _system_check_non_empty(a, b)

        # check the systems are dimensionally valid
        _system_check_dimensionally_valid(a, b)

        minmn = min(m, n)
        maxmn = max(m, n)

        # a is destroyed on exit, copy it
        acpy = _copy_to_fortran_order(a)

        # b is overwritten on exit with the solution, copy allocate
        bcpy = np.empty((nrhs, maxmn), dtype=np_dt).T
        # specialised copy in due to b being 1 or 2D
        _system_copy_in_b(bcpy, b, nrhs)

        # Allocate returns
        s = np.empty(minmn, dtype=real_dtype)
        rank_ptr = np.empty(1, dtype=np.int32)

        r = numba_ez_gelsd(
            kind,  # kind
            m,  # m
            n,  # n
            nrhs,  # nrhs
            acpy.ctypes,  # a
            m,  # lda
            bcpy.ctypes,  # a
            maxmn,  # ldb
            s.ctypes,  # s
            rcond,  # rcond
            rank_ptr.ctypes  # rank
        )
        _handle_err_maybe_convergence_problem(r)

        # set rank to that which was computed
        rank = rank_ptr[0]

        # compute residuals
        if rank < n or m <= n:
            res = np.empty((0), dtype=real_dtype)
        else:
            # this requires additional dispatch as there's a faster
            # impl if the result is in the real domain (no abs() required)
            res = _lstsq_residual(bcpy, n, nrhs)

        # extract 'x', the solution
        x = _lstsq_solution(b, bcpy, n)

        # help liveness analysis
        _dummy_liveness_func([acpy.size, bcpy.size, s.size, rank_ptr.size])
        return (x, res, rank, s[:minmn])

    return lstsq_impl


def _solve_compute_return(b, bcpy):
    """
    Extract 'x' (the solution) from the 'bcpy' scratch space.
    Note 'b' is only used to check the system input dimension...
    """
    raise NotImplementedError


@overload(_solve_compute_return)
def _solve_compute_return_impl(b, bcpy):
    if b.ndim == 1:
        def oneD_impl(b, bcpy):
            return bcpy.T.ravel()
        return oneD_impl
    else:
        def twoD_impl(b, bcpy):
            return bcpy
        return twoD_impl


@overload(np.linalg.solve)
def solve_impl(a, b):
    ensure_lapack()

    _check_linalg_matrix(a, "solve")
    _check_linalg_1_or_2d_matrix(b, "solve")

    _check_homogeneous_types("solve", a, b)

    np_dt = np_support.as_dtype(a.dtype)
    nb_dt = a.dtype

    # the lapack solver
    numba_xgesv = _LAPACK().numba_xgesv(a.dtype)

    kind = ord(get_blas_kind(nb_dt, "solve"))

    def solve_impl(a, b):
        n = a.shape[-1]
        nrhs = _system_compute_nrhs(b)

        # check the systems have no inf or NaN
        _check_finite_matrix(a)
        _check_finite_matrix(b)

        # check the systems are dimensionally valid
        _system_check_dimensionally_valid(a, b)

        # a is destroyed on exit, copy it
        acpy = _copy_to_fortran_order(a)

        # b is overwritten on exit with the solution, copy allocate
        bcpy = np.empty((nrhs, n), dtype=np_dt).T
        if n == 0:
            return _solve_compute_return(b, bcpy)

        # specialised copy in due to b being 1 or 2D
        _system_copy_in_b(bcpy, b, nrhs)

        # allocate pivot array (needs to be fortran int size)
        ipiv = np.empty(n, dtype=F_INT_nptype)

        r = numba_xgesv(
            kind,        # kind
            n,           # n
            nrhs,        # nhrs
            acpy.ctypes,  # a
            n,           # lda
            ipiv.ctypes,  # ipiv
            bcpy.ctypes,  # b
            n            # ldb
        )
        _inv_err_handler(r)

        # help liveness analysis
        _dummy_liveness_func([acpy.size, bcpy.size, ipiv.size])
        return _solve_compute_return(b, bcpy)

    return solve_impl


@overload(np.linalg.pinv)
def pinv_impl(a, rcond=1.e-15):
    ensure_lapack()

    _check_linalg_matrix(a, "pinv")

    # convert typing floats to numpy floats for use in the impl
    s_type = getattr(a.dtype, "underlying_float", a.dtype)
    s_dtype = np_support.as_dtype(s_type)

    numba_ez_gesdd = _LAPACK().numba_ez_gesdd(a.dtype)

    numba_xxgemm = _BLAS().numba_xxgemm(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "pinv"))
    JOB = ord('S')

    # need conjugate transposes
    TRANSA = ord('C')
    TRANSB = ord('C')

    # scalar constants
    dt = np_support.as_dtype(a.dtype)
    zero = np.array([0.], dtype=dt)
    one = np.array([1.], dtype=dt)

    def pinv_impl(a, rcond=1.e-15):

        # The idea is to build the pseudo-inverse via inverting the singular
        # value decomposition of a matrix `A`. Mathematically, this is roughly
        # A = U*S*V^H        [The SV decomposition of A]
        # A^+ = V*(S^+)*U^H  [The inverted SV decomposition of A]
        # where ^+ is pseudo inversion and ^H is Hermitian transpose.
        # As V and U are unitary, their inverses are simply their Hermitian
        # transpose. S has singular values on its diagonal and zero elsewhere,
        # it is inverted trivially by reciprocal of the diagonal values with
        # the exception that zero singular values remain as zero.
        #
        # The practical implementation can take advantage of a few things to
        # gain a few % performance increase:
        # * A is destroyed by the SVD algorithm from LAPACK so a copy is
        #   required, this memory is exactly the right size in which to return
        #   the pseudo-inverse and so can be reused for this purpose.
        # * The pseudo-inverse of S can be applied to either V or U^H, this
        #   then leaves a GEMM operation to compute the inverse via either:
        #   A^+ = (V*(S^+))*U^H
        #   or
        #   A^+ = V*((S^+)*U^H)
        #   however application of S^+ to V^H or U is more convenient as they
        #   are the result of the SVD algorithm. The application of the
        #   diagonal system is just a matrix multiplication which results in a
        #   row/column scaling (direction depending). To save effort, this
        #   "matrix multiplication" is applied to the smallest of U or V^H and
        #   only up to the point of "cut-off" (see next note) just as a direct
        #   scaling.
        # * The cut-off level for application of S^+ can be used to reduce
        #   total effort, this cut-off can come via rcond or may just naturally
        #   be present as a result of zeros in the singular values. Regardless
        #   there's no need to multiply by zeros in the application of S^+ to
        #   V^H or U as above. Further, the GEMM operation can be shrunk in
        #   effort by noting that the possible zero block generated by the
        #   presence of zeros in S^+ has no effect apart from wasting cycles as
        #   it is all fmadd()s where one operand is zero. The inner dimension
        #   of the GEMM operation can therefore be set as shrunk accordingly!

        n = a.shape[-1]
        m = a.shape[-2]

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        if m == 0 or n == 0:
            return acpy.T.ravel().reshape(a.shape).T

        minmn = min(m, n)

        u = np.empty((minmn, m), dtype=a.dtype)
        s = np.empty(minmn, dtype=s_dtype)
        vt = np.empty((n, minmn), dtype=a.dtype)

        r = numba_ez_gesdd(
            kind,         # kind
            JOB,          # job
            m,            # m
            n,            # n
            acpy.ctypes,  # a
            m,            # lda
            s.ctypes,     # s
            u.ctypes,     # u
            m,            # ldu
            vt.ctypes,    # vt
            minmn         # ldvt
        )
        _handle_err_maybe_convergence_problem(r)

        # Invert singular values under threshold. Also find the index of
        # the threshold value as this is the upper limit for the application
        # of the inverted singular values. Finding this value saves
        # multiplication by a block of zeros that would be created by the
        # application of these values to either U or V^H ahead of multiplying
        # them together. This is done by simply in BLAS parlance via
        # restricting the `k` dimension to `cut_idx` in `xgemm` whilst keeping
        # the leading dimensions correct.

        cut_at = s[0] * rcond
        cut_idx = 0
        for k in range(minmn):
            if s[k] > cut_at:
                s[k] = 1. / s[k]
                cut_idx = k
        cut_idx += 1

        # Use cut_idx so there's no scaling by 0.
        if m >= n:
            # U is largest so apply S^+ to V^H.
            for i in range(n):
                for j in range(cut_idx):
                    vt[i, j] = vt[i, j] * s[j]
        else:
            # V^H is largest so apply S^+ to U.
            for i in range(cut_idx):
                s_local = s[i]
                for j in range(minmn):
                    u[i, j] = u[i, j] * s_local

        # Do (v^H)^H*U^H (obviously one of the matrices includes the S^+
        # scaling) and write back to acpy. Note the innner dimension of cut_idx
        # taking account of the possible zero block.
        # We can store the result in acpy, given we had to create it
        # for use in the SVD, and it is now redundant and the right size
        # but wrong shape.

        r = numba_xxgemm(
            kind,
            TRANSA,       # TRANSA
            TRANSB,       # TRANSB
            n,            # M
            m,            # N
            cut_idx,      # K
            one.ctypes,   # ALPHA
            vt.ctypes,    # A
            minmn,        # LDA
            u.ctypes,     # B
            m,            # LDB
            zero.ctypes,  # BETA
            acpy.ctypes,  # C
            n             # LDC
        )

        # help liveness analysis
        #acpy.size
        #vt.size
        #u.size
        #s.size
        #one.size
        #zero.size
        _dummy_liveness_func([acpy.size, vt.size, u.size, s.size, one.size,
            zero.size])
        return acpy.T.ravel().reshape(a.shape).T

    return pinv_impl


def _get_slogdet_diag_walker(a):
    """
    Walks the diag of a LUP decomposed matrix
    uses that det(A) = prod(diag(lup(A)))
    and also that log(a)+log(b) = log(a*b)
    The return sign is adjusted based on the values found
    such that the log(value) stays in the real domain.
    """
    if isinstance(a.dtype, types.Complex):
        @register_jitable
        def cmplx_diag_walker(n, a, sgn):
            # walk diagonal
            csgn = sgn + 0.j
            acc = 0.
            for k in range(n):
                absel = np.abs(a[k, k])
                csgn = csgn * (a[k, k] / absel)
                acc = acc + np.log(absel)
            return (csgn, acc)
        return cmplx_diag_walker
    else:
        @register_jitable
        def real_diag_walker(n, a, sgn):
            # walk diagonal
            acc = 0.
            for k in range(n):
                v = a[k, k]
                if v < 0.:
                    sgn = -sgn
                    v = -v
                acc = acc + np.log(v)
            # sgn is a float dtype
            return (sgn + 0., acc)
        return real_diag_walker


@overload(np.linalg.slogdet)
def slogdet_impl(a):
    ensure_lapack()

    _check_linalg_matrix(a, "slogdet")

    numba_xxgetrf = _LAPACK().numba_xxgetrf(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "slogdet"))

    diag_walker = _get_slogdet_diag_walker(a)

    ONE = a.dtype(1)
    ZERO = getattr(a.dtype, "underlying_float", a.dtype)(0)

    def slogdet_impl(a):
        n = a.shape[-1]
        if a.shape[-2] != n:
            msg = "Last 2 dimensions of the array must be square."
            raise np.linalg.LinAlgError(msg)

        if n == 0:
            return (ONE, ZERO)

        _check_finite_matrix(a)

        acpy = _copy_to_fortran_order(a)

        ipiv = np.empty(n, dtype=F_INT_nptype)

        r = numba_xxgetrf(kind, n, n, acpy.ctypes, n, ipiv.ctypes)

        if r > 0:
            # factorisation failed, return same defaults as np
            return (0., -np.inf)
        _inv_err_handler(r)  # catch input-to-lapack problem

        # The following, prior to the call to diag_walker, is present
        # to account for the effect of possible permutations to the
        # sign of the determinant.
        # This is the same idea as in numpy:
        # File name `umath_linalg.c.src` e.g.
        # https://github.com/numpy/numpy/blob/master/numpy/linalg/umath_linalg.c.src
        # in function `@TYPE@_slogdet_single_element`.
        sgn = 1
        for k in range(n):
            sgn = sgn + (ipiv[k] != (k + 1))

        sgn = sgn & 1
        if sgn == 0:
            sgn = -1

        # help liveness analysis
        _dummy_liveness_func([ipiv.size])
        return diag_walker(n, acpy, sgn)

    return slogdet_impl


@overload(np.linalg.det)
def det_impl(a):

    ensure_lapack()

    _check_linalg_matrix(a, "det")

    def det_impl(a):
        (sgn, slogdet) = np.linalg.slogdet(a)
        return sgn * np.exp(slogdet)

    return det_impl


def _compute_singular_values(a):
    """
    Compute singular values of *a*.
    """
    raise NotImplementedError


@overload(_compute_singular_values)
def _compute_singular_values_impl(a):
    """
    Returns a function to compute singular values of `a`
    """
    numba_ez_gesdd = _LAPACK().numba_ez_gesdd(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "svd"))

    # Flag for "only compute `S`" to give to xgesdd
    JOBZ_N = ord('N')

    nb_ret_type = getattr(a.dtype, "underlying_float", a.dtype)
    np_ret_type = np_support.as_dtype(nb_ret_type)
    np_dtype = np_support.as_dtype(a.dtype)

    # These are not referenced in the computation but must be set
    # for MKL.
    u = np.empty((1, 1), dtype=np_dtype)
    vt = np.empty((1, 1), dtype=np_dtype)

    def sv_function(a):
        """
        Computes singular values.
        """
        # Don't use the np.linalg.svd impl instead
        # call LAPACK to shortcut doing the "reconstruct
        # singular vectors from reflectors" step and just
        # get back the singular values.
        n = a.shape[-1]
        m = a.shape[-2]
        if m == 0 or n == 0:
            raise np.linalg.LinAlgError('Arrays cannot be empty')
        _check_finite_matrix(a)

        ldu = m
        minmn = min(m, n)

        # need to be >=1 but aren't referenced
        ucol = 1
        ldvt = 1

        acpy = _copy_to_fortran_order(a)

        # u and vt are not referenced however need to be
        # allocated (as done above) for MKL as it
        # checks for ref is nullptr.
        s = np.empty(minmn, dtype=np_ret_type)

        r = numba_ez_gesdd(
            kind,        # kind
            JOBZ_N,      # jobz
            m,           # m
            n,           # n
            acpy.ctypes,  # a
            m,           # lda
            s.ctypes,    # s
            u.ctypes,    # u
            ldu,         # ldu
            vt.ctypes,   # vt
            ldvt         # ldvt
        )
        _handle_err_maybe_convergence_problem(r)

        # help liveness analysis
        _dummy_liveness_func([acpy.size, vt.size, u.size, s.size])
        return s

    return sv_function


def _oneD_norm_2(a):
    """
    Compute the L2-norm of 1D-array *a*.
    """
    raise NotImplementedError


@overload(_oneD_norm_2)
def _oneD_norm_2_impl(a):
    nb_ret_type = getattr(a.dtype, "underlying_float", a.dtype)
    np_ret_type = np_support.as_dtype(nb_ret_type)

    xxnrm2 = _BLAS().numba_xxnrm2(a.dtype)

    kind = ord(get_blas_kind(a.dtype, "norm"))

    def impl(a):
        # Just ignore order, calls are guarded to only come
        # from cases where order=None or order=2.
        n = len(a)
        # Call L2-norm routine from BLAS
        ret = np.empty((1,), dtype=np_ret_type)
        jmp = int(a.strides[0] / a.itemsize)
        r = xxnrm2(
            kind,      # kind
            n,         # n
            a.ctypes,  # x
            jmp,       # incx
            ret.ctypes  # result
        )
        if r < 0:
            fatal_error_func()
            assert 0   # unreachable

        # help liveness analysis
        #ret.size
        #a.size
        _dummy_liveness_func([ret.size, a.size])
        return ret[0]

    return impl


def _get_norm_impl(x, ord_flag):
    # This function is quite involved as norm supports a large
    # range of values to select different norm types via kwarg `ord`.
    # The implementation below branches on dimension of the input
    # (1D or 2D). The default for `ord` is `None` which requires
    # special handling in numba, this is dealt with first in each of
    # the dimension branches. Following this the various norms are
    # computed via code that is in most cases simply a loop version
    # of a ufunc based version as found in numpy.

    # The following is common to both 1D and 2D cases.
    # Convert typing floats to numpy floats for use in the impl.
    # The return type is always a float, numba differs from numpy in
    # that it returns an input precision specific value whereas numpy
    # always returns np.float64.
    nb_ret_type = getattr(x.dtype, "underlying_float", x.dtype)
    np_ret_type = np_support.as_dtype(nb_ret_type)

    np_dtype = np_support.as_dtype(x.dtype)

    xxnrm2 = _BLAS().numba_xxnrm2(x.dtype)

    kind = ord(get_blas_kind(x.dtype, "norm"))

    if x.ndim == 1:
        # 1D cases

        # handle "ord" being "None", must be done separately
        if ord_flag in (None, types.none):
            def oneD_impl(x, ord=None):
                return _oneD_norm_2(x)
        else:
            def oneD_impl(x, ord=None):
                n = len(x)

                # Shortcut to handle zero length arrays
                # this differs slightly to numpy in that
                # numpy raises a ValueError for kwarg ord=
                # +/-np.inf as the reduction operations like
                # max() and min() don't accept zero length
                # arrays
                if n == 0:
                    return 0.0

                # Note: on order == 2
                # This is the same as for ord=="None" but because
                # we have to handle "None" specially this condition
                # is separated
                if ord == 2:
                    return _oneD_norm_2(x)
                elif ord == np.inf:
                    # max(abs(x))
                    ret = abs(x[0])
                    for k in range(1, n):
                        val = abs(x[k])
                        if val > ret:
                            ret = val
                    return ret

                elif ord == -np.inf:
                    # min(abs(x))
                    ret = abs(x[0])
                    for k in range(1, n):
                        val = abs(x[k])
                        if val < ret:
                            ret = val
                    return ret

                elif ord == 0:
                    # sum(x != 0)
                    ret = 0.0
                    for k in range(n):
                        if x[k] != 0.:
                            ret += 1.
                    return ret

                elif ord == 1:
                    # sum(abs(x))
                    ret = 0.0
                    for k in range(n):
                        ret += abs(x[k])
                    return ret

                else:
                    # sum(abs(x)**ord)**(1./ord)
                    ret = 0.0
                    for k in range(n):
                        ret += abs(x[k])**ord
                    return ret**(1. / ord)
        return oneD_impl

    elif x.ndim == 2:
        # 2D cases

        # handle "ord" being "None"
        if ord_flag in (None, types.none):
            # Force `x` to be C-order, so that we can take a contiguous
            # 1D view.
            if x.layout == 'C':
                @register_jitable
                def array_prepare(x):
                    return x
            elif x.layout == 'F':
                @register_jitable
                def array_prepare(x):
                    # Legal since L2(x) == L2(x.T)
                    return x.T
            else:
                @register_jitable
                def array_prepare(x):
                    return x.copy()

            # Compute the Frobenius norm, this is the L2,2 induced norm of `x`
            # which is the L2-norm of x.ravel() and so can be computed via BLAS
            def twoD_impl(x, ord=None):
                n = x.size
                if n == 0:
                    # reshape() currently doesn't support zero-sized arrays
                    return 0.0
                x_c = array_prepare(x)
                return _oneD_norm_2(x_c.reshape(n))
        else:
            # max value for this dtype
            max_val = np.finfo(np_ret_type.type).max

            def twoD_impl(x, ord=None):
                n = x.shape[-1]
                m = x.shape[-2]

                # Shortcut to handle zero size arrays
                # this differs slightly to numpy in that
                # numpy raises errors for some ord values
                # and in other cases returns zero.
                if x.size == 0:
                    return 0.0

                if ord == np.inf:
                    # max of sum of abs across rows
                    # max(sum(abs(x)), axis=1)
                    global_max = 0.
                    for ii in range(m):
                        tmp = 0.
                        for jj in range(n):
                            tmp += abs(x[ii, jj])
                        if tmp > global_max:
                            global_max = tmp
                    return global_max

                elif ord == -np.inf:
                    # min of sum of abs across rows
                    # min(sum(abs(x)), axis=1)
                    global_min = max_val
                    for ii in range(m):
                        tmp = 0.
                        for jj in range(n):
                            tmp += abs(x[ii, jj])
                        if tmp < global_min:
                            global_min = tmp
                    return global_min
                elif ord == 1:
                    # max of sum of abs across cols
                    # max(sum(abs(x)), axis=0)
                    global_max = 0.
                    for ii in range(n):
                        tmp = 0.
                        for jj in range(m):
                            tmp += abs(x[jj, ii])
                        if tmp > global_max:
                            global_max = tmp
                    return global_max

                elif ord == -1:
                    # min of sum of abs across cols
                    # min(sum(abs(x)), axis=0)
                    global_min = max_val
                    for ii in range(n):
                        tmp = 0.
                        for jj in range(m):
                            tmp += abs(x[jj, ii])
                        if tmp < global_min:
                            global_min = tmp
                    return global_min

                # Results via SVD, singular values are sorted on return
                # by definition.
                elif ord == 2:
                    # max SV
                    return _compute_singular_values(x)[0]
                elif ord == -2:
                    # min SV
                    return _compute_singular_values(x)[-1]
                else:
                    # replicate numpy error
                    raise ValueError("Invalid norm order for matrices.")
        return twoD_impl
    else:
        assert 0  # unreachable


@overload(np.linalg.norm)
def norm_impl(x, ord=None):
    ensure_lapack()

    _check_linalg_1_or_2d_matrix(x, "norm")

    return _get_norm_impl(x, ord)


@overload(np.linalg.cond)
def cond_impl(x, p=None):
    ensure_lapack()

    _check_linalg_matrix(x, "cond")

    def impl(x, p=None):
        # This is extracted for performance, numpy does approximately:
        # `condition = norm(x) * norm(inv(x))`
        # in the cases of `p == 2` or `p ==-2` singular values are used
        # for computing norms. This costs numpy an svd of `x` then an
        # inversion of `x` and another svd of `x`.
        # Below is a different approach, which also gives a more
        # accurate answer as there is no inversion involved.
        # Recall that the singular values of an inverted matrix are the
        # reciprocal of singular values of the original matrix.
        # Therefore calling `svd(x)` once yields all the information
        # needed about both `x` and `inv(x)` without the cost or
        # potential loss of accuracy incurred through inversion.
        # For the case of `p == 2`, the result is just the ratio of
        # `largest singular value/smallest singular value`, and for the
        # case of `p==-2` the result is simply the
        # `smallest singular value/largest singular value`.
        # As a result of this, numba accepts non-square matrices as
        # input when p==+/-2 as well as when p==None.
        if p == 2 or p == -2 or p is None:
            s = _compute_singular_values(x)
            if p == 2 or p is None:
                r = np.divide(s[0], s[-1])
            else:
                r = np.divide(s[-1], s[0])
        else:  # cases np.inf, -np.inf, 1, -1
            norm_x = np.linalg.norm(x, p)
            norm_inv_x = np.linalg.norm(np.linalg.inv(x), p)
            r = norm_x * norm_inv_x
        # NumPy uses a NaN mask, if the input has a NaN, it will return NaN,
        # Numba calls ban NaN through the use of _check_finite_matrix but this
        # catches cases where NaN occurs through floating point use
        if np.isnan(r):
            return np.inf
        else:
            return r
    return impl


@register_jitable
def _get_rank_from_singular_values(sv, t):
    """
    Gets rank from singular values with cut-off at a given tolerance
    """
    rank = 0
    for k in range(len(sv)):
        if sv[k] > t:
            rank = rank + 1
        else:  # sv is ordered big->small so break on condition not met
            break
    return rank


@overload(np.linalg.matrix_rank)
def matrix_rank_impl(A, tol=None):
    """
    Computes rank for matrices and vectors.
    The only issue that may arise is that because numpy uses double
    precision lapack calls whereas numba uses type specific lapack
    calls, some singular values may differ and therefore counting the
    number of them above a tolerance may lead to different counts,
    and therefore rank, in some cases.
    """
    ensure_lapack()

    _check_linalg_1_or_2d_matrix(A, "matrix_rank")

    def _2d_matrix_rank_impl(A, tol):

        # handle the tol==None case separately for type inference to work
        if tol in (None, types.none):
            nb_type = getattr(A.dtype, "underlying_float", A.dtype)
            np_type = np_support.as_dtype(nb_type)
            eps_val = np.finfo(np_type).eps

            def _2d_tol_none_impl(A, tol=None):
                s = _compute_singular_values(A)
                # replicate numpy default tolerance calculation
                r = A.shape[0]
                c = A.shape[1]
                l = max(r, c)
                t = s[0] * l * eps_val
                return _get_rank_from_singular_values(s, t)
            return _2d_tol_none_impl
        else:
            def _2d_tol_not_none_impl(A, tol=None):
                s = _compute_singular_values(A)
                return _get_rank_from_singular_values(s, tol)
            return _2d_tol_not_none_impl

    def _get_matrix_rank_impl(A, tol):
        ndim = A.ndim
        if ndim == 1:
            # NOTE: Technically, the numpy implementation could be argued as
            # incorrect for the case of a vector (1D matrix). If a tolerance
            # is provided and a vector with a singular value below tolerance is
            # encountered this should report a rank of zero, the numpy
            # implementation does not do this and instead elects to report that
            # if any value in the vector is nonzero then the rank is 1.
            # An example would be [0, 1e-15, 0, 2e-15] which numpy reports as
            # rank 1 invariant of `tol`. The singular value for this vector is
            # obviously sqrt(5)*1e-15 and so a tol of e.g. sqrt(6)*1e-15 should
            # lead to a reported rank of 0 whereas a tol of 1e-15 should lead
            # to a reported rank of 1, numpy reports 1 regardless.
            # The code below replicates the numpy behaviour.
            def _1d_matrix_rank_impl(A, tol=None):
                for k in range(len(A)):
                    if A[k] != 0.:
                        return 1
                return 0
            return _1d_matrix_rank_impl
        elif ndim == 2:
            return _2d_matrix_rank_impl(A, tol)
        else:
            assert 0  # unreachable

    return _get_matrix_rank_impl(A, tol)


@overload(np.linalg.matrix_power)
def matrix_power_impl(a, n):
    """
    Computes matrix power. Only integer powers are supported in numpy.
    """

    _check_linalg_matrix(a, "matrix_power")
    np_dtype = np_support.as_dtype(a.dtype)

    nt = getattr(n, 'dtype', n)
    if not isinstance(nt, types.Integer):
        raise NumbaTypeError("Exponent must be an integer.")

    def matrix_power_impl(a, n):

        if n == 0:
            # this should be eye() but it doesn't support
            # the dtype kwarg yet so do it manually to save
            # the copy required by eye(a.shape[0]).asdtype()
            A = np.zeros(a.shape, dtype=np_dtype)
            for k in range(a.shape[0]):
                A[k, k] = 1.
            return A

        am, an = a.shape[-1], a.shape[-2]
        if am != an:
            raise ValueError('input must be a square array')

        # empty, return a copy
        if am == 0:
            return a.copy()

        # note: to be consistent over contiguousness, C order is
        # returned as that is what dot() produces and the most common
        # paths through matrix_power will involve that. Therefore
        # copies are made here to ensure the data ordering is
        # correct for paths not going via dot().

        if n < 0:
            A = np.linalg.inv(a).copy()
            if n == -1:  # return now
                return A
            n = -n
        else:
            if n == 1:  # return a copy now
                return a.copy()
            A = a  # this is safe, `a` is only read

        if n < 4:
            if n == 2:
                return np.dot(A, A)
            if n == 3:
                return np.dot(np.dot(A, A), A)
        else:

            acc = A
            exp = n

            # Initialise ret, SSA cannot see the loop will execute, without this
            # it appears as uninitialised.
            ret = acc
            # tried a loop split and branchless using identity matrix as
            # input but it seems like having a "first entry" flag is quicker
            flag = True
            while exp != 0:
                if exp & 1:
                    if flag:
                        ret = acc
                        flag = False
                    else:
                        ret = np.dot(ret, acc)
                acc = np.dot(acc, acc)
                exp = exp >> 1

            return ret

    return matrix_power_impl

# This is documented under linalg despite not being in the module


@overload(np.trace)
def matrix_trace_impl(a, offset=0):
    """
    Computes the trace of an array.
    """

    _check_linalg_matrix(a, "trace", la_prefix=False)

    if not isinstance(offset, (int, types.Integer)):
        raise NumbaTypeError("integer argument expected, got %s" % offset)

    def matrix_trace_impl(a, offset=0):
        rows, cols = a.shape
        k = offset
        if k < 0:
            rows = rows + k
        if k > 0:
            cols = cols - k
        n = max(min(rows, cols), 0)
        ret = 0
        if k >= 0:
            for i in range(n):
                ret += a[i, k + i]
        else:
            for i in range(n):
                ret += a[i - k, i]
        return ret

    return matrix_trace_impl


def _check_scalar_or_lt_2d_mat(a, func_name, la_prefix=True):
    prefix = "np.linalg" if la_prefix else "np"
    interp = (prefix, func_name)
    # checks that a matrix is 1 or 2D
    if isinstance(a, types.Array):
        if not a.ndim <= 2:
            raise TypingError("%s.%s() only supported on 1 and 2-D arrays "
                              % interp, highlighting=False)


@register_jitable
def outer_impl_none(a, b, out):
    aa = np.asarray(a)
    bb = np.asarray(b)
    return np.multiply(aa.ravel().reshape((aa.size, 1)),
                        bb.ravel().reshape((1, bb.size)))


@register_jitable
def outer_impl_arr(a, b, out):
    aa = np.asarray(a)
    bb = np.asarray(b)
    np.multiply(aa.ravel().reshape((aa.size, 1)),
                bb.ravel().reshape((1, bb.size)),
                out)
    return out


def _get_outer_impl(a, b, out):
    if out in (None, types.none):
        return outer_impl_none
    else:
        return outer_impl_arr


@overload(np.outer)
def outer_impl(a, b, out=None):

    _check_scalar_or_lt_2d_mat(a, "outer", la_prefix=False)
    _check_scalar_or_lt_2d_mat(b, "outer", la_prefix=False)

    impl = _get_outer_impl(a, b, out)

    def outer_impl(a, b, out=None):
        return impl(a, b, out)

    return outer_impl


def _kron_normaliser_impl(x):
    # makes x into a 2d array
    if isinstance(x, types.Array):
        if x.layout not in ('C', 'F'):
            raise TypingError("np.linalg.kron only supports 'C' or 'F' layout "
                              "input arrays. Received an input of "
                              "layout '{}'.".format(x.layout))
        elif x.ndim == 2:
            @register_jitable
            def nrm_shape(x):
                xn = x.shape[-1]
                xm = x.shape[-2]
                return x.reshape(xm, xn)
            return nrm_shape
        else:
            @register_jitable
            def nrm_shape(x):
                xn = x.shape[-1]
                return x.reshape(1, xn)
            return nrm_shape
    else:  # assume its a scalar
        @register_jitable
        def nrm_shape(x):
            a = np.empty((1, 1), type(x))
            a[0] = x
            return a
        return nrm_shape


def _kron_return(a, b):
    # transforms c into something that kron would return
    # based on the shapes of a and b
    a_is_arr = isinstance(a, types.Array)
    b_is_arr = isinstance(b, types.Array)
    if a_is_arr and b_is_arr:
        if a.ndim == 2 or b.ndim == 2:
            @register_jitable
            def ret(a, b, c):
                return c
            return ret
        else:
            @register_jitable
            def ret(a, b, c):
                return c.reshape(c.size)
            return ret
    else:  # at least one of (a, b) is a scalar
        if a_is_arr:
            @register_jitable
            def ret(a, b, c):
                return c.reshape(a.shape)
            return ret
        elif b_is_arr:
            @register_jitable
            def ret(a, b, c):
                return c.reshape(b.shape)
            return ret
        else:  # both scalars
            @register_jitable
            def ret(a, b, c):
                return c[0]
            return ret


@overload(np.kron)
def kron_impl(a, b):

    _check_scalar_or_lt_2d_mat(a, "kron", la_prefix=False)
    _check_scalar_or_lt_2d_mat(b, "kron", la_prefix=False)

    fix_a = _kron_normaliser_impl(a)
    fix_b = _kron_normaliser_impl(b)
    ret_c = _kron_return(a, b)

    # this is fine because the ufunc for the Hadamard product
    # will reject differing dtypes in a and b.
    dt = getattr(a, 'dtype', a)

    def kron_impl(a, b):

        aa = fix_a(a)
        bb = fix_b(b)

        am = aa.shape[-2]
        an = aa.shape[-1]
        bm = bb.shape[-2]
        bn = bb.shape[-1]

        cm = am * bm
        cn = an * bn

        # allocate c
        C = np.empty((cm, cn), dtype=dt)

        # In practice this is runs quicker than the more obvious
        # `each element of A multiplied by B and assigned to
        # a block in C` like alg.

        # loop over rows of A
        for i in range(am):
            # compute the column offset into C
            rjmp = i * bm
            # loop over rows of B
            for k in range(bm):
                # compute row the offset into C
                irjmp = rjmp + k
                # slice a given row of B
                slc = bb[k, :]
                # loop over columns of A
                for j in range(an):
                    # vectorized assignment of an element of A
                    # multiplied by the current row of B into
                    # a slice of a row of C
                    cjmp = j * bn
                    C[irjmp, cjmp:cjmp + bn] = aa[i, j] * slc

        return ret_c(a, b, C)

    return kron_impl