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