TABLE OF CONTENTS


ABINIT/complete_gamma [ Functions ]

[ Top ] [ Functions ]

NAME

 complete_gamma

FUNCTION

 Use the set of special q points calculated by the Monkhorst & Pack Technique.
 Check if all the informations for the q points are present in the input gamma matrices.
 Generate the gamma matrices (already summed over the FS) of the set of q points which
 samples homogeneously the entire Brillouin zone.

COPYRIGHT

 Copyright (C) 2009-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

 qpttoqpt = qpoint index mapping under symops

OUTPUT

 gamma_qpt = in/out: set of gamma matrix elements completed and symmetrized 
    gamma_qpt(2,nbranch**2,nsppol,nqpt_full)

PARENTS

      elphon,integrate_gamma_alt,m_phgamma

CHILDREN

      zgemm

SOURCE

 35 #if defined HAVE_CONFIG_H
 36 #include "config.h"
 37 #endif
 38 
 39 #include "abi_common.h"
 40 
 41 subroutine complete_gamma(Cryst,nbranch,nsppol,nqptirred,nqpt_full,ep_scalprod,qirredtofull,qpttoqpt,gamma_qpt)
 42 
 43  use defs_basis
 44  use defs_elphon
 45  use m_profiling_abi
 46  use m_errors
 47 
 48  use m_crystal,    only : crystal_t
 49 
 50 !This section has been created automatically by the script Abilint (TD).
 51 !Do not modify the following lines by hand.
 52 #undef ABI_FUNC
 53 #define ABI_FUNC 'complete_gamma'
 54 !End of the abilint section
 55 
 56  implicit none
 57 
 58 !Arguments ------------------------------------
 59 !scalars
 60  integer,intent(in) :: nsppol,nbranch,nqptirred,nqpt_full,ep_scalprod
 61  type(crystal_t),intent(in) :: Cryst
 62 !arrays
 63  integer,intent(in) :: qirredtofull(nqptirred)
 64  integer,intent(in) :: qpttoqpt(2,Cryst%nsym,nqpt_full)
 65  real(dp), intent(inout) :: gamma_qpt(2,nbranch**2,nsppol,nqpt_full)
 66 
 67 !Local variables-------------------------------
 68 !scalars
 69  integer :: ibranch,ieqqpt,ii,natom,nsym,iqpt,isppol,isym
 70  integer :: itim,jbranch,jj,kk,ll,neqqpt,iatom,ancestor_iatom,iqpt_fullbz
 71 !arrays
 72  integer :: symmetrized_qpt(nqpt_full)
 73  integer :: gkk_flag(nbranch,nbranch,nsppol,nqpt_full)
 74  real(dp) :: ss(3,3)
 75  real(dp) :: tmp_mat(2,nbranch,nbranch)
 76  real(dp) :: tmp_mat2(2,nbranch,nbranch)
 77  real(dp) :: ss_allatoms(2,nbranch,nbranch)
 78  real(dp) :: c_one(2), c_zero(2)
 79  real(dp),allocatable :: gkk_qpt_new(:,:,:),gkk_qpt_tmp(:,:,:)
 80 
 81 ! *********************************************************************
 82 
 83  c_one = (/one,zero/)
 84  c_zero = (/zero,zero/)
 85 
 86  natom = Cryst%natom
 87  nsym  = Cryst%nsym
 88 
 89 !Generation of the gkk matrices relative to the q points
 90 !of the set which samples the entire Brillouin zone
 91 
 92 !set up flags for gamma_qpt matrices we have
 93  gkk_flag = -1
 94  do iqpt=1,nqptirred
 95    iqpt_fullbz = qirredtofull(iqpt)
 96    gkk_flag(:,:,:,iqpt_fullbz) = 1
 97  end do
 98 
 99  symmetrized_qpt(:) = -1
