The Elk Code
 
Loading...
Searching...
No Matches
addbfsm.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine addbfsm
7use modmain
8implicit none
9! local variables
10integer idm,is,ias,npc
11real(8) t1
12! add the global fixed spin moment B-field to the Kohn-Sham field
13if (any(abs(fsmtype) == [1,3,4,6])) then
14 do idm=1,ndmag
15 t1=bfsmc(idm)
16 do ias=1,natmtot
17 is=idxis(ias)
18 npc=npcmt(is)
19 bsmt(1:npc,ias,idm)=bsmt(1:npc,ias,idm)+t1
20 end do
21 bsirc(1:ngtc,idm)=bsirc(1:ngtc,idm)+t1*cfrc(1:ngtc)
22 end do
23end if
24! add the muffin-tin fields
25if (any(abs(fsmtype) == [2,3,5,6])) then
26 do idm=1,ndmag
27 do ias=1,natmtot
28 is=idxis(ias)
29 npc=npcmt(is)
30 t1=bfsmcmt(idm,ias)
31 bsmt(1:npc,ias,idm)=bsmt(1:npc,ias,idm)+t1
32 end do
33 end do
34end if
35end subroutine
36
subroutine addbfsm
Definition addbfsm.f90:7
real(8), dimension(:), allocatable cfrc
Definition modmain.f90:438
integer fsmtype
Definition modmain.f90:251
real(8), dimension(:,:), allocatable bfsmcmt
Definition modmain.f90:263
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer ngtc
Definition modmain.f90:392
real(8), dimension(3) bfsmc
Definition modmain.f90:257
integer natmtot
Definition modmain.f90:40
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
real(8), dimension(:,:), pointer, contiguous bsirc
Definition modmain.f90:660
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition modmain.f90:656
integer ndmag
Definition modmain.f90:238