The Elk Code
vblocal.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine vblocal(vmt,vir,bmt,bir)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 real(8), intent(out) :: vmt(npcmtmax,natmtot),vir(ngtc)
12 real(8), intent(out) :: bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag)
13 ! local variables
14 integer idm,is,ias,nthd
15 integer nrc,nrci,npc
16 call 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)
22 do ias=1,natmtot
23  is=idxis(ias)
24  nrc=nrcmt(is)
25  nrci=nrcmti(is)
26 ! convert muffin-tin Kohn-Sham potential to coarse radial mesh
27  call rfmtftoc(nrc,nrci,vsmt(:,ias),vmt(:,ias))
28 ! convert to spherical coordinates
29  call rbshtip(nrc,nrci,vmt(:,ias))
30 ! multiply by radial integration weights
31  call rfcmtwr(nrc,nrci,wr2cmt(:,is),vmt(:,ias))
32 end do
33 !$OMP END DO NOWAIT
34 ! multiply interstitial Kohn-Sham potential by characteristic function and
35 ! convert to coarse grid
36 !$OMP SINGLE
37 call rfirftoc(vsir,vir)
38 !$OMP END SINGLE NOWAIT
39 ! repeat for the Kohn-Sham magnetic field
40 if (spinpol) then
41  do idm=1,ndmag
42 !$OMP DO SCHEDULE(DYNAMIC)
43  do ias=1,natmtot
44  is=idxis(ias)
45  nrc=nrcmt(is)
46  nrci=nrcmti(is)
47  npc=npcmt(is)
48  bmt(1:npc,ias,idm)=bsmt(1:npc,ias,idm)
49  call rfcmtwr(nrc,nrci,wr2cmt(:,is),bmt(:,ias,idm))
50  end do
51 !$OMP END DO NOWAIT
52  end do
53 ! convert Kohn-Sham magnetic field to coarse grid
54 !$OMP SINGLE
55  do idm=1,ndmag
56  call rfirftoc(bsir(:,idm),bir(:,idm))
57  end do
58 !$OMP END SINGLE
59 end if
60 !$OMP END PARALLEL
61 call freethd(nthd)
62 end subroutine
63 
subroutine rbshtip(nr, nri, rfmt)
Definition: rbshtip.f90:7
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
logical spinpol
Definition: modmain.f90:228
Definition: modomp.f90:6
real(8), dimension(:), allocatable vsir
Definition: modmain.f90:651
pure subroutine rfcmtwr(nr, nri, wr, rfmt)
Definition: rfcmtwr.f90:7
real(8), dimension(:,:), allocatable bsir
Definition: modmain.f90:658
real(8), dimension(:,:), allocatable wr2cmt
Definition: modmain.f90:189
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition: rfmtftoc.f90:7
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
subroutine vblocal(vmt, vir, bmt, bir)
Definition: vblocal.f90:7
subroutine rfirftoc(rfir, rfirc)
Definition: rfirftoc.f90:7
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition: modmain.f90:656