pwdft  0.1
PW-DFT code for education
fftw_wrapper Module Reference

Functions/Subroutines

subroutine init_fft ()
 
subroutine fft_r2g (VlocR, VlocG)
 v(r) e^{-iGr} -> V(G) More...
 
subroutine fft_g2r (VlocG, VlocR)
 V(G) e^{iGr} -> v(r) More...
 
subroutine fft_r2g_w (wfR, wfG)
 w(r) e^{-iGr} -> W(G) -> compress More...
 
subroutine fft_g2r_w (wfG, wfR)
 Uncompress -> W(G) e^{iGr} -> w(r) More...
 

Variables

integer, dimension(:), allocatable, save w2r
 (g_wfnpw) g_wfnpw -> g_rhnr More...
 

Function/Subroutine Documentation

◆ init_fft()

subroutine fftw_wrapper::init_fft

Definition at line 47 of file fftw_wrapper.F90.

47  !
48  use gvec, only : g_rh, g_wf
49  !
50  integer :: ipw, igv(3)
51  integer(c_int) :: fftwdim(3)
52  !
53  fftwdim(1:3) = (/g_rh%nft(3), g_rh%nft(2), g_rh%nft(1)/)
54  !
55  allocate(fft_in(g_rh%nr), fft_out(g_rh%nr))
56  !
57  plan_g2r = fftw_plan_dft(3, fftwdim, fft_in, fft_out, &
58  & fftw_backward, fftw_estimate)
59  !
60  plan_r2g = fftw_plan_dft(3, fftwdim, fft_in, fft_out, &
61  & fftw_forward, fftw_estimate)
62  !
63  allocate(w2r(g_wf%npw))
64  do ipw = 1, g_wf%npw
65  igv(1:3) = g_wf%mill(1:3, g_wf%map(ipw))
66  igv(1:3) = modulo(igv(1:3), g_rh%nft(1:3))
67  w2r(ipw) = 1 + igv(1) + g_rh%nft(1)*(igv(2) + g_rh%nft(2)*igv(3))
68  end do
69  !

References gvec::g_rh, gvec::g_wf, and w2r.

Referenced by pwdft().

Here is the caller graph for this function:

◆ fft_r2g()

subroutine fftw_wrapper::fft_r2g ( real(8), dimension(g_rh%nr), intent(in)  VlocR,
complex(8), dimension(g_rh%nr), intent(out)  VlocG 
)

v(r) e^{-iGr} -> V(G)

Definition at line 75 of file fftw_wrapper.F90.

75  !
76  use gvec, only : g_rh
77  !
78  real(8),intent(in) :: VlocR(g_rh%nr)
79  complex(8),intent(out) :: VlocG(g_rh%nr)
80  !
81  fft_in(1:g_rh%nr) = vlocr(1:g_rh%nr)
82  call fftw_execute_dft(plan_r2g, fft_in, fft_out)
83  vlocg(1:g_rh%nr) = fft_out(1:g_rh%nr) / dble(g_rh%nr)
84  !

References gvec::g_rh.

Referenced by diag_direct::direct(), energy::hartree(), and rho_v::hartree_pot().

Here is the caller graph for this function:

◆ fft_g2r()

subroutine fftw_wrapper::fft_g2r ( complex(8), dimension(g_rh%nr), intent(in)  VlocG,
real(8), dimension(g_rh%nr), intent(out)  VlocR 
)

V(G) e^{iGr} -> v(r)

Definition at line 90 of file fftw_wrapper.F90.

90  !
91  use gvec, only : g_rh
92  !
93  complex(8),intent(in) :: VlocG(g_rh%nr)
94  real(8),intent(out) :: VlocR(g_rh%nr)
95  !
96  fft_in(1:g_rh%nr) = vlocg(1:g_rh%nr)
97  call fftw_execute_dft(plan_g2r, fft_in, fft_out)
98  vlocr(1:g_rh%nr) = dble(fft_out(1:g_rh%nr))
99  !

References gvec::g_rh.

Referenced by rho_v::generate_vps(), and rho_v::hartree_pot().

Here is the caller graph for this function:

◆ fft_r2g_w()

subroutine fftw_wrapper::fft_r2g_w ( complex(8), dimension(g_rh%nr), intent(in)  wfR,
complex(8), dimension(g_wf%npw), intent(out)  wfG 
)

w(r) e^{-iGr} -> W(G) -> compress

Definition at line 105 of file fftw_wrapper.F90.

105  !
106  use atm_spec, only : vcell
107  use gvec, only : g_rh, g_wf
108  !
109  complex(8),intent(in) :: wfR(g_rh%nr)
110  complex(8),intent(out) :: wfG(g_wf%npw)
111  !
112  fft_in(1:g_rh%nr) = wfr(1:g_rh%nr)
113  call fftw_execute_dft(plan_r2g, fft_in, fft_out)
114  wfg(1:g_wf%npw) = fft_out(w2r(1:g_wf%npw)) * sqrt(vcell) / dble(g_rh%nr)
115  !

References gvec::g_rh, gvec::g_wf, atm_spec::vcell, and w2r.

Referenced by hamiltonian::h_psi().

Here is the caller graph for this function:

◆ fft_g2r_w()

subroutine fftw_wrapper::fft_g2r_w ( complex(8), dimension(g_wf%npw), intent(in)  wfG,
complex(8), dimension(g_rh%nr), intent(out)  wfR 
)

Uncompress -> W(G) e^{iGr} -> w(r)

Definition at line 121 of file fftw_wrapper.F90.

121  !
122  use atm_spec, only : vcell
123  use gvec, only : g_rh, g_wf
124  !
125  complex(8),intent(in) :: wfG(g_wf%npw)
126  complex(8),intent(out) :: wfR(g_rh%nr)
127  !
128  fft_in(1:g_rh%nr) = 0.0d0
129  fft_in(w2r(1:g_wf%npw)) = wfg(1:g_wf%npw)
130  call fftw_execute_dft(plan_g2r, fft_in, fft_out)
131  wfr(1:g_rh%nr) = fft_out(1:g_rh%nr) / sqrt(vcell)
132  !

References gvec::g_rh, gvec::g_wf, atm_spec::vcell, and w2r.

Referenced by hamiltonian::h_psi(), and k_point::ksum_rho().

Here is the caller graph for this function:

Variable Documentation

◆ w2r

integer, dimension(:), allocatable, save fftw_wrapper::w2r

(g_wfnpw) g_wfnpw -> g_rhnr

Definition at line 35 of file fftw_wrapper.F90.

35  integer,allocatable,save :: &
36  & w2r(:) !< (g_wf%npw) g_wf%npw -> g_rh%nr

Referenced by fft_g2r_w(), fft_r2g_w(), and init_fft().

atm_spec::vcell
real(8), save vcell
Unit cell volume [Bohr^3].
Definition: atm_spec.F90:50
gvec
Definition: gvec.F90:23
atm_spec
Definition: atm_spec.F90:23