fvnlib.f90 53 KB
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
module fvn
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! fvn : a f95 module replacement for some imsl routines
! it uses lapack for linear algebra
! it uses modified quadpack for integration
!
! William Daniau 2007
! william.daniau@femto-st.fr
!
! Routines naming scheme :
!
!           fvn_x_name
!           where x can be  s : real 
!                           d : real double precision
!                           c : complex
!                           z : double complex
!
!
! This piece of code is totally free! Do whatever you want with it. However
! if you find it usefull it would be kind to give credits ;-) Nevertheless, you
! may give credits to quadpack authors. 
!
! Version 1.1
!
! TO DO LIST :
! + Order eigenvalues and vectors in decreasing eigenvalue's modulus order -> atm 
! eigenvalues are given with no particular order.
! + Generic interface for fvn_x_name family  -> fvn_name
! + Make some parameters optional, status for example
! + use f95 kinds "double complex" -> complex(kind=8)
! + unify quadpack routines
! + ...
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

implicit none
! All quadpack routines are private to the module
private ::  d1mach,dqag,dqag_2d_inner,dqag_2d_outer,dqage,dqage_2d_inner, &
            dqage_2d_outer,dqk15,dqk15_2d_inner,dqk15_2d_outer,dqk21,dqk21_2d_inner,dqk21_2d_outer, &
            dqk31,dqk31_2d_inner,dqk31_2d_outer,dqk41,dqk41_2d_inner,dqk41_2d_outer, &
            dqk51,dqk51_2d_inner,dqk51_2d_outer,dqk61,dqk61_2d_inner,dqk61_2d_outer,dqpsrt
            
            
