The Elk Code
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(:),bmt(:,:,:),bir(:,:)
16 complex(8), allocatable :: se(:,:,:)
17 ! generate the momentum matrix elements
18 call genpmat
19 ! generate the inverse RPA response function
20 call epsinv
21 ! compute the matrix elements of -V_xc and -B_xc
22 allocate(vmt(npcmtmax,natmtot),vir(ngtc))
23 if (spinpol) then
24  allocate(bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag))
25 else
26  allocate(bmt(1,1,1),bir(1,1))
27 end if
28 call gwlocal(vmt,vir,bmt,bir)
29 ! synchronise MPI processes
30 call mpi_barrier(mpicom,ierror)
31 if (mp_mpi) write(*,*)
32 ! loop over reduced k-point set
33 call holdthd(nkpt/np_mpi,nthd)
34 !$OMP PARALLEL DEFAULT(SHARED) &
35 !$OMP PRIVATE(se) &
36 !$OMP NUM_THREADS(nthd)
37 allocate(se(nstsv,nstsv,0:nwfm))
38 !$OMP DO SCHEDULE(DYNAMIC)
39 do ik=1,nkpt
40 ! distribute among MPI processes
41  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
42 !$OMP CRITICAL(gwrhomag_)
43  write(*,'("Info(gwrhomag): ",I6," of ",I6," k-points")') ik,nkpt
44 !$OMP END CRITICAL(gwrhomag_)
45 ! determine the self-energy at the fermionic frequencies for current k-point
46  call gwsefmk(ik,vmt,vir,bmt,bir,se)
47 ! write the self-energy to file
48  call putgwsefm(ik,se)
49 end do
50 !$OMP END DO
51 deallocate(se)
52 !$OMP END PARALLEL
53 call freethd(nthd)
54 deallocate(vmt,vir,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 
integer npcmtmax
Definition: modmain.f90:216
character(256) filext
Definition: modmain.f90:1301
logical mp_mpi
Definition: modmpi.f90:17
integer ngtc
Definition: modmain.f90:392
subroutine epsinv
Definition: epsinv.f90:7
logical spinpol
Definition: modmain.f90:228
integer nkpt
Definition: modmain.f90:461
integer ndmag
Definition: modmain.f90:238
Definition: modomp.f90:6
subroutine rhomag
Definition: rhomag.f90:7
subroutine gwlocal(vmt, vir, bmt, bir)
Definition: gwlocal.f90:7
integer np_mpi
Definition: modmpi.f90:13
subroutine gwrhomag
Definition: gwrhomag.f90:7
integer nstsv
Definition: modmain.f90:889
subroutine gwsefmk(ikp, vmt, vir, bmt, bir, se)
Definition: gwsefmk.f90:7
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
subroutine gwefermi
Definition: gwefermi.f90:7
subroutine putgwsefm(ik, se)
Definition: putgwsefm.f90:7
Definition: modgw.f90:6
Definition: modmpi.f90:6
subroutine genpmat
Definition: genpmat.f90:7
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer nwfm
Definition: modgw.f90:19
subroutine putoccsv(fext, ik, occsvp)
Definition: putoccsv.f90:7
integer natmtot
Definition: modmain.f90:40
integer mpicom
Definition: modmpi.f90:11
subroutine gwdmatk(ik)
Definition: gwdmatk.f90:7
integer ierror
Definition: modmpi.f90:19