The Elk Code
 
Loading...
Searching...
No Matches
bfieldfsm.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
6!BOP
7! !ROUTINE: bfieldfsm
8! !INTERFACE:
9subroutine bfieldfsm
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Updates the effective magnetic field, ${\bf B}_{\rm FSM}$, required for
14! fixing the spin moment to a given value, ${\bf M}_{\rm FSM}$. This is done
15! by adding a vector to the field which is proportional to the difference
16! between the moment calculated in the $i$th self-consistent loop and the
17! required moment:
18! $$ {\bf B}_{\rm FSM}^{i+1}={\bf B}_{\rm FSM}^i+\lambda\left({\bf M}^i
19! -{\bf M}_{\rm FSM}\right), $$
20! where $\lambda$ is a scaling factor.
21!
22! !REVISION HISTORY:
23! Created March 2005 (JKD)
24!EOP
25!BOC
26implicit none
27! local variables
28integer is,ia,ias
29real(8) v1(3),v2(3),t1,t2
30if ((.not.spinpol).or.(fsmtype == 0)) return
31! fixed spin direction not valid for collinear magnetism
32if ((.not.ncmag).and.(fsmtype < 0)) return
33! determine the global effective field
34if ((abs(fsmtype) == 1).or.(abs(fsmtype) == 3)) then
35 if (ncmag) then
36 v1(:)=momtot(:)
37 else
38 v1(:)=0.d0
39 v1(3)=momtot(1)
40 end if
41 v2(:)=v1(:)-momfix(:)
42 if (ncmag) then
43 bfsmc(:)=bfsmc(:)+taufsm*v2(:)
44 else
45 bfsmc(1)=bfsmc(1)+taufsm*v2(3)
46 end if
47! make sure that the constraining field is perpendicular to the fixed moment
48! for fixed direction calculations (Y. Kvashnin and LN)
49 if (fsmtype < 0) call r3vo(momfix,bfsmc)
50end if
51! determine the muffin-tin fields for fixed local moments
52if ((abs(fsmtype) == 2).or.(abs(fsmtype) == 3)) then
53 do is=1,nspecies
54 do ia=1,natoms(is)
55 ias=idxas(ia,is)
56! if any component is >= 1000 then do not fix the moment
57 if (any(abs(mommtfix(1:3,ia,is)) >= 1000.d0)) cycle
58 if (ncmag) then
59 v1(:)=mommt(:,ias)
60 else
61 v1(:)=0.d0
62 v1(3)=mommt(1,ias)
63 end if
64 v2(:)=v1(:)-mommtfix(:,ia,is)
65 if (ncmag) then
66 bfsmcmt(:,ias)=bfsmcmt(:,ias)+taufsm*v2(:)
67 else
68 bfsmcmt(1,ias)=bfsmcmt(1,ias)+taufsm*v2(3)
69 end if
70! fixed spin direction
71 if (fsmtype < 0) call r3vo(mommtfix(:,ia,is),bfsmcmt(:,ias))
72 end do
73 end do
74end if
75! global fixed spin magnitude
76if ((fsmtype == 4).or.(fsmtype == 6)) then
78 bfsmc(1:ndmag)=bfsmc(1:ndmag)+t1*momtot(1:ndmag)
79end if
80! local fixed spin magnitude
81if ((fsmtype == 5).or.(fsmtype == 6)) then
82 do is=1,nspecies
83 do ia=1,natoms(is)
84 ias=idxas(ia,is)
85 t1=mommtfixm(ia,is)
86 if (t1 < 0.d0) cycle
87 t2=norm2(mommt(1:ndmag,ias))
88 t1=taufsm*(t2-t1)
89 bfsmcmt(1:ndmag,ias)=bfsmcmt(1:ndmag,ias)+t1*mommt(1:ndmag,ias)
90 end do
91 end do
92end if
93end subroutine
94!EOC
95
subroutine bfieldfsm
Definition bfieldfsm.f90:10
logical ncmag
Definition modmain.f90:240
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8) momfixm
Definition modmain.f90:255
real(8), dimension(3) momfix
Definition modmain.f90:253
real(8), dimension(:,:), allocatable mommt
Definition modmain.f90:744
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
logical spinpol
Definition modmain.f90:228
integer fsmtype
Definition modmain.f90:251
real(8), dimension(:,:), allocatable bfsmcmt
Definition modmain.f90:263
real(8) taufsm
Definition modmain.f90:265
real(8), dimension(maxatoms, maxspecies) mommtfixm
Definition modmain.f90:261
real(8), dimension(3) bfsmc
Definition modmain.f90:257
real(8), dimension(3) momtot
Definition modmain.f90:738
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition modmain.f90:259
real(8) momtotm
Definition modmain.f90:740
integer nspecies
Definition modmain.f90:34
integer ndmag
Definition modmain.f90:238
pure subroutine r3vo(x, y)
Definition r3vo.f90:7