The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine moke
7use modmain
8implicit none
9! local variables
10integer iw,ios
11complex(8) z1,z2,z3
12! allocatable arrays
13real(8), allocatable :: w(:),sig1(:,:),sig2(:,:)
14complex(8), allocatable :: kerr(:)
15! calculate dielectric function for the 11 and 12 components
17optcomp(1,1)=1
18optcomp(2,1)=1
19optcomp(1,2)=1
20optcomp(2,2)=2
21call dielectric
22! allocate local arrays
23allocate(w(nwplot))
24allocate(sig1(nwplot,2),sig2(nwplot,2))
25allocate(kerr(nwplot))
26! read diagonal contribution to optical conductivity
27open(50,file='SIGMA_11.OUT',status='OLD',form='FORMATTED',iostat=ios)
28if (ios /= 0) then
29 write(*,*)
30 write(*,'("Error(moke): error opening SIGMA_11.OUT")')
31 write(*,*)
32 stop
33end if
34do iw=1,nwplot
35 read(50,'(2G18.10)') w(iw),sig1(iw,1)
36end do
37read(50,*)
38do iw=1,nwplot
39 read(50,'(2G18.10)') w(iw),sig2(iw,1)
40end do
41close(50)
42! read off-diagonal contribution to optical conductivity
43open(50,file='SIGMA_12.OUT',status='OLD',form='FORMATTED',iostat=ios)
44if (ios /= 0) then
45 write(*,*)
46 write(*,'("Error(moke): error opening SIGMA_12.OUT")')
47 write(*,*)
48 stop
49end if
50do iw=1,nwplot
51 read(50,'(2G18.10)') w(iw),sig1(iw,2)
52end do
53read(50,*)
54do iw=1,nwplot
55 read(50,'(2G18.10)') w(iw),sig2(iw,2)
56end do
57close(50)
58! calculate the complex Kerr angle
59do 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*zi*z1/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
72end do
73open(50,file='KERR.OUT',form='FORMATTED')
74do iw=1,nwplot
75 write(50,'(2G18.10)') w(iw),dble(kerr(iw))*180.d0/pi
76end do
77write(50,'(" ")')
78do iw=1,nwplot
79 write(50,'(2G18.10)') w(iw),aimag(kerr(iw))*180.d0/pi
80end do
81close(50)
82write(*,*)
83write(*,'("Info(moke):")')
84write(*,'(" complex Kerr angle in degrees written to KERR.OUT")')
85deallocate(w,sig1,sig2,kerr)
86end subroutine
87
subroutine dielectric
subroutine moke
Definition moke.f90:7
real(8), parameter pi
Definition modmain.f90:1229
integer, dimension(3, 27) optcomp
Definition modmain.f90:1090
complex(8), parameter zi
Definition modmain.f90:1239
integer nwplot
Definition modmain.f90:1070
real(8), parameter fourpi
Definition modmain.f90:1231
integer noptcomp
Definition modmain.f90:1088