TABLE OF CONTENTS


ABINIT/smatrix_k_paw [ Functions ]

[ Top ] [ Functions ]

NAME

 smatrix_k_paw

FUNCTION

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

 smat_k_paw :: array of the on-site PAW parts of the overlaps between Bloch states at points
   k and k+b, for the various pairs of bands, that is, the on-site part of 
   <u_nk|u_mk+b>

SIDE EFFECTS

NOTES

 This routine assumes that the cprj are not explicitly ordered by 
 atom type.

PARENTS

      berryphase_new,cgwf,make_grad_berry

CHILDREN

SOURCE

 42 #if defined HAVE_CONFIG_H
 43 #include "config.h"
 44 #endif
 45 
 46 #include "abi_common.h"
 47 
 48  subroutine smatrix_k_paw(cprj_k,cprj_kb,dtefield,kdir,kfor,mband,natom,smat_k_paw,typat)
 49 
 50  use m_profiling_abi
 51 
 52  use defs_basis
 53  use m_errors
 54  use m_efield
 55  use m_pawcprj, only : pawcprj_type
 56 
 57 !This section has been created automatically by the script Abilint (TD).
 58 !Do not modify the following lines by hand.
 59 #undef ABI_FUNC
 60 #define ABI_FUNC 'smatrix_k_paw'
 61 !End of the abilint section
 62 
 63  implicit none
 64 
 65 !Arguments---------------------------
 66 !scalars
 67  integer,intent(in) :: kdir,kfor,mband,natom
 68  type(efield_type),intent(in) :: dtefield
 69  type(pawcprj_type),intent(in) :: cprj_k(natom,dtefield%nspinor*mband)
 70  type(pawcprj_type),intent(in) :: cprj_kb(natom,dtefield%nspinor*mband)
 71 
 72 !arrays
 73  integer,intent(in) :: typat(natom)
 74  real(dp),intent(out) :: smat_k_paw(2,dtefield%mband_occ,dtefield%mband_occ)
 75 
 76 !Local variables---------------------------
 77 !scalars
 78  integer :: iatom,iband,ibs,ilmn,ispinor,itypat
 79  integer :: jband,jbs,jlmn,klmn,nspinor
 80  complex(dpc) :: cpk,cpkb,cterm,paw_onsite
 81 
 82 ! *************************************************************************
 83 
 84 !initialize smat_k_paw
 85  smat_k_paw(:,:,:) = zero
 86 
 87  nspinor = dtefield%nspinor
 88 
 89  do iatom = 1, natom
 90    itypat = typat(iatom)
 91 
 92    do ilmn=1,dtefield%lmn_size(itypat)
 93      do jlmn=1,dtefield%lmn_size(itypat)
 94        klmn=max(ilmn,jlmn)*(max(ilmn,jlmn)-1)/2 + min(ilmn,jlmn)
 95        paw_onsite = cmplx(dtefield%qijb_kk(1,klmn,iatom,kdir),&
 96 &       dtefield%qijb_kk(2,klmn,iatom,kdir))
 97        if (kfor > 1) paw_onsite = conjg(paw_onsite)
 98        do iband = 1, dtefield%mband_occ
 99          do jband = 1, dtefield%mband_occ
100            do ispinor = 1, nspinor
101              ibs = nspinor*(iband-1) + ispinor
102              jbs = nspinor*(jband-1) + ispinor
103              cpk=cmplx(cprj_k(iatom,ibs)%cp(1,ilmn),cprj_k(iatom,ibs)%cp(2,ilmn))
104              cpkb=cmplx(cprj_kb(iatom,jbs)%cp(1,jlmn),cprj_kb(iatom,jbs)%cp(2,jlmn))
105              cterm = conjg(cpk)*paw_onsite*cpkb
106              smat_k_paw(1,iband,jband) = smat_k_paw(1,iband,jband)+dreal(cterm)
107              smat_k_paw(2,iband,jband) = smat_k_paw(2,iband,jband)+dimag(cterm)
108            end do ! end loop over ispinor
109          end do ! end loop over jband
110        end do ! end loop over iband
111      end do ! end loop over ilmn
112    end do ! end loop over jlmn
113 
114  end do ! end loop over atoms
115 
116  end subroutine    smatrix_k_paw