The Elk Code
writeefieldw.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2024 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 writeefieldw
7 use modmain
8 use modtddft
9 implicit none
10 ! local variables
11 integer iw,i
12 real(8) w1,w2,t0,t1
13 complex(8) z1
14 ! allocatable arrays
15 real(8), allocatable :: w(:),wt(:)
16 complex(8), allocatable :: ew(:,:)
17 ! read time-dependent A-field from file
18 call readafieldt
19 ! generate energy grid (always non-negative)
20 allocate(w(nwplot))
21 w1=max(wplot(1),0.d0)
22 w2=max(wplot(2),w1)
23 t1=(w2-w1)/dble(nwplot)
24 do iw=1,nwplot
25  w(iw)=w1+t1*dble(iw-1)
26 end do
27 ! determine the weights for the spline integration
28 allocate(wt(ntimes))
29 call wsplint(ntimes,times,wt)
30 ! compute the electric field from E = -1/c dA/dt and Fourier transform
31 allocate(ew(nwplot,3))
32 t0=-1.d0/solsc
33 do i=1,3
34 ! Fourier transform A(t) numerically to obtain A(ω)
35  call zftft(w,wt,3,afieldt(i,1),ew(:,i))
36 ! take the time derivative E(t)=-1/c dA(t)/dt analytically to get E(ω)
37  do iw=1,nwplot
38  z1=ew(iw,i)
39  ew(iw,i)=t0*w(iw)*cmplx(z1%im,-z1%re,8)
40  end do
41 ! filter the high-frequency components from E(ω) with a Lorentzian convolution
42  call zlrzncnv(nwplot,swidth,w,ew(:,i))
43 end do
44 ! write Fourier transform of electric field to file
45 open(50,file='EFIELDW.OUT',form='FORMATTED')
46 do i=1,3
47  do iw=1,nwplot
48  write(50,'(3G18.10)') w(iw),ew(iw,i)
49  end do
50  write(50,*)
51 end do
52 close(50)
53 write(*,*)
54 write(*,'("Info(writeefieldw):")')
55 write(*,'(" Fourier transform of electric field E(ω) written to EFIELDW.OUT")')
56 deallocate(w,wt,ew)
57 end subroutine
58 
subroutine zftft(w, wt, ld, ft, fw)
Definition: zftft.f90:7
subroutine writeefieldw
Definition: writeefieldw.f90:7
integer ntimes
Definition: modtddft.f90:42
real(8) swidth
Definition: modmain.f90:895
real(8), dimension(2) wplot
Definition: modmain.f90:1079
subroutine readafieldt
Definition: readafieldt.f90:7
real(8) solsc
Definition: modmain.f90:1253
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
subroutine zlrzncnv(n, s, w, zf)
Definition: zlrzncnv.f90:7
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
subroutine wsplint(n, x, w)
Definition: wsplint.f90:7
integer nwplot
Definition: modmain.f90:1073