TABLE OF CONTENTS
ABINIT/fourdp_6d [ Functions ]
NAME
fourdp_6d
FUNCTION
calculate a 6-dimensional Fast Fourier Transform isign=-1 : A(G1,G2) = Sum(r1,r2) A(r1,r2) exp(-iG1.r1) exp(+iG2.r2) ^ ^ isign=+1 : A(r1,r2) = Sum(G1,G2) A(G1,G2) exp(+iG1.r1) exp(-iG2.r2) ^ ^ isign=-1 and isign=1 form a transform/inverse-transform pair: calling one then the other will take you back to the original function, multiplied by a factor of (nl*nm*nn)**2. ------------------------------------------------------------------ input: a: A(r1,r2) [overwritten] output: a: A(G1,G2) in the format IGFFT ------------------------------------------------------------------
COPYRIGHT
Copyright (C) 2008-2018 ABINIT group (MG) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
INPUTS
OUTPUT
SIDE EFFECTS
NOTES
PARENTS
m_kxc
CHILDREN
fourdp
SOURCE
46 #if defined HAVE_CONFIG_H 47 #include "config.h" 48 #endif 49 50 #include "abi_common.h" 51 52 subroutine fourdp_6d(cplex,matrix,isign,MPI_enreg,nfft,ngfft,paral_kgb,tim_fourdp) 53 54 use defs_basis 55 use defs_abitypes 56 use m_profiling_abi 57 use m_errors 58 59 !This section has been created automatically by the script Abilint (TD). 60 !Do not modify the following lines by hand. 61 #undef ABI_FUNC 62 #define ABI_FUNC 'fourdp_6d' 63 use interfaces_53_ffts, except_this_one => fourdp_6d 64 !End of the abilint section 65 66 implicit none 67 68 !Arguments ------------------------------------ 69 !scalars 70 integer,intent(in) :: cplex,isign,nfft,paral_kgb,tim_fourdp 71 type(MPI_type),intent(in) :: MPI_enreg 72 !arrays 73 integer,intent(in) :: ngfft(18) 74 complex(gwpc),intent(inout) :: matrix(nfft,nfft) 75 76 !Local variables------------------------------- 77 !scalars 78 !integer,parameter :: cplex=2 79 integer :: i1,i2,i3,ifft 80 integer :: n1,n2,n3 81 !arrays 82 real(dp),allocatable :: fofg(:,:),fofr(:) 83 84 ! ************************************************************************* 85 86 !TODO check normalization factor, it is better if we use the GW conventions. 87 n1 = ngfft(1) 88 n2 = ngfft(2) 89 n3 = ngfft(3) 90 91 ABI_MALLOC(fofg,(2,nfft)) 92 ABI_MALLOC(fofr,(cplex*nfft)) 93 94 do i3=0,n3-1 95 do i2=0,n2-1 96 do i1=0,n1-1 97 98 ifft=1+i1+i2*n1+i3*n1*n2 99 if (isign==1) then ! G1 -> r1 transform for each G2 to form A(r1,G2) 100 fofg(1,:)=REAL (matrix(:,ifft)) 101 fofg(2,:)=AIMAG(matrix(:,ifft)) 102 else if (isign==-1) then ! r1 -> G1 transform for each r2 to form A(G1,r2) 103 fofr(1:nfft) =REAL (matrix(:,ifft)) 104 fofr(nfft+1:2*nfft)=AIMAG(matrix(:,ifft)) 105 else 106 MSG_ERROR("Wrong isign") 107 end if 108 109 call fourdp(cplex,fofg,fofr,isign,MPI_enreg,nfft,ngfft,paral_kgb,tim_fourdp) 110 111 if (isign==1) then ! Save A(r1,G2) 112 matrix(:,ifft)=CMPLX(fofr(1:nfft),fofr(nfft+1:2*nfft)) 113 else if (isign==-1) then ! Save A(G1,r2) 114 matrix(:,ifft)=CMPLX(fofg(1,:),fofg(2,:)) 115 end if 116 117 end do 118 end do 119 end do 120 121 do i3=0,n3-1 122 do i2=0,n2-1 123 do i1=0,n1-1 124 125 ifft=1+i1+i2*n1+i3*n1*n2 126 if (isign==1) then ! Do the G2 -> r2 transform of A(r1,G2) to get A(r1,r2) 127 fofr(1:nfft )=REAL (matrix(ifft,:)) 128 fofr(nfft+1:2*nfft)=AIMAG(matrix(ifft,:)) 129 else if (isign==-1) then 130 ! Do the r2 -> G2 transform of A(G1,r2) to get A(G1,G2) 131 fofg(1,:)=REAL (matrix(ifft,:)) 132 fofg(2,:)=AIMAG(matrix(ifft,:)) 133 end if 134 135 call fourdp(2,fofg,fofr,-isign,MPI_enreg,nfft,ngfft,paral_kgb,tim_fourdp) 136 137 if (isign==1) then 138 matrix(ifft,:)=CMPLX(fofg(1,:),fofg(2,:)) 139 else if (isign==-1) then 140 matrix(ifft,:)=CMPLX(fofr(1:nfft),fofr(nfft+1:2*nfft)) 141 end if 142 143 end do 144 end do 145 end do 146 147 ABI_FREE(fofg) 148 ABI_FREE(fofr) 149 150 end subroutine fourdp_6d