The Elk Code
zftft.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2025 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 zftft(w,wt,ld,ft,fw)
7 use modmain
8 use modtddft
9 use modomp
10 implicit none
11 ! arguments
12 real(8), intent(in) :: w(nwplot),wt(ntimes)
13 integer, intent(in) :: ld
14 real(8), intent(in) :: ft(ld,ntimes)
15 complex(8), intent(out) :: fw(nwplot)
16 ! local variables
17 integer iw,its,nthd
18 real(8) t1,t2
19 ! automatic arrays
20 real(8) f1(ntimes),f2(ntimes)
21 call holdthd(nwplot,nthd)
22 !$OMP PARALLEL DO DEFAULT(SHARED) &
23 !$OMP PRIVATE(f1,f2,its,t1,t2) &
24 !$OMP SCHEDULE(DYNAMIC) &
25 !$OMP NUM_THREADS(nthd)
26 do iw=1,nwplot
27  do its=1,ntimes
28  t1=ft(1,its)
29  t2=w(iw)*times(its)
30  f1(its)=t1*cos(t2)
31  f2(its)=t1*sin(t2)
32  end do
33  t1=dot_product(wt(:),f1(:))
34  t2=dot_product(wt(:),f2(:))
35  fw(iw)=cmplx(t1,t2,8)
36 end do
37 !$OMP END PARALLEL DO
38 call freethd(nthd)
39 end subroutine
40 
subroutine zftft(w, wt, ld, ft, fw)
Definition: zftft.f90:7
Definition: modomp.f90:6
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78