The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine genwgw
7use modmain
8use modgw
9implicit none
10! local variables
11integer ik,iw,jw,n
12real(8) de,t0,t1
14if (wmaxgw <= 0.d0) then
15! read the Fermi energy from file
16 call readfermi
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
23end if
24! number of Matsubara frequencies
25t1=pi*t0
26nwgw=2*nint(wmaxgw/t1)
27nwgw=max(nwgw,2)
28call nfftifc(npfftw,1,nwgw)
29! determine integer ranges for grid
30intwgw(1)=nwgw/2-nwgw+1
31intwgw(2)=nwgw/2
32if (allocated(iwfft)) deallocate(iwfft)
33allocate(iwfft(intwgw(1):intwgw(2)))
34if (allocated(wgw)) deallocate(wgw)
35allocate(wgw(intwgw(1):intwgw(2)))
36do 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
44end do
45n=minval(abs(intwgw(1:2)))
46if (n == 0) then
47 write(*,*)
48 write(*,'("Error(genwgw): not enough Matsubara frequencies")')
49 write(*,'("Increase wmaxgw")')
50 write(*,*)
51 stop
52end if
53if (mod(n,2) == 0) then
54 nwbs=n
55 nwfm=n-1
56else
57 nwfm=n
58 nwbs=n-1
59end if
60! store the complex fermionic frequencies
61if (allocated(wfm)) deallocate(wfm)
62allocate(wfm(0:nwfm))
63do iw=-nwfm,nwfm,2
64 jw=(iw+nwfm)/2
65 wfm(jw)=cmplx(0.d0,wgw(iw),8)
66end do
67! store the complex response function frequencies
68nwrf=nwbs+1
69if (allocated(wrf)) deallocate(wrf)
70allocate(wrf(nwrf))
71do iw=-nwbs,nwbs,2
72 jw=(iw+nwbs)/2+1
73 wrf(jw)=cmplx(0.d0,wgw(iw),8)
74end do
75end subroutine
76
subroutine genwgw
Definition genwgw.f90:7
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
Definition modgw.f90:6
real(8) wmaxgw
Definition modgw.f90:9
complex(8), dimension(:), allocatable wfm
Definition modgw.f90:25
integer nwgw
Definition modgw.f90:11
integer, dimension(2) intwgw
Definition modgw.f90:13
integer, dimension(:), allocatable iwfft
Definition modgw.f90:17
real(8), dimension(:), allocatable wgw
Definition modgw.f90:23
integer npfftw
Definition modgw.f90:15
integer nwfm
Definition modgw.f90:19
integer nwbs
Definition modgw.f90:21
real(8) efermi
Definition modmain.f90:904
real(8), parameter pi
Definition modmain.f90:1229
character(256) filext
Definition modmain.f90:1300
integer nwrf
Definition modmain.f90:1165
integer nkpt
Definition modmain.f90:461
real(8), parameter kboltz
Definition modmain.f90:1262
integer nstsv
Definition modmain.f90:886
complex(8), dimension(:), allocatable wrf
Definition modmain.f90:1167
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) tempk
Definition modmain.f90:684
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
subroutine nfftifc(np, nd, n)
Definition nfftifc.f90:10
subroutine readfermi
Definition readfermi.f90:10