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