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