The Elk Code
 
Loading...
Searching...
No Matches
gwlocal.f90
Go to the documentation of this file.
1
2! Copyright (C) 2017 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 gwlocal(vmt,vir,bmt,bir)
7use modmain
8use modomp
9implicit none
10! arguments
11real(8), intent(out) :: vmt(npcmtmax,natmtot),vir(ngtc)
12real(8), intent(out) :: bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag)
13! local variables
14integer idm,is,ias,nthd
15integer nrc,nrci,npc
16call holdthd(natmtot+1,nthd)
17!$OMP PARALLEL DEFAULT(SHARED) &
18!$OMP PRIVATE(ias,is) &
19!$OMP PRIVATE(nrc,nrci,npc,idm) &
20!$OMP NUM_THREADS(nthd)
21!$OMP DO SCHEDULE(DYNAMIC)
22do ias=1,natmtot
23 is=idxis(ias)
24 nrc=nrcmt(is)
25 nrci=nrcmti(is)
26 npc=npcmt(is)
27! convert exchange-correlation potential to a coarse radial mesh
28 call rfmtftoc(nrc,nrci,vxcmt(:,ias),vmt(:,ias))
29! negate because V_xc should be removed from the self-energy
30 vmt(1:npc,ias)=-vmt(1:npc,ias)
31! convert to spherical coordinates
32 call rbshtip(nrc,nrci,vmt(:,ias))
33! multiply by radial integration weights
34 call rfcmtwr(nrc,nrci,wr2cmt(:,is),vmt(:,ias))
35end do
36!$OMP END DO NOWAIT
37! multiply -V_xc by the characteristic function and convert to coarse grid
38!$OMP SINGLE
39call rfirftoc(vxcir,vir)
40vir(1:ngtc)=-vir(1:ngtc)
41!$OMP END SINGLE NOWAIT
42! do the same for B_xc in the spin-polarised case
43if (spinpol) then
44 do idm=1,ndmag
45!$OMP DO SCHEDULE(DYNAMIC)
46 do ias=1,natmtot
47 is=idxis(ias)
48 nrc=nrcmt(is)
49 nrci=nrcmti(is)
50 npc=npcmt(is)
51 call rfmtftoc(nrc,nrci,bxcmt(:,ias,idm),bmt(:,ias,idm))
52 bmt(1:npc,ias,idm)=-bmt(1:npc,ias,idm)
53 call rbshtip(nrc,nrci,bmt(:,ias,idm))
54 call rfcmtwr(nrc,nrci,wr2cmt(:,is),bmt(:,ias,idm))
55 end do
56!$OMP END DO NOWAIT
57 end do
58!$OMP SINGLE
59 do idm=1,ndmag
60 call rfirftoc(bxcir(:,idm),bir(:,idm))
61 bir(1:ngtc,idm)=-bir(1:ngtc,idm)
62 end do
63!$OMP END SINGLE
64end if
65!$OMP END PARALLEL
66call freethd(nthd)
67end subroutine
68
subroutine gwlocal(vmt, vir, bmt, bir)
Definition gwlocal.f90:7
real(8), dimension(:,:,:), allocatable bxcmt
Definition modmain.f90:636
logical spinpol
Definition modmain.f90:228
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
real(8), dimension(:,:), allocatable bxcir
Definition modmain.f90:636
real(8), dimension(:), allocatable vxcir
Definition modmain.f90:634
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
real(8), dimension(:,:), allocatable vxcmt
Definition modmain.f90:634
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rbshtip(nr, nri, rfmt)
Definition rbshtip.f90:7
pure subroutine rfcmtwr(nr, nri, wr, rfmt)
Definition rfcmtwr.f90:7
subroutine rfirftoc(rfir, rfirc)
Definition rfirftoc.f90:7
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition rfmtftoc.f90:7