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
|
! -*- f90 -*-
!
! Contains wrappers for the following LAPACK routines:
!
! Driver routines for standard eigenvalue and singular value problems:
! syev, heev (SEP symmetric/hermitian, eigenvalues/vectors)
! syevd, heevd (SEP symmetric/hermitian, eigenvalues/vectors, D&C)
! syevx, heevx (.., expert) - Not Implemented
! syevr, heevr (.., RRR)
! spev, hpev, spevd, hpevd, spevx, hpevx (..., packed storage) - Not Implemented
! sbev, hbev, sbevd, hbevd, sbevx, hbevx (..., band) - Not Implemented
! stev, stevd, stevx, stevr (..., tridiagonal) - Not Implemented
! gees (NEP, general, Schur factorization)
! geesx (NEP, general, Schur factorization, expert) - Not Implemented
! geev (NEP, general, eigenvalues/vectors)
! geevx (NEP, general, eigenvalues/vectors, expert) - Not Implemented
! gesvd (SVD, general, singular values/vectors) - Not Implemented
! gesdd (SVD, general, singular values/vectors, D&C)
!
!
! <sym=sy,\0,he,\2>
subroutine <prefix><sym>ev(compute_v,lower,n,w,a,work,lwork,<_1=,,rwork\,,\2>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.
! w,v,info = heev(a,compute_v=1,lower=0,lwork=2*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,<_2=,,rwork\,,\2>&info)
callprotoargument char*,char*,int*,<ctype>*,int*,<ctypereal>*,<ctype>*,int*,<_3=,,float*\,,double*\,>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)
<ftype> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
intent(in,copy,out,out=v) :: a
<ftypereal> dimension(n),intent(out),depend(n) :: w
! <_lwork=3*n-1,\0,2*n-1,\2>
integer optional,intent(in),depend(n) :: lwork=<_lwork>
check(lwork>=<_lwork>) :: lwork
<ftype> dimension(lwork),intent(hide,cache),depend(lwork) :: work
<ftypereal> dimension(3*n-1),intent(hide,cache),depend(n) :: rwork
integer intent(out) :: info
end subroutine <prefix><sym>ev
subroutine <prefix><sym>evd(compute_v,lower,n,w,a,work,lwork,iwork,liwork,<_1=,,rwork\,lrwork\,,\2>info)
! w,v,info = syevd(a,compute_v=1,lower=0,lwork=min_lwork,overwrite_a=0)
! Compute all eigenvalues and, optionally, eigenvectors of a
! real symmetric matrix A using D&C.
!
! Performance tip:
! If compute_v=0 then set also overwrite_a=1.
! w,v,info = heevd(a,compute_v=1,lower=0,lwork=min_lwork,overwrite_a=0)
! Compute all eigenvalues and, optionally, eigenvectors of a
! complex Hermitian matrix A using D&C.
!
! 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,<_2=,,rwork\,&lrwork\,,\2>iwork,&liwork,&info)
callprotoargument char*,char*,int*,<ctype>*,int*,<ctypereal>*,<ctype>*,int*,<_3=,,float*\,int*\,,double*\,int*\,>int*,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)
<ftype> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
intent(in,copy,out,out=v) :: a
<ftypereal> dimension(n),intent(out),depend(n) :: w
! <_lwork=(compute_v?1+6*n+2*n*n:2*n+1),\0,(compute_v?2*n+n*n:n+1),\2>
integer optional,intent(in),depend(n,compute_v) :: lwork=<_lwork>
check(lwork>=<_lwork>) :: lwork
<ftype> dimension(lwork),intent(hide,cache),depend(lwork) :: work
integer intent(hide),depend(n,compute_v) :: liwork = (compute_v?3+5*n:1)
integer dimension(liwork),intent(hide,cache),depend(liwork) :: iwork
! <_lrwork=,,(compute_v?1+5*n+2*n*n:n),\2>
integer intent(hide),depend(n,compute_v) :: lrwork = <_lrwork>
<ftypereal> dimension(lrwork),intent(hide,cache),depend(n,lrwork) :: rwork
integer intent(out) :: info
end subroutine <prefix><sym>evd
subroutine <prefix><sym>evr(n,a,compute_v,lower,vrange,irange,atol,w,z,m,ldz,isuppz,work,lwork,<,,rwork\,lrwork\,,\2>iwork,liwork,info)
! w,v,info = {sy|he}evr(a,compute_v=1,lower=0,vrange=None,irange=None,atol=-1,lwork=min_lwork,overwrite_a=0)
!
! Compute range of eigenvalues and, optionally, eigenvectors of a
! real symmetric matrix A using RRR.
!
! Performance tip:
! If compute_v=0 then set also overwrite_a=1.
! Warning:
! If compute_v=0 and overwrite_a=1, the contents of a is destroyed.
callstatement if(irange_capi==Py_None);else{irange[0]++;irange[1]++;}(*f2py_func)((compute_v?"V":"N"),(vrange_capi==Py_None?(irange_capi==Py_None?"A":"I"):"V"),(lower?"L":"U"),&n,a,&n,vrange,vrange+1,irange,irange+1,&atol,&m,w,z,&ldz,isuppz,work,&lwork,<_2=,,rwork\,&lrwork\,,\2>iwork,&liwork,&info);if(irange_capi==Py_None);else{irange[0]--;irange[1]--;}if(vrange_capi==Py_None);else {capi_w_tmp-\>dimensions[0]=capi_z_tmp-\>dimensions[1]=m;/*capi_z_tmp-\>strides[0]=m*capi_z_tmp-\>descr-\>elsize;*/}
callprotoargument char*,char*,char*,int*,<ctype>*,int*,<ctypereal>*,<ctypereal>*,int*,int*,<ctypereal>*,int*,<ctypereal>*,<ctype>*,int*,int*,<ctype>*,int*,<_3=,,float*\,int*\,,double*\,int*\,>int*,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)
<ftype> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
intent(in,copy) :: a
<ftypereal> optional,dimension(2),intent(in) :: vrange
integer optional,dimension(2),intent(in),depend(n) :: irange
check(irange_capi==Py_None || (irange[0]>=0 && irange[1]<n)) irange
<ftypereal> optional,intent(in) :: atol = -1.0
integer intent(hide),depend(vrange,irange,n) :: m = (irange_capi==Py_None?n:irange[1]-irange[0]+1)
<ftypereal> dimension(m),intent(out),depend(m) :: w
integer intent(hide),depend(compute_v,n) :: ldz = (compute_v?n:1)
<ftype> dimension(ldz,m),intent(out,out=v),depend(ldz,m) :: z
integer intent(hide),depend(m),dimension(2*m) :: isuppz
! <_lwork=26*n,\0,18*n,\2> Includes bug fix in `man zheevr`.
integer optional,intent(in),depend(n) :: lwork=<_lwork>
check(lwork>=<_lwork>) lwork
<ftype> dimension(lwork),intent(hide,cache),depend(lwork) :: work
integer intent(hide),depend(n,compute_v) :: liwork = 10*n
integer dimension(liwork),intent(hide,cache),depend(liwork) :: iwork
! <_lrwork=,,24*n,\2>
integer intent(hide),depend(n) :: lrwork = <_lrwork>
<ftypereal> dimension(lrwork),intent(hide,cache),depend(n,lrwork) :: rwork
integer intent(out) :: info
end subroutine <prefix><sym>evr
subroutine <prefix>gees(compute_v,sort_t,<prefix>select,n,a,nrows,sdim,<wr\,wi,\0,w,\2>,vs,ldvs,work,lwork,<,,rwork\,,\2>bwork,info)
! t,sdim,(wr,wi|w),vs,info = gees(zselect,a,compute_v=1,sort_t=0,lwork=3*n,zselect_extra_args=(),overwrite_a=0)
! 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
! t,sdim,wr,wi,vs,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_<prefix>select_in_gees__user__routines,&n,a,&nrows,&sdim,<wr\,wi,\0,w,\2>,vs,&ldvs,work,&lwork,<,,rwork\,,\2>bwork,&info,1,1)
callprotoargument char*,char*,int(*)(<float*\,float*,double*\,double*,complex_float*,complex_double*>),int*,<ctype>*,int*,int*,<ctype>*,<float*\,,double*\,,,><ctype>*,int*,<ctype>*,int*,<,,float*\,,double*\,>int*,int*,int,int
use 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 <prefix>select
integer intent(hide),depend(a) :: n = shape(a,1)
<ftype> 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
<ftype> intent(out),dimension(n) :: <wr\,wi,\0,w,\2>
<ftype> intent(out),depend(ldvs,n),dimension(ldvs,n) :: vs
integer intent(hide),depend(compute_v,n) :: ldvs=((compute_v==1)?n:1)
<ftype> intent(hide,cache),depend(lwork),dimension(lwork) :: work
integer optional,intent(in),check(lwork >= MAX(1,3*n)),depend(n) :: lwork = 3*n
<ftypereal> intent(hide,cache),depend(n),dimension(n) :: rwork
logical intent(hide,cache),depend(n),dimension(n) :: bwork
integer intent(out) :: info
end subroutine <prefix>gees
subroutine <prefix>geev(compute_vl,compute_vr,n,a,<wr\,wi,\0,w,\2>,vl,ldvl,vr,ldvr,work,lwork,<,,rwork\,,\2>info)
! wr,wi,vl,vr,info = geev(a,compute_vl=1,compute_vr=1,lwork=4*n,overwrite_a=0)
! 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,<wr\,wi,\0,w,\2>,vl,&ldvl,vr,&ldvr,work,&lwork,<,,rwork\,,\2>&info);}
callprotoargument char*,char*,int*,<ctype>*,int*,<ctype>*,<float*\,,double*\,,,><ctype>*,int*,<ctype>*,int*,<ctype>*,int*,<,,float*\,,double*\,>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)
<ftype> dimension(n,n),intent(in,copy) :: a
check(shape(a,0)==shape(a,1)) :: a
<ftype> dimension(n),intent(out),depend(n) :: <wr\,wi,\0,w,\2>
<ftype> dimension(ldvl,n),intent(out) :: vl
integer intent(hide),depend(n,compute_vl) :: ldvl=(compute_vl?n:1)
<ftype> dimension(ldvr,n),intent(out) :: vr
integer intent(hide),depend(n,compute_vr) :: ldvr=(compute_vr?n:1)
! <_lwork=(compute_vl||compute_vr)?4*n:3*n,\0,2*n,\2>
integer optional,intent(in),depend(n,compute_vl,compute_vr) :: lwork=<_lwork>
check(lwork>=<_lwork>) :: lwork
<ftype> dimension(lwork),intent(hide,cache),depend(lwork) :: work
<ftypereal> dimension(2*n),intent(hide,cache),depend(n) :: rwork
integer intent(out):: info
end subroutine <prefix>geev
subroutine <prefix>gesdd(m,n,minmn,du,dvt,a,compute_uv,u,s,vt,work,lwork,<,,rwork\,,\2>iwork,info)
! u,s,vh,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\,,\2>iwork,&info)
callprotoargument char*,int*,int*,<ctype>*,int*,<ctypereal>*,<ctype>*,int*,<ctype>*,int*,<ctype>*,int*,<,,float*\,,double*\,>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)
<ftype> dimension(m,n),intent(in,copy) :: a
<ftypereal> dimension(minmn),intent(out),depend(minmn) :: s
<ftype> dimension(du,du),intent(out),depend(du) :: u
<ftype> dimension(dvt,dvt),intent(out),depend(dvt) :: vt
<ftype> dimension(lwork),intent(hide,cache),depend(lwork) :: work
! <_lwork=(compute_uv?4*minmn*minmn+MAX(m\,n)+9*minmn:MAX(14*minmn+4\,10*minmn+827)+MAX(m\,n)),\0,(compute_uv?2*minmn*minmn+MAX(m\,n)+2*minmn:2*minmn+MAX(m\,n)),\2>
integer optional,intent(in),depend(minmn,compute_uv) &
:: lwork = <_lwork>
! 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
<ftypereal> dimension((compute_uv?5*minmn*minmn+7*minmn:5*minmn)),intent(hide,cache),depend(minmn,compute_uv) :: rwork
integer intent(hide,cache),dimension(8*minmn),depend(minmn) :: iwork
integer intent(out)::info
end subroutine <prefix>gesdd
|