The Elk Code
moke.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 S. Sharma, J. K. Dewhurst 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 subroutine moke
7 use modmain
8 implicit none
9 ! local variables
10 integer iw,ios
11 complex(8) z1,z2,z3
12 ! allocatable arrays
13 real(8), allocatable :: w(:),sig1(:,:),sig2(:,:)
14 complex(8), allocatable :: kerr(:)
15 ! calculate dielectric function for the 11 and 12 components
16 noptcomp=2
17 optcomp(1,1)=1
18 optcomp(2,1)=1
19 optcomp(1,2)=1
20 optcomp(2,2)=2
21 call dielectric
22 ! allocate local arrays
23 allocate(w(nwplot))
24 allocate(sig1(nwplot,2),sig2(nwplot,2))
25 allocate(kerr(nwplot))
26 ! read diagonal contribution to optical conductivity
27 open(50,file='SIGMA_11.OUT',status='OLD',form='FORMATTED',iostat=ios)
28 if (ios /= 0) then
29  write(*,*)
30  write(*,'("Error(moke): error opening SIGMA_11.OUT")')
31  write(*,*)
32  stop
33 end if
34 do iw=1,nwplot
35  read(50,'(2G18.10)') w(iw),sig1(iw,1)
36 end do
37 read(50,*)
38 do iw=1,nwplot
39  read(50,'(2G18.10)') w(iw),sig2(iw,1)
40 end do
41 close(50)
42 ! read off-diagonal contribution to optical conductivity
43 open(50,file='SIGMA_12.OUT',status='OLD',form='FORMATTED',iostat=ios)
44 if (ios /= 0) then
45  write(*,*)
46  write(*,'("Error(moke): error opening SIGMA_12.OUT")')
47  write(*,*)
48  stop
49 end if
50 do iw=1,nwplot
51  read(50,'(2G18.10)') w(iw),sig1(iw,2)
52 end do
53 read(50,*)
54 do iw=1,nwplot
55  read(50,'(2G18.10)') w(iw),sig2(iw,2)
56 end do
57 close(50)
58 ! calculate the complex Kerr angle
59 do iw=1,nwplot
60  if (w(iw) > 0.d0) then
61  z1=cmplx(sig1(iw,1),sig2(iw,1),8)
62  z2=cmplx(sig1(iw,2),sig2(iw,2),8)
63  z3=z1*sqrt(1.d0+fourpi*cmplx(-z1%im,z1%re,8)/w(iw))
64  if (abs(z3) > 1.d-8) then
65  kerr(iw)=-z2/z3
66  else
67  kerr(iw)=0.d0
68  end if
69  else
70  kerr(iw)=0.d0
71  end if
72 end do
73 open(50,file='KERR.OUT',form='FORMATTED')
74 do iw=1,nwplot
75  write(50,'(2G18.10)') w(iw),dble(kerr(iw))*180.d0/pi
76 end do
77 write(50,'(" ")')
78 do iw=1,nwplot
79  write(50,'(2G18.10)') w(iw),aimag(kerr(iw))*180.d0/pi
80 end do
81 close(50)
82 write(*,*)
83 write(*,'("Info(moke):")')
84 write(*,'(" complex Kerr angle in degrees written to KERR.OUT")')
85 deallocate(w,sig1,sig2,kerr)
86 end subroutine
87 
real(8), parameter pi
Definition: modmain.f90:1232
subroutine dielectric
Definition: dielectric.f90:10
integer, dimension(3, 27) optcomp
Definition: modmain.f90:1093
integer noptcomp
Definition: modmain.f90:1091
subroutine moke
Definition: moke.f90:7
real(8), parameter fourpi
Definition: modmain.f90:1234
integer nwplot
Definition: modmain.f90:1073