The Elk Code
Loading...
Searching...
No Matches
gwrhomag.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2018 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
gwrhomag
7
use
modmain
8
use
modgw
9
use
modmpi
10
use
modomp
11
implicit none
12
! local variables
13
integer
ik,lp,nthd
14
! allocatable arrays
15
real
(8),
allocatable
:: vmt(:,:),vir(:)
16
real
(8),
allocatable
:: bmt(:,:,:),bir(:,:)
17
complex(8)
,
allocatable
:: se(:,:,:)
18
! generate the momentum matrix elements
19
call
genpmat
20
! generate the inverse RPA response function
21
call
epsinv
22
! compute the matrix elements of -V_xc and -B_xc
23
allocate
(vmt(
npcmtmax
,
natmtot
),vir(
ngtc
))
24
if
(
spinpol
)
then
25
allocate
(bmt(
npcmtmax
,
natmtot
,
ndmag
),bir(
ngtc
,
ndmag
))
26
end if
27
call
gwlocal
(vmt,vir,bmt,bir)
28
! synchronise MPI processes
29
call
mpi_barrier(
mpicom
,
ierror
)
30
if
(
mp_mpi
)
write
(*,*)
31
! loop over reduced k-point set
32
call
holdthd
(
nkpt
/
np_mpi
,nthd)
33
!$OMP PARALLEL DEFAULT(SHARED) &
34
!$OMP PRIVATE(se) &
35
!$OMP NUM_THREADS(nthd)
36
allocate
(se(
nstsv
,
nstsv
,0:
nwfm
))
37
!$OMP DO SCHEDULE(DYNAMIC)
38
do
ik=1,
nkpt
39
! distribute among MPI processes
40
if
(mod(ik-1,
np_mpi
) /=
lp_mpi
) cycle
41
!$OMP CRITICAL(gwrhomag_)
42
write
(*,
'("Info(gwrhomag): ",I6," of ",I6," k-points")'
) ik,
nkpt
43
!$OMP END CRITICAL(gwrhomag_)
44
! determine the self-energy at the fermionic frequencies for current k-point
45
call
gwsefmk
(ik,vmt,vir,bmt,bir,se)
46
! write the self-energy to file
47
call
putgwsefm
(ik,se)
48
end do
49
!$OMP END DO
50
deallocate
(se)
51
!$OMP END PARALLEL
52
call
freethd
(nthd)
53
deallocate
(vmt,vir)
54
if
(
spinpol
)
deallocate
(bmt,bir)
55
! synchronise MPI processes
56
call
mpi_barrier(
mpicom
,
ierror
)
57
! determine the GW Fermi energy
58
call
gwefermi
59
! compute the GW density matrices and write the natural orbitals and occupation
60
! numbers to EVECSV.OUT and OCCSV.OUT, respectively
61
call
holdthd
(
nkpt
/
np_mpi
,nthd)
62
!$OMP PARALLEL DO DEFAULT(SHARED) &
63
!$OMP SCHEDULE(DYNAMIC) &
64
!$OMP NUM_THREADS(nthd)
65
do
ik=1,
nkpt
66
! distribute among MPI processes
67
if
(mod(ik-1,
np_mpi
) /=
lp_mpi
) cycle
68
call
gwdmatk
(ik)
69
end do
70
!$OMP END PARALLEL DO
71
call
freethd
(nthd)
72
! broadcast occupation number array to every MPI process
73
if
(
np_mpi
> 1)
then
74
do
ik=1,
nkpt
75
lp=mod(ik-1,
np_mpi
)
76
call
mpi_bcast(
occsv
(:,ik),
nstsv
,mpi_double_precision,lp,
mpicom
,
ierror
)
77
end do
78
end if
79
! write the occupation numbers to file
80
if
(
mp_mpi
)
then
81
do
ik=1,
nkpt
82
call
putoccsv
(
filext
,ik,
occsv
(:,ik))
83
end do
84
end if
85
! synchronise MPI processes
86
call
mpi_barrier(
mpicom
,
ierror
)
87
! determine the density and magnetisation in the usual way
88
call
rhomag
89
end subroutine
90
epsinv
subroutine epsinv
Definition
epsinv.f90:7
genpmat
subroutine genpmat
Definition
genpmat.f90:7
gwdmatk
subroutine gwdmatk(ik)
Definition
gwdmatk.f90:7
gwefermi
subroutine gwefermi
Definition
gwefermi.f90:7
gwlocal
subroutine gwlocal(vmt, vir, bmt, bir)
Definition
gwlocal.f90:7
gwrhomag
subroutine gwrhomag
Definition
gwrhomag.f90:7
gwsefmk
subroutine gwsefmk(ikp, vmt, vir, bmt, bir, se)
Definition
gwsefmk.f90:7
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::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::occsv
real(8), dimension(:,:), allocatable occsv
Definition
modmain.f90:902
modmain::ndmag
integer ndmag
Definition
modmain.f90:238
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
putoccsv
subroutine putoccsv(fext, ik, occsvp)
Definition
putoccsv.f90:7
rhomag
subroutine rhomag
Definition
rhomag.f90:7
gwrhomag.f90
Generated by
1.9.8