TABLE OF CONTENTS
ABINIT/overlap_k1k2_paw [ Functions ]
NAME
overlap_k1k2_paw
FUNCTION
compute PAW overlap between two k points, similar to smatrix_k_paw.F90 but more generic
COPYRIGHT
Copyright (C) 2005-2017 ABINIT group 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
cprj_k1 (pawcprj_type) :: cprj for occupied bands at point k1 cprj_k2 :: cprj for occupied bands at point k2 dk(3) :: vector k2 - k1 gprimd(3,3)=dimensioned primitive translations of reciprocal lattice lmn2max :: lmnmax*(lmnmax+1)/2 lmnsize(ntypat) :: lmnsize for each atom type mband :: number of bands natom=number of atoms in unit cell nspinor :: number of spinors (1 or 2) ntypat=number of types of atoms in unit cell pawang <type(pawang_type)>=paw angular mesh and related data pawrad(ntypat) <type(pawrad_type)>=paw radial mesh and related data pawtab(ntypat) <type(pawtab_type)>=paw tabulated starting data typat=typat(natom) list of atom types xred(natom,3) :: locations of atoms in cell
OUTPUT
k1k2_paw(2,mband,mband) :: array of the on-site PAW parts of the overlaps between Bloch states at points k1 and k2, for the various pairs of bands, that is, the on-site part of <u_nk1|u_mk2>
SIDE EFFECTS
NOTES
This routine assumes that the cprj are not explicitly ordered by atom type.
PARENTS
chern_number
CHILDREN
expibi,qijb_kk
SOURCE
52 #if defined HAVE_CONFIG_H 53 #include "config.h" 54 #endif 55 56 #include "abi_common.h" 57 58 subroutine overlap_k1k2_paw(cprj_k1,cprj_k2,dk,gprimd,k1k2_paw,lmn2max,lmnsize,mband,& 59 & natom,nspinor,ntypat,pawang,pawrad,pawtab,typat,xred) 60 61 use m_profiling_abi 62 63 use defs_basis 64 use m_errors 65 use m_pawang, only : pawang_type 66 use m_pawcprj, only : pawcprj_type 67 use m_pawrad, only : pawrad_type 68 use m_pawtab, only : pawtab_type 69 70 !This section has been created automatically by the script Abilint (TD). 71 !Do not modify the following lines by hand. 72 #undef ABI_FUNC 73 #define ABI_FUNC 'overlap_k1k2_paw' 74 use interfaces_65_paw, except_this_one => overlap_k1k2_paw 75 !End of the abilint section 76 77 implicit none 78 79 !Arguments--------------------------- 80 !scalars 81 integer,intent(in) :: lmn2max,mband,natom,nspinor,ntypat 82 type(pawang_type),intent(in) :: pawang 83 type(pawcprj_type),intent(in) :: cprj_k1(natom,mband),cprj_k2(natom,mband) 84 85 !arrays 86 integer,intent(in) :: lmnsize(ntypat),typat(natom) 87 real(dp),intent(in) :: dk(3),gprimd(3,3),xred(natom,3) 88 real(dp),intent(out) :: k1k2_paw(2,mband,mband) 89 type(pawrad_type),intent(in) :: pawrad(ntypat) 90 type(pawtab_type),intent(in) :: pawtab(ntypat) 91 92 !Local variables--------------------------- 93 !scalars 94 integer :: iatom,iband,ibs,ilmn,ispinor,itypat 95 integer :: jband,jbs,jlmn,klmn 96 complex(dpc) :: cpk1,cpk2,cterm,paw_onsite 97 98 ! arrays 99 real(dp),allocatable :: calc_expibi(:,:),calc_qijb(:,:,:) 100 101 ! ************************************************************************* 102 103 !initialize k1k2_paw output variable 104 k1k2_paw(:,:,:) = zero 105 106 ! obtain the atomic phase factors for the input k vector shift 107 ABI_ALLOCATE(calc_expibi,(2,natom)) 108 call expibi(calc_expibi,dk,natom,xred) 109 110 ! obtain the onsite PAW terms for the input k vector shift 111 ABI_ALLOCATE(calc_qijb,(2,lmn2max,natom)) 112 call qijb_kk(calc_qijb,dk,calc_expibi,gprimd,lmn2max,natom,ntypat,pawang,pawrad,pawtab,typat) 113 ABI_DEALLOCATE(calc_expibi) 114 115 do iatom = 1, natom 116 itypat = typat(iatom) 117 118 do ilmn=1,lmnsize(itypat) 119 do jlmn=1,lmnsize(itypat) 120 klmn=max(ilmn,jlmn)*(max(ilmn,jlmn)-1)/2 + min(ilmn,jlmn) 121 paw_onsite = cmplx(calc_qijb(1,klmn,iatom),calc_qijb(2,klmn,iatom)) 122 do iband = 1, mband 123 do jband = 1, mband 124 do ispinor = 1, nspinor 125 ibs = nspinor*(iband-1) + ispinor 126 jbs = nspinor*(jband-1) + ispinor 127 cpk1=cmplx(cprj_k1(iatom,ibs)%cp(1,ilmn),cprj_k1(iatom,ibs)%cp(2,ilmn)) 128 cpk2=cmplx(cprj_k2(iatom,jbs)%cp(1,jlmn),cprj_k2(iatom,jbs)%cp(2,jlmn)) 129 cterm = conjg(cpk1)*paw_onsite*cpk2 130 k1k2_paw(1,iband,jband) = k1k2_paw(1,iband,jband)+real(cterm) 131 k1k2_paw(2,iband,jband) = k1k2_paw(2,iband,jband)+aimag(cterm) 132 end do ! end loop over ispinor 133 end do ! end loop over jband 134 end do ! end loop over iband 135 end do ! end loop over ilmn 136 end do ! end loop over jlmn 137 138 end do ! end loop over atoms 139 140 ABI_DEALLOCATE(calc_qijb) 141 142 end subroutine overlap_k1k2_paw