contains 

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! Matrix inversion subroutines
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine fvn_s_matinv(d,a,inva,status)
    !
    ! Matrix inversion of a real matrix using BLAS and LAPACK
    !
    ! d (in) : matrix rank
    ! a (in) : input matrix
    ! inva (out) : inversed matrix
    ! status (ou) : =0 if something failed
    !
    integer, intent(in) :: d
    real, intent(in) :: a(d,d)
    real, intent(out) :: inva(d,d)
    integer, intent(out) :: status

    integer, allocatable :: ipiv(:)
    real, allocatable :: work(:)
    real twork(1)
    integer :: info
    integer :: lwork

    status=1

    allocate(ipiv(d))
    ! copy a into inva using BLAS
    !call scopy(d*d,a,1,inva,1)
    inva(:,:)=a(:,:)
    ! LU factorization using LAPACK
    call sgetrf(d,d,inva,d,ipiv,info)
    ! if info is not equal to 0, something went wrong we exit setting status to 0
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        return
    end if
    ! we use the query fonction of xxxtri to obtain the optimal workspace size
    call sgetri(d,inva,d,ipiv,twork,-1,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    ! Matrix inversion using LAPACK
    call sgetri(d,inva,d,ipiv,work,lwork,info)
    ! again if info is not equal to 0, we exit setting status to 0
    if (info /= 0) then
        status=0
    end if
    deallocate(work)
    deallocate(ipiv)
end subroutine

subroutine fvn_d_matinv(d,a,inva,status)
    !
    ! Matrix inversion of a double precision matrix using BLAS and LAPACK
    !
    ! d (in) : matrix rank
    ! a (in) : input matrix
    ! inva (out) : inversed matrix
    ! status (ou) : =0 if something failed
    !
    integer, intent(in) :: d
    double precision, intent(in) :: a(d,d)
    double precision, intent(out) :: inva(d,d)
    integer, intent(out) :: status

    integer, allocatable :: ipiv(:)
    double precision, allocatable :: work(:)
    double precision :: twork(1)
    integer :: info
    integer :: lwork

    status=1

    allocate(ipiv(d))
    ! copy a into inva using BLAS
    !call dcopy(d*d,a,1,inva,1)
    inva(:,:)=a(:,:)
    ! LU factorization using LAPACK
    call dgetrf(d,d,inva,d,ipiv,info)
    ! if info is not equal to 0, something went wrong we exit setting status to 0
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        return
    end if
    ! we use the query fonction of xxxtri to obtain the optimal workspace size
    call dgetri(d,inva,d,ipiv,twork,-1,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    ! Matrix inversion using LAPACK
    call dgetri(d,inva,d,ipiv,work,lwork,info)
    ! again if info is not equal to 0, we exit setting status to 0
    if (info /= 0) then
        status=0
    end if
    deallocate(work)
    deallocate(ipiv)
end subroutine

subroutine fvn_c_matinv(d,a,inva,status)
    !
    ! Matrix inversion of a complex matrix using BLAS and LAPACK
    !
    ! d (in) : matrix rank
    ! a (in) : input matrix
    ! inva (out) : inversed matrix
    ! status (ou) : =0 if something failed
    !
    integer, intent(in) :: d
    complex, intent(in) :: a(d,d)
    complex, intent(out) :: inva(d,d)
    integer, intent(out) :: status

    integer, allocatable :: ipiv(:)
    complex, allocatable :: work(:)
    complex :: twork(1)
    integer :: info
    integer :: lwork

    status=1

    allocate(ipiv(d))
    ! copy a into inva using BLAS
    !call ccopy(d*d,a,1,inva,1)
    inva(:,:)=a(:,:)
    
    ! LU factorization using LAPACK
    call cgetrf(d,d,inva,d,ipiv,info)
    ! if info is not equal to 0, something went wrong we exit setting status to 0
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        return
    end if
    ! we use the query fonction of xxxtri to obtain the optimal workspace size
    call cgetri(d,inva,d,ipiv,twork,-1,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    ! Matrix inversion using LAPACK
    call cgetri(d,inva,d,ipiv,work,lwork,info)
    ! again if info is not equal to 0, we exit setting status to 0
    if (info /= 0) then
        status=0
    end if
    deallocate(work)
    deallocate(ipiv)
end subroutine

subroutine fvn_z_matinv(d,a,inva,status)
    !
    ! Matrix inversion of a double complex matrix using BLAS and LAPACK
    !
    ! d (in) : matrix rank
    ! a (in) : input matrix
    ! inva (out) : inversed matrix
    ! status (ou) : =0 if something failed
    !
    integer, intent(in) :: d
    double complex, intent(in) :: a(d,d)
    double complex, intent(out) :: inva(d,d)
    integer, intent(out) :: status

    integer, allocatable :: ipiv(:)
    double complex, allocatable :: work(:)
    double complex :: twork(1)
    integer :: info
    integer :: lwork

    status=1

    allocate(ipiv(d))
    ! copy a into inva using BLAS
    !call zcopy(d*d,a,1,inva,1)
    inva(:,:)=a(:,:)
    
    ! LU factorization using LAPACK
    call zgetrf(d,d,inva,d,ipiv,info)
    ! if info is not equal to 0, something went wrong we exit setting status to 0
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        return
    end if
    ! we use the query fonction of xxxtri to obtain the optimal workspace size
    call zgetri(d,inva,d,ipiv,twork,-1,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    ! Matrix inversion using LAPACK
    call zgetri(d,inva,d,ipiv,work,lwork,info)
    ! again if info is not equal to 0, we exit setting status to 0
    if (info /= 0) then
        status=0
    end if
    deallocate(work)
    deallocate(ipiv)
end subroutine


!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! Determinants
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
function fvn_s_det(d,a,status)
    !
    ! Evaluate the determinant of a square matrix using lapack LU factorization
    !
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! status (out) : =0 if LU factorization failed
    !
    integer, intent(in) :: d
    real, intent(in) :: a(d,d)
    integer, intent(out) :: status
    real :: fvn_s_det
    
    real, allocatable :: wc_a(:,:)
    integer, allocatable :: ipiv(:)
    integer :: info,i
    
    status=1
    allocate(wc_a(d,d))
    allocate(ipiv(d))
    wc_a(:,:)=a(:,:)
    call sgetrf(d,d,wc_a,d,ipiv,info)
    if (info/= 0) then
        status=0
        fvn_s_det=0.e0
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if
    fvn_s_det=1.e0
    do i=1,d
        if (ipiv(i)==i) then
            fvn_s_det=fvn_s_det*wc_a(i,i)
        else
            fvn_s_det=-fvn_s_det*wc_a(i,i)
        end if
    end do
    deallocate(ipiv)
    deallocate(wc_a)

end function

function fvn_d_det(d,a,status)
    !
    ! Evaluate the determinant of a square matrix using lapack LU factorization
    !
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! status (out) : =0 if LU factorization failed
    !
    integer, intent(in) :: d
    double precision, intent(in) :: a(d,d)
    integer, intent(out) :: status
    double precision :: fvn_d_det
    
    double precision, allocatable :: wc_a(:,:)
    integer, allocatable :: ipiv(:)
    integer :: info,i
    
    status=1
    allocate(wc_a(d,d))
    allocate(ipiv(d))
    wc_a(:,:)=a(:,:)
    call dgetrf(d,d,wc_a,d,ipiv,info)
    if (info/= 0) then
        status=0
        fvn_d_det=0.d0
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if
    fvn_d_det=1.d0
    do i=1,d
        if (ipiv(i)==i) then
            fvn_d_det=fvn_d_det*wc_a(i,i)
        else
            fvn_d_det=-fvn_d_det*wc_a(i,i)
        end if
    end do
    deallocate(ipiv)
    deallocate(wc_a)

end function

function fvn_c_det(d,a,status)    !
    ! Evaluate the determinant of a square matrix using lapack LU factorization
    !
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! status (out) : =0 if LU factorization failed
    !
    integer, intent(in) :: d
    complex, intent(in) :: a(d,d)
    integer, intent(out) :: status
    complex :: fvn_c_det
    
    complex, allocatable :: wc_a(:,:)
    integer, allocatable :: ipiv(:)
    integer :: info,i
    
    status=1
    allocate(wc_a(d,d))
    allocate(ipiv(d))
    wc_a(:,:)=a(:,:)
    call cgetrf(d,d,wc_a,d,ipiv,info)
    if (info/= 0) then
        status=0
        fvn_c_det=(0.e0,0.e0)
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if
    fvn_c_det=(1.e0,0.e0)
    do i=1,d
        if (ipiv(i)==i) then
            fvn_c_det=fvn_c_det*wc_a(i,i)
        else
            fvn_c_det=-fvn_c_det*wc_a(i,i)
        end if
    end do
    deallocate(ipiv)
    deallocate(wc_a)

end function

function fvn_z_det(d,a,status)
    !
    ! Evaluate the determinant of a square matrix using lapack LU factorization
    !
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! det (out) : determinant
    ! status (out) : =0 if LU factorization failed
    !
    integer, intent(in) :: d
    double complex, intent(in) :: a(d,d)
    integer, intent(out) :: status
    double complex :: fvn_z_det
    
    double complex, allocatable :: wc_a(:,:)
    integer, allocatable :: ipiv(:)
    integer :: info,i
    
    status=1
    allocate(wc_a(d,d))
    allocate(ipiv(d))
    wc_a(:,:)=a(:,:)
    call zgetrf(d,d,wc_a,d,ipiv,info)
    if (info/= 0) then
        status=0
        fvn_z_det=(0.d0,0.d0)
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if
    fvn_z_det=(1.d0,0.d0)
    do i=1,d
        if (ipiv(i)==i) then
            fvn_z_det=fvn_z_det*wc_a(i,i)
        else
            fvn_z_det=-fvn_z_det*wc_a(i,i)
        end if
    end do
    deallocate(ipiv)
    deallocate(wc_a)

end function

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! Condition test
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! 1-norm 
! fonction lapack slange,dlange,clange,zlange pour obtenir la 1-norm
! fonction lapack sgecon,dgecon,cgecon,zgecon pour calculer la rcond
!
subroutine fvn_s_matcon(d,a,rcond,status)
    ! Matrix condition (reciprocal of condition number)
    ! 
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! rcond (out) : guess what
    ! status (out) : =0 if something went wrong
    !
    integer, intent(in) :: d
    real, intent(in) :: a(d,d)
    real, intent(out) :: rcond
    integer, intent(out) :: status
    
    real, allocatable :: work(:)
    integer, allocatable :: iwork(:)
    real :: anorm
    real, allocatable :: wc_a(:,:) ! working copy of a
    integer :: info
    integer, allocatable :: ipiv(:)
    
    real, external :: slange
    
    
    status=1
    
    anorm=slange('1',d,d,a,d,work) ! work is unallocated as it is only used when computing infinity norm
    
    allocate(wc_a(d,d))
    !call scopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    allocate(ipiv(d))
    call sgetrf(d,d,wc_a,d,ipiv,info)
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if
    allocate(work(4*d))
    allocate(iwork(d))
    call sgecon('1',d,wc_a,d,anorm,rcond,work,iwork,info)
    if (info /= 0) then
        status=0
    end if
    deallocate(iwork)
    deallocate(work)
    deallocate(ipiv)
    deallocate(wc_a)

end subroutine

subroutine fvn_d_matcon(d,a,rcond,status)
    ! Matrix condition (reciprocal of condition number)
    ! 
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! rcond (out) : guess what
    ! status (out) : =0 if something went wrong
    !
    integer, intent(in) :: d
    double precision, intent(in) :: a(d,d)
    double precision, intent(out) :: rcond
    integer, intent(out) :: status
    
    double precision, allocatable :: work(:)
    integer, allocatable :: iwork(:)
    double precision :: anorm
    double precision, allocatable :: wc_a(:,:) ! working copy of a
    integer :: info
    integer, allocatable :: ipiv(:)
    
    double precision, external :: dlange
    
    
    status=1
    
    anorm=dlange('1',d,d,a,d,work) ! work is unallocated as it is only used when computing infinity norm
    
    allocate(wc_a(d,d))
    !call dcopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    allocate(ipiv(d))
    call dgetrf(d,d,wc_a,d,ipiv,info)
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if

    allocate(work(4*d))
    allocate(iwork(d))
    call dgecon('1',d,wc_a,d,anorm,rcond,work,iwork,info)
    if (info /= 0) then
        status=0
    end if
    deallocate(iwork)
    deallocate(work)
    deallocate(ipiv)
    deallocate(wc_a)

end subroutine

subroutine fvn_c_matcon(d,a,rcond,status)
    ! Matrix condition (reciprocal of condition number)
    ! 
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! rcond (out) : guess what
    ! status (out) : =0 if something went wrong
    !
    integer, intent(in) :: d
    complex, intent(in) :: a(d,d)
    real, intent(out) :: rcond
    integer, intent(out) :: status
    
    real, allocatable :: rwork(:)
    complex, allocatable :: work(:)
    integer, allocatable :: iwork(:)
    real :: anorm
    complex, allocatable :: wc_a(:,:) ! working copy of a
    integer :: info
    integer, allocatable :: ipiv(:)
    
    real, external :: clange
    
    
    status=1
    
    anorm=clange('1',d,d,a,d,rwork) ! rwork is unallocated as it is only used when computing infinity norm
    
    allocate(wc_a(d,d))
    !call ccopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    allocate(ipiv(d))
    call cgetrf(d,d,wc_a,d,ipiv,info)
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if
    allocate(work(2*d))
    allocate(rwork(2*d))
    call cgecon('1',d,wc_a,d,anorm,rcond,work,rwork,info)
    if (info /= 0) then
        status=0
    end if
    deallocate(rwork)
    deallocate(work)
    deallocate(ipiv)
    deallocate(wc_a)
end subroutine

subroutine fvn_z_matcon(d,a,rcond,status)
    ! Matrix condition (reciprocal of condition number)
    ! 
    ! d (in) : matrix rank
    ! a (in) : The Matrix
    ! rcond (out) : guess what
    ! status (out) : =0 if something went wrong
    !
    integer, intent(in) :: d
    double complex, intent(in) :: a(d,d)
    double precision, intent(out) :: rcond
    integer, intent(out) :: status
    
    double complex, allocatable :: work(:)
    double precision, allocatable :: rwork(:)
    double precision :: anorm
    double complex, allocatable :: wc_a(:,:) ! working copy of a
    integer :: info
    integer, allocatable :: ipiv(:)
    
    double precision, external :: zlange
    
    
    status=1
    
    anorm=zlange('1',d,d,a,d,rwork) ! rwork is unallocated as it is only used when computing infinity norm
    
    allocate(wc_a(d,d))
    !call zcopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    allocate(ipiv(d))
    call zgetrf(d,d,wc_a,d,ipiv,info)
    if (info /= 0) then
        status=0
        deallocate(ipiv)
        deallocate(wc_a)
        return
    end if

    allocate(work(2*d))
    allocate(rwork(2*d))
    call zgecon('1',d,wc_a,d,anorm,rcond,work,rwork,info)
    if (info /= 0) then
        status=0
    end if
    deallocate(rwork)
    deallocate(work)
    deallocate(ipiv)
    deallocate(wc_a)
end subroutine

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! Valeurs propres/ Vecteurs propre
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

subroutine fvn_s_matev(d,a,evala,eveca,status)
    ! 
    ! integer d (in) : matrice rank
    ! real a(d,d) (in) : The Matrix
    ! complex evala(d) (out) : eigenvalues
    ! complex eveca(d,d) (out) : eveca(:,j) = jth eigenvector
    ! integer (out) : status =0 if something went wrong
    !
    ! interfacing Lapack routine SGEEV

    integer, intent(in) :: d
    real, intent(in) :: a(d,d)
    complex, intent(out) :: evala(d)
    complex, intent(out) :: eveca(d,d)
    integer, intent(out) :: status
    
    real, allocatable :: wc_a(:,:)  ! a working copy of a
    integer :: info
    integer :: lwork
    real, allocatable :: wr(:),wi(:)
    real :: vl      ! unused but necessary for the call
    real, allocatable :: vr(:,:)
    real, allocatable :: work(:)
    real :: twork(1)
    integer i
    integer j
    
    ! making a working copy of a
    allocate(wc_a(d,d))
    !call scopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    allocate(wr(d))
    allocate(wi(d))
    allocate(vr(d,d))
    ! query optimal work size
    call sgeev('N','V',d,wc_a,d,wr,wi,vl,1,vr,d,twork,-1,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    call sgeev('N','V',d,wc_a,d,wr,wi,vl,1,vr,d,work,lwork,info)

    if (info /= 0) then
        status=0
        deallocate(work)
        deallocate(vr)
        deallocate(wi)
        deallocate(wr)
        deallocate(wc_a)
        return
    end if

    ! now fill in the results
    i=1
    do while(i<=d)
        evala(i)=cmplx(wr(i),wi(i))
        if (wi(i) == 0.) then ! eigenvalue is real
            eveca(:,i)=cmplx(vr(:,i),0.)
        else ! eigenvalue is complex
            evala(i+1)=cmplx(wr(i+1),wi(i+1))
            eveca(:,i)=cmplx(vr(:,i),vr(:,i+1))
            eveca(:,i+1)=cmplx(vr(:,i),-vr(:,i+1))
            i=i+1
        end if
        i=i+1
    enddo
    deallocate(work)
    deallocate(vr)
    deallocate(wi)
    deallocate(wr)
    deallocate(wc_a)

end subroutine

subroutine fvn_d_matev(d,a,evala,eveca,status)
    ! 
    ! integer d (in) : matrice rank
    ! double precision a(d,d) (in) : The Matrix
    ! double complex evala(d) (out) : eigenvalues
    ! double complex eveca(d,d) (out) : eveca(:,j) = jth eigenvector
    ! integer (out) : status =0 if something went wrong
    !
    ! interfacing Lapack routine DGEEV
    integer, intent(in) :: d
    double precision, intent(in) :: a(d,d)
    double complex, intent(out) :: evala(d)
    double complex, intent(out) :: eveca(d,d)
    integer, intent(out) :: status
    
    double precision, allocatable :: wc_a(:,:)  ! a working copy of a
    integer :: info
    integer :: lwork
    double precision, allocatable :: wr(:),wi(:)
    double precision :: vl      ! unused but necessary for the call
    double precision, allocatable :: vr(:,:)
    double precision, allocatable :: work(:)
    double precision :: twork(1)
    integer i
    integer j
    
    ! making a working copy of a
    allocate(wc_a(d,d))
    !call dcopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    allocate(wr(d))
    allocate(wi(d))
    allocate(vr(d,d))
    ! query optimal work size
    call dgeev('N','V',d,wc_a,d,wr,wi,vl,1,vr,d,twork,-1,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    call dgeev('N','V',d,wc_a,d,wr,wi,vl,1,vr,d,work,lwork,info)

    if (info /= 0) then
        status=0
        deallocate(work)
        deallocate(vr)
        deallocate(wi)
        deallocate(wr)
        deallocate(wc_a)
        return
    end if

    ! now fill in the results
    i=1
    do while(i<=d)
        evala(i)=dcmplx(wr(i),wi(i))
        if (wi(i) == 0.) then ! eigenvalue is real
            eveca(:,i)=dcmplx(vr(:,i),0.)
        else ! eigenvalue is complex
            evala(i+1)=dcmplx(wr(i+1),wi(i+1))
            eveca(:,i)=dcmplx(vr(:,i),vr(:,i+1))
            eveca(:,i+1)=dcmplx(vr(:,i),-vr(:,i+1))
            i=i+1
        end if
        i=i+1
    enddo

    deallocate(work)
    deallocate(vr)
    deallocate(wi)
    deallocate(wr)
    deallocate(wc_a)

end subroutine

subroutine fvn_c_matev(d,a,evala,eveca,status)
    ! 
    ! integer d (in) : matrice rank
    ! complex a(d,d) (in) : The Matrix
    ! complex evala(d) (out) : eigenvalues
    ! complex eveca(d,d) (out) : eveca(:,j) = jth eigenvector
    ! integer (out) : status =0 if something went wrong
    !
    ! interfacing Lapack routine CGEEV

    integer, intent(in) :: d
    complex, intent(in) :: a(d,d)
    complex, intent(out) :: evala(d)
    complex, intent(out) :: eveca(d,d)
    integer, intent(out) :: status

    complex, allocatable :: wc_a(:,:) ! a working copy of a
    integer :: info
    integer :: lwork
    complex, allocatable :: work(:)
    complex :: twork(1)
    real, allocatable :: rwork(:)
    complex :: vl   ! unused but necessary for the call

    status=1

    ! making a working copy of a
    allocate(wc_a(d,d))
    !call ccopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    
    ! query optimal work size
    call cgeev('N','V',d,wc_a,d,evala,vl,1,eveca,d,twork,-1,rwork,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    allocate(rwork(2*d))
    call cgeev('N','V',d,wc_a,d,evala,vl,1,eveca,d,work,lwork,rwork,info)
    
    if (info /= 0) then
        status=0
    end if
    deallocate(rwork)
    deallocate(work)
    deallocate(wc_a)

end subroutine

subroutine fvn_z_matev(d,a,evala,eveca,status)
    ! 
    ! integer d (in) : matrice rank
    ! double complex a(d,d) (in) : The Matrix
    ! double complex evala(d) (out) : eigenvalues
    ! double complex eveca(d,d) (out) : eveca(:,j) = jth eigenvector
    ! integer (out) : status =0 if something went wrong
    !
    ! interfacing Lapack routine ZGEEV

    integer, intent(in) :: d
    double complex, intent(in) :: a(d,d)
    double complex, intent(out) :: evala(d)
    double complex, intent(out) :: eveca(d,d)
    integer, intent(out) :: status

    double complex, allocatable :: wc_a(:,:) ! a working copy of a
    integer :: info
    integer :: lwork
    double complex, allocatable :: work(:)
    double complex :: twork(1)
    double precision, allocatable :: rwork(:)
    double complex :: vl   ! unused but necessary for the call
    
    status=1

    ! making a working copy of a
    allocate(wc_a(d,d))
    !call zcopy(d*d,a,1,wc_a,1)
    wc_a(:,:)=a(:,:)
    
    ! query optimal work size
    call zgeev('N','V',d,wc_a,d,evala,vl,1,eveca,d,twork,-1,rwork,info)
    lwork=int(twork(1))
    allocate(work(lwork))
    allocate(rwork(2*d))
    call zgeev('N','V',d,wc_a,d,evala,vl,1,eveca,d,work,lwork,rwork,info)
    
    if (info /= 0) then
        status=0
    end if
    deallocate(rwork)
    deallocate(work)
    deallocate(wc_a)

end subroutine

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! Akima spline interpolation and spline evaluation
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

! Single precision
subroutine fvn_s_akima(n,x,y,br,co)
    implicit none
    integer, intent(in)  :: n
    real, intent(in) :: x(n)
    real, intent(in) :: y(n)
    real, intent(out) :: br(n)
    real, intent(out) :: co(4,n)
    
    real, allocatable :: var(:),z(:)
    real :: wi_1,wi
    integer :: i
    real :: dx,a,b

    ! br is just a copy of x
    br(:)=x(:)
    
    allocate(var(n))
    allocate(z(n))
    ! evaluate the variations
    do i=1, n-1
        var(i+2)=(y(i+1)-y(i))/(x(i+1)-x(i))
    end do
    var(n+2)=2.e0*var(n+1)-var(n)
    var(n+3)=2.e0*var(n+2)-var(n+1)
    var(2)=2.e0*var(3)-var(4)
    var(1)=2.e0*var(2)-var(3)
  
    do i = 1, n
    wi_1=abs(var(i+3)-var(i+2))
    wi=abs(var(i+1)-var(i))
    if ((wi_1+wi).eq.0.e0) then
        z(i)=(var(i+2)+var(i+1))/2.e0
    else
        z(i)=(wi_1*var(i+1)+wi*var(i+2))/(wi_1+wi)
    end if
    end do
    
    do i=1, n-1
        dx=x(i+1)-x(i)
        a=(z(i+1)-z(i))*dx                      ! coeff intermediaires pour calcul wd
        b=y(i+1)-y(i)-z(i)*dx                   ! coeff intermediaires pour calcul wd
        co(1,i)=y(i)
        co(2,i)=z(i)
        !co(3,i)=-(a-3.*b)/dx**2                ! méthode wd
        !co(4,i)=(a-2.*b)/dx**3                 ! méthode wd
        co(3,i)=(3.e0*var(i+2)-2.e0*z(i)-z(i+1))/dx   ! méthode JP Moreau
        co(4,i)=(z(i)+z(i+1)-2.e0*var(i+2))/dx**2  !
        ! les coefficients donnés par imsl sont co(3,i)*2 et co(4,i)*6
        ! etrangement la fonction csval corrige et donne la bonne valeur ...
    end do
    co(1,n)=y(n)
    co(2,n)=z(n)
    co(3,n)=0.e0
    co(4,n)=0.e0

    deallocate(z)
    deallocate(var)

end subroutine

! Double precision
subroutine fvn_d_akima(n,x,y,br,co)

    implicit none
    integer, intent(in)  :: n
    double precision, intent(in) :: x(n)
    double precision, intent(in) :: y(n)
    double precision, intent(out) :: br(n)
    double precision, intent(out) :: co(4,n)
    
    double precision, allocatable :: var(:),z(:)
    double precision :: wi_1,wi
    integer :: i
    double precision :: dx,a,b
    
    ! br is just a copy of x
    br(:)=x(:)

    allocate(var(n))
    allocate(z(n))
    ! evaluate the variations
    do i=1, n-1
        var(i+2)=(y(i+1)-y(i))/(x(i+1)-x(i))
    end do
    var(n+2)=2.d0*var(n+1)-var(n)
    var(n+3)=2.d0*var(n+2)-var(n+1)
    var(2)=2.d0*var(3)-var(4)
    var(1)=2.d0*var(2)-var(3)
  
    do i = 1, n
    wi_1=dabs(var(i+3)-var(i+2))
    wi=dabs(var(i+1)-var(i))
    if ((wi_1+wi).eq.0.d0) then
        z(i)=(var(i+2)+var(i+1))/2.d0
    else
        z(i)=(wi_1*var(i+1)+wi*var(i+2))/(wi_1+wi)
    end if
    end do
    
    do i=1, n-1
        dx=x(i+1)-x(i)
        a=(z(i+1)-z(i))*dx                      ! coeff intermediaires pour calcul wd
        b=y(i+1)-y(i)-z(i)*dx                   ! coeff intermediaires pour calcul wd
        co(1,i)=y(i)
        co(2,i)=z(i)
        !co(3,i)=-(a-3.*b)/dx**2                ! méthode wd
        !co(4,i)=(a-2.*b)/dx**3                 ! méthode wd
        co(3,i)=(3.d0*var(i+2)-2.d0*z(i)-z(i+1))/dx   ! méthode JP Moreau
        co(4,i)=(z(i)+z(i+1)-2.d0*var(i+2))/dx**2  !
        ! les coefficients donnés par imsl sont co(3,i)*2 et co(4,i)*6
        ! etrangement la fonction csval corrige et donne la bonne valeur ...
    end do
    co(1,n)=y(n)
    co(2,n)=z(n)
    co(3,n)=0.d0
    co(4,n)=0.d0

    deallocate(z)
    deallocate(var)

end subroutine

!
! Single precision spline evaluation
!
function fvn_s_spline_eval(x,n,br,co)
    implicit none
    real, intent(in) :: x           ! x must be br(1)<= x <= br(n+1) otherwise value is extrapolated
    integer, intent(in) :: n        ! number of intervals
    real, intent(in) :: br(n+1)     ! breakpoints
    real, intent(in) :: co(4,n+1)   ! spline coeeficients
    real :: fvn_s_spline_eval
    
    integer :: i
    real :: dx
    
    if (x<=br(1)) then
        i=1
    else if (x>=br(n+1)) then
        i=n
    else
    i=1
    do while(x>=br(i))
        i=i+1
    end do
    i=i-1
    end if
    dx=x-br(i)
    fvn_s_spline_eval=co(1,i)+co(2,i)*dx+co(3,i)*dx**2+co(4,i)*dx**3

end function

! Double precision spline evaluation
function fvn_d_spline_eval(x,n,br,co)
    implicit none
    double precision, intent(in) :: x           ! x must be br(1)<= x <= br(n+1) otherwise value is extrapolated
    integer, intent(in) :: n        ! number of intervals
    double precision, intent(in) :: br(n+1)     ! breakpoints
    double precision, intent(in) :: co(4,n+1)   ! spline coeeficients
    double precision :: fvn_d_spline_eval
    
    integer :: i
    double precision :: dx
    
    
    if (x<=br(1)) then
        i=1
    else if (x>=br(n+1)) then
        i=n
    else
    i=1
    do while(x>=br(i))
        i=i+1
    end do
    i=i-1
    end if
    
    dx=x-br(i)
    fvn_d_spline_eval=co(1,i)+co(2,i)*dx+co(3,i)*dx**2+co(4,i)*dx**3

end function


!
! Muller
!
!
!
! William Daniau 2007
!
! This routine is a fortran 90 port of Hans D. Mittelmann's routine muller.f
! http://plato.asu.edu/ftp/other_software/muller.f
!
! it can be used as a replacement for imsl routine dzanly with minor changes
!
!-----------------------------------------------------------------------
!
!   purpose             - zeros of an analytic complex function
!                           using the muller method with deflation
!
!   usage               - call fvn_z_muller (f,eps,eps1,kn,n,nguess,x,itmax,
!                           infer,ier)
!
!   arguments    f      - a complex function subprogram, f(z), written
!                           by the user specifying the equation whose
!                           roots are to be found.  f must appear in
!                           an external statement in the calling pro-
!                           gram.
!                eps    - 1st stopping criterion.  let fp(z)=f(z)/p
!                           where p = (z-z(1))*(z-z(2))*,,,*(z-z(k-1))
!                           and z(1),...,z(k-1) are previously found
!                           roots.  if ((cdabs(f(z)).le.eps) .and.
!                           (cdabs(fp(z)).le.eps)), then z is accepted
!                           as a root. (input)
!                eps1   - 2nd stopping criterion.  a root is accepted
!                           if two successive approximations to a given
!                           root agree within eps1. (input)
!                             note. if either or both of the stopping
!                             criteria are fulfilled, the root is
!                             accepted.
!                kn     - the number of known roots which must be stored
!                           in x(1),...,x(kn), prior to entry to muller
!                nguess - the number of initial guesses provided. these
!                           guesses must be stored in x(kn+1),...,
!                           x(kn+nguess).  nguess must be set equal
!                           to zero if no guesses are provided. (input)
!                n      - the number of new roots to be found by
!                           muller (input)
!                x      - a complex vector of length kn+n.  x(1),...,
!                           x(kn) on input must contain any known
!                           roots.  x(kn+1),..., x(kn+n) on input may,
!                           on user option, contain initial guesses for
!                           the n new roots which are to be computed.
!                           if the user does not provide an initial
!                           guess, zero is used.
!                           on output, x(kn+1),...,x(kn+n) contain the
!                           approximate roots found by muller.
!                itmax  - the maximum allowable number of iterations
!                           per root (input)
!                infer  - an integer vector of length kn+n.  on
!                           output infer(j) contains the number of
!                           iterations used in finding the j-th root
!                           when convergence was achieved.  if
!                           convergence was not obtained in itmax
!                           iterations, infer(j) will be greater than
!                           itmax (output).
!                ier    - error parameter (output)
!                         warning error
!                           ier = 33 indicates failure to converge with-
!                             in itmax iterations for at least one of
!                             the (n) new roots.
!
!
!   remarks      muller always returns the last approximation for root j
!                in x(j). if the convergence criterion is satisfied,
!                then infer(j) is less than or equal to itmax. if the
!                convergence criterion is not satisified, then infer(j)
!                is set to either itmax+1 or itmax+k, with k greater
!                than 1. infer(j) = itmax+1 indicates that muller did
!                not obtain convergence in the allowed number of iter-
!                ations. in this case, the user may wish to set itmax
!                to a larger value. infer(j) = itmax+k means that con-
!                vergence was obtained (on iteration k) for the defla-
!                ted function
!                              fp(z) = f(z)/((z-z(1)...(z-z(j-1)))
!
!                but failed for f(z). in this case, better initial
!                guesses might help or, it might be necessary to relax
!                the convergence criterion.
!
!-----------------------------------------------------------------------
!
subroutine fvn_z_muller (f,eps,eps1,kn,nguess,n,x,itmax,infer,ier)
     implicit none
      double precision :: rzero,rten,rhun,rp01,ax,eps1,qz,eps,tpq
      double complex ::   d,dd,den,fprt,frt,h,rt,t1,t2,t3, &
                          tem,z0,z1,z2,bi,xx,xl,y0,y1,y2,x0, &
                          zero,p1,one,four,p5
      
      double complex, external :: f
      integer :: ickmax,kn,nguess,n,itmax,ier,knp1,knpn,i,l,ic, &
                    knpng,jk,ick,nn,lm1,errcode
      double complex :: x(kn+n)
      integer :: infer(kn+n)
      
      
      data                zero/(0.0,0.0)/,p1/(0.1,0.0)/, &
                          one/(1.0,0.0)/,four/(4.0,0.0)/, &
                          p5/(0.5,0.0)/, &
                          rzero/0.0/,rten/10.0/,rhun/100.0/, &
                          ax/0.1/,ickmax/3/,rp01/0.01/

            ier = 0
            if (n .lt. 1) then ! What the hell are doing here then ...
                return
            end if
            !eps1 = rten **(-nsig)
            eps1 = min(eps1,rp01)
            
            knp1 = kn+1
            knpn = kn+n
            knpng = kn+nguess
            do i=1,knpn
                infer(i) = 0
                if (i .gt. knpng) x(i) = zero
            end do
            l= knp1

            ic=0
rloop:      do while (l<=knpn)   ! Main loop over new roots
                jk = 0
                ick = 0
                xl = x(l)
icloop:         do
                    ic = 0
                    h = ax
                    h = p1*h
                    if (cdabs(xl) .gt. ax) h = p1*xl
!                                  first three points are
!                                    xl+h,  xl-h,  xl
                    rt = xl+h
                    call deflated_work(errcode)
                    if (errcode == 1) then
                        exit icloop
                    end if

                    z0 = fprt
                    y0 = frt
                    x0 = rt
                    rt = xl-h
                    call deflated_work(errcode)
                    if (errcode == 1) then
                        exit icloop
                    end if

                    z1 = fprt
                    y1 = frt
                    h = xl-rt
                    d = h/(rt-x0)
                    rt = xl

                    call deflated_work(errcode)
                    if (errcode == 1) then
                        exit icloop
                    end if

   
                    z2 = fprt
                    y2 = frt
!                                  begin main algorithm
 iloop:             do
                        dd = one + d
                        t1 = z0*d*d
                        t2 = z1*dd*dd
                        xx = z2*dd
                        t3 = z2*d
                        bi = t1-t2+xx+t3
                        den = bi*bi-four*(xx*t1-t3*(t2-xx))
!                                  use denominator of maximum amplitude 
                        t1 = cdsqrt(den)
                        qz = rhun*max(cdabs(bi),cdabs(t1))
                        t2 = bi + t1
                        tpq = cdabs(t2)+qz
                        if (tpq .eq. qz) t2 = zero
                        t3 = bi - t1
                        tpq = cdabs(t3) + qz
                        if (tpq .eq. qz) t3 = zero
                        den = t2
                        qz = cdabs(t3)-cdabs(t2)
                        if (qz .gt. rzero) den = t3
!                                  test for zero denominator            
                        if (cdabs(den) .eq. rzero) then
                            call trans_rt()
                            call deflated_work(errcode)
                            if (errcode == 1) then
                                exit icloop
                            end if
                            z2 = fprt
                            y2 = frt
                            cycle iloop
                        end if


                        d = -xx/den
                        d = d+d
                        h = d*h
                        rt = rt + h
!                                  check convergence of the first kind  
                        if (cdabs(h) .le. eps1*max(cdabs(rt),ax)) then
                            if (ic .ne. 0) then
                                exit icloop
                            end if
                            ic = 1
                            z0 = y1
                            z1 = y2
                            z2 = f(rt)
                            xl = rt
                            ick = ick+1
                            if (ick .le. ickmax) then
                                cycle iloop 
                            end if
!                                  warning error, itmax = maximum
                            jk = itmax + jk
                            ier = 33
                        end if
                        if (ic .ne. 0) then
                            cycle icloop
                        end if
                        call deflated_work(errcode)
                        if (errcode == 1) then
                            exit icloop
                        end if

                        do while ( (cdabs(fprt)-cdabs(z2)*rten) .ge. rzero)
                            !   take remedial action to induce
                            !   convergence
                            d = d*p5
                            h = h*p5
                            rt = rt-h
                            call deflated_work(errcode)
                            if (errcode == 1) then
                                exit icloop
                            end if
                        end do
                        z0 = z1
                        z1 = z2
                        z2 = fprt
                        y0 = y1
                        y1 = y2
                        y2 = frt
                    end do iloop
                end do icloop
        x(l) = rt
        infer(l) = jk
        l = l+1
      end do rloop
      
      contains
        subroutine trans_rt()
           tem = rten*eps1
           if (cdabs(rt) .gt. ax) tem = tem*rt
           rt = rt+tem
           d = (h+tem)*d/h
           h = h+tem
        end subroutine trans_rt
        
        subroutine deflated_work(errcode)
            ! errcode=0 => no errors
            ! errcode=1 => jk>itmax or convergence of second kind achieved
            integer :: errcode,flag
            
            flag=1
    loop1:  do while(flag==1)
                errcode=0
                jk = jk+1
                if (jk .gt. itmax) then
                    ier=33
                    errcode=1
                    return
                end if
                frt = f(rt)
                fprt = frt
                if (l /= 1) then
                    lm1 = l-1
                    do i=1,lm1
                        tem = rt - x(i)
                        if (cdabs(tem) .eq. rzero) then
                        !if (ic .ne. 0) go to 15 !! ?? possible?
                            call trans_rt()
                            cycle loop1
                        end if
                        fprt = fprt/tem
                    end do
                end if
                flag=0
            end do loop1
 
            if (cdabs(fprt) .le. eps .and. cdabs(frt) .le. eps) then
                errcode=1
                return
            end if
            
        end subroutine deflated_work
      
      end subroutine


!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
!   Integration
!
!   Only double precision coded atm
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!


subroutine fvn_d_gauss_legendre(n,qx,qw)
!
! This routine compute the n Gauss Legendre abscissas and weights
! Adapted from Numerical Recipes routine gauleg
!
! n (in) : number of points
! qx(out) : abscissas
! qw(out) : weights
!
implicit none
double precision,parameter :: pi=3.141592653589793d0
integer, intent(in) :: n
double precision, intent(out) :: qx(n),qw(n)

integer :: m,i,j
double precision :: z,z1,p1,p2,p3,pp
m=(n+1)/2
do i=1,m
    z=cos(pi*(dble(i)-0.25d0)/(dble(n)+0.5d0))
iloop:  do 
            p1=1.d0
            p2=0.d0
            do j=1,n
                p3=p2
                p2=p1
                p1=((2.d0*dble(j)-1.d0)*z*p2-(dble(j)-1.d0)*p3)/dble(j)
            end do
            pp=dble(n)*(z*p1-p2)/(z*z-1.d0)
            z1=z
            z=z1-p1/pp
            if (dabs(z-z1)<=epsilon(z)) then
                exit iloop
            end if
        end do iloop
    qx(i)=-z
    qx(n+1-i)=z
    qw(i)=2.d0/((1.d0-z*z)*pp*pp)
    qw(n+1-i)=qw(i)
end do
end subroutine



subroutine fvn_d_gl_integ(f,a,b,n,res)
!
! This is a simple non adaptative integration routine 
! using n gauss legendre abscissas and weights
!
!   f(in)   : the function to integrate
!   a(in)   : lower bound
!   b(in)   : higher bound
!   n(in)   : number of gauss legendre pairs
!   res(out): the evaluation of the integral
!
double precision,external :: f
double precision, intent(in) :: a,b
integer, intent(in):: n
double precision, intent(out) :: res

double precision, allocatable :: qx(:),qw(:)
double precision :: xm,xr
integer :: i

! First compute n gauss legendre abs and weight
allocate(qx(n))
allocate(qw(n))
call fvn_d_gauss_legendre(n,qx,qw)

xm=0.5d0*(b+a)
xr=0.5d0*(b-a)

res=0.d0

do i=1,n
    res=res+qw(i)*f(xm+xr*qx(i))
end do

res=xr*res

deallocate(qw)
deallocate(qx)

end subroutine

!!!!!!!!!!!!!!!!!!!!!!!!
!
! Simple and double adaptative Gauss Kronrod integration based on
! a modified version of quadpack ( http://www.netlib.org/quadpack
!
! Common parameters :
!
!       key (in)
!       epsabs
!       epsrel
!
!
!!!!!!!!!!!!!!!!!!!!!!!!

subroutine fvn_d_integ_1_gk(f,a,b,epsabs,epsrel,key,res,abserr,ier,limit)
!
! Evaluate the integral of function f(x) between a and b
!
! f(in) : the function
! a(in) : lower bound
! b(in) : higher bound
! epsabs(in) : desired absolute error
! epsrel(in) : desired relative error
! key(in) : gauss kronrod rule
!                     1:   7 - 15 points
!                     2:  10 - 21 points
!                     3:  15 - 31 points
!                     4:  20 - 41 points
!                     5:  25 - 51 points
!                     6:  30 - 61 points
!
! limit(in) : maximum number of subintervals in the partition of the 
!               given integration interval (a,b). A value of 500 will give the same
!               behaviour as the imsl routine dqdag
!
! res(out) : estimated integral value
! abserr(out) : estimated absolute error
! ier(out) : error flag from quadpack routines
!               0 : no error
!               1 : maximum number of subdivisions allowed
!                   has been achieved. one can allow more
!                   subdivisions by increasing the value of
!                   limit (and taking the according dimension
!                   adjustments into account). however, if
!                   this yield no improvement it is advised
!                   to analyze the integrand in order to
!                   determine the integration difficulaties.
!                   if the position of a local difficulty can
!                   be determined (i.e.singularity,
!                   discontinuity within the interval) one
!                   will probably gain from splitting up the
!                   interval at this point and calling the
!                   integrator on the subranges. if possible,
!                   an appropriate special-purpose integrator
!                   should be used which is designed for
!                   handling the type of difficulty involved.
!               2 : the occurrence of roundoff error is
!                   detected, which prevents the requested
!                   tolerance from being achieved.
!               3 : extremely bad integrand behaviour occurs
!                   at some points of the integration
!                   interval.
!               6 : the input is invalid, because
!                   (epsabs.le.0 and
!                   epsrel.lt.max(50*rel.mach.acc.,0.5d-28))
!                   or limit.lt.1 or lenw.lt.limit*4.
!                   result, abserr, neval, last are set
!                   to zero.
!                   except when lenw is invalid, iwork(1),
!                   work(limit*2+1) and work(limit*3+1) are
!                   set to zero, work(1) is set to a and
!                   work(limit+1) to b.

implicit none
double precision, external :: f
double precision, intent(in) :: a,b,epsabs,epsrel
integer, intent(in) :: key
integer, intent(in) :: limit
double precision, intent(out) :: res,abserr
integer, intent(out) :: ier

double precision, allocatable :: work(:)
integer, allocatable :: iwork(:)
integer :: lenw,neval,last

! imsl value for limit is 500
lenw=limit*4

allocate(iwork(limit))
allocate(work(lenw))

call dqag(f,a,b,epsabs,epsrel,key,res,abserr,neval,ier,limit,lenw,last,iwork,work)

deallocate(work)
deallocate(iwork)

end subroutine



subroutine fvn_d_integ_2_gk(f,a,b,g,h,epsabs,epsrel,key,res,abserr,ier,limit)
!
! Evaluate the double integral of function f(x,y) for x between a and b and y between g(x) and h(x)
!
! f(in) : the function
! a(in) : lower bound
! b(in) : higher bound
! g(in) : external function describing lower bound for y
! h(in) : external function describing higher bound for y
! epsabs(in) : desired absolute error
! epsrel(in) : desired relative error
! key(in) : gauss kronrod rule
!                     1:   7 - 15 points
!                     2:  10 - 21 points
!                     3:  15 - 31 points
!                     4:  20 - 41 points
!                     5:  25 - 51 points
!                     6:  30 - 61 points
!
! limit(in) : maximum number of subintervals in the partition of the 
!               given integration interval (a,b). A value of 500 will give the same
!               behaviour as the imsl routine dqdag
!
! res(out) : estimated integral value
! abserr(out) : estimated absolute error
! ier(out) : error flag from quadpack routines
!               0 : no error
!               1 : maximum number of subdivisions allowed
!                   has been achieved. one can allow more
!                   subdivisions by increasing the value of
!                   limit (and taking the according dimension
!                   adjustments into account). however, if
!                   this yield no improvement it is advised
!                   to analyze the integrand in order to
!                   determine the integration difficulaties.
!                   if the position of a local difficulty can
!                   be determined (i.e.singularity,
!                   discontinuity within the interval) one
!                   will probably gain from splitting up the
!                   interval at this point and calling the
!                   integrator on the subranges. if possible,
!                   an appropriate special-purpose integrator
!                   should be used which is designed for
!                   handling the type of difficulty involved.
!               2 : the occurrence of roundoff error is
!                   detected, which prevents the requested
!                   tolerance from being achieved.
!               3 : extremely bad integrand behaviour occurs
!                   at some points of the integration
!                   interval.
!               6 : the input is invalid, because
!                   (epsabs.le.0 and
!                   epsrel.lt.max(50*rel.mach.acc.,0.5d-28))
!                   or limit.lt.1 or lenw.lt.limit*4.
!                   result, abserr, neval, last are set
!                   to zero.
!                   except when lenw is invalid, iwork(1),
!                   work(limit*2+1) and work(limit*3+1) are
!                   set to zero, work(1) is set to a and
!                   work(limit+1) to b.

implicit none
double precision, external:: f,g,h
double precision, intent(in) :: a,b,epsabs,epsrel
integer, intent(in) :: key,limit
integer, intent(out) :: ier
double precision, intent(out) :: res,abserr


double precision, allocatable :: work(:)
integer, allocatable :: iwork(:)
integer :: lenw,neval,last

! imsl value for limit is 500
lenw=limit*4
allocate(work(lenw))
allocate(iwork(limit))

call dqag_2d_outer(f,a,b,g,h,epsabs,epsrel,key,res,abserr,neval,ier,limit,lenw,last,iwork,work)

deallocate(iwork)
deallocate(work)
end subroutine



subroutine fvn_d_integ_2_inner_gk(f,x,a,b,epsabs,epsrel,key,res,abserr,ier,limit)
!
! Evaluate the single integral of function f(x,y) for y between a and b with a
! given x value
!
! This function is used for the evaluation of the double integral fvn_d_integ_2_gk
!
! f(in) : the function
! x(in) : x
! a(in) : lower bound
! b(in) : higher bound
! epsabs(in) : desired absolute error
! epsrel(in) : desired relative error
! key(in) : gauss kronrod rule
!                     1:   7 - 15 points
!                     2:  10 - 21 points
!                     3:  15 - 31 points
!                     4:  20 - 41 points
!                     5:  25 - 51 points
!                     6:  30 - 61 points
!
! limit(in) : maximum number of subintervals in the partition of the 
!               given integration interval (a,b). A value of 500 will give the same
!               behaviour as the imsl routine dqdag
!
! res(out) : estimated integral value
! abserr(out) : estimated absolute error
! ier(out) : error flag from quadpack routines
!               0 : no error
!               1 : maximum number of subdivisions allowed
!                   has been achieved. one can allow more
!                   subdivisions by increasing the value of
!                   limit (and taking the according dimension
!                   adjustments into account). however, if
!                   this yield no improvement it is advised
!                   to analyze the integrand in order to
!                   determine the integration difficulaties.
!                   if the position of a local difficulty can
!                   be determined (i.e.singularity,
!                   discontinuity within the interval) one
!                   will probably gain from splitting up the
!                   interval at this point and calling the
!                   integrator on the subranges. if possible,
!                   an appropriate special-purpose integrator
!                   should be used which is designed for
!                   handling the type of difficulty involved.
!               2 : the occurrence of roundoff error is
!                   detected, which prevents the requested
!                   tolerance from being achieved.
!               3 : extremely bad integrand behaviour occurs
!                   at some points of the integration
!                   interval.
!               6 : the input is invalid, because
!                   (epsabs.le.0 and
!                   epsrel.lt.max(50*rel.mach.acc.,0.5d-28))
!                   or limit.lt.1 or lenw.lt.limit*4.
!                   result, abserr, neval, last are set
!                   to zero.
!                   except when lenw is invalid, iwork(1),
!                   work(limit*2+1) and work(limit*3+1) are
!                   set to zero, work(1) is set to a and
!                   work(limit+1) to b.

implicit none
double precision, external:: f
double precision, intent(in) :: x,a,b,epsabs,epsrel
integer, intent(in) :: key,limit
integer, intent(out) :: ier
double precision, intent(out) :: res,abserr


double precision, allocatable :: work(:)
integer, allocatable :: iwork(:)
integer :: lenw,neval,last

! imsl value for limit is 500
lenw=limit*4
allocate(work(lenw))
allocate(iwork(limit))

call dqag_2d_inner(f,x,a,b,epsabs,epsrel,key,res,abserr,neval,ier,limit,lenw,last,iwork,work)

deallocate(iwork)
deallocate(work)
end subroutine

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Include the modified quadpack files
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
include "fvn_quadpack/dqag_2d_inner.f"
include "fvn_quadpack/dqk15_2d_inner.f"
include "fvn_quadpack/dqk31_2d_outer.f"
include "fvn_quadpack/d1mach.f"
include "fvn_quadpack/dqk31_2d_inner.f"
include "fvn_quadpack/dqage.f"
include "fvn_quadpack/dqk15.f"
include "fvn_quadpack/dqk21.f"
include "fvn_quadpack/dqk31.f"
include "fvn_quadpack/dqk41.f"
include "fvn_quadpack/dqk51.f"
include "fvn_quadpack/dqk61.f"
include "fvn_quadpack/dqk41_2d_outer.f"
include "fvn_quadpack/dqk41_2d_inner.f"
include "fvn_quadpack/dqag_2d_outer.f"
include "fvn_quadpack/dqpsrt.f"
include "fvn_quadpack/dqag.f"
include "fvn_quadpack/dqage_2d_outer.f"
include "fvn_quadpack/dqage_2d_inner.f"
include "fvn_quadpack/dqk51_2d_outer.f"
include "fvn_quadpack/dqk51_2d_inner.f"
include "fvn_quadpack/dqk61_2d_outer.f"
include "fvn_quadpack/dqk21_2d_outer.f"
include "fvn_quadpack/dqk61_2d_inner.f"
include "fvn_quadpack/dqk21_2d_inner.f"
include "fvn_quadpack/dqk15_2d_outer.f"





end module fvn