The Elk Code
getgwsefm.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2017 A. Davydov, A. Sanna, J. K. Dewhurst, S. Sharma and
3 ! E. K. U. Gross. This file is distributed under the terms of the GNU General
4 ! Public License. See the file COPYING for license details.
5 
6 subroutine getgwsefm(ik,se)
7 use modmain
8 use modgw
9 use modramdisk
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik
13 complex(8), intent(out) :: se(nstsv,nstsv,0:nwfm)
14 ! local variables
15 logical tgs
16 integer recl,nstsv_,nwfm_
17 real(8) vkl_(3),t1
18 !$OMP CRITICAL(u280)
19 ! read from RAM disk if required
20 if (ramdisk) then
21  call getrd('GWSEFM.OUT',ik,tgs,v1=vkl_,n1=nstsv_,n2=nwfm_, &
22  nzv=nstsv*nstsv*(nwfm+1),zva=se)
23  if (tgs) goto 10
24 end if
25 ! find the record length
26 inquire(iolength=recl) vkl_,nstsv_,nwfm_,se
27 open(280,file='GWSEFM.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
28 read(280,rec=ik) vkl_,nstsv_,nwfm_,se
29 close(280)
30 10 continue
31 !$OMP END CRITICAL(u280)
32 t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
33 if (t1 > epslat) then
34  write(*,*)
35  write(*,'("Error(getgwsefm): differing vectors for k-point ",I8)') ik
36  write(*,'(" current : ",3G18.10)') vkl(:,ik)
37  write(*,'(" GWSEFM.OUT : ",3G18.10)') vkl_
38  write(*,*)
39  stop
40 end if
41 if (nstsv /= nstsv_) then
42  write(*,*)
43  write(*,'("Error(getgwsefm): differing nstsv for k-point ",I8)') ik
44  write(*,'(" current : ",I8)') nstsv
45  write(*,'(" GWSEFM.OUT : ",I8)') nstsv_
46  write(*,*)
47  stop
48 end if
49 if (nwfm /= nwfm_) then
50  write(*,*)
51  write(*,'("Error(getgwsefm): differing nwfm for k-point ",I8)') ik
52  write(*,'(" current : ",I8)') nwfm
53  write(*,'(" GWSEFM.OUT : ",I8)') nwfm_
54  write(*,*)
55  stop
56 end if
57 end subroutine
58 
subroutine getgwsefm(ik, se)
Definition: getgwsefm.f90:7
logical ramdisk
Definition: modramdisk.f90:9
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
Definition: modramdisk.f90:214
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
Definition: modgw.f90:6
real(8) epslat
Definition: modmain.f90:24