The Elk Code
 
Loading...
Searching...
No Matches
cftwfir.f90
Go to the documentation of this file.
1
2! Copyright (C) 2017 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
6subroutine cftwfir(ngp,igpig,wfir)
7use modmain
8use modomp
9implicit none
10! arguments
11integer, intent(in) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
12complex(4), intent(inout) :: wfir(ngtc,nspinor,nstsv)
13! local variables
14integer ist,ispn,jspn
15integer n,igp,nthd
16real(4) t0
17! automatic arrays
18complex(4) c(ngkmax)
19t0=1.e0/sqrt(omega)
20call holdthd(nstsv,nthd)
21!$OMP PARALLEL DO DEFAULT(SHARED) &
22!$OMP PRIVATE(c,ispn,jspn,n,igp) &
23!$OMP SCHEDULE(DYNAMIC) &
24!$OMP NUM_THREADS(nthd)
25do ist=1,nstsv
26 do ispn=1,nspinor
27 jspn=jspnfv(ispn)
28 n=ngp(jspn)
29 c(1:n)=wfir(1:n,ispn,ist)
30 wfir(1:ngtc,ispn,ist)=0.e0
31 do igp=1,n
32 wfir(igfc(igpig(igp,jspn)),ispn,ist)=t0*c(igp)
33 end do
34 call cfftifc(3,ngdgc,1,wfir(:,ispn,ist))
35 end do
36end do
37!$OMP END PARALLEL DO
38call freethd(nthd)
39end subroutine
40
subroutine cfftifc(nd, n, sgn, c)
subroutine cftwfir(ngp, igpig, wfir)
Definition cftwfir.f90:7
integer, dimension(3) ngdgc
Definition modmain.f90:388
real(8) omega
Definition modmain.f90:20
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(2) jspnfv
Definition modmain.f90:291
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106