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