TABLE OF CONTENTS
ABINIT/ph1d3d [ Functions ]
NAME
ph1d3d
FUNCTION
Compute the three-dimensional phase factor $e^{i 2 \pi (k+G) cdot xred}$ from the three one-dimensional factors, the k point coordinates, and the atom positions, for all planewaves which fit in the fft box.
COPYRIGHT
Copyright (C) 1998-2017 ABINIT group (XG,DR) 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
iatom, jatom= bounds of atom indices in ph1d for which ph3d has to be computed kg_k(3,npw_k)=reduced planewave coordinates. matblk= dimension of ph3d natom= dimension of ph1d npw=number of plane waves n1,n2,n3=dimensions of fft box (ngfft(3)). phkxred(2,natom)=phase factors exp(2 pi k.xred) ph1d(2,(2*n1+1)*natom+(2*n2+1)*natom+(2*n3+1)*natom)=exp(2Pi i G xred) for vectors (Gx,0,0), (0,Gy,0) and (0,0,Gz) with components ranging from -nj <= Gj <= nj
OUTPUT
ph3d(2,npw_k,matblk)=$e^{2 i \pi (k+G) cdot xred}$ for vectors (Gx,Gy,Gz), and for atoms in the range iatom to jatom with respect to ph1d
PARENTS
cg_rotate,ctocprj,getcprj,m_cut3d,m_epjdos,m_fock,m_hamiltonian,m_wfd nonlop_pl,nonlop_ylm,partial_dos_fractions,suscep_stat,wfconv
CHILDREN
SOURCE
42 #if defined HAVE_CONFIG_H 43 #include "config.h" 44 #endif 45 46 #include "abi_common.h" 47 48 subroutine ph1d3d(iatom,jatom,kg_k,matblk,natom,npw_k,n1,n2,n3,phkxred,ph1d,ph3d) 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 'ph1d3d' 58 !End of the abilint section 59 60 implicit none 61 62 !Arguments ------------------------------------ 63 !scalars 64 integer,intent(in) :: iatom,jatom,matblk,n1,n2,n3,natom,npw_k 65 !arrays 66 integer,intent(in) :: kg_k(3,npw_k) 67 real(dp),intent(in) :: ph1d(2,(2*n1+1+2*n2+1+2*n3+1)*natom) 68 real(dp),intent(in) :: phkxred(2,natom) 69 real(dp),intent(out) :: ph3d(2,npw_k,matblk) 70 71 !Local variables------------------------------- 72 !scalars 73 integer :: i1,ia,iatblk,ig,shift1,shift2,shift3 74 real(dp) :: ph12i,ph12r,ph1i,ph1r,ph2i,ph2r,ph3i,ph3r,phkxi,phkxr 75 character(len=500) :: message 76 !arrays 77 real(dp),allocatable :: ph1kxred(:,:) 78 79 ! ************************************************************************* 80 81 if(matblk-1 < jatom-iatom)then 82 write(message,'(a,a,a,a,a,i0,a,a,i0,a,i0,a)')& 83 & 'Input natom-1 must be larger or equal to jatom-iatom,',ch10,& 84 & 'while their value is : ',ch10,& 85 & 'natom-1 = ',natom-1,ch10,& 86 & 'jatom=',jatom,', iatom=',iatom,'.' 87 MSG_BUG(message) 88 end if 89 90 ABI_ALLOCATE(ph1kxred,(2,-n1:n1)) 91 92 !ia runs from iatom to jatom 93 do ia=iatom,jatom 94 95 ! iatblk runs from 1 to matblk 96 iatblk=ia-iatom+1 97 !write(87,*) iatblk 98 shift1=1+n1+(ia-1)*(2*n1+1) 99 shift2=1+n2+(ia-1)*(2*n2+1)+natom*(2*n1+1) 100 shift3=1+n3+(ia-1)*(2*n3+1)+natom*(2*n1+1+2*n2+1) 101 ! Compute product of phkxred by phase for the first component of G vector 102 phkxr=phkxred(1,ia) 103 phkxi=phkxred(2,ia) 104 ! DEBUG (needed to compare with version prior to 2.0) 105 ! phkxr=1.0d0 106 ! phkxi=0.0d0 107 ! ENDDEBUG 108 do i1=-n1,n1 109 ph1kxred(1,i1)=ph1d(1,i1+shift1)*phkxr-ph1d(2,i1+shift1)*phkxi 110 ph1kxred(2,i1)=ph1d(2,i1+shift1)*phkxr+ph1d(1,i1+shift1)*phkxi 111 end do 112 113 ! Compute tri-dimensional phase factor 114 !$OMP PARALLEL DO PRIVATE(ig,ph1r,ph1i,ph2r,ph2i,ph3r,ph3i,ph12r,ph12i) 115 do ig=1,npw_k 116 ph1r=ph1kxred(1,kg_k(1,ig)) 117 ph1i=ph1kxred(2,kg_k(1,ig)) 118 ph2r=ph1d(1,kg_k(2,ig)+shift2) 119 ph2i=ph1d(2,kg_k(2,ig)+shift2) 120 ph3r=ph1d(1,kg_k(3,ig)+shift3) 121 ph3i=ph1d(2,kg_k(3,ig)+shift3) 122 ph12r=ph1r*ph2r-ph1i*ph2i 123 ph12i=ph1r*ph2i+ph1i*ph2r 124 !if(ig==487) then 125 !write(87,*)iatblk,ph3d(1,ig,iatblk),ph12r,ph3r,ph12i,ph3i 126 !endif 127 ph3d(1,ig,iatblk)=ph12r*ph3r-ph12i*ph3i 128 ph3d(2,ig,iatblk)=ph12r*ph3i+ph12i*ph3r 129 end do 130 !$OMP END PARALLEL DO 131 end do 132 !write(87,*)ph3d(1,487,8) 133 ABI_DEALLOCATE(ph1kxred) 134 135 end subroutine ph1d3d