The Elk Code
 
Loading...
Searching...
No Matches
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:
9subroutine bdipole
10! !USES:
11use 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
33implicit none
34! local variables
35integer idm,is,ias
36integer nrc,nrci,np,npc
37real(8) cb,t1,t2
38! automatic arrays
39real(8) rfmt(npcmtmax)
40! allocatable arrays
41real(8), allocatable :: rvfmt(:,:,:),rvfir(:,:)
42complex(8), allocatable :: zrhomt(:,:),zrhoir(:)
43complex(8), allocatable :: zvclmt(:,:),zvclir(:)
44! external functions
45real(8), external :: rfmtint
46! prefactor for the spin dipole magnetic field
47cb=gfacte/(4.d0*solsc)
48! allocate local arrays
49allocate(rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3))
50allocate(zrhomt(npmtmax,natmtot),zrhoir(ngtot))
51allocate(zvclmt(npmtmax,natmtot),zvclir(ngtot))
52! compute the curl of the magnetisation density, i.e. the magnetisation current
53call curlrvf(magmt,magir,rvfmt,rvfir)
54! multiply by prefactor
55do idm=1,3
56 do ias=1,natmtot
57 is=idxis(ias)
58 np=npmt(is)
59 rvfmt(1:np,ias,idm)=cb*rvfmt(1:np,ias,idm)
60 end do
61end do
62rvfir(1:ngtot,1:3)=cb*rvfir(1:ngtot,1:3)
63! add the gauge-invariant current density if required
64if (tjr.and.(iscl >= 1)) then
65 call genjr
66 t1=1.d0/solsc
67 do idm=1,3
68 do ias=1,natmtot
69 is=idxis(ias)
70 np=npmt(is)
71 rvfmt(1:np,ias,idm)=rvfmt(1:np,ias,idm)+t1*jrmt(1:np,ias,idm)
72 end do
73 end do
74 rvfir(1:ngtot,1:3)=rvfir(1:ngtot,1:3)+t1*jrir(1:ngtot,1:3)
75end if
76do idm=1,3
77! transform to complex spherical harmonics
78 do ias=1,natmtot
79 is=idxis(ias)
80 call rtozfmt(nrmt(is),nrmti(is),rvfmt(:,ias,idm),zrhomt(:,ias))
81 end do
82! solve Poisson's equation in the muffin-tin to find the A-field
83 call genzvclmt(nrmt,nrmti,nrmtmax,rlmt,wprmt,npmtmax,zrhomt,zvclmt)
84 zrhoir(1:ngtot)=rvfir(1:ngtot,idm)
85! solve in the entire unit cell
87 ngvec,jlgrmt,ylmg,sfacg,zrhoir,npmtmax,zvclmt,zvclir)
88! convert muffin-tin A-field to real spherical harmonics
89 do ias=1,natmtot
90 is=idxis(ias)
91 call ztorfmt(nrmt(is),nrmti(is),zvclmt(:,ias),rvfmt(:,ias,idm))
92 end do
93! store the real part of the interstitial A-field
94 rvfir(1:ngtot,idm)=dble(zvclir(1:ngtot))
95end do
96! compute the curl of A to obtain the dipole B-field
97call curlrvf(rvfmt,rvfir,bdmt,bdir)
98! scale dipole B-field if required (by scaling the prefactor)
99cb=cb*bdipscf
100! add to the Kohn-Sham field
101do idm=1,3
102 do ias=1,natmtot
103 is=idxis(ias)
104 nrc=nrcmt(is)
105 nrci=nrcmti(is)
106 npc=npcmt(is)
107! convert to coarse radial mesh
108 call rfmtftoc(nrc,nrci,bdmt(:,ias,idm),rfmt)
109! convert to spherical coordinates
110 call rbshtip(nrc,nrci,rfmt)
111 bsmt(1:npc,ias,idm)=bsmt(1:npc,ias,idm)+cb*rfmt(1:npc)
112! store the average dipole field in each muffin-tin
113 t1=rfmtint(nrmt(is),nrmti(is),wr2mt(:,is),bdmt(:,ias,idm))
114 t2=(4.d0/3.d0)*pi*rmt(is)**3
115 bdmta(idm,ias)=t1/t2
116 end do
117end do
118do idm=1,3
119 bsir(1:ngtot,idm)=bsir(1:ngtot,idm)+cb*bdir(1:ngtot,idm)
120end do
121deallocate(rvfmt,rvfir)
122deallocate(zrhomt,zrhoir,zvclmt,zvclir)
123end subroutine
124!EOC
125
subroutine bdipole
Definition bdipole.f90:10
subroutine curlrvf(rvfmt, rvfir, curlmt, curlir)
Definition curlrvf.f90:7
subroutine genjr
Definition genjr.f90:7
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zrhomt, zvclmt)
Definition genzvclmt.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter gfacte
Definition modmain.f90:1276
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngridg
Definition modmain.f90:386
real(8), parameter pi
Definition modmain.f90:1229
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
real(8), dimension(:,:), allocatable bdir
Definition modmain.f90:638
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(maxspecies) rmt
Definition modmain.f90:162
real(8), dimension(:,:,:), allocatable jlgrmt
Definition modmain.f90:426
integer ngvec
Definition modmain.f90:396
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
complex(8), dimension(:,:), allocatable sfacg
Definition modmain.f90:430
logical tjr
Definition modmain.f90:620
integer natmtot
Definition modmain.f90:40
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
real(8), dimension(:,:), pointer, contiguous magir
Definition modmain.f90:616
real(8), dimension(:,:), allocatable bdmta
Definition modmain.f90:640
real(8), dimension(:,:,:), allocatable wprmt
Definition modmain.f90:185
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:,:), allocatable jrmt
Definition modmain.f90:622
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable bsir
Definition modmain.f90:658
real(8) solsc
Definition modmain.f90:1252
real(8), dimension(:,:,:), allocatable bdmt
Definition modmain.f90:638
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
real(8), dimension(:), allocatable gclg
Definition modmain.f90:424
real(8) bdipscf
Definition modmain.f90:645
integer nrmtmax
Definition modmain.f90:152
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition modmain.f90:656
complex(8), dimension(:,:), allocatable ylmg
Definition modmain.f90:428
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
integer iscl
Definition modmain.f90:1048
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
real(8), dimension(:,:), allocatable jrir
Definition modmain.f90:622
subroutine rbshtip(nr, nri, rfmt)
Definition rbshtip.f90:7
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition rfmtftoc.f90:7
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition rtozfmt.f90:7
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, zrhoir, ld3, zvclmt, zvclir)
Definition zpotcoul.f90:11
pure subroutine ztorfmt(nr, nri, zfmt, rfmt)
Definition ztorfmt.f90:7