!%f90 -*- f90 -*- ! Signatures for f2py-wrappers of FORTRAN LEVEL 1 BLAS functions. ! ! Author: Pearu Peterson ! Created: Jan-Feb 2002 ! Modified: Fabian Pedregosa, 2011 ! ! Implemented: ! ! rotg, rotmg, rot, rotm, swap, scal, copy, axpy, dot, dotu, dotc ! nrm2, asum, amax, iamax ! NOTE: Avoiding wrappers hack does not work under 64-bit Gentoo system ! with single precision routines, so they are removed. ! ! Shorthand notations ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! Level 1 BLAS ! subroutine rotg(a,b,c,s) ! Computes the parameters for a Givens rotation. ! ! Given the Cartesian coordinates (a, b) of a point, these routines return ! the parameters c, s, r, and z associated with the Givens rotation. The ! parameters c and s define a unitary matrix such that: ! ! ( c s ) ( a ) ( r ) ! ( ) ( ) = ( ) ! (-s c ) ( b ) = ( 0 ) ! ! The parameter z is defined such that if |a| > |b|, z is s; otherwise if c ! is not 0 z is 1/c; otherwise z is 1. ! ! ! Parameters ! ---------- ! a : float or complex number ! Provides the x-coordinate for the point. ! ! b : float or complex number ! Provides the y-coordinate. ! ! Returns ! ------- ! c, s : ! Parameter c associated with the Givens rotation. ! ! Notes ! ----- ! Unlike the FORTRAN implementation, this function will not return ! parameters r and z, as these can easily be computed from the ! returned parameters. ! callprotoargument *,*,*,* intent(in) :: a, b intent(out,out=c) :: c intent(out,out=s) :: s end subroutine rotg subroutine rotmg(d1,d2,x1,y1,param) ! Computes the parameters for a modified Givens rotation. ! ! Given Cartesian coordinates (x1, y1) of an input vector, this ! routine compute the components of a modified Givens ! transformation matrix H that zeros the y-component of the ! resulting vector: ! ! [x] [sqrt(d1) x1] ! [ ] = H [ ] ! [0] [sqrt(d2) y1] ! callstatement (*f2py_func)(&d1,&d2,&x1,&y1,param) callprotoargument *,*,*,*,* intent(in) :: d1, d2, x1, y1 intent(out), dimension(5) :: param end subroutine rotmg subroutine rot(n,x,offx,incx,y,offy,incy,c,s) ! Applies a plane rotation with real cosine and complex sine to a ! pair of complex vectors and returns the modified vectors. ! ! x, y are input vectors and c, s are values that define a rotation: ! ! [ c s] ! [ ] ! [-conj(s) c] ! ! where c*c + s*conjg(s) = 1.0. ! callstatement (*f2py_func)(&n,x+offx,&incx,y+offy,&incy,&c,&s) callprotoargument F_INT*,*,F_INT*,*,F_INT*,*,* dimension(*),intent(in,out,copy) :: x,y intent(in) :: c, s integer optional, intent(in), check(incx>0||incx<0) :: incx = 1 integer optional, intent(in), check(incy>0||incy<0) :: incy = 1 integer optional, intent(in), depend(x) :: offx=0 integer optional, intent(in), depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end subroutine rot subroutine rotm(n,x,offx,incx,y,offy,incy,param) ! Performs rotation of points in the modified plane ! ! Given two complex vectors x and y, each vector element of these vectors is ! replaced as follows: ! ! x(i) = H*x(i) + H*y(i) ! y(i) = H*y(i) - H*x(i) ! ! where H is a modified Givens transformation matrix whose values are stored ! in the param(2) through param(5) array. callstatement (*f2py_func)(&n,x+offx,&incx,y+offy,&incy,param) callprotoargument F_INT*,*,F_INT*,*,F_INT*,* dimension(*), intent(in,out,copy) :: x, y dimension(5), intent(in) :: param integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional, intent(in),depend(x) :: offx=0 integer optional, intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end subroutine rotm subroutine swap(n,x,offx,incx,y,offy,incy) ! Swap two arrays, `x` and `y` callstatement (*f2py_func)(&n,x+offx,&incx,y+offy,&incy) callprotoargument F_INT*,*,F_INT*,*,F_INT* dimension(*), intent(in,out) :: x, y integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional, intent(in),depend(x) :: offx=0 integer optional, intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end subroutine swap subroutine scal(n,a,x,offx,incx) ! Computes the product of a vector by a scalar: y = a*x callstatement (*f2py_func)(&n,&a,x+offx,&incx) callprotoargument F_INT*,*,*,F_INT* intent(in):: a dimension(*), intent(in,out) :: x integer optional, intent(in), check(incx>0||incx<0) :: incx = 1 integer optional, intent(in), depend(x) :: offx=0 check(offx>=0 && offx(n-1)*abs(incx)) :: n end subroutine scal subroutine scal(n,a,x,offx,incx) ! Computes the product of a vector by a scalar, scales a complex ! vector by a real constant ! y = a*x callstatement (*f2py_func)(&n,&a,x+offx,&incx) callprotoargument F_INT*,*,*,F_INT* intent(in) :: a dimension(*), intent(in,out,copy) :: x integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),depend(x) :: offx=0 check(offx>=0 && offx(n-1)*abs(incx)) :: n end subroutine scal subroutine copy(n,x,offx,incx,y,offy,incy) ! Copy y <- x callstatement (*f2py_func)(&n,x+offx,&incx,y+offy,&incy) callprotoargument F_INT*,*,F_INT*,*,F_INT* dimension(*), intent(in) :: x dimension(*), intent(in,out) :: y integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional, intent(in),depend(x) :: offx=0 integer optional, intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end subroutine copy subroutine axpy(n,a,x,offx,incx,y,offy,incy) ! Calculate z = a*x+y, where a is scalar. callstatement (*f2py_func)(&n,&a,x+offx,&incx,y+offy,&incy) callprotoargument F_INT*,*,*,F_INT*,*,F_INT* dimension(*), intent(in) :: x dimension(*), intent(in,out,out=z) :: y optional, intent(in):: a=<1.0,\0,(1.0\,0.0),\2> integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional, intent(in),depend(x) :: offx=0 integer optional, intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end subroutine axpy function sdot(n,x,offx,incx,y,offy,incy) result (xy) ! Computes a vector-vector dot product. fortranname sdot callstatement (*f2py_func)(&sdot,&n,x+offx,&incx,y+offy,&incy) callprotoargument float*,F_INT*,float*,F_INT*,float*,F_INT* real dimension(*), intent(in) :: x real dimension(*), intent(in) :: y real sdot,xy integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional, intent(in),depend(x) :: offx=0 integer optional, intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end function sdot function ddot(n,x,offx,incx,y,offy,incy) result (xy) ! Computes a vector-vector dot product. callstatement (*f2py_func)(&ddot,&n,x+offx,&incx,y+offy,&incy) callprotoargument double*,F_INT*,double*,F_INT*,double*,F_INT* double precision dimension(*), intent(in) :: x double precision dimension(*), intent(in) :: y double precision ddot,xy integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional, intent(in),depend(x) :: offx=0 integer optional, intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end function ddot function dotu(n,x,offx,incx,y,offy,incy) result(xy) :: dotu, xy fortranname wdotu callstatement (*f2py_func)(&dotu,&n,x+offx,&incx,y+offy,&incy) callprotoargument *,F_INT*,*,F_INT*,*,F_INT* dimension(*),intent(in) :: x dimension(*),intent(in) :: y integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional,intent(in),depend(x) :: offx=0 integer optional,intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end function dotu function dotc(n,x,offx,incx,y,offy,incy) result(xy) :: dotc, xy fortranname wdotc callstatement (*f2py_func)(&dotc,&n,x+offx,&incx,y+offy,&incy) callprotoargument *,F_INT*,*,F_INT*,*,F_INT* dimension(*),intent(in) :: x dimension(*),intent(in) :: y integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional, intent(in),check(incy>0||incy<0) :: incy = 1 integer optional,intent(in),depend(x) :: offx=0 integer optional,intent(in),depend(y) :: offy=0 check(offx>=0 && offx=0 && offy(n-1)*abs(incx)) :: n check(len(y)-offy>(n-1)*abs(incy)) :: n end function dotc function nrm2(n,x,offx,incx) result(n2) nrm2, n2 callstatement (*f2py_func)(&nrm2, &n,x+offx,&incx) callprotoargument *,F_INT*,*,F_INT* dimension(*),intent(in) :: x integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional,intent(in),depend(x) :: offx=0 check(offx>=0 && offx(n-1)*abs(incx)) :: n end function nrm2 function nrm2(n,x,offx,incx) result(n2) nrm2, n2 callstatement (*f2py_func)(&nrm2, &n,x+offx,&incx) callprotoargument *,F_INT*,*,F_INT* dimension(*),intent(in) :: x integer optional, intent(in),check(incx>0||incx<0) :: incx = 1 integer optional,intent(in),depend(x) :: offx=0 check(offx>=0 && offx(n-1)*abs(incx)) :: n end function nrm2 function asum(n,x,offx,incx) result (s) ! Computes the sum of magnitudes of the vector elements callstatement (*f2py_func)(&asum,&n,x+offx,&incx) callprotoargument *,F_INT*,*,F_INT* dimension(*), intent(in) :: x asum,s integer optional, intent(in), check(incx>0||incx<0) :: incx = 1 integer optional, intent(in), depend(x) :: offx=0 check(offx>=0 && offx(n-1)*abs(incx)) :: n end function asum function asum(n,x,offx,incx) result (s) ! Computes the sum of magnitudes of the vector elements callstatement (*f2py_func)(&asum,&n,x+offx,&incx) callprotoargument *,F_INT*,*,F_INT* dimension(*), intent(in) :: x asum,s integer optional, intent(in), check(incx>0||incx<0) :: incx = 1 integer optional, intent(in), depend(x) :: offx=0 check(offx>=0 && offx(n-1)*abs(incx)) :: n end function asum function iamax(n,x,offx,incx) result(k) ! Finds the index of the element with maximum absolute value. callstatement iamax_return_value = (*f2py_func)(&n,x+offx,&incx) - 1 callprotoargument F_INT*,*,F_INT* ! This is to avoid Fortran wrappers. integer iamax,k fortranname F_FUNC(iamax,IAMAX) intent(c) iamax dimension(*), intent(in) :: x integer optional, intent(in), check(incx>0||incx<0) :: incx = 1 integer optional, intent(in), depend(x) :: offx=0 check(offx>=0 && offx(n-1)*abs(incx)) :: n end function iamax