The Elk Code
Loading...
Searching...
No Matches
bdipoleu.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2024 Wenhan Chen, 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
6
subroutine
bdipoleu
7
use
modmain
8
use
modulr
9
implicit none
10
! local variables
11
integer
idm,iq,ifq
12
real
(8) cb,t1
13
complex(8)
zv(3),z1
14
! automatic arrays
15
real
(8) rfft(nqpt)
16
complex(8)
zvq(nfqrz,3)
17
if
(.not.
ncmag
)
then
18
write
(*,*)
19
write
(*,
'("Error(bdipoleu): non-collinear magnetism required for inclusion &
20
&of the dipole field")'
)
21
write
(*,*)
22
stop
23
end if
24
! prefactor for the electron spin dipole magnetic field
25
cb=
gfacte
/(4.d0*
solsc
)
26
! Fourier transform the R-dependent magnetisation to Q-space
27
do
idm=1,3
28
rfft(1:nqpt)=
momtotru
(idm,1:nqpt)/
omega
29
call
rzfftifc
(3,
ngridq
,-1,rfft,zvq(:,idm))
30
end do
31
do
ifq=1,nfqrz
32
iq=
iqrzf
(ifq)
33
! compute B(Q) = 4π (gₛ/4c)² [ m(Q) - Q (Q⋅m(Q))/Q² ]
34
zv(1:3)=zvq(ifq,1:3)
35
if
(iq > 1)
then
36
z1=
vqc
(1,iq)*zv(1)+
vqc
(2,iq)*zv(2)+
vqc
(3,iq)*zv(3)
37
t1=
vqc
(1,iq)**2+
vqc
(2,iq)**2+
vqc
(3,iq)**2
38
z1=z1/t1
39
zv(1:3)=zv(1:3)-z1*
vqc
(:,iq)
40
end if
41
t1=
bdipscf
*
fourpi
*cb**2
42
bdipq
(1:3,ifq)=t1*zv(1:3)
43
end do
44
end subroutine
45
bdipoleu
subroutine bdipoleu
Definition
bdipoleu.f90:7
modmain
Definition
modmain.f90:6
modmain::gfacte
real(8), parameter gfacte
Definition
modmain.f90:1276
modmain::ncmag
logical ncmag
Definition
modmain.f90:240
modmain::omega
real(8) omega
Definition
modmain.f90:20
modmain::fourpi
real(8), parameter fourpi
Definition
modmain.f90:1231
modmain::vqc
real(8), dimension(:,:), allocatable vqc
Definition
modmain.f90:547
modmain::iqrzf
integer, dimension(:), allocatable iqrzf
Definition
modmain.f90:543
modmain::solsc
real(8) solsc
Definition
modmain.f90:1252
modmain::bdipscf
real(8) bdipscf
Definition
modmain.f90:645
modmain::ngridq
integer, dimension(3) ngridq
Definition
modmain.f90:515
modulr
Definition
modulr.f90:6
modulr::bdipq
complex(8), dimension(:,:), allocatable bdipq
Definition
modulr.f90:80
modulr::momtotru
real(8), dimension(:,:), allocatable momtotru
Definition
modulr.f90:57
rzfftifc
subroutine rzfftifc(nd, n, sgn, r, z)
Definition
zfftifc_fftw.f90:32
bdipoleu.f90
Generated by
1.9.8