! -*- f90 -*- ! Signatures for f2py-wrappers of FORTRAN LEVEL 2 BLAS functions. ! ! Author: Pearu Peterson ! Created: Jan-Feb 2002 ! Modified: Fabian Pedregosa, 2011; Eric Moore, 2014 ! ! Implemented: ! gemv, hemv, symv, trmv, ger, geru, gerc, her, syr, hpr, spr, her2, syr2, ! spmv, hpmv, hbmv, sbmv, hpr2, spr2, gbmv, trsv, tbmv, tpmv, tbsv, tpsv ! ! ! Level 2 BLAS ! subroutine gemv(m,n,alpha,a,x,beta,y,offx,incx,offy,incy,trans,rows,cols,ly) ! Computes a matrix-vector product using a general matrix ! ! y = gemv(alpha,a,x,beta=0,y=0,offx=0,incx=1,offy=0,incy=0,trans=0) ! Calculate y <- alpha * op(A) * x + beta * y callstatement (*f2py_func)((trans?(trans==2?"C":"T"):"N"),&m,&n,&alpha,a,&m, & x+offx,&incx,&beta,y+offy,&incy) callprotoargument char*,F_INT*,F_INT*,*,*,F_INT*,*,F_INT*,*, & *,F_INT* integer optional, intent(in), check(trans>=0 && trans <=2) :: trans = 0 integer optional, intent(in), check(incx>0||incx<0) :: incx = 1 integer optional, intent(in), check(incy>0||incy<0) :: incy = 1 intent(in) :: alpha intent(in), optional :: beta = <0.0,\0,(0.0\,0.0),\2> dimension(*), intent(in) :: x dimension(ly), intent(in,copy,out), depend(ly),optional :: y integer intent(hide), depend(incy,rows,offy) :: ly = & (y_capi==Py_None?1+offy+(rows-1)*abs(incy):-1) dimension(m,n), intent(in) :: a integer depend(a), intent(hide):: m = shape(a,0) integer depend(a), intent(hide):: n = shape(a,1) integer optional, intent(in) :: offx=0 integer optional, intent(in) :: offy=0 check(offx>=0 && offxoffx+(cols-1)*abs(incx)) :: x depend(offx,cols,incx) :: x check(offy>=0 && offyoffy+(rows-1)*abs(incy)) :: y depend(offy,rows,incy) :: y integer depend(m,n,trans), intent(hide) :: rows = (trans?n:m) integer depend(m,n,trans), intent(hide) :: cols = (trans?m:n) end subroutine gemv subroutine gbmv(m,n,kl,ku,alpha,a,lda,x,incx,offx,beta,y,incy,offy,trans,ly) ! Performs one of the matrix-vector operations ! ! y := alpha*A*x + beta*y, or y := alpha*A**T*x + beta*y, ! or y := alpha*A**H*x + beta*y, ! ! where alpha and beta are scalars, x and y are vectors and A is an ! m by n band matrix, with kl sub-diagonals and ku super-diagonals. callstatement (*f2py_func)((trans?(trans==2?"C":"T"):"N"),&m,&n,&kl,&ku,&alpha,a,&lda,x+offx,&incx,&beta,y+offy,&incy) callprotoargument char*,F_INT*,F_INT*,F_INT*,F_INT*,*,*,F_INT*,*,F_INT*,*,*,F_INT* integer optional,intent(in),check(trans>=0 && trans <=2) :: trans = 0 integer intent(in), depend(ku,kl),check(m>=ku+kl+1) :: m integer intent(in),check(n>=0&&n==shape(a,1)),depend(a) :: n integer intent(in),check(kl>=0) :: kl integer intent(in),check(ku>=0) :: ku integer intent(hide),depend(a) :: lda = MAX(shape(a,0),1) integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer intent(hide),depend(m,n,incy,offy,trans) :: ly = & (y_capi==Py_None?1+offy+(trans==0?m-1:n-1)*abs(incy):-1) integer optional, intent(in) :: offx=0 integer optional, intent(in) :: offy=0 intent(in) :: alpha intent(in),optional :: beta = <0.0,\0,(0.0\,0.0),\2> dimension(lda,n),intent(in) :: a dimension(ly), intent(in,out,copy,out=yout),depend(ly),optional :: y check(offy>=0 && offyoffy+(trans==0?m-1:n-1)*abs(incy)) :: y depend(offy,n,incy) :: y dimension(*), intent(in) :: x check(offx>=0 && offxoffx+(trans==0?m-1:n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine gbmv subroutine bmv(n,k,alpha,a,lda,x,incx,offx,beta,y,incy,offy,lower,ly) ! Performs the matrix-vector operation ! ! y := alpha*A*x + beta*y, ! ! where alpha and beta are scalars, x and y are n element vectors and ! A is an n by n symmetric band matrix, with k super-diagonals. callstatement (*f2py_func)((lower?"L":"U"),&n,&k,&alpha,a,&lda,x+offx,&incx,&beta,y+offy,&incy) callprotoargument char*,F_INT*,F_INT*,*,*,F_INT*,*,F_INT*,*,*,F_INT* integer optional, intent(in),check(lower==0||lower==1) :: lower = 0 integer intent(hide),depend(a),check(n>=0) :: n = shape(a,1) integer intent(hide),depend(a) :: lda = MAX(shape(a,0),1) integer intent(in),depend(lda),check(k>=0&&k<=lda-1) :: k integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer intent(hide),depend(incy,n,offy) :: ly = (y_capi==Py_None?1+offy+(n-1)*abs(incy):-1) integer optional, intent(in) :: offx=0 integer optional, intent(in) :: offy=0 intent(in) :: alpha intent(in),optional :: beta = <0.0,\0,(0.0\,0.0),\2> dimension(lda,n),intent(in) :: a dimension(ly), intent(in,out,copy,out=yout),depend(ly),optional :: y check(offy>=0 && offyoffy+(n-1)*abs(incy)) :: y depend(offy,n,incy) :: y dimension(*), intent(in) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine bmv subroutine pmv(n,alpha,ap,x,incx,offx,beta,y,incy,offy,lower,ly) ! Computes a matrix-vector product for a symmetric matrix ! ! Calculate y <- alpha * A * x + beta * y, A is symmmetric in packed form. callstatement (*f2py_func)((lower?"L":"U"),&n,&alpha,ap,x+offx,&incx,&beta,y+offy,&incy) callprotoargument char*,F_INT*,*,*,*,F_INT*,*,*,F_INT* integer optional, intent(in),check(lower==0||lower==1) :: lower = 0 integer intent(in),check(n>=0) :: n integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer intent(hide),depend(incy,n,offy) :: ly = (y_capi==Py_None?1+offy+(n-1)*abs(incy):-1) integer optional, intent(in) :: offx=0 integer optional, intent(in) :: offy=0 intent(in) :: alpha intent(in),optional :: beta = <0.0,\0,(0.0\,0.0),\2,\2,\2> dimension(*),depend(n),intent(in),check(len(ap)>=(n*(n+1)/2)) :: ap dimension(ly), intent(in,out,copy,out=yout),depend(ly),optional :: y check(offy>=0 && offyoffy+(n-1)*abs(incy)) :: y depend(offy,n,incy) :: y dimension(*), intent(in) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine pmv subroutine (n,alpha,a,x,beta,y,offx,incx,offy,incy,lower,ly) ! Computes a matrix-vector product for a symmetric/hermitian matrix ! ! Calculate y <- alpha * A * x + beta * y, A is symmmetric/hermitian callstatement (*f2py_func)((lower?"L":"U"),&n,&alpha,a,&n,x+offx,&incx,&beta, & y+offy,&incy) callprotoargument char*,F_INT*,*,*,F_INT*,*,F_INT*,*, & *,F_INT* integer optional, intent(in),check(lower==0||lower==1) :: lower = 0 integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 intent(in) :: alpha intent(in),optional :: beta = <0.0,\0,(0.0\,0.0),\2> dimension(*), intent(in) :: x dimension(ly), intent(in,copy,out),depend(ly),optional :: y integer intent(hide),depend(incy,n,offy) :: ly = & (y_capi==Py_None?1+offy+(n-1)*abs(incy):-1) dimension(n,n), intent(in),check(shape(a,0)==shape(a,1)) :: a integer depend(a), intent(hide):: n = shape(a,0) integer optional, intent(in) :: offx=0 integer optional, intent(in) :: offy=0 check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x check(offy>=0 && offyoffy+(n-1)*abs(incy)) :: y depend(offy,n,incy) :: y end subroutine subroutine ger<,,u,u,c,c>(m,n,alpha,x,incx,y,incy,a,lda) ! Performs a rank-1 update of a general matrix. ! ! Calculate a <- alpha*x*y^T + a ! Calculate a <- alpha*x*y^H + a ! integer intent(hide),depend(x) :: m = len(x) integer intent(hide),depend(y) :: n = len(y) intent(in) :: alpha dimension(m), intent(in,overwrite) :: x integer optional, intent(in),check(incx==1||incx==-1) :: incx = 1 dimension(n), intent(in,overwrite) :: y integer optional, intent(in),check(incy==1||incy==-1) :: incy = 1 dimension(m,n), intent(in,out,copy),optional :: & a = <0.0,\0,(0.0\,0.0),\2,\2,\2> integer intent(hide), depend(m) :: lda=m end subroutine ger<,,u,u,c,c> subroutine r(alpha,x,lower,incx,offx,n,a) ! Performs a rank-1 update of a symmetric/hermitian matrix. ! ! Calculate a <- alpha*x*x^T + a ! Calculate a <- alpha*x*x^H + a ! callstatement (*f2py_func)((lower?"L":"U"),&n,&alpha,x+offx,&incx,a,&n) callprotoargument char*, F_INT*, *, *, F_INT*, *, F_INT* integer, optional, intent(in), check(lower == 0 || lower == 1) :: lower = 0 intent(in) :: alpha dimension(*), depend(offx) :: x check(offx >= 0 && offx < len(x)) :: x integer, intent(in), optional :: offx = 0 integer, intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer, intent(in), optional :: n = (len(x)-1-offx)/abs(incx)+1 check(n >= 0) :: n check(n <= (len(x)-1-offx)/abs(incx)+1) :: n depend(x, offx, incx) :: n dimension(n,n), intent(in,copy,out), optional :: a depend(x, offx, incx, n) :: a end subroutine r subroutine r2(alpha,x,y,lower,incx,offx,incy,offy,n,a) ! Performs a rank-2 update of a symmetric/hermitian matrix. ! ! Calculate a <- alpha*x*y^T + alpha*y*x^T + a ! Calculate a <- alpha*x*y^H + alpha*y*x^H + a ! callstatement (*f2py_func)((lower?"L":"U"),&n,&alpha,x+offx,&incx,y+offy,&incy,a,&n) callprotoargument char*, F_INT*, *, *, F_INT*, *, F_INT*, *, F_INT* integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 intent(in) :: alpha dimension(*), intent(in), check(offx >= 0 && offx < len(x)), depend(offx) :: x dimension(*), intent(in), check(offy >= 0 && offy < len(y)), depend(offy) :: y integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(in), optional :: offx = 0 integer intent(in), optional, check(incy>0||incy<0) :: incy = 1 integer intent(in), optional :: offy = 0 integer intent(in), optional :: n = ((len(x)-1-offx)/abs(incx)+1 <= (len(y)-1-offy)/abs(incy)+1 ? (len(x)-1-offx)/abs(incx)+1 : (len(y)-1-offy)/abs(incy)+1) depend(x,incx,offx,y,incy,offy) :: n check(n>=0) :: n check(n <= (len(x)-1-offx)/abs(incx)+1) :: n check(n <= (len(y)-1-offy)/abs(incy)+1) :: n dimension(n,n), intent(in,copy,out), optional :: a depend(incx, offx, x, incy, offy, y, n) :: a end subroutine r2 subroutine pr(n,alpha,x,incx,offx,ap,lower) ! Performs the symmetric rank 1 operation ! ! A := alpha*x*x**T(H) + A, ! ! where alpha is a scalar, x is an n element vector and A is an n x n ! symmetric/hermitian matrix, supplied in packed form. callstatement (*f2py_func)((lower?"L":"U"),&n,&alpha,x+offx,&incx,ap) callprotoargument char*,F_INT*,*,*,F_INT*,* integer intent(in),check(n>=0) :: n integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(in), optional :: offx = 0 intent(in) :: alpha dimension(*), intent(in) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x dimension(*),depend(n),intent(in,out,copy,out=apu) :: ap check(len(ap)>=(n*(n+1)/2)) :: ap end subroutine pr subroutine pr2(n,alpha,x,incx,offx,y,incy,offy,ap,lower) ! Performs the symmetric rank 2 operation ! ! A := alpha*x*y**T + alpha*y*x**T + A, ! ! where alpha is a scalar, x and y are n element vectors and A is an ! n by n symmetric/hermitian matrix, supplied in packed form. callstatement (*f2py_func)((lower?"L":"U"),&n,&alpha,x+offx,&incx,y+offy,&incy,ap) callprotoargument char*,F_INT*,*,*,F_INT*,*,F_INT*,* integer intent(in),check(n>=0) :: n integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(in), optional, check(incy>0||incy<0) :: incy = 1 integer intent(in), optional :: offx = 0 integer intent(in), optional :: offy = 0 intent(in) :: alpha dimension(*), intent(in) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x dimension(*), intent(in) :: y check(offy>=0 && offyoffy+(n-1)*abs(incy)) :: y depend(offy,n,incy) :: y dimension(*),depend(n),intent(in,out,copy,out=apu) :: ap check(len(ap)>=(n*(n+1)/2)) :: ap end subroutinepr2 subroutine tbsv(n,k,a,lda,x,incx,offx,lower,trans,diag) ! Solves one of the systems of equations ! ! A*xout = x, or A**T*xout = x, or A**H*xout = x, ! ! where xout and x are n element vectors and A is an n by n unit, or ! non-unit, upper or lower triangular matrix. ! ! No test for singularity or near-singularity is included in this ! routine. Such tests must be performed before calling this routine. callstatement (*f2py_func)((lower?"L":"U"),(trans?(trans==2?"C":"T"):"N"),(diag?"U":"N"),&n,&k,a,&lda,x+offx,&incx) callprotoargument char*,char*,char*,F_INT*,F_INT*,*,F_INT*,*,F_INT* integer intent(hide),check(n>=0),depend(a) :: n = shape(a,1) integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 integer intent(in), optional, check(trans >= 0 || trans <= 2) :: trans = 0 integer intent(in), optional, check(diag == 0 || diag == 1) :: diag = 0 integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(hide), depend(a) :: lda = MAX(shape(a,0),1) integer intent(in),depend(lda),check(k>=0&&k<=lda-1) :: k integer intent(in), optional :: offx = 0 intent(in), dimension(lda,n) :: a dimension(*), intent(in,out,copy,out=xout) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine tbsv subroutine tpsv(n,ap,x,incx,offx,lower,trans,diag) ! Solves one of the systems of equations ! ! A*xout = x, or A**T*xout = x, or A**H*xout = x, ! ! where xout and x are n element vectors and A is an n by n unit, or ! non-unit, upper or lower triangular matrix, supplied in packed form. ! ! No test for singularity or near-singularity is included in this ! routine. Such tests must be performed before calling this routine. callstatement (*f2py_func)((lower?"L":"U"),(trans?(trans==2?"C":"T"):"N"),(diag?"U":"N"),&n,ap,x+offx,&incx) callprotoargument char*,char*,char*,F_INT*,*,*,F_INT* integer intent(in),check(n>=0) :: n integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 integer intent(in), optional, check(trans >= 0 || trans <= 2) :: trans = 0 integer intent(in), optional, check(diag == 0 || diag == 1) :: diag = 0 integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(in), optional :: offx = 0 dimension(*),depend(n),intent(in),check(len(ap)>=(n*(n+1)/2)) :: ap dimension(*), intent(in,out,copy,out=xout) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine tpsv subroutine trmv(n,a,x,offx,incx,lower,trans,diag) ! Computes a matrix-vector product using a triangular matrix ! ! x <- op(A) * x, A is triangular ! callstatement (*f2py_func)((lower?"L":"U"), (trans?(trans==2?"C":"T"):"N"), & (diag?"U":"N"),&n,a,&n,x+offx,&incx) callprotoargument char*,char*,char*,F_INT*,*,F_INT*,*,F_INT* integer optional, intent(in), check(trans>=0 && trans <=2) :: trans = 0 integer optional, intent(in), check(lower==0||lower==1) :: lower = 0 integer optional, intent(in), check(diag==0||diag==1) :: diag = 0 integer optional, intent(in), check(incx>0||incx<0) :: incx = 1 dimension(*), intent(in,out,copy) :: x dimension(n,n), intent(in),check(shape(a,0)==shape(a,1)) :: a integer depend(a), intent(hide):: n = shape(a,0) integer optional, intent(in), depend(x) :: offx=0 check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: n depend(x,offx,incx) :: n end subroutine trmv subroutine trsv(n,a,lda,x,incx,offx,lower,trans,diag) ! Solves one of the systems of equations ! ! A*xout = x, or A**T*xout = x, or A**H*xout = x, ! ! where xout and x are n element vectors and A is an n by n unit, or ! non-unit, upper or lower triangular matrix. ! ! No test for singularity or near-singularity is included in this ! routine. Such tests must be performed before calling this routine. callstatement (*f2py_func)((lower?"L":"U"),(trans?(trans==2?"C":"T"):"N"),(diag?"U":"N"),&n,a,&lda,x+offx,&incx) callprotoargument char*,char*,char*,F_INT*,*,F_INT*,*,F_INT* integer intent(hide),check(n>=0),depend(a) :: n = shape(a,0) integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 integer intent(in), optional, check(trans >= 0 || trans <= 2) :: trans = 0 integer intent(in), optional, check(diag == 0 || diag == 1) :: diag = 0 integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(hide), depend(a) :: lda = MAX(shape(a,0),1) integer intent(in), optional :: offx = 0 intent(in), dimension(n,n), check(shape(a,0)==shape(a,1)) :: a dimension(*), intent(in,out,copy,out=xout) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine trsv subroutine tbmv(n,k,a,lda,x,incx,offx,lower,trans,diag) ! Solves one of the systems of equations ! ! A*xout = x, or A**T*xout = x, or A**H*xout = x, ! ! where xout and x are an n element vector and A is an n by n unit, or non-unit, ! upper or lower triangular band matrix, with ( k + 1 ) diagonals. callstatement (*f2py_func)((lower?"L":"U"),(trans?(trans==2?"C":"T"):"N"),(diag?"U":"N"),&n,&k,a,&lda,x+offx,&incx) callprotoargument char*,char*,char*,F_INT*,F_INT*,*,F_INT*,*,F_INT* integer intent(hide),check(n>=0),depend(a) :: n = shape(a,1) integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 integer intent(in), optional, check(trans >= 0 || trans <= 2) :: trans = 0 integer intent(in), optional, check(diag == 0 || diag == 1) :: diag = 0 integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(hide), depend(a) :: lda = MAX(shape(a,0),1) integer intent(in),depend(lda),check(k>=0&&k<=lda-1) :: k integer intent(in), optional :: offx = 0 intent(in), dimension(lda,n) :: a dimension(*), intent(in,out,copy,out=xout) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine tbmv subroutine tpmv(n,ap,x,incx,offx,lower,trans,diag) ! performs one of the matrix-vector operations ! ! x := A*x, or x := A**T*x, ! ! where x is n element vectors and A is an n by n unit, or ! non-unit, upper or lower triangular matrix, supplied in packed form. callstatement (*f2py_func)((lower?"L":"U"),(trans?(trans==2?"C":"T"):"N"),(diag?"U":"N"),&n,ap,x+offx,&incx) callprotoargument char*,char*,char*,F_INT*,*,*,F_INT* integer intent(in),check(n>=0) :: n integer intent(in), optional, check(lower == 0 || lower == 1) :: lower = 0 integer intent(in), optional, check(trans >= 0 || trans <= 2) :: trans = 0 integer intent(in), optional, check(diag == 0 || diag == 1) :: diag = 0 integer intent(in), optional, check(incx>0||incx<0) :: incx = 1 integer intent(in), optional :: offx = 0 dimension(*),depend(n),intent(in),check(len(ap)>=(n*(n+1)/2)) :: ap dimension(*), intent(in,out,copy,out=xout) :: x check(offx>=0 && offxoffx+(n-1)*abs(incx)) :: x depend(offx,n,incx) :: x end subroutine tpmv