Skip to content

Commit

Permalink
BUG: Prevent crash in LAPACK Hermitian eigendecomposition wrapper.
Browse files Browse the repository at this point in the history
In the LAPACK wrapper for ?heevr, enlarge the size of the isuppz array
parameter in order to prevent out-of-bound access by the underlying
library. For both ?heevr and ?syevr, more robust checks and default
values are also put in place.

Reference: scipy#11709
  • Loading branch information
congma committed Mar 28, 2020
1 parent a09b5fe commit 0fdce89
Showing 1 changed file with 22 additions and 17 deletions.
39 changes: 22 additions & 17 deletions scipy/linalg/flapack_sym_herm.pyf.src
Original file line number Diff line number Diff line change
Expand Up @@ -757,22 +757,22 @@ subroutine <prefix2>syevr(compute_v,range,lower,n,a,lda,vl,vu,il,iu,abstol,w,z,m
integer optional,intent(in),check(compute_v==1||compute_v==0):: compute_v = 1
character optional,intent(in),check(*range=='A'||*range=='V' ||*range=='I') :: range='A'
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer optional,intent(in),check(il>=1) :: il=1
integer optional,intent(in),check(iu>=il||iu<=n),depend(il,n) :: iu=n
integer optional,intent(in),check(il>=1&&il<=n),depend(n) :: il=1
integer optional,intent(in),check(iu>=il&&iu<=n),depend(il,n) :: iu=n
<ftype2> optional,intent(in) :: vl=0.0
<ftype2> optional,intent(in),check(vu>=vl),depend(vl) :: vu=1.0
<ftype2> optional,intent(in),check(vu>vl),depend(vl) :: vu=1.0
<ftype2> intent(in) :: abstol=0.0
integer optional,intent(in),depend(n),check(lwork>=1||lwork==-1) :: lwork=max(26*n,1)
integer optional,intent(in),depend(n),check(liwork>=1||liwork==-1):: liwork= max(1,10*n)
integer optional,intent(in),depend(n),check(lwork>=max(1,26*n)||lwork==-1) :: lwork=max(26*n,1)
integer optional,intent(in),depend(n),check(liwork>=max(1,10*n)||liwork==-1):: liwork= max(1,10*n)

integer intent(hide),depend(a) :: n=shape(a,0)
integer intent(hide),depend(n) :: lda=max(1,n)
integer intent(hide),depend(z) :: ldz=max(1,shape(z,0))
integer intent(hide),depend(compute_v,n) :: ldz=(compute_v?max(1,n):1)
<ftype2> intent(hide),dimension(lwork),depend(lwork) :: work
integer intent(hide),dimension(liwork),depend(liwork) :: iwork

<ftype2> intent(out),dimension(n),depend(n) :: w
<ftype2> intent(out),dimension((compute_v?MAX(0,n):0),(compute_v?(*range=='I'?iu-il+1:MAX(1,n)):0)),depend(n,compute_v,range,iu,il) :: z
<ftype2> intent(out),dimension((compute_v?ldz:0),(compute_v?(*range=='I'?iu-il+1:MAX(1,n)):0)),depend(n,compute_v,ldz,range,iu,il) :: z
integer intent(out) :: m
! Only returned if range=='A' or range=='I' and il, iu = 1, n
integer intent(out),dimension((compute_v?(2*(*range=='A'||(*range=='I' && iu-il+1==n)?n:0)):0)),depend(n,iu,il,compute_v,range) :: isuppz
Expand Down Expand Up @@ -827,27 +827,32 @@ subroutine <prefix2c>heevr(compute_v,range,lower,n,a,lda,vl,vu,il,iu,abstol,w,z,
integer optional,intent(in),check(compute_v==1||compute_v==0):: compute_v = 1
character optional,intent(in),check(*range=='A'||*range=='V' ||*range=='I') :: range='A'
integer optional,intent(in),check(lower==0||lower==1) :: lower = 0
integer optional,intent(in) :: il=1
integer optional,intent(in),depend(n) :: iu=n
integer optional,intent(in),check(il>=1&&il<=n),depend(n) :: il=1
integer optional,intent(in),check(iu>=il&&iu<=n),depend(il,n) :: iu=n
<ftype2> optional,intent(in) :: vl=0.0
<ftype2> optional,intent(in),check(vu>vl),depend(vl) :: vu=1.0
<ftype2> intent(in) :: abstol=0.0
integer optional,intent(in),depend(n),check(lwork>=1||lwork==-1) :: lwork=max(2*n,1)
integer optional,intent(in),depend(n),check(lwork>=1||lwork==-1) :: lrwork=max(24*n,1)
integer optional,intent(in),depend(n),check(liwork>=1||liwork==-1):: liwork= max(1,10*n)
integer optional,intent(in),depend(n),check(lwork>=max(2*n,1)||lwork==-1) :: lwork=max(2*n,1)
integer optional,intent(in),depend(n),check(lrwork>=max(24*n,1)||lrwork==-1) :: lrwork=max(24*n,1)
integer optional,intent(in),depend(n),check(liwork>=max(1,10*n)||liwork==-1):: liwork= max(1,10*n)

integer intent(hide),depend(a) :: n=shape(a,0)
integer intent(hide),depend(n) :: lda=max(1,n)
integer intent(hide),depend(z) :: ldz=max(1,shape(z,0))
integer intent(hide),depend(compute_v,n) :: ldz=(compute_v?max(1,n):1)
<ftype2c> intent(hide),dimension(lwork),depend(lwork) :: work
<ftype2> intent(hide),dimension(lrwork) :: rwork
<ftype2> intent(hide),dimension(lrwork),depend(lrwork) :: rwork
integer intent(hide),dimension(liwork),depend(liwork) :: iwork

<ftype2> intent(out),dimension(n),depend(n) :: w
<ftype2c> intent(out),dimension((compute_v?MAX(0,n):0),(compute_v?(*range=='I'?iu-il+1:MAX(1,n)):0)),depend(n,compute_v,range,iu,il) :: z
<ftype2c> intent(out),dimension((compute_v?ldz:0),(compute_v?(*range=='I'?iu-il+1:MAX(1,n)):0)),depend(n,ldz,compute_v,range,iu,il) :: z
integer intent(out) :: m
! Only returned if range=='A' or range=='I' and il, iu = 1, n
integer intent(out),dimension((compute_v?(2*(*range=='A'||(*range=='I' && iu-il+1==n)?n:0)):0)),depend(n,iu,il,range,compute_v) :: isuppz
! According to the NetLib LAPACK reference-implementation the arrray isuppz
! is only "implemented" if range=='A' or range=='I' and il, iu = 1, n.
! The array is nevertheless allocated with a size that is large enough
! as precaution against out-of-bound access. The underlying LAPACK
! implementation may as well use the buffer in any case, even if the
! eigenvectors are not requested.
integer intent(out),dimension(2*max(1,n)),depend(n) :: isuppz
integer intent(out) :: info

end subroutine <prefix2c>heevr
Expand Down

0 comments on commit 0fdce89

Please sign in to comment.