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