TABLE OF CONTENTS
ABINIT/getph [ Functions ]
NAME
getph
FUNCTION
Compute three factors of one-dimensional structure factor phase for input atomic coordinates, for all planewaves which fit in fft box. The storage of these atomic factors is made according to the values provided by the index table atindx. This will save time in nonlop.
COPYRIGHT
Copyright (C) 1998-2017 ABINIT group (DCA, 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 . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
atindx(natom)=index table for atoms (see gstate.f) natom=number of atoms in cell. n1,n2,n3=dimensions of fft box (ngfft(3)). xred(3,natom)=reduced atomic coordinates.
OUTPUT
ph1d(2,(2*n1+1)*natom+(2*n2+1)*natom+(2*n3+1)*natom)=exp(2Pi i G.xred) for integer vector G with components ranging from -nj <= G <= nj. Real and imag given in usual Fortran convention.
PARENTS
afterscfloop,bethe_salpeter,cg_rotate,dfpt_looppert,dfptnl_loop extrapwf,gstate,m_cut3d,m_fock,m_gkk,m_hamiltonian,m_phgamma,m_phpi m_sigmaph,m_wfd,partial_dos_fractions,prcref,prcref_PMA,respfn,scfcv screening,sigma,update_e_field_vars,wfconv
CHILDREN
SOURCE
41 #if defined HAVE_CONFIG_H 42 #include "config.h" 43 #endif 44 45 #include "abi_common.h" 46 47 48 subroutine getph(atindx,natom,n1,n2,n3,ph1d,xred) 49 50 use defs_basis 51 use m_profiling_abi 52 use m_errors 53 54 !This section has been created automatically by the script Abilint (TD). 55 !Do not modify the following lines by hand. 56 #undef ABI_FUNC 57 #define ABI_FUNC 'getph' 58 !End of the abilint section 59 60 implicit none 61 62 !Arguments ------------------------------------ 63 !scalars 64 integer,intent(in) :: n1,n2,n3,natom 65 !arrays 66 integer,intent(in) :: atindx(natom) 67 real(dp),intent(in) :: xred(3,natom) 68 real(dp),intent(out) :: ph1d(:,:) 69 70 !Local variables------------------------------- 71 !scalars 72 integer,parameter :: im=2,re=1 73 integer :: i1,i2,i3,ia,ii,ph1d_size1,ph1d_size2,ph1d_sizemin 74 !character(len=500) :: msg 75 real(dp) :: arg 76 77 ! ************************************************************************* 78 79 ph1d_size1=size(ph1d,1);ph1d_size2=size(ph1d,2) 80 ph1d_sizemin=(2*n1+1+2*n2+1+2*n3+1)*natom 81 if (ph1d_size1/=2.or.ph1d_size2<ph1d_sizemin) then 82 MSG_BUG('Wrong ph1d sizes!') 83 end if 84 85 do ia=1,natom 86 87 ! Store the phase factor of atom number ia in place atindx(ia) 88 i1=(atindx(ia)-1)*(2*n1+1) 89 i2=(atindx(ia)-1)*(2*n2+1)+natom*(2*n1+1) 90 i3=(atindx(ia)-1)*(2*n3+1)+natom*(2*n1+1+2*n2+1) 91 92 do ii=1,2*n1+1 93 arg=two_pi*dble(ii-1-n1)*xred(1,ia) 94 ph1d(re,ii+i1)=dcos(arg) 95 ph1d(im,ii+i1)=dsin(arg) 96 end do 97 98 do ii=1,2*n2+1 99 arg=two_pi*dble(ii-1-n2)*xred(2,ia) 100 ph1d(re,ii+i2)=dcos(arg) 101 ph1d(im,ii+i2)=dsin(arg) 102 end do 103 104 do ii=1,2*n3+1 105 arg=two_pi*dble(ii-1-n3)*xred(3,ia) 106 ph1d(re,ii+i3)=dcos(arg) 107 ph1d(im,ii+i3)=dsin(arg) 108 end do 109 110 end do 111 112 !This is to avoid uninitialized ph1d values 113 if (ph1d_sizemin<ph1d_size2) then 114 ph1d(:,ph1d_sizemin+1:ph1d_size2)=zero 115 end if 116 117 end subroutine getph