The Elk Code
Loading...
Searching...
No Matches
magnetoelt.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2023 J. K. Dewhurst and S. Sharma.
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
magnetoelt
7
use
modmain
8
use
modmpi
9
use
modtest
10
implicit none
11
! local variables
12
integer
i,j
13
real
(8) pvl1(3),pvl2(3)
14
real
(8) felt(3,3),vc(3),t1
15
! initialise universal variables
16
call
init0
17
call
init1
18
! store original parameters
19
bfieldc00
(:)=
bfieldc0
(:)
20
ngridk0
(:)=
ngridk
(:)
21
maxscl0
=
maxscl
22
spinpol0
=
spinpol
23
reducebf0
=
reducebf
24
tshift0
=
tshift
25
! no shifting of the atomic basis
26
tshift
=.false.
27
! enable spin-polarisation
28
spinpol
=.true.
29
! no magnetic field reduction
30
reducebf
=1.d0
31
! initial ground-state run should start from atomic densities
32
trdstate
=.false.
33
do
j=1,3
34
if
(
mp_mpi
)
then
35
write
(*,
'("Info(magnetoelt): working on magnetic field component ",I1)'
) j
36
end if
37
! apply negative change to magnetic field
38
bfieldc0
(j)=
bfieldc00
(j)-0.5d0*
deltabf
39
! run the ground-state calculation
40
call
gndstate
41
! subsequent calculations will read in the previous potential
42
trdstate
=.true.
43
! compute the first polarisation in lattice coordinates
44
call
polar
(pvl1)
45
! apply positive change to magnetic field
46
bfieldc0
(j)=
bfieldc00
(j)+0.5d0*
deltabf
47
! run the ground-state calculation again
48
call
gndstate
49
! compute the second polarisation
50
call
polar
(pvl2)
51
do
i=1,3
52
! add multiple of 2*pi to bring polarisation vectors into coincidence
53
pvl1(i)=modulo(pvl1(i),
twopi
)
54
pvl2(i)=modulo(pvl2(i),
twopi
)
55
t1=pvl1(i)-pvl2(i)
56
if
(abs(t1-
twopi
) < abs(t1))
then
57
pvl1(i)=pvl1(i)-
twopi
58
else
if
(abs(t1+
twopi
) < abs(t1))
then
59
pvl1(i)=pvl1(i)+
twopi
60
end if
61
! calculate the magnetoelectric tensor element from difference in polarisations
62
t1=
wkptnr
*
occmax
*dble(
nkspolar
*
ngridk
(i))/(
twopi
*
deltabf
)
63
felt(i,j)=t1*(pvl2(i)-pvl1(i))
64
end do
65
end do
66
if
(
mp_mpi
)
then
67
open
(50,file=
'MAGNETOELT.OUT'
,form=
'FORMATTED'
)
68
write
(50,*)
69
write
(50,
'("The magnetoelectric tensor is the change in the polarisation")'
)
70
write
(50,ᵢⱼ
'("with respect to the external magnetic field dP/dB, for")'
)
71
write
(50,
'("components i,j=1...3")'
)
72
do
j=1,3
73
write
(50,*)
74
write
(50,
'("Magnetic field Cartesian component j : ",I1)'
) j
75
write
(50,ᵢⱼ
'("Magnetoelectric tensor components dP/dB, i=1...3")'
)
76
write
(50,
'(" lattice coordinates : ",3G18.10)'
) felt(:,j)
77
call
r3mv
(
avec
,felt(:,j),vc)
78
write
(50,
'(" Cartesian coordinates : ",3G18.10)'
) vc
79
write
(50,
'(" length : ",G18.10)'
) norm2(vc(1:3))
80
end do
81
close
(50)
82
write
(*,*)
83
write
(*,
'("Info(magnetoelt):")'
)
84
write
(*,
'(" Magnetoelectric tensor writtent to MAGNETOELT.OUT")'
)
85
end if
86
! write test file if required
87
call
writetest
(390,
'Magnetoelectric tensor'
,nv=9,tol=1.d-5,rva=felt)
88
! restore original parameters
89
tshift
=
tshift0
90
spinpol
=
spinpol0
91
bfieldc0
(:)=
bfieldc00
(:)
92
reducebf
=
reducebf0
93
end subroutine
94
gndstate
subroutine gndstate
Definition
gndstate.f90:10
init0
subroutine init0
Definition
init0.f90:10
init1
subroutine init1
Definition
init1.f90:10
magnetoelt
subroutine magnetoelt
Definition
magnetoelt.f90:7
modmain
Definition
modmain.f90:6
modmain::bfieldc00
real(8), dimension(3) bfieldc00
Definition
modmain.f90:271
modmain::tshift
logical tshift
Definition
modmain.f90:352
modmain::ngridk0
integer, dimension(3) ngridk0
Definition
modmain.f90:448
modmain::deltabf
real(8) deltabf
Definition
modmain.f90:281
modmain::reducebf0
real(8) reducebf0
Definition
modmain.f90:279
modmain::spinpol
logical spinpol
Definition
modmain.f90:228
modmain::twopi
real(8), parameter twopi
Definition
modmain.f90:1230
modmain::avec
real(8), dimension(3, 3) avec
Definition
modmain.f90:12
modmain::tshift0
logical tshift0
Definition
modmain.f90:352
modmain::nkspolar
integer nkspolar
Definition
modmain.f90:485
modmain::spinpol0
logical spinpol0
Definition
modmain.f90:228
modmain::ngridk
integer, dimension(3) ngridk
Definition
modmain.f90:448
modmain::maxscl0
integer maxscl0
Definition
modmain.f90:1046
modmain::trdstate
logical trdstate
Definition
modmain.f90:682
modmain::bfieldc0
real(8), dimension(3) bfieldc0
Definition
modmain.f90:271
modmain::maxscl
integer maxscl
Definition
modmain.f90:1046
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmain::reducebf
real(8) reducebf
Definition
modmain.f90:279
modmain::wkptnr
real(8) wkptnr
Definition
modmain.f90:477
modmpi
Definition
modmpi.f90:6
modmpi::mp_mpi
logical mp_mpi
Definition
modmpi.f90:17
modtest
Definition
modtest.f90:6
modtest::writetest
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition
modtest.f90:16
polar
subroutine polar(pvl)
Definition
polar.f90:10
r3mv
pure subroutine r3mv(a, x, y)
Definition
r3mv.f90:10
magnetoelt.f90
Generated by
1.9.8