The Elk Code
potuinit.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 T. Mueller, 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 potuinit
7 use modmain
8 use modulr
9 use modomp
10 implicit none
11 ! local variables
12 integer ifq,idm,is,ias
13 integer nrc,nrci,npc
14 ! automatic arrays
15 real(8) rfmt(npcmtmax)
16 ! set the Q=0 muffin-tin potential equal to that of the normal ground-state
17 do 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),rfmt)
23  call rbshtip(nrc,nrci,rfmt)
24  vsqmt(1:npc,ias,1)=rfmt(1:npc)
25 end do
26 ! zero the muffin-tin potential for non-zero Q
27 do ifq=2,nfqrz
28  do ias=1,natmtot
29  is=idxis(ias)
30  vsqmt(1:npcmt(is),ias,ifq)=0.d0
31  end do
32 end do
33 ! repeat for the interstitial potential
34 vsqir(1:ngtot,1)=vsir(1:ngtot)
35 vsqir(1:ngtot,2:nfqrz)=0.d0
36 if (.not.spinpol) return
37 ! set the Q=0 muffin-tin magnetic field equal to that of the normal ground-state
38 do idm=1,ndmag
39  do ias=1,natmtot
40  is=idxis(ias)
41  npc=npcmt(is)
42  bsqmt(1:npc,ias,idm,1)=bsmt(1:npc,ias,idm)
43  end do
44 end do
45 ! zero the magnetic field for non-zero Q
46 do ifq=2,nfqrz
47  do idm=1,ndmag
48  do ias=1,natmtot
49  is=idxis(ias)
50  bsqmt(1:npcmt(is),ias,idm,ifq)=0.d0
51  end do
52  end do
53 end do
54 ! repeat for the interstitial magnetic field
55 bsqir(1:ngtot,1:ndmag,1)=bsir(1:ngtot,1:ndmag)
56 bsqir(1:ngtot,1:ndmag,2:nfqrz)=0.d0
57 end subroutine
58 
subroutine potuinit
Definition: potuinit.f90:7
complex(8), dimension(:,:,:), pointer, contiguous vsqmt
Definition: modulr.f90:84
subroutine rbshtip(nr, nri, rfmt)
Definition: rbshtip.f90:7
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer ngtot
Definition: modmain.f90:390
logical spinpol
Definition: modmain.f90:228
complex(8), dimension(:,:), pointer, contiguous vsqir
Definition: modulr.f90:84
integer ndmag
Definition: modmain.f90:238
Definition: modomp.f90:6
complex(8), dimension(:,:,:), pointer, contiguous bsqir
Definition: modulr.f90:85
real(8), dimension(:), allocatable vsir
Definition: modmain.f90:651
real(8), dimension(:,:), allocatable bsir
Definition: modmain.f90:658
complex(8), dimension(:,:,:,:), pointer, contiguous bsqmt
Definition: modulr.f90:85
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer nfqrz
Definition: modmain.f90:539
integer natmtot
Definition: modmain.f90:40
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
Definition: modulr.f90:6
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition: modmain.f90:656