The Elk Code
mae.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2013 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 mae
7 use modmain
8 use modmpi
9 use moddelf
10 use modtest
11 implicit none
12 ! local variables
13 integer i,j,i0,i1
14 real(8) e0,e1,de
15 real(8) v1(3),v2(3),th
16 real(8) a(3,3),b(3,3)
17 ! initialise global variables
18 call init0
19 ! store original parameters
20 avec0(:,:)=avec(:,:)
24 bfieldc00(:)=bfieldc0(:)
27 vkloff0(:)=vkloff(:)
28 ! enable spin-orbit coupling
29 spinorb=.true.
30 ! enforce collinear magnetism in the z-direction
31 cmagz=.true.
32 ! no fixed spin moment calculation: the crystal is rotated instead
33 fsmtype=0
34 ! if task=28 then start from atomic densities; if task=29 read STATE.OUT
35 trdstate=(task == 29)
36 ! zero k-point offset
37 vkloff(:)=0.d0
38 ! start with large magnetic field
39 bfieldc0(1:2)=0.d0
40 bfieldc0(3)=-2.d0
41 ! reduce the external magnetic field after each s.c. loop
42 reducebf=0.75d0
43 ! generate the spin moment directions in (theta,phi) coordinates
44 call gentpmae
45 ! open MAE_INFO.OUT
46 if (mp_mpi) then
47  open(71,file='MAE_INFO.OUT',form='FORMATTED')
48  write(71,*)
49  write(71,'("Scale factor of spin-orbit coupling term : ",G18.10)') socscf
50 end if
51 i0=1; i1=1
52 e0=1.d8; e1=-1.d8
53 ! loop over points on sphere
54 do i=1,npmae
55  if (mp_mpi) then
56  write(*,'("Info(mae): fixed spin moment direction ",I6," of ",I6)') i,npmae
57  end if
58 ! rotate lattice vectors instead of moment (thanks to J. Glasbrenner,
59 ! K. Bussmann and I. Mazin)
60 ! first by -phi around the z-axis
61  v1(:)=0.d0
62  v1(3)=1.d0
63  th=-tpmae(2,i)
64  call axangrot(v1,th,a)
65 ! then by -theta around the y-axis
66  v1(:)=0.d0
67  v1(2)=1.d0
68  th=-tpmae(1,i)
69  call axangrot(v1,th,b)
70  call r3mm(b,a,rotsht)
71  call r3mm(rotsht,avec0,avec)
72 ! find the corresponding moment direction vector
73  call r3minv(rotsht,a)
74  v1(:)=0.d0
75  v1(3)=1.d0
76  call r3mv(a,v1,v2)
77  do j=1,3
78  if (abs(v2(j)) < epslat) v2(j)=0.d0
79  end do
80 ! rotate the spherical cover used for the spherical harmonic transform
81  trotsht=.true.
82 ! run the ground-state calculation
83  call gndstate
84 ! subsequent calculations should read the previous density
85  trdstate=.true.
86 ! make external magnetic field small
87  bfieldc0(3)=-0.01d0
88  if (mp_mpi) then
89  write(71,*)
90  write(71,'("Fixed spin moment direction point ",I6," of ",I6)') i,npmae
91  write(71,'("Spherical coordinates of direction : ",2G18.10)') tpmae(:,i)
92  write(71,'("Direction vector (Cartesian coordinates) : ",3G18.10)') v2
93  write(71,'("Calculated total moment magnitude : ",G18.10)') momtotm
94  write(71,'("Total energy : ",G24.14)') engytot
95  flush(71)
96  end if
97 ! check for minimum and maximum total energy
98  if (engytot < e0) then
99  e0=engytot
100  i0=i
101  end if
102  if (engytot > e1) then
103  e1=engytot
104  i1=i
105  end if
106 ! delete the eigenvector files
107  call delfiles(evec=.true.)
108 ! synchronise MPI processes
109  call mpi_barrier(mpicom,ierror)
110 end do
111 ! magnetic anisotropy energy
112 de=e1-e0
113 if (mp_mpi) then
114  write(71,*)
115  write(71,'("Minimum energy point : ",I6)') i0
116  write(71,'("Maximum energy point : ",I6)') i1
117  write(71,*)
118  write(71,'("Estimated magnetic anisotropy energy (MAE) : ",G18.10)') de
119  write(71,*)
120  write(71,'("MAE per unit volume : ",G18.10)') de/omega
121  close(71)
122  open(50,file='MAE.OUT',form='FORMATTED')
123  write(50,'(G18.10)') de
124  close(50)
125  open(50,file='MAEPUV.OUT',form='FORMATTED')
126  write(50,'(G18.10)') de/omega
127  close(50)
128  write(*,*)
129  write(*,'("Info(mae):")')
130  write(*,'(" Estimated magnetic anisotropy energy written to MAE.OUT")')
131  write(*,'(" MAE per unit volume written to MAEPUV.OUT")')
132  write(*,*)
133  write(*,'(" Number of fixed spin moment directions used : ",I6)') npmae
134  write(*,*)
135  write(*,'(" Additional information written to MAE_INFO.OUT")')
136 end if
137 ! write the MAE to test file
138 call writetest(28,'magnetic anisotropy energy',tol=1.d-5,rv=de)
139 ! restore original input parameters
140 avec(:,:)=avec0(:,:)
145 bfieldc0(:)=bfieldc00(:)
147 vkloff(:)=vkloff0(:)
148 trotsht=.false.
149 end subroutine
150 
real(8) socscf
Definition: modmain.f90:232
real(8), dimension(:,:), allocatable tpmae
Definition: modmain.f90:304
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
subroutine gndstate
Definition: gndstate.f90:10
integer task
Definition: modmain.f90:1299
logical mp_mpi
Definition: modmpi.f90:17
integer fsmtype0
Definition: modmain.f90:251
logical spinpol
Definition: modmain.f90:228
real(8) reducebf
Definition: modmain.f90:279
real(8) momtotm
Definition: modmain.f90:740
real(8) omega
Definition: modmain.f90:20
subroutine r3minv(a, b)
Definition: r3minv.f90:10
real(8), dimension(3) vkloff0
Definition: modmain.f90:450
real(8), dimension(3) vkloff
Definition: modmain.f90:450
subroutine gentpmae
Definition: gentpmae.f90:7
real(8) reducebf0
Definition: modmain.f90:279
real(8), dimension(3, 3) rotsht
Definition: modmain.f90:563
real(8) engytot
Definition: modmain.f90:983
pure subroutine axangrot(v, th, rot)
Definition: axangrot.f90:10
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
logical cmagz
Definition: modmain.f90:242
logical spinorb0
Definition: modmain.f90:230
logical spinpol0
Definition: modmain.f90:228
real(8), dimension(3, 3) avec0
Definition: modmain.f90:12
subroutine mae
Definition: mae.f90:7
subroutine delfiles(evec, devec, eval, occ, pmat, epsi)
Definition: moddelf.f90:25
real(8) epslat
Definition: modmain.f90:24
Definition: modmpi.f90:6
logical trotsht
Definition: modmain.f90:561
logical cmagz0
Definition: modmain.f90:242
logical spinorb
Definition: modmain.f90:230
real(8), dimension(3) bfieldc0
Definition: modmain.f90:271
real(8), dimension(3) bfieldc00
Definition: modmain.f90:271
logical trdstate
Definition: modmain.f90:682
subroutine init0
Definition: init0.f90:10
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
pure subroutine r3mm(a, b, c)
Definition: r3mm.f90:10
integer mpicom
Definition: modmpi.f90:11
integer fsmtype
Definition: modmain.f90:251
integer npmae
Definition: modmain.f90:302
integer ierror
Definition: modmpi.f90:19