The Elk Code
dsygvxi.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2026 J. K. Dewhurst and S. Sharma.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine dsygvxi(n,m,ld1,a,b,w,ld2,z,lwork,work)
7 use modomp
8 implicit none
9 ! arguments
10 integer, intent(in) :: n,m,ld1
11 real(8), intent(in) :: a(ld1,*),b(ld1,*)
12 real(8), intent(out) :: w(m)
13 integer, intent(in) :: ld2
14 real(8), intent(out) :: z(ld2,m)
15 integer, intent(in) :: lwork
16 real(8), intent(out) :: work(lwork)
17 ! local variables
18 integer nts,p,info,nthd
19 real(8) vl,vu
20 ! enable MKL parallelism
21 call holdthd(maxthdmkl,nthd)
22 nts=mkl_set_num_threads_local(nthd)
23 block
24 integer iwork(5*n),ifail(n)
25 real(8) wn(n)
26 ! find the first m eigenvalues and eigenvectors
27 call dsygvx(1,'V','I','U',n,a,ld1,b,ld1,vl,vu,1,m,-1.d0,p,wn,z,ld2,work,lwork, &
28  iwork,ifail,info)
29 w(1:m)=wn(1:m)
30 end block
31 nts=mkl_set_num_threads_local(0)
32 call freethd(nthd)
33 if (info /= 0) then
34  write(*,*)
35  write(*,'("Error(dsygvxi): diagonalisation failed")')
36  write(*,'(" DSYGVX returned INFO = ",I0)') info
37  write(*,*)
38  stop
39 end if
40 end subroutine
41 
Definition: modomp.f90:6
integer maxthdmkl
Definition: modomp.f90:15
subroutine dsygvxi(n, m, ld1, a, b, w, ld2, z, lwork, work)
Definition: dsygvxi.f90:7
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78