The Elk Code
 
Loading...
Searching...
No Matches
dielectric_tdrt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2017 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
7use modmain
8use modtddft
9implicit none
10! local variables
11integer iw,i,j
12real(8) w1,w2,t0,t1
13complex(8) eta,z1,z2
14character(256) fname
15! allocatable arrays
16real(8), allocatable :: w(:),wt(:),jt(:,:)
17complex(8), allocatable :: ew(:,:),jw(:,:),eps(:)
18! initialise global variables
19call init0
20call init1
21! generate energy grid (always non-negative)
22allocate(w(nwplot))
23w1=max(wplot(1),0.d0)
24w2=max(wplot(2),w1)
25t1=(w2-w1)/dble(nwplot)
26do iw=1,nwplot
27 w(iw)=w1+t1*dble(iw-1)
28end do
29! i divided by the complex relaxation time
30eta=cmplx(0.d0,swidth,8)
31! determine the weights for the spline integration
32allocate(wt(ntimes))
33call wsplint(ntimes,times,wt)
34! compute the electric field from E = -1/c dA/dt and Fourier transform
35allocate(ew(nwplot,3))
36t0=-1.d0/solsc
37do i=1,3
38 if (task == 480) then
39! Fourier transform A(t) numerically to obtain A(ω)
40 call zftft(w,wt,3,afieldt(i,1),ew(:,i))
41! take the time derivative E(t)=-1/c dA(t)/dt analytically to get E(ω)
42 ew(:,i)=t0*zmi*w(:)*ew(:,i)
43! filter the high-frequency components from E(ω) with a Lorentzian convolution
44 call zlrzncnv(nwplot,swidth,w,ew(:,i))
45 else
46! analytic Fourier transform of E(t) assumed to be a delta function at t=0
47 t1=t0*afieldt(i,1)
48 ew(:,i)=t1
49 end if
50end do
51! read in the total current from file
52allocate(jt(3,ntimes))
53call readjtot(jt)
54! divide by the unit cell volume
55jt(:,:)=jt(:,:)/omega
56! set the constant part of J(t) to zero if required; this effectively removes
57! the Drude term
58if (jtconst0) then
59 do i=1,3
60 t1=dot_product(wt(:),jt(i,:))
61 t1=t1/times(ntimes)
62 jt(i,1:ntimes)=jt(i,1:ntimes)-t1
63 end do
64end if
65! Fourier transform the current
66allocate(jw(nwplot,3))
67do i=1,3
68 call zftft(w,wt,3,jt(i,1),jw(:,i))
69! filter the high-frequency components from J(ω) with a Lorentzian convolution
70 call zlrzncnv(nwplot,swidth,w,jw(:,i))
71end do
72deallocate(wt,jt)
73! compute the dielectric function and write to file
74allocate(eps(nwplot))
75do i=1,3
76 do j=1,3
77 do iw=1,nwplot
78 z1=jw(iw,i)
79 z2=ew(iw,j)
80 t1=abs(z2%re)+abs(z2%im)
81 if (t1 > 1.d-8) then
82 z1=z1/z2
83 else
84 z1=0.d0
85 end if
86 z1=fourpi*zi*z1/(w(iw)+eta)
87 if (i == j) z1=z1+1.d0
88 eps(iw)=z1
89 end do
90! filter the high-frequency components of ϵ
91 call zlrzncnv(nwplot,2.d0*swidth,w,eps)
92 write(fname,'("EPSILON_TDRT_",2I1,".OUT")') i,j
93 open(50,file=trim(fname),form='FORMATTED')
94 do iw=1,nwplot
95 write(50,'(2G18.10)') w(iw),dble(eps(iw))
96 end do
97 write(50,*)
98 do iw=1,nwplot
99 write(50,'(2G18.10)') w(iw),aimag(eps(iw))
100 end do
101 close(50)
102 end do
103end do
104! write Fourier transform of total current to file
105open(50,file='JTOTW.OUT',form='FORMATTED')
106do i=1,3
107 do iw=1,nwplot
108 write(50,'(3G18.10)') w(iw),jw(iw,i)
109 end do
110 write(50,*)
111end do
112close(50)
113write(*,*)
114write(*,'("Info(dielectric_tdrt):")')
115write(*,'(" dielectric tensor determined from real-time evolution")')
116write(*,'(" written to EPSILON_TDRT_ij.OUT for components i,j = 1,2,3")')
117write(*,*)
118write(*,'(" (Note that only those components which are not orthogonal to the")')
119write(*,'(" applied A-field will be calculated correctly)")')
120write(*,*)
121write(*,ω'(" Fourier transform of total current J() written to JTOTW.OUT")')
122deallocate(w,ew,jw,eps)
123end subroutine
124
subroutine dielectric_tdrt
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
complex(8), parameter zmi
Definition modmain.f90:1239
real(8) omega
Definition modmain.f90:20
complex(8), parameter zi
Definition modmain.f90:1239
integer nwplot
Definition modmain.f90:1070
real(8), parameter fourpi
Definition modmain.f90:1231
real(8) swidth
Definition modmain.f90:892
real(8), dimension(2) wplot
Definition modmain.f90:1076
integer task
Definition modmain.f90:1298
real(8) solsc
Definition modmain.f90:1252
real(8), dimension(:), allocatable times
Definition modtddft.f90:48
real(8), dimension(:,:), allocatable afieldt
Definition modtddft.f90:58
logical jtconst0
Definition modtddft.f90:108
integer ntimes
Definition modtddft.f90:42
subroutine readjtot(jt)
Definition readjtot.f90:7
subroutine wsplint(n, x, w)
Definition wsplint.f90:7
subroutine zftft(w, wt, ld, ft, fw)
Definition zftft.f90:7
subroutine zlrzncnv(n, s, w, zf)
Definition zlrzncnv.f90:7