TABLE OF CONTENTS


ABINIT/nmsq_pure_gkk [ Functions ]

[ Top ] [ Functions ]

NAME

 nmsq_pure_gkk

FUNCTION

  Calculate gamma matrices for pure gkk case, ie when the
  scalar product with the displacement vector is done later
  Sum over bands is carried out later.

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

   displ_red = phonon displacement in reduced coordinates (used to calculate the ph linewidth)
   elph_ds = datastructure with gkk matrix elements
   FSfullpqtofull = mapping of k+q to k
   kpt_phon = coordinates of kpoints near to FS
   h1_mat_el_sq = matrix elements $<psi_{k+q,m} | H^{1} | psi_{k,n}>$ squared
   iqptirred = index of present qpoint

OUTPUT

   elph_ds%gkq filled
   accum_mat = matrix for accumulating FS average of gkk (gamma matrix -> linewidths)
   accum_mat2 = complex array whose real part contains the phonon linewidth

PARENTS

      normsq_gkq

CHILDREN

      gam_mult_displ

SOURCE

 40 #if defined HAVE_CONFIG_H
 41 #include "config.h"
 42 #endif
 43 
 44 #include "abi_common.h"
 45 
 46 
 47 subroutine nmsq_pure_gkk(accum_mat,accum_mat2,displ_red,elph_ds,FSfullpqtofull,&
 48 &   h1_mat_el_sq,iqptirred)
 49 
 50  use defs_basis
 51  use defs_elphon
 52  use m_profiling_abi
 53  use m_errors
 54 
 55 !This section has been created automatically by the script Abilint (TD).
 56 !Do not modify the following lines by hand.
 57 #undef ABI_FUNC
 58 #define ABI_FUNC 'nmsq_pure_gkk'
 59  use interfaces_77_ddb, except_this_one => nmsq_pure_gkk
 60 !End of the abilint section
 61 
 62  implicit none
 63 
 64 !Arguments ------------------------------------
 65 !scalars
 66  integer,intent(in) :: iqptirred
 67  type(elph_type),intent(inout) :: elph_ds
 68 !arrays
 69  integer,intent(in) :: FSfullpqtofull(elph_ds%k_phon%nkpt,elph_ds%nqpt_full)
 70  real(dp),intent(in) :: displ_red(2,elph_ds%nbranch,elph_ds%nbranch)
 71  real(dp),intent(inout) :: &
 72 & h1_mat_el_sq(2,elph_ds%nFSband*elph_ds%nFSband,elph_ds%nbranch*elph_ds%nbranch,elph_ds%k_phon%my_nkpt,elph_ds%nsppol)
 73  real(dp),intent(inout) :: accum_mat(2,elph_ds%nbranch,elph_ds%nbranch,elph_ds%nsppol)
 74  real(dp),intent(inout) :: accum_mat2(2,elph_ds%nbranch,elph_ds%nbranch,elph_ds%nsppol)
 75 
 76 !Local variables-------------------------------
 77 !scalars
 78  integer :: ikpt_phon,ikpt_phonq,ib1,ib2,ibeff,ipert1,isppol
 79  integer :: iqpt_fullbz
 80  integer :: ik_this_proc
 81  real(dp) :: sd1,sd2
 82  character(len=500) :: message
 83 !arrays
 84  real(dp) :: gkq_sum_bands(2,elph_ds%nbranch,elph_ds%nbranch)
 85  real(dp) :: tmp_mat2(2,elph_ds%nbranch,elph_ds%nbranch)
 86  real(dp) :: zgemm_tmp_mat(2,elph_ds%nbranch,elph_ds%nbranch)
 87 
 88 ! *************************************************************************
 89 
 90  if (elph_ds%ep_keepbands /= 1) then
 91    message = ' elph_ds%ep_keepbands should be 1 to keep bands!'
 92    MSG_ERROR(message)
 93  end if
 94 
 95  iqpt_fullbz = elph_ds%qirredtofull(iqptirred)
 96 
 97 !h1_mat_el_sq is already fine here - nothing to do
 98 
 99 
100 !MG20060603 NOTE:
101 !accum_mat and accum_mat2 are real, the imaginary part is used for debugging purpose
102 !accum_mat2 is used to store the phonon-linewidhts before interpolation
103 
104 !MJV 20070525 NOTE:
105 !in some of the nmsq routines, in particular this one, the work done to
106 !calculate accum_mat,accum_mat2 is completely superfluous and will be re-done
107 !on the interpolated values.
108 !MG uses them for the QPT output, however, so keep it for consistency for the
109 !moment.
110 
111  do isppol=1,elph_ds%nsppol
112    do ik_this_proc =1, elph_ds%k_phon%my_nkpt
113      ikpt_phon = elph_ds%k_phon%my_ikpt(ik_this_proc)
114 
115      ikpt_phonq = FSfullpqtofull(ikpt_phon,iqpt_fullbz)
116 
117      gkq_sum_bands(:,:,:) = zero
118 
119 !    gkq_sum_bands = \sum_{ib1,ib2} \langle k+q \mid H^{(1)}_{q,\tau_i,\alpha_i} \mid k   \rangle
120 !    \cdot \langle k   \mid H^{(1)}_{q,\tau_j,\alpha_j} \mid k+q \rangle
121 !    where ibranch -> \tau_i,\alpha_i  and  jbranch -> \tau_j,\alpha_j
122 
123      do ib1=1,elph_ds%nFSband
124 
125        sd1 = elph_ds%k_phon%wtk(ib1,ikpt_phon,isppol)      !  weights for distance from the fermi surface
126 
127        do ib2=1,elph_ds%nFSband
128 
129          sd2 = elph_ds%k_phon%wtk(ib2,ikpt_phonq,isppol)  !  weights for distance from the fermi surface
130          ibeff = ib2+(ib1-1)*elph_ds%nFSband
131 
132          gkq_sum_bands = gkq_sum_bands + &
133 &         sd1*sd2*reshape(h1_mat_el_sq(:,ibeff,:,ik_this_proc,isppol),(/2,elph_ds%nbranch,elph_ds%nbranch/))
134 
135        end do !ib2
136      end do !ib1
137 !    END loops over bands
138 
139 
140      accum_mat(:,:,:,isppol) = accum_mat(:,:,:,isppol) + gkq_sum_bands(:,:,:)
141    end do
142 !  END loop over kpt_phon
143 
144 !  MG20060603
145 !  do scalar product with the displ_red to calculate the ph lwdth before interpolation (stored in accum_mat2)
146 
147    zgemm_tmp_mat = accum_mat(:,:,:,isppol)
148 
149    call gam_mult_displ(elph_ds%nbranch, displ_red, zgemm_tmp_mat, tmp_mat2)
150 
151    do ipert1=1,elph_ds%nbranch
152      accum_mat2(1,ipert1,ipert1,isppol) = accum_mat2(1,ipert1,ipert1,isppol) + tmp_mat2(1,ipert1,ipert1)
153    end do
154 
155 !  ENDMG
156 
157  end do ! isppol
158 
159 end subroutine nmsq_pure_gkk