The Elk Code
rfzfftq.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 T. Mueller, 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 
6 subroutine rfzfftq(sgn,nf,ngt,rfmt,rfir,zfmt,zfir)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 integer, intent(in) :: sgn,nf,ngt
12 real(8), intent(inout) :: rfmt(npcmtmax,natmtot,nf,nqpt)
13 real(8), intent(inout) :: rfir(ngt,nf,nqpt)
14 complex(8), intent(inout) :: zfmt(npcmtmax,natmtot,nf,nfqrz)
15 complex(8), intent(inout) :: zfir(ngt,nf,nfqrz)
16 ! local variables
17 integer jf,is,ias,npc,i,nthd
18 ! automatic arrays
19 real(8) r(nqpt)
20 complex(8) z(nfqrz)
21 call holdthd(npcmtmax+ngt,nthd)
22 !$OMP PARALLEL DEFAULT(SHARED) &
23 !$OMP PRIVATE(r,z,jf,ias,is,npc,i) &
24 !$OMP NUM_THREADS(nthd)
25 if (sgn == -1) then
26 ! loop over the number of functions
27  do jf=1,nf
28 ! Fourier transform the muffin-tin function
29  do ias=1,natmtot
30  is=idxis(ias)
31  npc=npcmt(is)
32 !$OMP DO SCHEDULE(DYNAMIC)
33  do i=1,npc
34  r(1:nqpt)=rfmt(i,ias,jf,1:nqpt)
35  call rzfftifc(3,ngridq,-1,r,z)
36  zfmt(i,ias,jf,1:nfqrz)=z(1:nfqrz)
37  end do
38 !$OMP END DO NOWAIT
39  end do
40 ! Fourier transform the interstitial function
41 !$OMP DO SCHEDULE(DYNAMIC)
42  do i=1,ngt
43  r(1:nqpt)=rfir(i,jf,1:nqpt)
44  call rzfftifc(3,ngridq,-1,r,z)
45  zfir(i,jf,1:nfqrz)=z(1:nfqrz)
46  end do
47 !$OMP END DO NOWAIT
48 ! end loop over number of functions
49  end do
50 else
51 ! loop over the number of functions
52  do jf=1,nf
53  do ias=1,natmtot
54  is=idxis(ias)
55  npc=npcmt(is)
56 !$OMP DO SCHEDULE(DYNAMIC)
57  do i=1,npc
58  z(1:nfqrz)=zfmt(i,ias,jf,1:nfqrz)
59  call rzfftifc(3,ngridq,1,r,z)
60  rfmt(i,ias,jf,1:nqpt)=r(1:nqpt)
61  end do
62 !$OMP END DO NOWAIT
63  end do
64 !$OMP DO SCHEDULE(DYNAMIC)
65  do i=1,ngt
66  z(1:nfqrz)=zfir(i,jf,1:nfqrz)
67  call rzfftifc(3,ngridq,1,r,z)
68  rfir(i,jf,1:nqpt)=r(1:nqpt)
69  end do
70 !$OMP END DO NOWAIT
71 ! end loop over number of functions
72  end do
73 end if
74 !$OMP END PARALLEL
75 call freethd(nthd)
76 end subroutine
77 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
Definition: modomp.f90:6
integer, dimension(3) ngridq
Definition: modmain.f90:515
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine rzfftifc(nd, n, sgn, r, z)
subroutine rfzfftq(sgn, nf, ngt, rfmt, rfir, zfmt, zfir)
Definition: rfzfftq.f90:7