TABLE OF CONTENTS
ABINIT/mkfsqgrid [ Functions ]
NAME
mkfsqgrid
FUNCTION
This routine sets up the qpoints between the full FS kpt grid points
COPYRIGHT
Copyright (C) 2004-2018 ABINIT group (MVer) This file is distributed under the terms of the GNU General Public Licence, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
kpt_phon = kpoints on the FS nkpt_phon = number of kpts on the FS
OUTPUT
nFSqpt = full number of qpoints between points on the FS FStoqpt = qpoint index for each pair of kpt on the FS tmpFSqpt = temp array with coordinates of the qpoints between points on the FS
NOTES
might need the inverse indexing qpt -> (kpt1,kpt2) (kpt3,kpt4) ...
PARENTS
CHILDREN
wrap2_pmhalf
SOURCE
37 #if defined HAVE_CONFIG_H 38 #include "config.h" 39 #endif 40 41 #include "abi_common.h" 42 43 44 subroutine mkfsqgrid(kpt_phon,FStoqpt,nkpt_phon,nFSqpt,tmpFSqpt) 45 46 use defs_basis 47 use m_profiling_abi 48 use m_numeric_tools, only : wrap2_pmhalf 49 50 !This section has been created automatically by the script Abilint (TD). 51 !Do not modify the following lines by hand. 52 #undef ABI_FUNC 53 #define ABI_FUNC 'mkfsqgrid' 54 !End of the abilint section 55 56 implicit none 57 58 !Arguments ------------------------------------ 59 !scalars 60 integer,intent(in) :: nkpt_phon 61 integer,intent(out) :: nFSqpt 62 !arrays 63 integer,intent(out) :: FStoqpt(nkpt_phon,nkpt_phon) 64 real(dp),intent(in) :: kpt_phon(3,nkpt_phon) 65 real(dp),intent(out) :: tmpFSqpt(3,nkpt_phon*nkpt_phon) 66 67 !Local variables------------------------------- 68 !scalars 69 integer :: iFSqpt,ikpt1,ikpt2,new 70 real(dp) :: shift,ss 71 !arrays 72 real(dp) :: k1(3),kpt(3) 73 74 ! ************************************************************************* 75 76 nFSqpt=0 77 tmpFSqpt(:,:)=zero 78 do ikpt1=1,nkpt_phon 79 do ikpt2=1,nkpt_phon 80 k1(:) = kpt_phon(:,ikpt1) - kpt_phon(:,ikpt2) 81 call wrap2_pmhalf(k1(1),kpt(1),shift) 82 call wrap2_pmhalf(k1(2),kpt(2),shift) 83 call wrap2_pmhalf(k1(3),kpt(3),shift) 84 85 new=1 86 ! is kpt among the FS qpts found already? 87 do iFSqpt=1,nFSqpt 88 ss=(kpt(1)-tmpFSqpt(1,iFSqpt))**2 + & 89 & (kpt(2)-tmpFSqpt(2,iFSqpt))**2 + & 90 & (kpt(3)-tmpFSqpt(3,iFSqpt))**2 91 if (ss < tol6) then 92 FStoqpt(ikpt1,ikpt2) = iFSqpt 93 new=0 94 exit 95 end if 96 end do 97 if (new == 1) then 98 nFSqpt=nFSqpt+1 99 tmpFSqpt(:,nFSqpt) = kpt(:) 100 FStoqpt(ikpt1,ikpt2) = nFSqpt 101 end if 102 103 end do 104 end do 105 106 !got nFSqpt,tmpFSqpt,FStoqpt 107 108 end subroutine mkfsqgrid