TABLE OF CONTENTS
ABINIT/test_ftgkk [ Functions ]
NAME
test_ftgkk
FUNCTION
Test the fourier transform routine ftgkk for the el-phon matrix elements
COPYRIGHT
Copyright (C) 2004-2018 ABINIT group (MVerstra) 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
elph_ds = elphon datastructure with matrix elements gprim = reciprocal lattice vectors natom = number of atoms nrpt = number of real space points for FT interpolation rpt = coordinates of real space points for FT interpolation qpt_full = qpoint coordinates wghatm = weights for pairs of atoms in FT interpolation
OUTPUT
SIDE EFFECTS
NOTES
MJV 18/5/2008 reverted to old syntax/use for ftgkk, with all ft being done in a batch. Might come back to 5.5 version with atomic FT in ftgkk, but later.
PARENTS
CHILDREN
ftgkk
SOURCE
40 #if defined HAVE_CONFIG_H 41 #include "config.h" 42 #endif 43 44 #include "abi_common.h" 45 46 47 subroutine test_ftgkk(elph_ds,gprim,natom,nrpt,rpt,qpt_full,wghatm) 48 49 use defs_basis 50 use defs_elphon 51 use m_profiling_abi 52 53 !This section has been created automatically by the script Abilint (TD). 54 !Do not modify the following lines by hand. 55 #undef ABI_FUNC 56 #define ABI_FUNC 'test_ftgkk' 57 use interfaces_77_ddb, except_this_one => test_ftgkk 58 !End of the abilint section 59 60 implicit none 61 62 !Arguments ------------------------------------ 63 !scalars 64 integer,intent(in) :: natom,nrpt 65 type(elph_type),intent(inout) :: elph_ds 66 !arrays 67 real(dp),intent(in) :: gprim(3,3),rpt(3,nrpt),qpt_full(3,elph_ds%nqpt_full) 68 real(dp),intent(in) :: wghatm(natom,natom,nrpt) 69 70 !Local variables------------------------------- 71 !scalars 72 integer :: ikpt_phon,iqpt,isppol,qtor,sz1,sz2 73 !arrays 74 real(dp),allocatable :: gkq_disk(:,:,:,:,:),tmp_gkq(:,:,:,:,:) 75 76 ! ************************************************************************* 77 78 !for each qpt do FT to recuperate original values 79 80 isppol = 1 81 qtor = 0 82 sz1=elph_ds%ngkkband*elph_ds%ngkkband 83 sz2=elph_ds%nbranch*elph_ds%nbranch 84 ABI_ALLOCATE(gkq_disk,(2,sz1,sz2,elph_ds%k_phon%nkpt,elph_ds%nsppol)) 85 ABI_ALLOCATE(tmp_gkq,(2,sz1,sz2,elph_ds%k_phon%nkpt,elph_ds%nsppol)) 86 87 do iqpt=1,elph_ds%nqpt_full 88 tmp_gkq(:,:,:,:,:) = zero 89 90 call ftgkk (wghatm,tmp_gkq,elph_ds%gkk_rpt,elph_ds%gkqwrite,& 91 & elph_ds%gkk_rptwrite,gprim,1,natom,& 92 & elph_ds%k_phon%nkpt,elph_ds%ngkkband,elph_ds%k_phon%nkpt,1,& 93 & nrpt,elph_ds%nsppol,qtor,rpt,qpt_full,elph_ds%unit_gkk_rpt,elph_ds%unitgkq) 94 95 if (elph_ds%gkqwrite == 0) then 96 do ikpt_phon=1,10 97 write (93,*) tmp_gkq(:,:,:,ikpt_phon,isppol)-elph_ds%gkk_qpt(:,:,:,ikpt_phon,isppol,iqpt) 98 end do 99 else 100 do ikpt_phon=1, elph_ds%k_phon%nkpt 101 read (elph_ds%unitgkq,REC=((iqpt-1)*elph_ds%k_phon%nkpt+ikpt_phon)) gkq_disk(:,:,:,ikpt_phon,:) 102 end do 103 do ikpt_phon=1,10 104 write (93,*) tmp_gkq(:,:,:,ikpt_phon,isppol)-gkq_disk(:,:,:,ikpt_phon,isppol) 105 end do 106 end if 107 end do 108 109 end subroutine test_ftgkk