The Elk Code
hflocal.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2012 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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 hflocal(vmt,vir,bmt,bir)
7 use modmain
8 implicit none
9 ! arguments
10 real(8), intent(out) :: vmt(npcmtmax,natmtot),vir(ngtc)
11 real(8), intent(out) :: bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag)
12 ! local variables
13 integer idm,is,ias
14 integer np,nrc,nrci
15 ! automatic arrays
16 real(8) rfmt(npmtmax),rfir(ngtot)
17 ! compute the Coulomb potential
18 call potcoul
19 ! convert to spherical coordinates and store in output arrays
20 if (hybrid) then
21 ! hybrid functional case
24  do ias=1,natmtot
25  is=idxis(ias)
26  nrc=nrcmt(is)
27  nrci=nrcmti(is)
28  np=npmt(is)
29  rfmt(1:np)=vclmt(1:np,ias)+vxcmt(1:np,ias)
30  call rfmtftoc(nrc,nrci,rfmt,vmt(:,ias))
31  call rbshtip(nrc,nrci,vmt(:,ias))
32  call rfcmtwr(nrc,nrci,wr2cmt(:,is),vmt(:,ias))
33  end do
34  rfir(1:ngtot)=vclir(1:ngtot)+vxcir(1:ngtot)
35 ! multiply interstitial potential by characteristic function and convert to
36 ! coarse grid
37  call rfirftoc(rfir,vir)
38  if (spinpol) then
39  do idm=1,ndmag
40  do ias=1,natmtot
41  is=idxis(ias)
42  nrc=nrcmt(is)
43  nrci=nrcmti(is)
44  call rfmtftoc(nrc,nrci,bxcmt(:,ias,idm),bmt(:,ias,idm))
45  call rbshtip(nrc,nrci,bmt(:,ias,idm))
46  call rfcmtwr(nrc,nrci,wr2cmt(:,is),bmt(:,ias,idm))
47  end do
48  call rfirftoc(bxcir(:,idm),bir(:,idm))
49  end do
50  end if
51 else
52 ! normal Hartree-Fock case
53  do ias=1,natmtot
54  is=idxis(ias)
55  nrc=nrcmt(is)
56  nrci=nrcmti(is)
57  call rfmtftoc(nrc,nrci,vclmt(:,ias),vmt(:,ias))
58  call rbshtip(nrc,nrci,vmt(:,ias))
59  call rfcmtwr(nrc,nrci,wr2cmt(:,is),vmt(:,ias))
60  end do
61  call rfirftoc(vclir,vir)
62 end if
63 end subroutine
64 
real(8), dimension(:), allocatable wxcir
Definition: modmain.f90:676
subroutine rbshtip(nr, nri, rfmt)
Definition: rbshtip.f90:7
real(8), dimension(:,:), allocatable vxcmt
Definition: modmain.f90:634
integer, dimension(3) xctype
Definition: modmain.f90:588
real(8), dimension(:), allocatable ecir
Definition: modmain.f90:632
logical spinpol
Definition: modmain.f90:228
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
logical hybrid
Definition: modmain.f90:1152
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
real(8), dimension(:,:), allocatable ecmt
Definition: modmain.f90:632
real(8), dimension(:,:), allocatable exmt
Definition: modmain.f90:630
pure subroutine rfcmtwr(nr, nri, wr, rfmt)
Definition: rfcmtwr.f90:7
subroutine potxc(tsh, xctype_, rhomt_, rhoir_, magmt_, magir_, taumt_, tauir_, exmt_, exir_, ecmt_, ecir_, vxcmt_, vxcir_, bxcmt_, bxcir_, wxcmt_, wxcir_)
Definition: potxc.f90:11
real(8), dimension(:), allocatable vxcir
Definition: modmain.f90:634
real(8), dimension(:,:), allocatable wxcmt
Definition: modmain.f90:676
real(8), dimension(:,:,:), allocatable bxcmt
Definition: modmain.f90:636
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
real(8), dimension(:), allocatable vclir
Definition: modmain.f90:624
real(8), dimension(:,:), allocatable bxcir
Definition: modmain.f90:636
subroutine hflocal(vmt, vir, bmt, bir)
Definition: hflocal.f90:7
real(8), dimension(:,:), allocatable wr2cmt
Definition: modmain.f90:189
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:), allocatable tauir
Definition: modmain.f90:672
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
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
subroutine potcoul
Definition: potcoul.f90:10
real(8), dimension(:,:,:), allocatable taumt
Definition: modmain.f90:672
subroutine rfirftoc(rfir, rfirc)
Definition: rfirftoc.f90:7
real(8), dimension(:), allocatable exir
Definition: modmain.f90:630