TABLE OF CONTENTS
ABINIT/kgindex [ Functions ]
NAME
kgindex
FUNCTION
Compute the index of each plane wave on a FFT grid.
COPYRIGHT
Copyright (C) 2001-2018 ABINIT group (XG) 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
kg_k(3,npw_k)=dimensionless coords of G vecs (integer) mpi_enreg=information about MPI parallelization ngfft(18)=contain all needed information about 3D FFT, see ~abinit/doc/variables/vargs.htm#ngfft npw_k=number of planewaves
OUTPUT
indpw_k(npw_k)=linear list number (in fft box) of given G vector for the current processor (local adress) =0 if kg_k(ipw) is not treated by this procesor mask(npw_k)=True if kg_k(ipw) belongs to this procesor, false otherwise.
NOTES
mpi_enreg is not necessary in this case (the info is also in ngfft), but much more easy to read...
PARENTS
m_fft_prof,m_gsphere,m_screening,m_shirley,m_wfd,prcref,prcref_PMA
CHILDREN
SOURCE
36 #if defined HAVE_CONFIG_H 37 #include "config.h" 38 #endif 39 40 #include "abi_common.h" 41 42 subroutine kgindex(indpw_k,kg_k,mask,mpi_enreg,ngfft,npw_k) 43 44 use defs_basis 45 use defs_abitypes 46 use m_profiling_abi 47 use m_errors 48 49 !This section has been created automatically by the script Abilint (TD). 50 !Do not modify the following lines by hand. 51 #undef ABI_FUNC 52 #define ABI_FUNC 'kgindex' 53 !End of the abilint section 54 55 implicit none 56 57 !Arguments ------------------------------------ 58 !scalars 59 integer,intent(in) :: npw_k 60 type(MPI_type),intent(in) :: mpi_enreg 61 !arrays 62 integer,intent(in) :: kg_k(3,npw_k),ngfft(18) 63 integer,intent(out) :: indpw_k(npw_k) 64 logical,intent(out) :: mask(npw_k) 65 !Local variables------------------------------- 66 !scalars 67 integer :: ig,ig1,ig2,ig3,me_fft,n1,n2,n3,nd2 68 character(len=500) :: msg 69 !arrays 70 integer, ABI_CONTIGUOUS pointer :: fftn2_distrib(:),ffti2_local(:) 71 !integer, ABI_CONTIGUOUS pointer :: fftn3_distrib(:),ffti3_local(:) 72 73 ! ************************************************************************* 74 75 n1=ngfft(1); n2=ngfft(2); n3=ngfft(3) 76 77 !Use the following indexing (N means ngfft of the adequate direction) 78 !0 1 2 3 ... N/2 -(N-1)/2 ... -1 <= kg 79 !1 2 3 4 ....N/2+1 N/2+2 ... N <= index 80 81 me_fft=mpi_enreg%me_fft 82 nd2=(n2-1)/mpi_enreg%nproc_fft+1 83 84 if (n2== mpi_enreg%distribfft%n2_coarse) then 85 fftn2_distrib => mpi_enreg%distribfft%tab_fftdp2_distrib 86 ffti2_local => mpi_enreg%distribfft%tab_fftdp2_local 87 else if (n2 == mpi_enreg%distribfft%n2_fine) then 88 fftn2_distrib => mpi_enreg%distribfft%tab_fftdp2dg_distrib 89 ffti2_local => mpi_enreg%distribfft%tab_fftdp2dg_local 90 else 91 MSG_BUG("Unable to find an allocated distrib for this fft grid") 92 end if 93 94 !call ptabs_fourwf(mpi_enreg,n2,n3,fftn2_distrib,ffti2_local,fftn3_distrib,ffti3_local) 95 96 do ig=1,npw_k 97 ig1=modulo(kg_k(1,ig),n1) 98 ig2=modulo(kg_k(2,ig),n2) 99 ig3=modulo(kg_k(3,ig),n3) 100 if(me_fft==fftn2_distrib(ig2+1)) then 101 ig2=ffti2_local(ig2+1) - 1 102 indpw_k(ig)=ig1+1+n1*(ig2+nd2*ig3) 103 mask(ig)=.true. 104 else 105 indpw_k(ig)=0 106 mask(ig)=.false. 107 end if 108 if ( ANY(kg_k(:,ig)>ngfft(1:3)/2) .or. ANY(kg_k(:,ig)<-(ngfft(1:3)-1)/2) ) then 109 write(msg,'(a,i0,a)')" The G-vector with ig: ",ig," falls outside the FFT box." 110 MSG_ERROR(msg) 111 end if 112 end do 113 114 end subroutine kgindex