The Elk Code
getoccsv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 getoccsv(fext,ikp,vpl,occsvp)
7 use modmain
8 implicit none
9 ! arguments
10 character(*), intent(in) :: fext
11 integer, intent(in) :: ikp
12 real(8), intent(in) :: vpl(3)
13 real(8), intent(out) :: occsvp(nstsv)
14 ! local variables
15 integer isym,ik
16 integer recl,nstsv_
17 real(8) vkl_(3),t1
18 if (ikp > 0) then
19  ik=ikp
20 else
21 ! find the k-point number
22  call findkpt(vpl,isym,ik)
23 end if
24 ! find the record length
25 inquire(iolength=recl) vkl_,nstsv_,occsvp
26 !$OMP CRITICAL(u208)
27 open(208,file='OCCSV'//trim(fext),form='UNFORMATTED',access='DIRECT',recl=recl)
28 read(208,rec=ik) vkl_,nstsv_,occsvp
29 close(208)
30 !$OMP END CRITICAL(u208)
31 t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
32 if (t1 > epslat) then
33  write(*,*)
34  write(*,'("Error(getoccsv): differing vectors for k-point ",I8)') ik
35  write(*,'(" current : ",3G18.10)') vkl(:,ik)
36  write(*,'(" OCCSV.OUT : ",3G18.10)') vkl_
37  write(*,*)
38  stop
39 end if
40 if (nstsv /= nstsv_) then
41  write(*,*)
42  write(*,'("Error(getoccsv): differing nstsv for k-point ",I8)') ik
43  write(*,'(" current : ",I8)') nstsv
44  write(*,'(" OCCSV.OUT : ",I8)') nstsv_
45  write(*,*)
46  stop
47 end if
48 end subroutine
49 
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition: getoccsv.f90:7
subroutine findkpt(vpl, isym, ik)
Definition: findkpt.f90:7