The Elk Code
 
Loading...
Searching...
No Matches
gwbandstr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2017 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine gwbandstr
7use modmain
8use modgw
9use modmpi
10use modomp
11implicit none
12! local variables
13integer ip,iw
14real(8) dw
15! allocatable arrays
16real(8), allocatable :: vmt(:,:),vir(:)
17real(8), allocatable :: bmt(:,:,:),bir(:,:)
18real(8), allocatable :: wr(:),sf(:)
19complex(8), allocatable :: se(:,:,:)
20! store original parameters
21vkloff0(:)=vkloff(:)
22! initialise universal variables
23call init0
24call init1
25call init2
26call init3
27! read density and potentials from file
28call readstate
29! generate the core wavefunctions and densities
30call gencore
31! Fourier transform Kohn-Sham potential to G-space
32call genvsig
33! generate k-points along a path for band structure plots
35! compute the matrix elements of -V_xc and -B_xc
36allocate(vmt(npcmtmax,natmtot),vir(ngtc))
37if (spinpol) then
38 allocate(bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag))
39end if
40call gwlocal(vmt,vir,bmt,bir)
41! real axis frequencies
42allocate(wr(nwplot))
43dw=(wplot(2)-wplot(1))/dble(nwplot)
44do iw=1,nwplot
45 wr(iw)=dw*dble(iw-1)+wplot(1)
46end do
47allocate(sf(nwplot),se(nstsv,nstsv,0:nwfm))
48if (mp_mpi) then
49 if (ip01d == 1) then
50 open(85,file='GWBAND.OUT',form='FORMATTED',action='WRITE')
51 write(85,'(2I6," : grid size")') nwplot,npp1d
52 else
53 open(85,file='GWBAND.OUT',form='FORMATTED',action='WRITE',position='APPEND')
54 end if
55end if
56! loop over plot points along path
57do ip=ip01d,npp1d
58 if (mp_mpi) then
59 write(*,'("Info(gwbandstr): ",I6," of ",I6," plot points")') ip,npp1d
60 end if
61! reset the OpenMP thread variables
62 call omp_reset
63! change the k-point offset
64 vkloff(:)=vplp1d(:,ip)*ngridk(:)
65! generate the new k-point set
66 call init1
67! determine the Kohn-Sham ground-state for this k-point offset
68 call linengy
69 call genapwlofr
70 call gensocfr
71 call genevfsv
72 call occupy
73! write the momentum matrix elements to file
74 call genpmat
75! generate the inverse dielectric function and write to file
76 call epsinv
77! determine the self-energy for the first k-point
78 if (mp_mpi) then
79 write(*,'("Info(gwbandstr): calculating self-energy for first k-point")')
80 end if
81 call gwsefmk(1,vmt,vir,bmt,bir,se)
82! solve the Dyson equation on the real axis
83 call dysonr(1,wr,se,sf)
84 if (mp_mpi) then
85 do iw=1,nwplot
86 write(85,'(3G18.10)') dpp1d(ip),wr(iw),sf(iw)
87 end do
88 flush(85)
89 end if
90! synchronise MPI processes
91 call mpi_barrier(mpicom,ierror)
92end do
93deallocate(vmt,vir,wr,sf,se)
94if (spinpol) deallocate(bmt,bir)
95if (mp_mpi) then
96 close(85)
97 write(*,*)
98 write(*,'("Info(gwbandstr):")')
99 write(*,'(" GW spectral function band structure written to GWBAND.OUT")')
100end if
101! restore original input parameters
102vkloff(:)=vkloff0(:)
103! synchronise MPI processes
104call mpi_barrier(mpicom,ierror)
105end subroutine
106
subroutine dysonr(ik, wr, sem, sf)
Definition dysonr.f90:7
subroutine epsinv
Definition epsinv.f90:7
subroutine genapwlofr
Definition genapwlofr.f90:7
subroutine gencore
Definition gencore.f90:10
subroutine genevfsv
Definition genevfsv.f90:7
subroutine genpmat
Definition genpmat.f90:7
subroutine gensocfr
Definition gensocfr.f90:7
subroutine genvsig
Definition genvsig.f90:10
subroutine gwbandstr
Definition gwbandstr.f90:7
subroutine gwlocal(vmt, vir, bmt, bir)
Definition gwlocal.f90:7
subroutine gwsefmk(ikp, vmt, vir, bmt, bir, se)
Definition gwsefmk.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init2
Definition init2.f90:7
subroutine init3
Definition init3.f90:7
subroutine linengy
Definition linengy.f90:10
Definition modgw.f90:6
integer nwfm
Definition modgw.f90:19
integer npp1d
Definition modmain.f90:1113
real(8), dimension(3) vkloff0
Definition modmain.f90:450
integer nvp1d
Definition modmain.f90:1111
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
logical spinpol
Definition modmain.f90:228
integer nwplot
Definition modmain.f90:1070
real(8), dimension(:,:), allocatable vvlp1d
Definition modmain.f90:1117
integer ip01d
Definition modmain.f90:1115
integer ngtc
Definition modmain.f90:392
integer natmtot
Definition modmain.f90:40
real(8), dimension(2) wplot
Definition modmain.f90:1076
integer npcmtmax
Definition modmain.f90:216
integer, dimension(3) ngridk
Definition modmain.f90:448
real(8), dimension(:,:), allocatable vplp1d
Definition modmain.f90:1121
real(8), dimension(:), allocatable dvp1d
Definition modmain.f90:1119
integer nstsv
Definition modmain.f90:886
real(8), dimension(:), allocatable dpp1d
Definition modmain.f90:1123
real(8), dimension(3) vkloff
Definition modmain.f90:450
integer ndmag
Definition modmain.f90:238
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
logical mp_mpi
Definition modmpi.f90:17
subroutine omp_reset
Definition modomp.f90:71
subroutine occupy
Definition occupy.f90:10
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition plotpt1d.f90:10
subroutine readstate
Definition readstate.f90:10