The Elk Code
bdipole.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 !BOP
7 ! !ROUTINE: bdipole
8 ! !INTERFACE:
9 subroutine bdipole
10 ! !USES:
11 use modmain
12 ! !DESCRIPTION:
13 ! Calculates the magnetic dipole field arising from the spin and orbital
14 ! current. The total current density is
15 ! $$ {\bf j}({\bf r}) = {\rm Im}\sum_{i{\bf k}}^{\rm occ}
16 ! \varphi_{i{\bf k}}^{\dag}({\bf r})\nabla\varphi_{i{\bf k}}({\bf r})
17 ! -\frac{1}{c}{\bf A}\,\rho({\bf r})
18 ! +\frac{g_s}{4}\nabla\times{\bf m}({\bf r}), $$
19 ! where $g_s$ is the electron spin $g$-factor. The vector potential arising
20 ! from ${\bf j}({\bf r})$ is calculated by
21 ! $$ {\bf A}({\bf r})=\frac{1}{c}\int d^3r'\,\frac{{\bf j}({\bf r}')}
22 ! {|{\bf r}-{\bf r}'|}, $$
23 ! using the Poisson equation solver {\tt zpotcoul}. Finally, the magnetic
24 ! field is determined from ${\bf B}({\bf r})=\nabla\times{\bf A}({\bf r})$.
25 ! This field is included as a Zeeman term in the second-variational
26 ! Hamiltonian:
27 ! $$ \hat{H}\rightarrow\hat{H}+\frac{g_s}{4c}{\bf B}\cdot\boldsymbol\sigma. $$
28 !
29 ! !REVISION HISTORY:
30 ! Created April 2018 (T. Mueller)
31 !EOP
32 !BOC
33 implicit none
34 ! local variables
35 integer idm,is,ias
36 integer nrc,nrci,np,npc
37 real(8) cb,t1,t2
38 ! automatic arrays
39 real(8) rfmt(npcmtmax)
40 ! allocatable arrays
41 real(8), allocatable :: rvfmt(:,:,:),rvfir(:,:)
42 complex(8), allocatable :: zvclmt(:,:),zvclir(:)
43 ! external functions
44 real(8), external :: rfmtint
45 ! prefactor for the spin dipole magnetic field
46 cb=gfacte/(4.d0*solsc)
47 ! allocate local arrays
48 allocate(rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3))
49 allocate(zvclmt(npmtmax,natmtot),zvclir(ngtot))
50 ! compute the curl of the magnetisation density, i.e. the magnetisation current
51 call curlrvf(magmt,magir,rvfmt,rvfir)
52 ! multiply by prefactor
53 do idm=1,3
54  do ias=1,natmtot
55  is=idxis(ias)
56  np=npmt(is)
57  rvfmt(1:np,ias,idm)=cb*rvfmt(1:np,ias,idm)
58  end do
59 end do
60 rvfir(1:ngtot,1:3)=cb*rvfir(1:ngtot,1:3)
61 ! add the gauge-invariant current density if required
62 if (tjr.and.(iscl >= 1)) then
63  call genjr
64  t1=1.d0/solsc
65  do idm=1,3
66  do ias=1,natmtot
67  is=idxis(ias)
68  np=npmt(is)
69  rvfmt(1:np,ias,idm)=rvfmt(1:np,ias,idm)+t1*jrmt(1:np,ias,idm)
70  end do
71  end do
72  rvfir(1:ngtot,1:3)=rvfir(1:ngtot,1:3)+t1*jrir(1:ngtot,1:3)
73 end if
74 do idm=1,3
75 ! transform to complex spherical harmonics
76  do ias=1,natmtot
77  is=idxis(ias)
78  call rtozfmt(nrmt(is),nrmti(is),rvfmt(:,ias,idm),zvclmt(:,ias))
79  end do
80 ! solve Poisson's equation in the muffin-tin to find the A-field
82  zvclir(1:ngtot)=rvfir(1:ngtot,idm)
83 ! solve in the entire unit cell
85  ngvec,jlgrmt,ylmg,sfacg,npmtmax,zvclmt,zvclir)
86 ! convert muffin-tin A-field to real spherical harmonics
87  do ias=1,natmtot
88  is=idxis(ias)
89  call ztorfmt(nrmt(is),nrmti(is),zvclmt(:,ias),rvfmt(:,ias,idm))
90  end do
91 ! store the real part of the interstitial A-field
92  rvfir(1:ngtot,idm)=dble(zvclir(1:ngtot))
93 end do
94 ! compute the curl of A to obtain the dipole B-field
95 call curlrvf(rvfmt,rvfir,bdmt,bdir)
96 ! scale dipole B-field if required (by scaling the prefactor)
97 cb=cb*bdipscf
98 ! add to the Kohn-Sham field
99 do idm=1,3
100  do ias=1,natmtot
101  is=idxis(ias)
102  nrc=nrcmt(is)
103  nrci=nrcmti(is)
104  npc=npcmt(is)
105 ! convert to coarse radial mesh
106  call rfmtftoc(nrc,nrci,bdmt(:,ias,idm),rfmt)
107 ! convert to spherical coordinates
108  call rbshtip(nrc,nrci,rfmt)
109  bsmt(1:npc,ias,idm)=bsmt(1:npc,ias,idm)+cb*rfmt(1:npc)
110 ! store the average dipole field in each muffin-tin
111  t1=rfmtint(nrmt(is),nrmti(is),wr2mt(:,is),bdmt(:,ias,idm))
112  t2=(4.d0/3.d0)*pi*rmt(is)**3
113  bdmta(idm,ias)=t1/t2
114  end do
115 end do
116 do idm=1,3
117  bsir(1:ngtot,idm)=bsir(1:ngtot,idm)+cb*bdir(1:ngtot,idm)
118 end do
119 deallocate(rvfmt,rvfir,zvclmt,zvclir)
120 end subroutine
121 !EOC
122 
complex(8), dimension(:,:), allocatable sfacg
Definition: modmain.f90:430
subroutine rbshtip(nr, nri, rfmt)
Definition: rbshtip.f90:7
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
logical tjr
Definition: modmain.f90:620
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer ngtot
Definition: modmain.f90:390
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
integer iscl
Definition: modmain.f90:1045
real(8), dimension(:,:,:), allocatable bdmt
Definition: modmain.f90:638
subroutine genjr
Definition: genjr.f90:7
real(8), dimension(:,:), allocatable bdmta
Definition: modmain.f90:640
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition: rtozfmt.f90:7
real(8), parameter pi
Definition: modmain.f90:1226
complex(8), dimension(:,:), allocatable ylmg
Definition: modmain.f90:428
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
real(8), dimension(:,:), allocatable bsir
Definition: modmain.f90:658
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
subroutine curlrvf(rvfmt, rvfir, curlmt, curlir)
Definition: curlrvf.f90:7
subroutine bdipole
Definition: bdipole.f90:10
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zvclmt)
Definition: genzvclmt.f90:7
real(8), dimension(:,:), allocatable bdir
Definition: modmain.f90:638
integer ngvec
Definition: modmain.f90:396
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, ld3, zvclmt, zvclir)
Definition: zpotcoul.f90:11
real(8), parameter gfacte
Definition: modmain.f90:1271
real(8), dimension(maxspecies) rmt
Definition: modmain.f90:162
real(8) solsc
Definition: modmain.f90:1247
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:,:), allocatable jlgrmt
Definition: modmain.f90:426
pure subroutine ztorfmt(nr, nri, zfmt, rfmt)
Definition: ztorfmt.f90:7
real(8), dimension(:), allocatable gclg
Definition: modmain.f90:424
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
real(8) bdipscf
Definition: modmain.f90:645
integer npmtmax
Definition: modmain.f90:216
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
real(8), dimension(:,:,:), allocatable jrmt
Definition: modmain.f90:622
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:,:), allocatable wprmt
Definition: modmain.f90:185
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
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), dimension(:,:), allocatable jrir
Definition: modmain.f90:622
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition: modmain.f90:656
real(8), dimension(:,:), allocatable wr2mt
Definition: modmain.f90:183
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150