TABLE OF CONTENTS


ABINIT/mkpwind_k [ Functions ]

[ Top ] [ Functions ]

NAME

 mkpwind_k

FUNCTION

 Make plane wave index at k point for basis at second k point,
 needed to compute overlaps $\langle u_{k,n}|u_{k+b,n}\rangle$
 as appear in Berry phase derived quantities

COPYRIGHT

 Copyright (C) 2003-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 .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt.

INPUTS

 dk(3)=real vector difference of ket kpt - bra kpt
 dtset <type(dataset_type)>=all input variables in this dataset
 fnkpt=number of kpts in full BZ
 fkptns=kpts in full BZ
 gmet(3,3)=metric in reciprocal space
 indkk_f2ibz(fnkpt,6)=information on folding from FBZ to IBZ (see initberry or initorbmag)
 ikpt=index of bra k pt
 ikpt1=index of neighbour ket k pt
 kg(3,dtset%mpw*dtset%mkmem)=planewave basis data
 kgindex(dtset%nkpt)= index of kg per kpt
 mpi_enreg=information about MPI parallelization
 npw_k=number of planewaves at k
 symrec(3,3,nsym) = symmetries in reciprocal space in terms of
   reciprocal space primitive translations

OUTPUT

 pwind_k1(dtset%mpw)=output index of ikpt1 basis states refered to ikpt

SIDE EFFECTS

TODO

NOTES

PARENTS

      chern_number

CHILDREN

      kpgsph

SOURCE

 51 #if defined HAVE_CONFIG_H
 52 #include "config.h"
 53 #endif
 54 
 55 #include "abi_common.h"
 56 
 57 subroutine mkpwind_k(dk,dtset,fnkpt,fkptns,gmet,indkk_f2ibz,ikpt,ikpt1,&
 58 & kg,kgindex,mpi_enreg,npw_k,pwind_k1,symrec)
 59 
 60  use defs_basis
 61  use defs_abitypes
 62  use defs_datatypes
 63  use m_xmpi
 64  use m_errors
 65  use m_profiling_abi
 66 
 67  use m_fftcore, only : kpgsph
 68 
 69 !This section has been created automatically by the script Abilint (TD).
 70 !Do not modify the following lines by hand.
 71 #undef ABI_FUNC
 72 #define ABI_FUNC 'mkpwind_k'
 73 !End of the abilint section
 74 
 75  implicit none
 76 
 77 !Arguments ------------------------------------
 78 !scalars
 79  integer,intent(in) :: fnkpt,ikpt,ikpt1,npw_k
 80  type(dataset_type),intent(in) :: dtset
 81  type(MPI_type), intent(inout) :: mpi_enreg
 82 
 83 !arrays
 84  integer,intent(in) :: indkk_f2ibz(fnkpt,6),kg(3,dtset%mpw*dtset%mkmem),kgindex(dtset%nkpt)
 85  integer,intent(in) :: symrec(3,3,dtset%nsym)
 86  integer,intent(out) :: pwind_k1(dtset%mpw)
 87  real(dp),intent(in) :: dk(3),fkptns(3,fnkpt),gmet(3,3)
 88  
 89 !Local variables -------------------------
 90 !scalars
 91  integer :: exchn2n3d,idum1,ikg1,ipw,istwf_k,isym,isym1,jpw,npw_k1
 92  real(dp) :: ecut_eff
 93  
 94 !arrays
 95  integer,allocatable :: kg1_k(:,:)
 96  real(dp) :: dg(3),dum33(3,3),kpt1(3),iadum(3),iadum1(3)
 97 
 98 ! ***********************************************************************
 99 
100  ABI_ALLOCATE(kg1_k,(3,dtset%mpw))
101  
102  ecut_eff = dtset%ecut*(dtset%dilatmx)**2
103  exchn2n3d = 0 ; istwf_k = 1 ; ikg1 = 0
104 
105  ! Build basis sphere of plane waves for the nearest neighbour of the k-point 
106 
107  kg1_k(:,:) = 0
108  kpt1(:) = dtset%kptns(:,ikpt1)
109  call kpgsph(ecut_eff,exchn2n3d,gmet,ikg1,ikpt1,istwf_k,kg1_k,kpt1,1,mpi_enreg,dtset%mpw,npw_k1)
110 
111  !        
112  !        Deal with symmetry transformations
113  !        
114 
115  !        bra k-point k(b) and IBZ k-point kIBZ(b) related by
116  !        k(b) = alpha(b) S(b)^t kIBZ(b) + G(b)
117  !        where alpha(b), S(b) and G(b) are given by indkk_f2ibz
118  !        
119  !        For the ket k-point:
120  !        k(k) = alpha(k) S(k)^t kIBZ(k) + G(k) - GBZ(k)
121  !        where GBZ(k) takes k(k) to the BZ
122  !        
123  
124  isym  = indkk_f2ibz(ikpt,2)
125  isym1 = indkk_f2ibz(ikpt1,2)
126 
127  !        Construct transformed G vector that enters the matching condition:
128  !        alpha(k) S(k)^{t,-1} ( -G(b) - GBZ(k) + G(k) )
129 
130  dg(:) = -indkk_f2ibz(ikpt,3:5) &
131 & -nint(-fkptns(:,ikpt) - dk(:) - tol10 + &
132 & fkptns(:,ikpt1)) &
133 & +indkk_f2ibz(ikpt1,3:5)
134 
135  iadum(:) = MATMUL(TRANSPOSE(dtset%symrel(:,:,isym1)),dg(:))
136 
137  dg(:) = iadum(:)
138  
139  !        Construct S(k)^{t,-1} S(b)^{t}
140 
141  dum33(:,:) = MATMUL(TRANSPOSE(dtset%symrel(:,:,isym1)),symrec(:,:,isym))
142 
143  !        Construct alpha(k) alpha(b)
144  
145  pwind_k1(:) = 0
146  do ipw = 1, npw_k
147 
148     !          NOTE: the bra G vector is taken for the sym-related IBZ k point,
149     !          not for the FBZ k point
150    iadum(:) = kg(:,kgindex(ikpt) + ipw)
151    
152     !          to determine r.l.v. matchings, we transformed the bra vector
153     !          Rotation
154    iadum1(:)=0
155    do idum1=1,3
156      iadum1(:)=iadum1(:)+dum33(:,idum1)*iadum(idum1)
157    end do
158    iadum(:)=iadum1(:)
159    iadum(:) = iadum(:) + dg(:)
160 
161    do jpw = 1, npw_k1
162      iadum1(1:3) = kg1_k(1:3,jpw)
163      if ( (iadum(1) == iadum1(1)).and. &
164 &     (iadum(2) == iadum1(2)).and. &
165 &     (iadum(3) == iadum1(3)) ) then
166        pwind_k1(ipw) = jpw
167           ! write(std_out,'(a,2i4)')'JWZ debug : bg ipw == jpw ',ipw,jpw
168        exit
169      end if
170    end do
171  end do
172 
173  ABI_DEALLOCATE(kg1_k)
174 
175 end subroutine mkpwind_k