TABLE OF CONTENTS
ABINIT/dsdr_k_paw [ Functions ]
NAME
dsdr_k_paw
FUNCTION
compute on-site terms for forces and stresses for finite electric fields with PAW
COPYRIGHT
Copyright (C) 2005-2018 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_k (pawcprj_type) :: cprj for occupied bands at point k cprj_kb :: cprj for occupied bands at point k+b dtefield :: structure referring to all efield and berry's phase variables kdir :: integer giving direction along which overlap is computed for ket kfor :: integer indicating whether to compute forward (1) or backward (2) along kpt string natom :: number of atoms in cell typat :: typat(natom) type of each atom
OUTPUT
SIDE EFFECTS
dsdr :: array of the on-site PAW parts of the derivatives with respect to atm positions and/or strains of the overlaps between Bloch states at points k and k+b, for the various pairs of bands
NOTES
This routine assumes that the cprj are not explicitly ordered by atom type.
PARENTS
berryphase_new
CHILDREN
SOURCE
43 #if defined HAVE_CONFIG_H 44 #include "config.h" 45 #endif 46 47 #include "abi_common.h" 48 49 subroutine dsdr_k_paw(cprj_k,cprj_kb,dsdr,dtefield,kdir,kfor,mband,natom,ncpgr,typat) 50 51 use m_profiling_abi 52 53 use defs_basis 54 use m_errors 55 use m_efield 56 use m_pawcprj, only : pawcprj_type 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 'dsdr_k_paw' 62 !End of the abilint section 63 64 implicit none 65 66 !Arguments--------------------------- 67 !scalars 68 integer,intent(in) :: kdir,kfor,mband,natom,ncpgr 69 character(len=500) :: message 70 type(efield_type),intent(in) :: dtefield 71 type(pawcprj_type),intent(in) :: cprj_k(natom,dtefield%nspinor*mband) 72 type(pawcprj_type),intent(in) :: cprj_kb(natom,dtefield%nspinor*mband) 73 74 !arrays 75 integer,intent(in) :: typat(natom) 76 real(dp),intent(inout) :: dsdr(2,natom,ncpgr,dtefield%mband_occ,dtefield%mband_occ) 77 78 !Local variables--------------------------- 79 !scalars 80 integer :: iatom,iband,ibs,icpgr,ilmn,ispinor,itypat 81 integer :: jband,jbs,jlmn,klmn,nspinor 82 complex(dpc) :: cpk,cpkb,dcpk,dcpkb,cterm,paw_onsite 83 84 ! ************************************************************************* 85 86 !initialize dsdr 87 dsdr(:,:,:,:,:) = zero 88 89 ! if 3 gradients we are in the ctocprj choice 2 case 90 ! and the 3 gradients are due to the atomic displacements 91 ! if 6 gradients we are in the ctocprj choice 3 case 92 ! and the 6 gradients are due to the strains 93 ! if 9 gradients we are in the ctocprj choice 23 case 94 ! and the first six are due to strain, last three due to displacements 95 if (ncpgr /= 3 .and. ncpgr /= 6 .and. ncpgr /= 9) then 96 message = ' dsdr_k_paw called with ncpgr /= 3, 6, or 9 (no gradients) ' 97 MSG_BUG(message) 98 end if 99 100 nspinor = dtefield%nspinor 101 102 do iatom = 1, natom 103 itypat = typat(iatom) 104 105 do ilmn=1,dtefield%lmn_size(itypat) 106 do jlmn=1,dtefield%lmn_size(itypat) 107 klmn=max(ilmn,jlmn)*(max(ilmn,jlmn)-1)/2 + min(ilmn,jlmn) 108 paw_onsite = cmplx(dtefield%qijb_kk(1,klmn,iatom,kdir),& 109 & dtefield%qijb_kk(2,klmn,iatom,kdir)) 110 if (kfor > 1) paw_onsite = conjg(paw_onsite) 111 do iband = 1, dtefield%mband_occ 112 do jband = 1, dtefield%mband_occ 113 do ispinor = 1, nspinor 114 do icpgr = 1, ncpgr 115 ibs = nspinor*(iband-1) + ispinor 116 jbs = nspinor*(jband-1) + ispinor 117 cpk=cmplx(cprj_k(iatom,ibs)%cp(1,ilmn),cprj_k(iatom,ibs)%cp(2,ilmn)) 118 dcpk=cmplx(cprj_k(iatom,ibs)%dcp(1,icpgr,ilmn),cprj_k(iatom,ibs)%dcp(2,icpgr,ilmn)) 119 cpkb=cmplx(cprj_kb(iatom,jbs)%cp(1,jlmn),cprj_kb(iatom,jbs)%cp(2,jlmn)) 120 dcpkb=cmplx(cprj_kb(iatom,jbs)%dcp(1,icpgr,jlmn),cprj_kb(iatom,jbs)%dcp(2,icpgr,jlmn)) 121 cterm=paw_onsite*(conjg(dcpk)*cpkb+conjg(cpk)*dcpkb) 122 dsdr(1,iatom,icpgr,iband,jband) = dsdr(1,iatom,icpgr,iband,jband)+real(cterm) 123 dsdr(2,iatom,icpgr,iband,jband) = dsdr(2,iatom,icpgr,iband,jband)+aimag(cterm) 124 end do ! end loop over icpgr 125 end do ! end loop over ispinor 126 end do ! end loop over jband 127 end do ! end loop over iband 128 end do ! end loop over ilmn 129 end do ! end loop over jlmn 130 131 end do ! end loop over atoms 132 133 end subroutine dsdr_k_paw