The Elk Code
 
Loading...
Searching...
No Matches
vblocalu.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
6subroutine vblocalu
7use modmain
8use modulr
9use modomp
10implicit none
11! local variables
12integer ifq,idm,is,ias
13integer nrc,nrci,npc,nthd
14! automatic arrays
15real(8) rfmt1(npcmtmax),rfmt2(npcmtmax)
16! subtract the normal Kohn-Sham potential for Q=0
17do ias=1,natmtot
18 is=idxis(ias)
19 nrc=nrcmt(is)
20 nrci=nrcmti(is)
21 npc=npcmt(is)
22 call rfmtftoc(nrc,nrci,vsmt(:,ias),rfmt1)
23 call rbsht(nrc,nrci,rfmt1,rfmt2)
24 vsqmt(1:npc,ias,1)=vsqmt(1:npc,ias,1)-rfmt2(1:npc)
25end do
26vsqir(1:ngtot,1)=vsqir(1:ngtot,1)-vsir(1:ngtot)
27call holdthd(2*nfqrz,nthd)
28!$OMP PARALLEL DO DEFAULT(SHARED) &
29!$OMP PRIVATE(ias,is,nrc,nrci) &
30!$OMP SCHEDULE(DYNAMIC) &
31!$OMP NUM_THREADS(nthd)
32do ifq=1,nfqrz
33 do ias=1,natmtot
34 is=idxis(ias)
35 nrc=nrcmt(is)
36 nrci=nrcmti(is)
37! multiply the muffin-tin potential by the radial integration weights
38 call zfcmtwr(nrc,nrci,wr2cmt(:,is),vsqmt(:,ias,ifq))
39 end do
40end do
41!$OMP END PARALLEL DO
42call freethd(nthd)
43if (.not.spinpol) return
44! subtract the normal Kohn-Sham magnetic field for Q=0 in the muffin-tins
45do idm=1,ndmag
46 do ias=1,natmtot
47 is=idxis(ias)
48 npc=npcmt(is)
49 bsqmt(1:npc,ias,idm,1)=bsqmt(1:npc,ias,idm,1)-bsmt(1:npc,ias,idm)
50 end do
51end do
52call holdthd(2*nfqrz,nthd)
53!$OMP PARALLEL DO DEFAULT(SHARED) &
54!$OMP PRIVATE(idm,ias,is,nrc,nrci) &
55!$OMP SCHEDULE(DYNAMIC) &
56!$OMP NUM_THREADS(nthd)
57do ifq=1,nfqrz
58 do idm=1,ndmag
59 do ias=1,natmtot
60 is=idxis(ias)
61 nrc=nrcmt(is)
62 nrci=nrcmti(is)
63! multiply the muffin-tin field by the radial integration weights
64 call zfcmtwr(nrc,nrci,wr2cmt(:,is),bsqmt(:,ias,idm,ifq))
65 end do
66 end do
67end do
68!$OMP END PARALLEL DO
69call freethd(nthd)
70! subtract the normal Kohn-Sham magnetic field for Q=0 in the interstitial
71do idm=1,ndmag
72 bsqir(1:ngtot,idm,1)=bsqir(1:ngtot,idm,1)-bsir(1:ngtot,idm)
73end do
74end subroutine
75
integer ngtot
Definition modmain.f90:390
integer nfqrz
Definition modmain.f90:539
logical spinpol
Definition modmain.f90:228
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer natmtot
Definition modmain.f90:40
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
real(8), dimension(:,:), pointer, contiguous vsmt
Definition modmain.f90:649
real(8), dimension(:,:), allocatable bsir
Definition modmain.f90:658
real(8), dimension(:), allocatable vsir
Definition modmain.f90:651
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition modmain.f90:656
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer ndmag
Definition modmain.f90:238
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
complex(8), dimension(:,:,:,:), pointer, contiguous bsqmt
Definition modulr.f90:84
complex(8), dimension(:,:,:), pointer, contiguous vsqmt
Definition modulr.f90:83
complex(8), dimension(:,:,:), pointer, contiguous bsqir
Definition modulr.f90:84
complex(8), dimension(:,:), pointer, contiguous vsqir
Definition modulr.f90:83
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition rbsht.f90:7
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition rfmtftoc.f90:7
subroutine vblocalu
Definition vblocalu.f90:7
pure subroutine zfcmtwr(nr, nri, wr, zfmt)
Definition zfcmtwr.f90:7