TABLE OF CONTENTS
ABINIT/mkqptequiv [ Functions ]
NAME
mkqptequiv
FUNCTION
This routine determines the equivalence between 1) qpoints and fermi surface kpoints 2) qpoints under symmetry operations
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
Cryst<crystal_t>=Info on unit cell and symmetries. kpt_phon = fermi surface kpoints nkpt_phon = number of kpoints in the full FS set nqpt = number of qpoints qpt_full = qpoint coordinates
OUTPUT
FSfullpqtofull = mapping of k + q onto k' for k and k' in full BZ qpttoqpt(itim,isym,iqpt) = qpoint index which transforms to iqpt under isym and with time reversal itim.
NOTES
REMOVED 3/6/2008: much too large matrix, and not used at present FStoqpt = mapping of kpoint pairs (1 irreducible and 1 full) to qpoints
PARENTS
elphon,get_tau_k
CHILDREN
destroy_kptrank,get_rank_1kpt,mkkptrank,wrtout
SOURCE
42 #if defined HAVE_CONFIG_H 43 #include "config.h" 44 #endif 45 46 #include "abi_common.h" 47 48 subroutine mkqptequiv(FSfullpqtofull,Cryst,kpt_phon,nkpt_phon,nqpt,qpttoqpt,qpt_full,mqtofull) 49 50 use defs_basis 51 use defs_elphon 52 use m_kptrank 53 use m_errors 54 use m_profiling_abi 55 56 use m_crystal, only : crystal_t 57 58 !This section has been created automatically by the script Abilint (TD). 59 !Do not modify the following lines by hand. 60 #undef ABI_FUNC 61 #define ABI_FUNC 'mkqptequiv' 62 use interfaces_14_hidewrite 63 !End of the abilint section 64 65 implicit none 66 67 !Arguments ------------------------------------ 68 !scalars 69 integer,intent(in) :: nkpt_phon,nqpt 70 type(crystal_t),intent(in) :: Cryst 71 !arrays 72 integer,intent(out) :: FSfullpqtofull(nkpt_phon,nqpt),qpttoqpt(2,Cryst%nsym,nqpt) 73 integer,intent(out),optional :: mqtofull(nqpt) 74 real(dp),intent(in) :: kpt_phon(3,nkpt_phon),qpt_full(3,nqpt) 75 76 !Local variables------------------------------- 77 !scalars 78 integer :: ikpt_phon,iFSqpt,iqpt,isym,symrankkpt_phon 79 character(len=500) :: message 80 type(kptrank_type) :: kptrank_t 81 !arrays 82 real(dp) :: tmpkpt(3),gamma_kpt(3) 83 84 ! ************************************************************************* 85 86 call wrtout(std_out,' mkqptequiv : making rankkpt_phon and invrankkpt_phon',"COLL") 87 88 call mkkptrank (kpt_phon,nkpt_phon,kptrank_t) 89 90 FSfullpqtofull = -999 91 gamma_kpt(:) = zero 92 93 do ikpt_phon=1,nkpt_phon 94 do iqpt=1,nqpt 95 ! tmpkpt = jkpt = ikpt + qpt 96 tmpkpt(:) = kpt_phon(:,ikpt_phon) + qpt_full(:,iqpt) 97 98 ! which kpt is it among the full FS kpts? 99 call get_rank_1kpt (tmpkpt,symrankkpt_phon,kptrank_t) 100 101 FSfullpqtofull(ikpt_phon,iqpt) = kptrank_t%invrank(symrankkpt_phon) 102 if (FSfullpqtofull(ikpt_phon,iqpt) == -1) then 103 MSG_ERROR("looks like no kpoint equiv to k+q !!!") 104 end if 105 106 end do 107 end do 108 109 if (present(mqtofull)) then 110 do iqpt=1,nqpt 111 tmpkpt(:) = gamma_kpt(:) - qpt_full(:,iqpt) 112 113 ! which kpt is it among the full FS kpts? 114 call get_rank_1kpt (tmpkpt,symrankkpt_phon,kptrank_t) 115 116 mqtofull(iqpt) = kptrank_t%invrank(symrankkpt_phon) 117 if (mqtofull(iqpt) == -1) then 118 MSG_ERROR("looks like no kpoint equiv to -q !!!") 119 end if 120 end do 121 end if 122 123 call destroy_kptrank (kptrank_t) 124 125 !start over with q grid 126 call wrtout(std_out,' mkqptequiv : FSfullpqtofull made. Do qpttoqpt',"COLL") 127 128 call mkkptrank (qpt_full,nqpt,kptrank_t) 129 130 qpttoqpt(:,:,:) = -1 131 do iFSqpt=1,nqpt 132 do isym=1,Cryst%nsym 133 tmpkpt(:) = Cryst%symrec(:,1,isym)*qpt_full(1,iFSqpt) & 134 & + Cryst%symrec(:,2,isym)*qpt_full(2,iFSqpt) & 135 & + Cryst%symrec(:,3,isym)*qpt_full(3,iFSqpt) 136 137 call get_rank_1kpt (tmpkpt,symrankkpt_phon,kptrank_t) 138 if (kptrank_t%invrank(symrankkpt_phon) == -1) then 139 message = "looks like no kpoint equiv to q by symmetry without time reversal!!!" 140 MSG_ERROR(message) 141 end if 142 qpttoqpt(1,isym,kptrank_t%invrank(symrankkpt_phon)) = iFSqpt 143 144 tmpkpt = -tmpkpt 145 call get_rank_1kpt (tmpkpt,symrankkpt_phon,kptrank_t) 146 if (kptrank_t%invrank(symrankkpt_phon) == -1) then 147 message = ' mkqptequiv : Error : looks like no kpoint equiv to q by symmetry with time reversal!!!' 148 MSG_ERROR(message) 149 end if 150 qpttoqpt(2,isym,kptrank_t%invrank(symrankkpt_phon)) = iFSqpt 151 end do 152 end do 153 154 call destroy_kptrank (kptrank_t) 155 156 end subroutine mkqptequiv