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
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
! allocatable arrays
14
real
(8),
allocatable
:: w(:),wt(:)
15
complex(8)
,
allocatable
:: ew(:,:)
16
! read time-dependent A-field from file
17
call
readafieldt
18
! generate energy grid (always non-negative)
19
allocate
(w(
nwplot
))
20
w1=max(
wplot
(1),0.d0)
21
w2=max(
wplot
(2),w1)
22
t1=(w2-w1)/dble(
nwplot
)
23
do
iw=1,
nwplot
24
w(iw)=w1+t1*dble(iw-1)
25
end do
26
! determine the weights for the spline integration
27
allocate
(wt(
ntimes
))
28
call
wsplint
(
ntimes
,
times
,wt)
29
! compute the electric field from E = -1/c dA/dt and Fourier transform
30
allocate
(ew(
nwplot
,3))
31
t0=-1.d0/
solsc
32
do
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))
39
end do
40
! write Fourier transform of electric field to file
41
open
(50,file=
'EFIELDW.OUT'
,form=
'FORMATTED'
)
42
do
i=1,3
43
do
iw=1,
nwplot
44
write
(50,
'(3G18.10)'
) w(iw),ew(iw,i)
45
end do
46
write
(50,*)
47
end do
48
close
(50)
49
write
(*,*)
50
write
(*,
'("Info(writeefieldw):")'
)
51
write
(*,ω
'(" Fourier transform of electric field E() written to EFIELDW.OUT")'
)
52
deallocate
(w,wt,ew)
53
end subroutine
54
modmain
Definition
modmain.f90:6
modmain::zmi
complex(8), parameter zmi
Definition
modmain.f90:1239
modmain::nwplot
integer nwplot
Definition
modmain.f90:1070
modmain::swidth
real(8) swidth
Definition
modmain.f90:892
modmain::wplot
real(8), dimension(2) wplot
Definition
modmain.f90:1076
modmain::solsc
real(8) solsc
Definition
modmain.f90:1252
modtddft
Definition
modtddft.f90:6
modtddft::times
real(8), dimension(:), allocatable times
Definition
modtddft.f90:48
modtddft::afieldt
real(8), dimension(:,:), allocatable afieldt
Definition
modtddft.f90:58
modtddft::ntimes
integer ntimes
Definition
modtddft.f90:42
readafieldt
subroutine readafieldt
Definition
readafieldt.f90:7
writeefieldw
subroutine writeefieldw
Definition
writeefieldw.f90:7
wsplint
subroutine wsplint(n, x, w)
Definition
wsplint.f90:7
zftft
subroutine zftft(w, wt, ld, ft, fw)
Definition
zftft.f90:7
zlrzncnv
subroutine zlrzncnv(n, s, w, zf)
Definition
zlrzncnv.f90:7
writeefieldw.f90
Generated by
1.9.8