The Elk Code
genwgw.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2016 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 genwgw
7 use modmain
8 use modgw
9 implicit none
10 ! local variables
11 integer ik,iw,jw,n
12 real(8) de,t0,t1
13 t0=kboltz*tempk
14 if (wmaxgw <= 0.d0) then
15 ! read the Fermi energy from file
16  call readefm
17 ! find the maximum eigenvalue range over all k-points
18  do ik=1,nkpt
19  call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
20  end do
21  de=maxval(abs(evalsv(1:nstsv,1:nkpt)-efermi))
22  wmaxgw=abs(wmaxgw)*de
23 end if
24 ! number of Matsubara frequencies
25 t1=pi*t0
26 nwgw=2*nint(wmaxgw/t1)
27 nwgw=max(nwgw,2)
28 call nfftifc(npfftw,1,nwgw)
29 ! determine integer ranges for grid
30 intwgw(1)=nwgw/2-nwgw+1
31 intwgw(2)=nwgw/2
32 if (allocated(iwfft)) deallocate(iwfft)
33 allocate(iwfft(intwgw(1):intwgw(2)))
34 if (allocated(wgw)) deallocate(wgw)
35 allocate(wgw(intwgw(1):intwgw(2)))
36 do iw=intwgw(1),intwgw(2)
37  if (iw >= 0) then
38  jw=iw
39  else
40  jw=nwgw+iw
41  end if
42  iwfft(iw)=jw+1
43  wgw(iw)=dble(iw)*t1
44 end do
45 n=minval(abs(intwgw(1:2)))
46 if (n == 0) then
47  write(*,*)
48  write(*,'("Error(genwgw): not enough Matsubara frequencies")')
49  write(*,'("Increase wmaxgw")')
50  write(*,*)
51  stop
52 end if
53 if (mod(n,2) == 0) then
54  nwbs=n
55  nwfm=n-1
56 else
57  nwfm=n
58  nwbs=n-1
59 end if
60 ! store the complex fermionic frequencies
61 if (allocated(wfm)) deallocate(wfm)
62 allocate(wfm(0:nwfm))
63 do iw=-nwfm,nwfm,2
64  jw=(iw+nwfm)/2
65  wfm(jw)=cmplx(0.d0,wgw(iw),8)
66 end do
67 ! store the complex response function frequencies
68 nwrf=nwbs+1
69 if (allocated(wrf)) deallocate(wrf)
70 allocate(wrf(nwrf))
71 do iw=-nwbs,nwbs,2
72  jw=(iw+nwbs)/2+1
73  wrf(jw)=cmplx(0.d0,wgw(iw),8)
74 end do
75 end subroutine
76 
real(8) efermi
Definition: modmain.f90:907
subroutine genwgw
Definition: genwgw.f90:7
integer, dimension(:), allocatable iwfft
Definition: modgw.f90:17
integer nwrf
Definition: modmain.f90:1168
character(256) filext
Definition: modmain.f90:1301
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
integer nkpt
Definition: modmain.f90:461
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition: getevalsv.f90:7
real(8), parameter kboltz
Definition: modmain.f90:1263
complex(8), dimension(:), allocatable wrf
Definition: modmain.f90:1170
integer nwbs
Definition: modgw.f90:21
real(8), parameter pi
Definition: modmain.f90:1232
integer nstsv
Definition: modmain.f90:889
subroutine readefm
Definition: readefm.f90:10
integer, dimension(2) intwgw
Definition: modgw.f90:13
real(8) wmaxgw
Definition: modgw.f90:9
real(8) tempk
Definition: modmain.f90:684
integer nwgw
Definition: modgw.f90:11
real(8), dimension(:), allocatable wgw
Definition: modgw.f90:23
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
Definition: modgw.f90:6
integer npfftw
Definition: modgw.f90:15
integer nwfm
Definition: modgw.f90:19
complex(8), dimension(:), allocatable wfm
Definition: modgw.f90:25
subroutine nfftifc(np, nd, n)
Definition: nfftifc.f90:10