100 
101  ABI_ALLOCATE(gkk_qpt_new,(2,nbranch**2,nsppol))
102  ABI_ALLOCATE(gkk_qpt_tmp,(2,nbranch**2,nsppol))
103 
104  do iqpt=1,nqpt_full
105 !  
106 !  Already symmetrized?
107    if (symmetrized_qpt(iqpt) == 1) cycle
108 
109    gkk_qpt_new(:,:,:) = zero
110 
111 !  loop over qpoints equivalent to iqpt
112    neqqpt=0
113 !  do not use time reversal symmetry to complete the qpoints:
114 !  do not know what happens to the gamma matrices
115 !  11/2011: MJV: time reversal is needed here if inversion is absent
116 !  - used in read_gkk and all reductions of q-points by symmetry.
117 
118    do itim=1,2
119      do isym=1,nsym
120 !      ieqqpt is sent onto iqpt by itim/isym
121        ieqqpt = qpttoqpt(itim,isym,iqpt)
122 
123 
124        if (gkk_flag(1,1,1,ieqqpt) == -1) cycle
125 !      if we have information on this qpt
126 !      iqpt is equivalent to ieqqpt: get it from file or memory
127        gkk_qpt_tmp(:,:,:) = gamma_qpt(:,:,:,ieqqpt)
128 
129        neqqpt=neqqpt+1
130 
131 !      
132 !      MJV note 02/2010:
133 !      the correspondence of symrel and symrec in the different cases, symmetrizing there
134 !      and back, has been fixed in the cases with and without scalprod (ie cartesian
135 !      and reduced real space coordinates) with respect to a calculation with no symmetries
136 !      I believe everything is settled, but still do not know why the 2 versions of the ss
137 !      matrices here use different rel/rec, instead of just being multiplied by the rprim gprim...
138 !      
139        if (ep_scalprod==1) then
140          do ii=1,3
141            do jj=1,3
142              ss(ii,jj)=zero
143              do kk=1,3
144                do ll=1,3
145                  ss(ii,jj)=ss(ii,jj)+Cryst%rprimd(ii,kk)*Cryst%symrel(kk,ll,isym)*Cryst%gprimd(ll,jj)
146                end do
147              end do
148            end do
149          end do
150        else
151          do ii=1,3
152            do jj=1,3
153              ss(ii,jj) = Cryst%symrec(ii,jj,isym)
154            end do
155          end do
156        end if
157 
158        ss_allatoms(:,:,:) = zero
159        do iatom=1,natom
160          ancestor_iatom = Cryst%indsym(4,isym,iatom)
161          ss_allatoms(1, (ancestor_iatom-1)*3+1:(ancestor_iatom-1)*3+3, (iatom-1)*3+1:(iatom-1)*3+3) = ss(1:3,1:3)
162        end do
163 
164 
165 !      NOTE   ssinv(ii,jj)=ssinv(ii,jj)+Cryst%gprimd(ii,kk)*rprimd(jj,ll)*Cryst%symrec(ll,kk,isym)
166 
167        do isppol=1,nsppol
168 !        multiply by the ss matrices
169          tmp_mat2(:,:,:) = zero
170          tmp_mat(:,:,:) = reshape(gkk_qpt_tmp(:,:,isppol),(/2,nbranch,nbranch/))
171 
172          call ZGEMM ('N','N',nbranch,nbranch,nbranch,&
173 &         c_one,ss_allatoms,nbranch,tmp_mat,nbranch,c_zero,tmp_mat2,nbranch)
174 
175          call ZGEMM ('N','T',nbranch,nbranch,nbranch,&
176 &         c_one,tmp_mat2,nbranch,ss_allatoms,nbranch,c_zero,tmp_mat,nbranch)
177 
178 !        add to gkk_qpt_new
179          do ibranch =1,nbranch
180            do jbranch =1,nbranch
181              gkk_qpt_new(:,(jbranch-1)*nbranch+ibranch,isppol) = &
182 &             gkk_qpt_new(:,(jbranch-1)*nbranch+ibranch,isppol) + tmp_mat(:,jbranch,ibranch)
183            end do
184          end do
185        end do ! isppol
186 !      
187      end do ! isym
188    end do ! itim
189 !  
190    ABI_CHECK(neqqpt>0,'no q-points found equivalent to iqpt ')
191 !  Divide by number of equivalent qpts found.
192    gkk_qpt_new(:,:,:) = gkk_qpt_new(:,:,:)/neqqpt
193 
194 !  copy the symmetrized version into all the equivalent qpoints, appropriately transformed
195    do itim=1,2
196      do isym=1,nsym
197 !      ieqqpt is sent onto iqpt by itim/isym
198        ieqqpt = qpttoqpt(itim,isym,iqpt)
199 
200        if (symmetrized_qpt(ieqqpt) /= -1) cycle
201        gkk_qpt_tmp(:,:,:) = zero
202 
203 !      use symrec matrices to get inverse transform from isym^{-1}
204        if (ep_scalprod==1) then
205          do ii=1,3
206            do jj=1,3
207              ss(ii,jj)=zero
208              do kk=1,3
209                do ll=1,3
210 !                Use inverse of symop matrix here to get back to ieqqpt (inv+transpose is in symrec and in gprimd)
211                  ss(ii,jj)=ss(ii,jj)+Cryst%rprimd(ii,kk)*Cryst%symrec(ll,kk,isym)*Cryst%gprimd(ll,jj)
212                end do
213              end do
214            end do
215          end do
216        else
217          do ii=1,3
218            do jj=1,3
219              ss(ii,jj) = Cryst%symrel(jj,ii,isym)
220            end do
221          end do
222        end if
223 
224        ss_allatoms(:,:,:) = zero
225        do iatom=1,natom
226          ancestor_iatom = Cryst%indsym(4,isym,iatom)
227          ss_allatoms(1, (ancestor_iatom-1)*3+1:(ancestor_iatom-1)*3+3, (iatom-1)*3+1:(iatom-1)*3+3) = ss(1:3,1:3)
228        end do
229 
230 !      ! Use inverse of symop matrix here to get back to ieqqpt
231 !      ssinv(ii,jj)=ssinv(ii,jj)+gprimd(ii,kk)*rprimd(jj,ll)*Cryst%symrel(kk,ll,isym)
232 
233        do isppol=1,nsppol
234 !        multiply by the ss^{-1} matrices
235          tmp_mat2(:,:,:) = zero
236          tmp_mat(:,:,:) = reshape(gkk_qpt_new(:,:,isppol),(/2,nbranch,nbranch/))
237 
238          call ZGEMM ('N','N',nbranch,nbranch,nbranch,&
239 &         c_one,ss_allatoms,nbranch,tmp_mat,nbranch,c_zero,tmp_mat2,nbranch)
240 
241          call ZGEMM ('N','T',nbranch,nbranch,nbranch,&
242 &         c_one,tmp_mat2,nbranch,ss_allatoms,nbranch,c_zero,tmp_mat,nbranch)
243 
244 !        FIXME: the following could just be a reshape
245          do ibranch =1,nbranch
246            do jbranch =1,nbranch
247              gkk_qpt_tmp(:,(jbranch-1)*nbranch+ibranch,isppol) =&
248 &             tmp_mat(:,jbranch,ibranch)
249            end do
250          end do
251          if (gkk_flag (1,1,isppol,ieqqpt) == -1) gkk_flag (:,:,isppol,ieqqpt) = 0
252        end do ! end isppol do
253 
254 !      save symmetrized matrices for qpt ieqqpt
255        gamma_qpt(:,:,:,ieqqpt) = gkk_qpt_tmp(:,:,:)
256 
257        symmetrized_qpt(ieqqpt) = 1
258 
259      end do !isym
260    end do !itim
261  end do !iqpt
262 
263  ABI_DEALLOCATE(gkk_qpt_new)
264  ABI_DEALLOCATE(gkk_qpt_tmp)
265 
266 end subroutine complete_gamma