The Elk Code
maginit.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2016 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 maginit
7 use modmain
8 implicit none
9 ! local variables
10 integer idm,is,ia,ias,np
11 ! magnetisation as fraction of density
12 real(8), parameter :: fmr=0.15d0
13 real(8) v(3),t1
14 ! initialise muffin-tin magnetisation
15 do ias=1,natmtot
16  is=idxis(ias)
17  ia=idxia(ias)
18  np=npmt(is)
19  v(1:3)=bfcmt(1:3,ia,is)+bfieldc(1:3)
20  t1=sqrt(v(1)**2+v(2)**2+v(3)**2)
21  if (t1 > 1.d-8) then
22  t1=-fmr/t1
23  v(1:3)=t1*v(1:3)
24  if (.not.ncmag) v(1)=v(3)
25  do idm=1,ndmag
26  t1=v(idm)
27  magmt(1:np,ias,idm)=t1*rhomt(1:np,ias)
28  end do
29  else
30  magmt(1:np,ias,1:ndmag)=0.d0
31  end if
32 end do
33 ! initialise interstitial magnetisation
34 v(1:3)=bfieldc(1:3)
35 t1=sqrt(v(1)**2+v(2)**2+v(3)**2)
36 if (t1 > 1.d-8) then
37  t1=-fmr/t1
38  v(1:3)=t1*v(1:3)
39  if (.not.ncmag) v(1)=v(3)
40  do idm=1,ndmag
41  t1=v(idm)
42  magir(1:ngtot,idm)=t1*rhoir(1:ngtot)
43  end do
44 else
45  magir(1:ngtot,1:ndmag)=0.d0
46 end if
47 end subroutine
48 
integer ngtot
Definition: modmain.f90:390
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
integer ndmag
Definition: modmain.f90:238
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
real(8), dimension(3, maxatoms, maxspecies) bfcmt
Definition: modmain.f90:273
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine maginit
Definition: maginit.f90:7
integer, dimension(maxatoms *maxspecies) idxia
Definition: modmain.f90:45
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
integer natmtot
Definition: modmain.f90:40
logical ncmag
Definition: modmain.f90:240
real(8), dimension(3) bfieldc
Definition: modmain.f90:269