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
|
!%f90 -*- f90 -*-
! Signatures for f2py wrappers of FORTRAN LAPACK functions.
!
! Author: Pearu Peterson
! Created: Jan-Feb 2002
! $Revision: 1.10 $ $Date: 2004/10/04 00:49:04 $
!
! Additions by Travis Oliphant
!
! Usage:
! f2py -c flapack.pyf -L/usr/local/lib/atlas -llapack -lf77blas -lcblas -latlas -lg2c
python module generic_flapack
interface
subroutine <tchar=s,d,c,z>gebal(scale,permute,n,a,m,lo,hi,pivscale,info)
!
! ba,lo,hi,pivscale,info = gebal(a,scale=0,permute=0,overwrite_a=0)
! Balance general matrix a.
! hi,lo are such that ba[i][j]==0 if i>j and j=0...lo-1 or i=hi+1..n-1
! pivscale([0:lo], [lo:hi+1], [hi:n+1]) = (p1,d,p2) where (p1,p2)[j] is
! the index of the row and column interchanged with row and column j.
! d[j] is the scaling factor applied to row and column j.
! The order in which the interchanges are made is n-1 to hi+1, then 0 to lo-1.
!
! P * A * P = [[T1,X,Y],[0,B,Z],[0,0,T2]]
! BA = [[T1,X*D,Y],[0,inv(D)*B*D,ind(D)*Z],[0,0,T2]]
! where D = diag(d), T1,T2 are upper triangular matrices.
! lo,hi mark the starting and ending columns of submatrix B.
callstatement { (*f2py_func)((permute?(scale?"B":"P"):(scale?"S":"N")),&n,a,&m,&lo,&hi,pivscale,&info); hi--; lo--; }
callprotoargument char*,int*,<type_in_c>*,int*,int*,int*,<rtype_in_c>*,int*
integer intent(in),optional :: permute = 0
integer intent(in),optional :: scale = 0
integer intent(hide),depend(a,n) :: m = shape(a,0)
integer intent(hide),depend(a) :: n = shape(a,1)
check(m>=n) m
integer intent(out) :: hi,lo
<rtype_in> dimension(n),intent(out),depend(n) :: pivscale
<type_in> dimension(m,n),intent(in,out,copy,out=ba) :: a
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>gebal
subroutine <tchar=s,d,c,z>gehrd(n,lo,hi,a,tau,work,lwork,info)
!
! hq,tau,info = gehrd(a,lo=0,hi=n-1,lwork=n,overwrite_a=0)
! Reduce general matrix A to upper Hessenberg form H by unitary similarity
! transform Q^H * A * Q = H
!
! Q = H(lo) * H(lo+1) * ... * H(hi-1)
! H(i) = I - tau * v * v^H
! v[0:i+1] = 0, v[i+1]=1, v[hi+1:n] = 0
! v[i+2:hi+1] is stored in hq[i+2:hi+i,i]
! tau is tau[i]
!
! hq for n=7,lo=1,hi=5:
! [a a h h h h a
! a h h h h a
! h h h h h h
! v2h h h h h
! v2v3h h h h
! v2v3v4h h h
! a]
!
callstatement { hi++; lo++; (*f2py_func)(&n,&lo,&hi,a,&n,tau,work,&lwork,&info); }
callprotoargument int*,int*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,int*,int*
integer intent(hide),depend(a) :: n = shape(a,0)
<type_in> dimension(n,n),intent(in,out,copy,out=ht),check(shape(a,0)==shape(a,1)) :: a
integer intent(in),optional :: lo = 0
integer intent(in),optional,depend(n) :: hi = n-1
<type_in> dimension(n-1),intent(out),depend(n) :: tau
<type_in> dimension(lwork),intent(cahce,hide),depend(lwork) :: work
integer intent(in),optional,depend(n),check(lwork>=MAX(n,1)) :: lwork = MAX(n,1)
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>gehrd
subroutine <tchar=s,d,c,z>gbsv(n,kl,ku,nrhs,ab,piv,b,info)
!
! lub,piv,x,info = gbsv(kl,ku,ab,b,overwrite_ab=0,overwrite_b=0)
! Solve A * X = B
! A = P * L * U
! A is a band matrix of order n with kl subdiagonals and ku superdiagonals
! starting at kl-th row.
! X, B are n-by-nrhs matrices
!
callstatement {int i=2*kl+ku+1;(*f2py_func)(&n,&kl,&ku,&nrhs,ab,&i,piv,b,&n,&info);for(i=0;i<n;--piv[i++]);}
callprotoargument int*,int*,int*,int*,<type_in_c>*,int*,int*,<type_in_c>*,int*,int*
integer depend(ab),intent(hide):: n = shape(ab,1)
integer intent(in) :: kl
integer intent(in) :: ku
integer depend(b),intent(hide) :: nrhs = shape(b,1)
<type_in> dimension(2*kl+ku+1,n),depend(kl,ku), check(2*kl+ku+1==shape(ab,0)) :: ab
integer dimension(n),depend(n),intent(out) :: piv
<type_in> dimension(n,nrhs),depend(n),check(shape(ab,1)==shape(b,0)) :: b
integer intent(out) :: info
intent(in,out,copy,out=x) b
intent(in,out,copy,out=lub) ab
end subroutine <tchar=s,d,c,z>gbsv
subroutine <tchar=s,d,c,z>gesv(n,nrhs,a,piv,b,info)
! lu,piv,x,info = gesv(a,b,overwrite_a=0,overwrite_b=0)
! Solve A * X = B.
! A = P * L * U
! U is upper diagonal triangular, L is unit lower triangular,
! piv pivots columns.
callstatement {int i;(*f2py_func)(&n,&nrhs,a,&n,piv,b,&n,&info);for(i=0;i<n;--piv[i++]);}
callprotoargument int*,int*,<type_in_c>*,int*,int*,<type_in_c>*,int*,int*
integer depend(a),intent(hide):: n = shape(a,0)
integer depend(b),intent(hide):: nrhs = shape(b,1)
<type_in> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
integer dimension(n),depend(n),intent(out) :: piv
<type_in> dimension(n,nrhs),check(shape(a,0)==shape(b,0)),depend(n) :: b
integer intent(out)::info
intent(in,out,copy,out=x) b
intent(in,out,copy,out=lu) a
end subroutine <tchar=s,d,c,z>gesv
subroutine <tchar=s,d,c,z>getrf(m,n,a,piv,info)
! lu,piv,info = getrf(a,overwrite_a=0)
! Compute an LU factorization of a general M-by-N matrix A.
! A = P * L * U
threadsafe
callstatement {int i;(*f2py_func)(&m,&n,a,&m,piv,&info);for(i=0,n=MIN(m,n);i<n;--piv[i++]);}
callprotoargument int*,int*,<type_in_c>*,int*,int*,int*
integer depend(a),intent(hide):: m = shape(a,0)
integer depend(a),intent(hide):: n = shape(a,1)
<type_in> dimension(m,n),intent(in,out,copy,out=lu) :: a
integer dimension(MIN(m,n)),depend(m,n),intent(out) :: piv
integer intent(out):: info
end subroutine <tchar=s,d,c,z>getrf
subroutine <tchar=s,d,c,z>getrs(n,nrhs,lu,piv,b,info,trans)
! x,info = getrs(lu,piv,b,trans=0,overwrite_b=0)
! Solve A * X = B if trans=0
! Solve A^T * X = B if trans=1
! Solve A^H * X = B if trans=2
! A = P * L * U
threadsafe
callstatement {int i;for(i=0;i<n;++piv[i++]);(*f2py_func)((trans?(trans==2?"C":"T"):"N"),&n,&nrhs,lu,&n,piv,b,&n,&info);for(i=0;i<n;--piv[i++]);}
callprotoargument char*,int*,int*,<type_in_c>*,int*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(trans>=0 && trans <=2) :: trans = 0
integer depend(lu),intent(hide):: n = shape(lu,0)
integer depend(b),intent(hide):: nrhs = shape(b,1)
<type_in> dimension(n,n),intent(in) :: lu
check(shape(lu,0)==shape(lu,1)) :: lu
integer dimension(n),intent(in),depend(n) :: piv
<type_in> dimension(n,nrhs),intent(in,out,copy,out=x),depend(n),check(shape(lu,0)==shape(b,0)) :: b
integer intent(out):: info
end subroutine <tchar=s,d,c,z>getrs
subroutine <tchar=s,d,c,z>getri(n,lu,piv,work,lwork,info)
! inv_a,info = getri(lu,piv,lwork=3*n,overwrite_lu=0)
! Find A inverse A^-1.
! A = P * L * U
callstatement {int i;for(i=0;i<n;++piv[i++]);(*f2py_func)(&n,lu,&n,piv,work,&lwork,&info);for(i=0;i<n;--piv[i++]);}
callprotoargument int*,<type_in_c>*,int*,int*,<type_in_c>*,int*,int*
integer depend(lu),intent(hide):: n = shape(lu,0)
<type_in> dimension(n,n),intent(in,out,copy,out=inv_a) :: lu
check(shape(lu,0)==shape(lu,1)) :: lu
integer dimension(n),intent(in),depend(n) :: piv
integer intent(out):: info
integer optional,intent(in),depend(n),check(lwork>=n) :: lwork=3*n
<type_in> dimension(lwork),intent(hide,cache),depend(lwork) :: work
end subroutine <tchar=s,d,c,z>getri
subroutine <tchar=s,d>gesdd(m,n,minmn,du,dvt,a,compute_uv,u,s,vt,work,lwork,iwork,info)
! u,s,vt,info = gesdd(a,compute_uv=1,lwork=..,overwrite_a=0)
! Compute the singular value decomposition (SVD):
! A = U * SIGMA * transpose(V)
! A - M x N matrix
! U - M x M matrix
! SIGMA - M x N zero matrix with a main diagonal filled with min(M,N)
! singular values
! transpose(V) - N x N matrix
callstatement (*f2py_func)((compute_uv?"A":"N"),&m,&n,a,&m,s,u,&du,vt,&dvt,work,&lwork,iwork,&info)
callprotoargument char*,int*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,int*,int*
integer intent(in),optional,check(compute_uv==0||compute_uv==1):: compute_uv = 1
integer intent(hide),depend(a):: m = shape(a,0)
integer intent(hide),depend(a):: n = shape(a,1)
integer intent(hide),depend(m,n):: minmn = MIN(m,n)
integer intent(hide),depend(compute_uv,minmn) :: du = (compute_uv?m:1)
integer intent(hide),depend(compute_uv,n) :: dvt = (compute_uv?n:1)
<type_in> dimension(m,n),intent(in,copy) :: a
<type_in> dimension(minmn),intent(out),depend(minmn) :: s
<type_in> dimension(du,du),intent(out),depend(du) :: u
<type_in> dimension(dvt,dvt),intent(out),depend(dvt) :: vt
<type_in> dimension(lwork),intent(hide,cache),depend(lwork) :: work
integer optional,intent(in),depend(minmn,compute_uv) &
:: lwork = (compute_uv?4*minmn*minmn+MAX(m,n)+9*minmn:MAX(14*minmn+4,10*minmn+2+25*(25+8))+MAX(m,n))
! gesdd docs are mess: optimal turns out to be less than minimal in docs
! check(lwork>=(compute_uv?3*minmn*minmn+MAX(MAX(m,n),4*minmn*(minmn+1)):MAX(14*minmn+4,10*minmn+2+25*(25+8))+MAX(m,n))) :: lwork
integer intent(hide,cache),dimension(8*minmn),depend(minmn) :: iwork
integer intent(out)::info
end subroutine <tchar=s,d>gesdd
subroutine <tchar=c,z>gesdd(m,n,minmn,du,dvt,a,compute_uv,u,s,vt,work,rwork,lwork,iwork,info)
! u,s,vt,info = gesdd(a,compute_uv=1,lwork=..,overwrite_a=0)
! Compute the singular value decomposition (SVD):
! A = U * SIGMA * conjugate-transpose(V)
! A - M x N matrix
! U - M x M matrix
! SIGMA - M x N zero matrix with a main diagonal filled with min(M,N)
! singular values
! conjugate-transpose(V) - N x N matrix
callstatement (*f2py_func)((compute_uv?"A":"N"),&m,&n,a,&m,s,u,&du,vt,&dvt,work,&lwork,rwork,iwork,&info)
callprotoargument char*,int*,int*,<type_in_c>*,int*,<rtype_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<rtype_in_c>*,int*,int*
integer intent(in),optional,check(compute_uv==0||compute_uv==1):: compute_uv = 1
integer intent(hide),depend(a):: m = shape(a,0)
integer intent(hide),depend(a):: n = shape(a,1)
integer intent(hide),depend(m,n):: minmn = MIN(m,n)
integer intent(hide),depend(compute_uv,minmn) :: du = (compute_uv?m:1)
integer intent(hide),depend(compute_uv,n) :: dvt = (compute_uv?n:1)
<type_in> dimension(m,n),intent(in,copy) :: a
<rtype_in> dimension(minmn),intent(out),depend(minmn) :: s
<type_in> dimension(du,du),intent(out),depend(du) :: u
<type_in> dimension(dvt,dvt),intent(out),depend(dvt) :: vt
<type_in> dimension(lwork),intent(hide,cache),depend(lwork) :: work
<rtype_in> dimension((compute_uv?5*minmn*minmn+7*minmn:5*minmn)),intent(hide,cache),depend(minmn,compute_uv) :: rwork
integer optional,intent(in),depend(minmn,compute_uv) &
:: lwork = (compute_uv?2*minmn*minmn+MAX(m,n)+2*minmn:2*minmn+MAX(m,n))
check(lwork>=(compute_uv?2*minmn*minmn+MAX(m,n)+2*minmn:2*minmn+MAX(m,n))) :: lwork
integer intent(hide,cache),dimension(8*minmn),depend(minmn) :: iwork
integer intent(out)::info
end subroutine <tchar=c,z>gesdd
subroutine <tchar=s,d>gelss(m,n,minmn,maxmn,nrhs,a,b,s,cond,r,work,lwork,info)
! v,x,s,rank,info = gelss(a,b,cond=-1.0,overwrite_a=0,overwrite_b=0)
! Solve Minimize 2-norm(A * X - B).
callstatement (*f2py_func)(&m,&n,&nrhs,a,&m,b,&maxmn,s,&cond,&r,work,&lwork,&info)
callprotoargument int*,int*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,int*
integer intent(hide),depend(a):: m = shape(a,0)
integer intent(hide),depend(a):: n = shape(a,1)
integer intent(hide),depend(m,n):: minmn = MIN(m,n)
integer intent(hide),depend(m,n):: maxmn = MAX(m,n)
<type_in> dimension(m,n),intent(in,out,copy,out=v) :: a
integer depend(b),intent(hide):: nrhs = shape(b,1)
<type_in> dimension(maxmn,nrhs),check(maxmn==shape(b,0)),depend(maxmn) :: b
intent(in,out,copy,out=x) b
integer intent(out)::info
integer optional,intent(in),depend(nrhs,minmn,maxmn), &
check(lwork>=1) &
:: lwork=3*minmn+MAX(2*minmn,MAX(maxmn,nrhs))
!check(lwork>=3*minmn+MAX(2*minmn,MAX(maxmn,nrhs)))
<type_in> dimension(lwork),intent(hide),depend(lwork) :: work
<type_in> intent(in),optional :: cond = -1.0
integer intent(out,out=rank) :: r
<type_in> intent(out),dimension(minmn),depend(minmn) :: s
end subroutine <tchar=s,d>gelss
subroutine <tchar=c,z>gelss(m,n,minmn,maxmn,nrhs,a,b,s,cond,r,work,rwork,lwork,info)
! v,x,s,rank,info = gelss(a,b,cond=-1.0,overwrite_a=0,overwrite_b=0)
! Solve Minimize 2-norm(A * X - B).
callstatement (*f2py_func)(&m,&n,&nrhs,a,&m,b,&maxmn,s,&cond,&r,work,&lwork,rwork,&info)
callprotoargument int*,int*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<rtype_in_c>*,<rtype_in_c>*,int*,<type_in_c>*,int*,<rtype_in_c>*,int*
integer intent(hide),depend(a):: m = shape(a,0)
integer intent(hide),depend(a):: n = shape(a,1)
integer intent(hide),depend(m,n):: minmn = MIN(m,n)
integer intent(hide),depend(m,n):: maxmn = MAX(m,n)
<type_in> dimension(m,n),intent(in,out,copy,out=v) :: a
integer depend(b),intent(hide):: nrhs = shape(b,1)
<type_in> dimension(maxmn,nrhs),check(maxmn==shape(b,0)),depend(maxmn) :: b
intent(in,out,copy,out=x) b
integer intent(out)::info
integer optional,intent(in),depend(nrhs,minmn,maxmn), &
check(lwork>=1) &
:: lwork=2*minmn+MAX(maxmn,nrhs)
! check(lwork>=2*minmn+MAX(maxmn,nrhs))
<type_in> dimension(lwork),intent(hide),depend(lwork) :: work
<rtype_in> dimension(5*minmn-1),intent(hide),depend(lwork) :: rwork
<rtype_in> intent(in),optional :: cond = -1.0
integer intent(out,out=rank) :: r
<rtype_in> intent(out),dimension(minmn),depend(minmn) :: s
end subroutine <tchar=c,z>gelss
subroutine <tchar=s,d,c,z>geqrf(m,n,a,tau,work,lwork,info)
! qr_a,tau,work,info = geqrf(a,lwork=3*n,overwrite_a=0)
! Compute a QR factorization of a real M-by-N matrix A:
! A = Q * R.
callstatement (*f2py_func)(&m,&n,a,&m,tau,work,&lwork,&info)
callprotoargument int*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,int*,int*
integer intent(hide),depend(a):: m = shape(a,0)
integer intent(hide),depend(a):: n = shape(a,1)
<type_in> dimension(m,n),intent(in,out,copy,out=qr) :: a
<type_in> dimension(MIN(m,n)),intent(out) :: tau
integer optional,intent(in),depend(n),check(lwork>=n||lwork==-1) :: lwork=3*n
<type_in> dimension(MAX(lwork,1)),intent(out),depend(lwork) :: work
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>geqrf
subroutine <tchar=s,d>geev(compute_vl,compute_vr,n,a,wr,wi,vl,ldvl,vr,ldvr,work,lwork,info)
! wr,wi,vl,vr,info = geev(a,compute_vl=1,compute_vr=1,lwork=4*n,overwrite_a=0)
callstatement {(*f2py_func)((compute_vl?"V":"N"),(compute_vr?"V":"N"),&n,a,&n,wr,wi,vl,&ldvl,vr,&ldvr,work,&lwork,&info);}
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,int*
integer optional,intent(in):: compute_vl = 1
check(compute_vl==1||compute_vl==0) compute_vl
integer optional,intent(in):: compute_vr = 1
check(compute_vr==1||compute_vr==0) compute_vr
integer intent(hide),depend(a) :: n = shape(a,0)
<type_in> dimension(n,n),intent(in,copy) :: a
check(shape(a,0)==shape(a,1)) :: a
<type_in> dimension(n),intent(out),depend(n) :: wr
<type_in> dimension(n),intent(out),depend(n) :: wi
<type_in> dimension(ldvl,n),intent(out) :: vl
integer intent(hide),depend(n,compute_vl) :: ldvl=(compute_vl?n:1)
<type_in> dimension(ldvr,n),intent(out) :: vr
integer intent(hide),depend(n,compute_vr) :: ldvr=(compute_vr?n:1)
integer optional,intent(in),depend(n,compute_vl,compute_vr) :: lwork=4*n
check(lwork>=((compute_vl||compute_vr)?4*n:3*n)) :: lwork
<type_in> dimension(lwork),intent(hide,cache),depend(lwork) :: work
integer intent(out):: info
end subroutine <tchar=s,d>geev
subroutine <tchar=c,z>geev(compute_vl,compute_vr,n,a,w,vl,ldvl,vr,ldvr,work,lwork,rwork,info)
! w,vl,vr,info = geev(a,compute_vl=1,compute_vr=1,lwork=2*n,overwrite_a=0)
callstatement (*f2py_func)((compute_vl?"V":"N"),(compute_vr?"V":"N"),&n,a,&n,w,vl,&ldvl,vr,&ldvr,work,&lwork,rwork,&info)
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<rtype_in_c>*,int*
integer optional,intent(in):: compute_vl = 1
check(compute_vl==1||compute_vl==0) compute_vl
integer optional,intent(in):: compute_vr = 1
check(compute_vr==1||compute_vr==0) compute_vr
integer intent(hide),depend(a) :: n = shape(a,0)
<type_in> dimension(n,n),intent(in,copy) :: a
check(shape(a,0)==shape(a,1)) :: a
<type_in> dimension(n),intent(out),depend(n) :: w
<type_in> dimension(ldvl,n),depend(ldvl),intent(out) :: vl
integer intent(hide),depend(compute_vl,n) :: ldvl=(compute_vl?n:1)
<type_in> dimension(ldvr,n),depend(ldvr),intent(out) :: vr
integer intent(hide),depend(compute_vr,n) :: ldvr=(compute_vr?n:1)
integer optional,intent(in),depend(n) :: lwork=2*n
check(lwork>=2*n) :: lwork
<type_in> dimension(lwork),intent(hide),depend(lwork) :: work
<rtype_in> dimension(2*n),intent(hide,cache),depend(n) :: rwork
integer intent(out):: info
end subroutine <tchar=c,z>geev
subroutine <tchar=s,d>gegv(compute_vl,compute_vr,n,a,b,alphar,alphai,beta,vl,ldvl,vr,ldvr,work,lwork,info)
! Compute the generalized eigenvalues (alphar +/- alphai*i, beta)
! of the real nonsymmetric matrices A and B: det(A-w*B)=0 where w=alpha/beta.
! Optionally, compute the left and/or right generalized eigenvectors:
! (A - w B) r = 0, l^H * (A - w B) = 0
!
! alphar,alphai,beta,vl,vr,info = gegv(a,b,compute_vl=1,compute_vr=1,lwork=8*n,overwrite_a=0,overwrite_b=0)
callstatement (*f2py_func)((compute_vl?"V":"N"),(compute_vr?"V":"N"),&n,a,&n,b,&n,alphar,alphai,beta,vl,&ldvl,vr,&ldvr,work,&lwork,&info)
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,int*
integer optional,intent(in):: compute_vl = 1
check(compute_vl==1||compute_vl==0) compute_vl
integer optional,intent(in):: compute_vr = 1
check(compute_vr==1||compute_vr==0) compute_vr
integer intent(hide),depend(a) :: n = shape(a,0)
<type_in> dimension(n,n),intent(in,copy) :: a
check(shape(a,0)==shape(a,1)) :: a
<type_in> dimension(n,n),depend(n),intent(in,copy) :: b
check(shape(b,0)==shape(b,1) && shape(b,0)==n) :: b
<type_in> dimension(n),depend(n),intent(out) :: alphar
<type_in> dimension(n),depend(n),intent(out) :: alphai
<type_in> dimension(n),depend(n),intent(out) :: beta
<type_in> dimension(ldvl,n),intent(out),depend(ldvl) :: vl
integer intent(hide),depend(compute_vl,n) :: ldvl=(compute_vl?n:1)
<type_in> dimension(ldvr,n),intent(out),depend(ldvr) :: vr
integer intent(hide),depend(compute_vr,n) :: ldvr=(compute_vr?n:1)
integer optional,intent(in),depend(n) :: lwork=8*n
check(lwork>=8*n) :: lwork
<type_in> dimension(lwork),intent(hide),depend(lwork) :: work
integer intent(out):: info
end subroutine <tchar=s,d>gegv
subroutine <tchar=c,z>gegv(compute_vl,compute_vr,n,a,b,alpha,beta,vl,ldvl,vr,ldvr,work,lwork,rwork,info)
! Compute the generalized eigenvalues (alpha, beta)
! of the comples nonsymmetric matrices A and B: det(A-w*B)=0 where w=alpha/beta.
! Optionally, compute the left and/or right generalized eigenvectors:
! (A - w B) r = 0, l^H * (A - w B) = 0
!
! alpha,beta,vl,vr,info = gegv(a,b,compute_vl=1,compute_vr=1,lwork=2*n,overwrite_a=0,overwrite_b=0)
callstatement (*f2py_func)((compute_vl?"V":"N"),(compute_vr?"V":"N"),&n,a,&n,b,&n,alpha,beta,vl,&ldvl,vr,&ldvr,work,&lwork,rwork,&info)
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<rtype_in_c>*,int*
integer optional,intent(in):: compute_vl = 1
check(compute_vl==1||compute_vl==0) compute_vl
integer optional,intent(in):: compute_vr = 1
check(compute_vr==1||compute_vr==0) compute_vr
integer intent(hide),depend(a) :: n = shape(a,0)
<type_in> dimension(n,n),intent(in,copy) :: a
check(shape(a,0)==shape(a,1)) :: a
<type_in> dimension(n,n),depend(n),intent(in,copy) :: b
check(shape(b,0)==shape(b,1) && shape(b,0)==n) :: b
<type_in> dimension(n),depend(n),intent(out) :: alpha
<type_in> dimension(n),depend(n),intent(out) :: beta
<type_in> dimension(ldvl,n),intent(out),depend(ldvl) :: vl
integer intent(hide),depend(compute_vl,n) :: ldvl=(compute_vl?n:1)
<type_in> dimension(ldvr,n),intent(out),depend(ldvr) :: vr
integer intent(hide),depend(compute_vr,n) :: ldvr=(compute_vr?n:1)
integer optional,intent(in),depend(n) :: lwork=2*n
check(lwork>=2*n) :: lwork
<type_in> dimension(lwork),intent(hide),depend(lwork) :: work
<rtype_in> dimension(8*n),intent(hide),depend(n) :: rwork
integer intent(out):: info
end subroutine <tchar=c,z>gegv
subroutine <tchar=s,d>syev(compute_v,lower,n,w,a,work,lwork,info)
! w,v,info = syev(a,compute_v=1,lower=0,lwork=3*n-1,overwrite_a=0)
! Compute all eigenvalues and, optionally, eigenvectors of a
! real symmetric matrix A.
!
! Performance tip:
! If compute_v=0 then set also overwrite_a=1.
callstatement (*f2py_func)((compute_v?"V":"N"),(lower?"L":"U"),&n,a,&n,w,work,&lwork,&info)
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,int*,int*
integer optional,intent(in):: compute_v = 1
check(compute_v==1||compute_v==0) compute_v
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer intent(hide),depend(a):: n = shape(a,0)
<type_in> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
intent(in,copy,out,out=v) :: a
<type_in> dimension(n),intent(out),depend(n) :: w
integer optional,intent(in),depend(n) :: lwork=3*n-1
check(lwork>=3*n-1) :: lwork
<type_in> dimension(lwork),intent(hide),depend(lwork) :: work
integer intent(out) :: info
end subroutine <tchar=s,d>syev
subroutine <tchar=c,z>heev(compute_v,lower,n,w,a,work,lwork,rwork,info)
! w,v,info = syev(a,compute_v=1,lower=0,lwork=3*n-1,overwrite_a=0)
! Compute all eigenvalues and, optionally, eigenvectors of a
! complex Hermitian matrix A.
!
! Warning:
! If compute_v=0 and overwrite_a=1, the contents of a is destroyed.
callstatement (*f2py_func)((compute_v?"V":"N"),(lower?"L":"U"),&n,a,&n,w,work,&lwork,rwork,&info)
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,int*,<rtype_in_c>*,int*
integer optional,intent(in):: compute_v = 1
check(compute_v==1||compute_v==0) compute_v
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer intent(hide),depend(a):: n = shape(a,0)
<type_in> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
intent(in,copy,out,out=v) :: a
<type_in> dimension(n),intent(out),depend(n) :: w
integer optional,intent(in),depend(n) :: lwork=2*n-1
check(lwork>=2*n-1) :: lwork
<type_in> dimension(lwork),intent(hide),depend(lwork) :: work
<rtype_in> dimension(3*n-1),intent(hide),depend(n) :: rwork
integer intent(out) :: info
end subroutine <tchar=c,z>heev
subroutine <tchar=s,d,c,z>posv(n,nrhs,a,b,info,lower)
! c,x,info = posv(a,b,lower=0,overwrite_a=0,overwrite_b=0)
! Solve A * X = B.
! A is symmetric positive defined
! A = U^T * U, C = U if lower = 0
! A = L * L^T, C = L if lower = 1
! C is triangular matrix of the corresponding Cholesky decomposition.
callstatement (*f2py_func)((lower?"L":"U"),&n,&nrhs,a,&n,b,&n,&info)
callprotoargument char*,int*,int*,<type_in_c>*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer depend(a),intent(hide):: n = shape(a,0)
integer depend(b),intent(hide):: nrhs = shape(b,1)
<type_in> dimension(n,n),intent(in,out,copy,out=c) :: a
check(shape(a,0)==shape(a,1)) :: a
<type_in> dimension(n,nrhs),intent(in,out,copy,out=x),depend(n):: b
check(shape(a,0)==shape(b,0)) :: b
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>posv
subroutine <tchar=s,d>potrf(n,a,info,lower,clean)
! c,info = potrf(a,lower=0,clean=1,overwrite_a=0)
! Compute Cholesky decomposition of symmetric positive defined matrix:
! A = U^T * U, C = U if lower = 0
! A = L * L^T, C = L if lower = 1
! C is triangular matrix of the corresponding Cholesky decomposition.
! clean==1 zeros strictly lower or upper parts of U or L, respectively
callstatement (*f2py_func)((lower?"L":"U"),&n,a,&n,&info); if(clean){int i,j;if(lower){for(i=0;i<n;++i) for(j=i+1;j<n;++j) *(a+j*n+i)=0.0;} else {for(i=0;i<n;++i) for(j=i+1;j<n;++j) *(a+i*n+j)=0.0;}}
callprotoargument char*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer optional,intent(in),check(clean==0||clean==1) :: clean = 1
integer depend(a),intent(hide):: n = shape(a,0)
<type_in> dimension(n,n),intent(in,out,copy,out=c) :: a
check(shape(a,0)==shape(a,1)) :: a
integer intent(out) :: info
end subroutine <tchar=s,d>potrf
subroutine <tchar=c,z>potrf(n,a,info,lower,clean)
! c,info = potrf(a,lower=0,clean=1,overwrite_a=0)
! Compute Cholesky decomposition of symmetric positive defined matrix:
! A = U^H * U, C = U if lower = 0
! A = L * L^H, C = L if lower = 1
! C is triangular matrix of the corresponding Cholesky decomposition.
! clean==1 zeros strictly lower or upper parts of U or L, respectively
callstatement (*f2py_func)((lower?"L":"U"),&n,a,&n,&info); if(clean){int i,j,k;if(lower){for(i=0;i<n;++i) for(j=i+1;j<n;++j) {k=j*n+i;(a+k)->r=(a+k)->i=0.0;}} else {for(i=0;i<n;++i) for(j=i+1;j<n;++j) {k=i*n+j;(a+k)->r=(a+k)->i=0.0;}}}
callprotoargument char*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer optional,intent(in),check(clean==0||clean==1) :: clean = 1
integer depend(a),intent(hide):: n = shape(a,0)
<type_in> dimension(n,n),intent(in,out,copy,out=c) :: a
check(shape(a,0)==shape(a,1)) :: a
integer intent(out) :: info
end subroutine <tchar=c,z>potrf
subroutine <tchar=s,d,c,z>potrs(n,nrhs,c,b,info,lower)
! x,info = potrs(c,b,lower=0=1,overwrite_b=0)
! Solve A * X = B.
! A is symmetric positive defined
! A = U^T * U, C = U if lower = 0
! A = L * L^T, C = L if lower = 1
! C is triangular matrix of the corresponding Cholesky decomposition.
callstatement (*f2py_func)((lower?"L":"U"),&n,&nrhs,c,&n,b,&n,&info)
callprotoargument char*,int*,int*,<type_in_c>*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer depend(c),intent(hide):: n = shape(c,0)
integer depend(b),intent(hide):: nrhs = shape(b,1)
<type_in> dimension(n,n),intent(in) :: c
check(shape(c,0)==shape(c,1)) :: c
<type_in> dimension(n,nrhs),intent(in,out,copy,out=x),depend(n):: b
check(shape(c,0)==shape(b,0)) :: b
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>potrs
subroutine <tchar=s,d,c,z>potri(n,c,info,lower)
! inv_a,info = potri(c,lower=0,overwrite_c=0)
! Compute A inverse A^-1.
! A = U^T * U, C = U if lower = 0
! A = L * L^T, C = L if lower = 1
! C is triangular matrix of the corresponding Cholesky decomposition.
callstatement (*f2py_func)((lower?"L":"U"),&n,c,&n,&info)
callprotoargument char*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer depend(c),intent(hide):: n = shape(c,0)
<type_in> dimension(n,n),intent(c,in,out,copy,out=inv_a) :: c
check(shape(c,0)==shape(c,1)) :: c
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>potri
subroutine <tchar=s,d,c,z>lauum(n,c,info,lower)
! a,info = lauum(c,lower=0,overwrite_c=0)
! Compute product
! U^T * U, C = U if lower = 0
! L * L^T, C = L if lower = 1
! C is triangular matrix of the corresponding Cholesky decomposition.
callstatement (*f2py_func)((lower?"L":"U"),&n,c,&n,&info)
callprotoargument char*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer depend(c),intent(hide):: n = shape(c,0)
<type_in> dimension(n,n),intent(in,out,copy,out=a) :: c
check(shape(c,0)==shape(c,1)) :: c
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>lauum
subroutine <tchar=s,d,c,z>trtri(n,c,info,lower,unitdiag)
! inv_c,info = trtri(c,lower=0,unitdiag=1,overwrite_c=0)
! Compute C inverse C^-1 where
! C = U if lower = 0
! C = L if lower = 1
! C is non-unit triangular matrix if unitdiag = 0
! C is unit triangular matrix if unitdiag = 1
callstatement (*f2py_func)((lower?"L":"U"),(unitdiag?"U":"N"),&n,c,&n,&info)
callprotoargument char*,char*,int*,<type_in_c>*,int*,int*
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer optional,intent(in),check(unitdiag==0||unitdiag==1) :: unitdiag = 0
integer depend(c),intent(hide):: n = shape(c,0)
<type_in> dimension(n,n),intent(in,out,copy,out=inv_c) :: c
check(shape(c,0)==shape(c,1)) :: c
integer intent(out) :: info
end subroutine <tchar=s,d,c,z>trtri
subroutine <tchar=s,d,c,z>laswp(n,a,nrows,k1,k2,piv,off,inc,m)
! a = laswp(a,piv,k1=0,k2=len(piv)-1,off=0,inc=1,overwrite_a=0)
! Perform row interchanges on the matrix A for each of row k1 through k2
!
! piv pivots rows.
callstatement {int i;m=len(piv);for(i=0;i<m;++piv[i++]);++k1;++k2; (*f2py_func)(&n,a,&nrows,&k1,&k2,piv+off,&inc); for(i=0;i<m;--piv[i++]);}
callprotoargument int*,<type_in_c>*,int*,int*,int*,int*,int*
integer depend(a),intent(hide):: nrows = shape(a,0)
integer depend(a),intent(hide):: n = shape(a,1)
<type_in> dimension(nrows,n),intent(in,out,copy) :: a
integer dimension(*),intent(in),depend(nrows) :: piv
check(len(piv)<=nrows) :: piv
!XXX: how to check that all elements in piv are < n?
integer optional,intent(in) :: k1 = 0
check(0<=k1) :: k1
integer optional,intent(in),depend(k1,piv,off) :: k2 = len(piv)-1
check(k1<=k2 && k2<len(piv)-off) :: k2
integer optional, intent(in),check(inc>0||inc<0) :: inc = 1
integer optional,intent(in),depend(piv) :: off=0
check(off>=0 && off<len(piv)) :: off
integer intent(hide),depend(piv,inc,off) :: m = (len(piv)-off)/abs(inc)
check(len(piv)-off>(m-1)*abs(inc)) :: m
end subroutine <tchar=s,d,c,z>laswp
subroutine <tchar=c,z>gees(compute_v,sort_t,<tchar>select,n,a,nrows,sdim,w,vs,ldvs,work,lwork,rwork,bwork,info)
! t,sdim,w,vs,work,info=gees(compute_v=1,sort_t=0,select,a,lwork=3*n)
! For an NxN matrix compute the eigenvalues, the schur form T, and optionally
! the matrix of Schur vectors Z. This gives the Schur factorization
! A = Z * T * Z^H -- a complex matrix is in Schur form if it is upper
! triangular
callstatement (*f2py_func)((compute_v?"V":"N"),(sort_t?"S":"N"),cb_<tchar>select_in_<tchar>gees__user__routines,&n,a,&nrows,&sdim,w,vs,&ldvs,work,&lwork,rwork,bwork,&info,1,1)
callprotoargument char*,char*,int(*)(<type_in_c>*),int*,<type_in_c>*,int*,int*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<rtype_in_c>*,int*,int*,int,int
use <tchar>gees__user__routines
integer optional,intent(in),check(compute_v==0||compute_v==1) :: compute_v = 1
integer optional,intent(in),check(sort_t==0||sort_t==1) :: sort_t = 0
external <tchar>select
integer intent(hide),depend(a) :: n = shape(a,1)
<type_in> intent(in,out,copy,out=t),check(shape(a,0)==shape(a,1)),dimension(n,n) :: a
integer intent(hide),depend(a) :: nrows=shape(a,0)
integer intent(out) :: sdim=0
<type_in> intent(out),dimension(n) :: w
<type_in> intent(out),depend(ldvs,n),dimension(ldvs,n) :: vs
integer intent(hide),depend(compute_v,n) :: ldvs=((compute_v==1)?n:1)
<type_in> intent(out),depend(lwork),dimension(MAX(lwork,1)) :: work
integer optional,intent(in),check((lwork==-1)||(lwork >= MAX(1,2*n))),depend(n) :: lwork = 3*n
<rtype_in> optional,intent(hide),depend(n),dimension(n) :: rwork
logical optional,intent(hide),depend(n),dimension(n) :: bwork
integer intent(out) :: info
end subroutine <tchar=c,z>gees
subroutine <tchar=d,s>gees(compute_v,sort_t,<tchar>select,n,a,nrows,sdim,wr,wi,vs,ldvs,work,lwork,bwork,info)
! t,sdim,w,vs,work,info=gees(compute_v=1,sort_t=0,select,a,lwork=3*n)
! For an NxN matrix compute the eigenvalues, the schur form T, and optionally
! the matrix of Schur vectors Z. This gives the Schur factorization
! A = Z * T * Z^H -- a real matrix is in Schur form if it is upper quasi-
! triangular with 1x1 and 2x2 blocks.
callstatement (*f2py_func)((compute_v?"V":"N"),(sort_t?"S":"N"),cb_<tchar>select_in_<tchar>gees__user__routines,&n,a,&nrows,&sdim,wr,wi,vs,&ldvs,work,&lwork,bwork,&info,1,1)
callprotoargument char*,char*,int(*)(<type_in_c>*,<type_in_c>*),int*,<type_in_c>*,int*,int*,<type_in_c>*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,int*,int*,int,int
use <tchar>gees__user__routines
integer optional,intent(in),check(compute_v==0||compute_v==1) :: compute_v = 1
integer optional,intent(in),check(sort_t==0||sort_t==1) :: sort_t = 0
external <tchar>select
integer intent(hide),depend(a) :: n = shape(a,1)
<type_in> intent(in,out,copy,out=t),check(shape(a,0)==shape(a,1)),dimension(n,n) :: a
integer intent(hide),depend(a) :: nrows=shape(a,0)
integer intent(out) :: sdim=0
<type_in> intent(out),dimension(n) :: wr
<type_in> intent(out),dimension(n) :: wi
<type_in> intent(out),depend(ldvs,n),dimension(ldvs,n) :: vs
integer intent(hide),depend(compute_v,n) :: ldvs=((compute_v==1)?n:1)
<type_in> intent(out),depend(lwork),dimension(MAX(lwork,1)) :: work
integer optional,intent(in),check((lwork==-1)||(lwork >= MAX(1,2*n))),depend(n) :: lwork = 3*n
<rtype_in> optional,intent(hide),depend(n),dimension(n) :: rwork
logical optional,intent(hide),depend(n),dimension(n) :: bwork
integer intent(out) :: info
end subroutine <tchar=s,d>gees
subroutine <tchar=s,d>ggev(compute_vl,compute_vr,n,a,b,alphar,alphai,beta,vl,ldvl,vr,ldvr,work,lwork,info)
callstatement {(*f2py_func)((compute_vl?"V":"N"),(compute_vr?"V":"N"),&n,a,&n,b,&n,alphar,alphai,beta,vl,&ldvl,vr,&ldvr,work,&lwork,&info);}
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,int*
integer optional,intent(in):: compute_vl = 1
check(compute_vl==1||compute_vl==0) compute_vl
integer optional,intent(in):: compute_vr = 1
check(compute_vr==1||compute_vr==0) compute_vr
integer intent(hide),depend(a) :: n = shape(a,0)
<type_in> dimension(n,n),intent(in,copy) :: a
check(shape(a,0)==shape(a,1)) :: a
<type_in> intent(in,copy), dimension(n,n) :: b
check(shape(b,0)==shape(b,1)) :: b
<type_in> intent(out), dimension(n), depend(n) :: alphar
<type_in> intent(out), dimension(n), depend(n) :: alphai
<type_in> intent(out), dimension(n), depend(n) :: beta
<type_in> depend(ldvl,n), dimension(ldvl,n),intent(out) :: vl
integer intent(hide),depend(n,compute_vl) :: ldvl=(compute_vl?n:1)
<type_in> depend(ldvr,n), dimension(ldvr,n),intent(out) :: vr
integer intent(hide),depend(n,compute_vr) :: ldvr=(compute_vr?n:1)
integer optional,intent(in),depend(n,compute_vl,compute_vr) :: lwork=8*n
check((lwork==-1) || (lwork>=MAX(1,8*n))) :: lwork
<type_in> intent(out), dimension(MAX(lwork,1)), depend(lwork) :: work
integer intent(out):: info
end subroutine <tchar>ggev
subroutine <tchar=c,z>ggev(compute_vl,compute_vr,n,a,b,alpha,beta,vl,ldvl,vr,ldvr,work,lwork,rwork,info)
callstatement {(*f2py_func)((compute_vl?"V":"N"),(compute_vr?"V":"N"),&n,a,&n,b,&n,alpha,beta,vl,&ldvl,vr,&ldvr,work,&lwork,rwork,&info);}
callprotoargument char*,char*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,<type_in_c>*,<type_in_c>*,int*,<type_in_c>*,int*,<type_in_c>*,int*,<rtype_in_c>*,int*
integer optional,intent(in):: compute_vl = 1
check(compute_vl==1||compute_vl==0) compute_vl
integer optional,intent(in):: compute_vr = 1
check(compute_vr==1||compute_vr==0) compute_vr
integer intent(hide),depend(a) :: n = shape(a,0)
<type_in> dimension(n,n),intent(in,copy) :: a
check(shape(a,0)==shape(a,1)) :: a
<type_in> intent(in,copy), dimension(n,n) :: b
check(shape(b,0)==shape(b,1)) :: b
<type_in> intent(out), dimension(n), depend(n) :: alpha
<type_in> intent(out), dimension(n), depend(n) :: beta
<type_in> depend(ldvl,n), dimension(ldvl,n),intent(out) :: vl
integer intent(hide),depend(n,compute_vl) :: ldvl=(compute_vl?n:1)
<type_in> depend(ldvr,n), dimension(ldvr,n),intent(out) :: vr
integer intent(hide),depend(n,compute_vr) :: ldvr=(compute_vr?n:1)
integer optional,intent(in),depend(n,compute_vl,compute_vr) :: lwork=8*n
check((lwork==-1) || (lwork>=MAX(1,8*n))) :: lwork
<type_in> intent(out), dimension(MAX(lwork,1)), depend(lwork) :: work
<rtype_in> intent(hide), dimension(8*n), depend(n) :: rwork
integer intent(out):: info
end subroutine <tchar>ggev
end interface
end python module generic_flapack
! This file was auto-generated with f2py (version:2.10.173).
! See http://cens.ioc.ee/projects/f2py2e/
|