The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine getgwsefm(ik,se)
7use modmain
8use modgw
10implicit none
11! arguments
12integer, intent(in) :: ik
13complex(8), intent(out) :: se(nstsv,nstsv,0:nwfm)
14! local variables
15logical tgs
16integer recl,nstsv_,nwfm_
17real(8) vkl_(3),t1
18!$OMP CRITICAL(u280)
19! read from RAM disk if required
20if (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
24end if
25! find the record length
26inquire(iolength=recl) vkl_,nstsv_,nwfm_,se
27open(280,file='GWSEFM.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
28read(280,rec=ik) vkl_,nstsv_,nwfm_,se
29close(280)
3010 continue
31!$OMP END CRITICAL(u280)
32t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
33if (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
40end if
41if (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
48end if
49if (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
56end if
57end subroutine
58
subroutine getgwsefm(ik, se)
Definition getgwsefm.f90:7
Definition modgw.f90:6
real(8) epslat
Definition modmain.f90:24
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
type(file_t), dimension(:), allocatable, private file
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
logical ramdisk
Definition modramdisk.f90:9