The Elk Code
Loading...
Searching...
No Matches
addbfsmu.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2019 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
addbfsmu
7
use
modmain
8
use
modulr
9
implicit none
10
! local variables
11
integer
idm,is,ias,npc
12
real
(8) t1
13
! add the global fixed spin moment B-field to the Kohn-Sham field
14
if
(any(abs(
fsmtype
) == [1,3,4,6]))
then
15
do
idm=1,
ndmag
16
t1=
bfsmc
(idm)
17
do
ias=1,
natmtot
18
is=
idxis
(ias)
19
npc=
npcmt
(is)
20
bsqmt
(1:npc,ias,idm,1)=
bsqmt
(1:npc,ias,idm,1)+t1
21
end do
22
bsqir
(1:
ngtot
,idm,1)=
bsqir
(1:
ngtot
,idm,1)+t1
23
end do
24
end if
25
! add the muffin-tin fields
26
if
(any(abs(
fsmtype
) == [2,3,5,6]))
then
27
do
idm=1,
ndmag
28
do
ias=1,
natmtot
29
is=
idxis
(ias)
30
npc=
npcmt
(is)
31
t1=
bfsmcmt
(idm,ias)
32
bsqmt
(1:npc,ias,idm,1)=
bsqmt
(1:npc,ias,idm,1)+t1
33
end do
34
end do
35
end if
36
end subroutine
37
addbfsmu
subroutine addbfsmu
Definition
addbfsmu.f90:7
modmain
Definition
modmain.f90:6
modmain::ngtot
integer ngtot
Definition
modmain.f90:390
modmain::fsmtype
integer fsmtype
Definition
modmain.f90:251
modmain::bfsmcmt
real(8), dimension(:,:), allocatable bfsmcmt
Definition
modmain.f90:263
modmain::idxis
integer, dimension(maxatoms *maxspecies) idxis
Definition
modmain.f90:44
modmain::bfsmc
real(8), dimension(3) bfsmc
Definition
modmain.f90:257
modmain::natmtot
integer natmtot
Definition
modmain.f90:40
modmain::npcmt
integer, dimension(maxspecies) npcmt
Definition
modmain.f90:214
modmain::ndmag
integer ndmag
Definition
modmain.f90:238
modulr
Definition
modulr.f90:6
modulr::bsqmt
complex(8), dimension(:,:,:,:), pointer, contiguous bsqmt
Definition
modulr.f90:84
modulr::bsqir
complex(8), dimension(:,:,:), pointer, contiguous bsqir
Definition
modulr.f90:84
addbfsmu.f90
Generated by
1.9.8