The Elk Code
Loading...
Searching...
No Matches
gwsefm.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
gwsefm
7
use
modmain
8
use
modgw
9
use
modmpi
10
use
modomp
11
implicit none
12
! local variables
13
integer
ik,nthd
14
! allocatable arrays
15
real
(8),
allocatable
:: vmt(:,:),vir(:)
16
real
(8),
allocatable
:: bmt(:,:,:),bir(:,:)
17
complex(8)
,
allocatable
:: se(:,:,:)
18
! initialise universal variables
19
call
init0
20
call
init1
21
call
init2
22
call
init3
23
! read density and potentials from file
24
call
readstate
25
! generate the core wavefunctions and densities
26
call
gencore
27
! find the new linearisation energies
28
call
linengy
29
! generate the APW radial functions
30
call
genapwfr
31
! generate the local-orbital radial functions
32
call
genlofr
33
! get the eigenvalues and occupation numbers from file
34
do
ik=1,
nkpt
35
call
getevalsv
(
filext
,ik,
vkl
(:,ik),
evalsv
(:,ik))
36
call
getoccsv
(
filext
,ik,
vkl
(:,ik),
occsv
(:,ik))
37
end do
38
! write the momentum matrix elements in the second-variational basis to file
39
call
genpmat
40
! generate the inverse dielectric function and write to file if required
41
if
(
task
/= 601)
call
epsinv
42
! compute the matrix elements of -V_xc and -B_xc
43
allocate
(vmt(
npcmtmax
,
natmtot
),vir(
ngtc
))
44
if
(
spinpol
)
then
45
allocate
(bmt(
npcmtmax
,
natmtot
,
ndmag
),bir(
ngtc
,
ndmag
))
46
end if
47
call
gwlocal
(vmt,vir,bmt,bir)
48
if
(
mp_mpi
)
write
(*,*)
49
! synchronise MPI processes
50
call
mpi_barrier(
mpicom
,
ierror
)
51
! loop over reduced k-point set
52
call
holdthd
(
nkpt
/
np_mpi
,nthd)
53
!$OMP PARALLEL DEFAULT(SHARED) &
54
!$OMP PRIVATE(se) &
55
!$OMP NUM_THREADS(nthd)
56
allocate
(se(
nstsv
,
nstsv
,0:
nwfm
))
57
!$OMP DO SCHEDULE(DYNAMIC)
58
do
ik=1,
nkpt
59
! distribute among MPI processes
60
if
(mod(ik-1,
np_mpi
) /=
lp_mpi
) cycle
61
!$OMP CRITICAL(gwsefm_)
62
write
(*,
'("Info(gwsefm): ",I6," of ",I6," k-points")'
) ik,
nkpt
63
!$OMP END CRITICAL(gwsefm_)
64
! determine the self-energy at the fermionic frequencies for current k-point
65
call
gwsefmk
(ik,vmt,vir,bmt,bir,se)
66
! write the self-energy to file
67
call
putgwsefm
(ik,se)
68
end do
69
!$OMP END DO
70
deallocate
(se)
71
!$OMP END PARALLEL
72
call
freethd
(nthd)
73
if
(
mp_mpi
)
then
74
write
(*,*)
75
write
(*,
'("Info(gwsefm): GW self-energy at the fermionic frequencies &
76
&written to GWSEFM.OUT")'
)
77
end if
78
deallocate
(vmt,vir)
79
if
(
spinpol
)
deallocate
(bmt,bir)
80
end subroutine
81
epsinv
subroutine epsinv
Definition
epsinv.f90:7
genapwfr
subroutine genapwfr
Definition
genapwfr.f90:10
gencore
subroutine gencore
Definition
gencore.f90:10
genlofr
subroutine genlofr
Definition
genlofr.f90:10
genpmat
subroutine genpmat
Definition
genpmat.f90:7
getevalsv
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition
getevalsv.f90:7
getoccsv
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition
getoccsv.f90:7
gwlocal
subroutine gwlocal(vmt, vir, bmt, bir)
Definition
gwlocal.f90:7
gwsefm
subroutine gwsefm
Definition
gwsefm.f90:7
gwsefmk
subroutine gwsefmk(ikp, vmt, vir, bmt, bir, se)
Definition
gwsefmk.f90:7
init0
subroutine init0
Definition
init0.f90:10
init1
subroutine init1
Definition
init1.f90:10
init2
subroutine init2
Definition
init2.f90:7
init3
subroutine init3
Definition
init3.f90:7
linengy
subroutine linengy
Definition
linengy.f90:10
modgw
Definition
modgw.f90:6
modgw::nwfm
integer nwfm
Definition
modgw.f90:19
modmain
Definition
modmain.f90:6
modmain::filext
character(256) filext
Definition
modmain.f90:1300
modmain::spinpol
logical spinpol
Definition
modmain.f90:228
modmain::ngtc
integer ngtc
Definition
modmain.f90:392
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::natmtot
integer natmtot
Definition
modmain.f90:40
modmain::npcmtmax
integer npcmtmax
Definition
modmain.f90:216
modmain::task
integer task
Definition
modmain.f90:1298
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::vkl
real(8), dimension(:,:), allocatable vkl
Definition
modmain.f90:471
modmain::occsv
real(8), dimension(:,:), allocatable occsv
Definition
modmain.f90:902
modmain::ndmag
integer ndmag
Definition
modmain.f90:238
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition
modmain.f90:918
modmpi
Definition
modmpi.f90:6
modmpi::lp_mpi
integer lp_mpi
Definition
modmpi.f90:15
modmpi::ierror
integer ierror
Definition
modmpi.f90:19
modmpi::mpicom
integer mpicom
Definition
modmpi.f90:11
modmpi::np_mpi
integer np_mpi
Definition
modmpi.f90:13
modmpi::mp_mpi
logical mp_mpi
Definition
modmpi.f90:17
modomp
Definition
modomp.f90:6
modomp::holdthd
subroutine holdthd(nloop, nthd)
Definition
modomp.f90:78
modomp::freethd
subroutine freethd(nthd)
Definition
modomp.f90:106
putgwsefm
subroutine putgwsefm(ik, se)
Definition
putgwsefm.f90:7
readstate
subroutine readstate
Definition
readstate.f90:10
gwsefm.f90
Generated by
1.9.8