TABLE OF CONTENTS


ABINIT/get_kpt_fullbz [ Functions ]

[ Top ] [ Functions ]

NAME

 get_kpt_fullbz

FUNCTION

 create full grid of kpoints from kptrlatt and shiftk

COPYRIGHT

 Copyright (C) 2002-2018 ABINIT group (MVer,XG)
 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

  kptrlatt(3,3)=lattice vectors for full kpoint grid
  nkpt_fullbz=number of kpoints in full brillouin zone
  nshiftk=number of kpoint grid shifts
  shiftk(3,nshiftk)=kpoint shifts

OUTPUT

  kpt_fullbz(3,nkpt_fullbz)=kpoints in full brillouin zone

NOTES

PARENTS

      get_full_kgrid,invars2

CHILDREN

      mati3det,matr3inv,wrap2_pmhalf

SOURCE

 35 #if defined HAVE_CONFIG_H
 36 #include "config.h"
 37 #endif
 38 
 39 #include "abi_common.h"
 40 
 41 subroutine get_kpt_fullbz(kpt_fullbz,kptrlatt,nkpt_fullbz,nshiftk,shiftk)
 42 
 43  use defs_basis
 44  use m_profiling_abi
 45  use m_errors
 46 
 47  use m_numeric_tools,   only : wrap2_pmhalf
 48 
 49 !This section has been created automatically by the script Abilint (TD).
 50 !Do not modify the following lines by hand.
 51 #undef ABI_FUNC
 52 #define ABI_FUNC 'get_kpt_fullbz'
 53  use interfaces_32_util
 54 !End of the abilint section
 55 
 56  implicit none
 57 
 58 !Arguments ------------------------------------
 59 !scalars
 60  integer,intent(in) :: nkpt_fullbz,nshiftk
 61 !arrays
 62  integer,intent(in) :: kptrlatt(3,3)
 63  real(dp),intent(in) :: shiftk(3,nshiftk)
 64  real(dp),intent(out) :: kpt_fullbz(3,nkpt_fullbz)
 65 
 66 !Local variables-------------------------------
 67 !scalars
 68  integer, parameter :: max_number_of_prime=47,mshiftk=210
 69  integer :: det,ii,ikshft,iprim,jj,kk,nn
 70  character(len=500) :: message
 71 
 72 !arrays
 73  integer :: boundmax(3),boundmin(3),common_factor(3)
 74  integer, parameter :: prime_factor(max_number_of_prime)=(/2,3,5,7,9, 11,13,17,19,23,&
 75 &  29,31,37,41,43, 47,53,59,61,67,&
 76 &  71,73,79,83,89, 97,101,103,107,109,&
 77 &  113,127,131,137,139, 149,151,157,163,167,&
 78 &  173,179,181,191,193, 197,199/)
 79  real(dp) :: k1(3),k2(3),klatt(3,3),rlatt(3,3),shift(3),test_rlatt(3,3)
 80 
 81 ! *********************************************************************
 82 
 83 !Identify first factors that can be used to rescale the three kptrlatt vectors
 84 !Only test a large set of prime factors, though ...
 85  do jj=1,3
 86    common_factor(jj)=1
 87    rlatt(:,jj)=kptrlatt(:,jj)
 88    do iprim=1,max_number_of_prime
 89      test_rlatt(:,jj)=rlatt(:,jj)/dble(prime_factor(iprim))
 90 !    If one of the components is lower than 1 in absolute value, then it is not worth to continue the search.
 91      if(minval(abs(abs(test_rlatt(:,jj))-half))<half-tol8)exit
 92      do
 93        if(sum(abs(test_rlatt(:,jj)-nint(test_rlatt(:,jj)) ))<tol8)then
 94          common_factor(jj)=prime_factor(iprim)*common_factor(jj)
 95          rlatt(:,jj)=rlatt(:,jj)/dble(prime_factor(iprim))
 96          test_rlatt(:,jj)=test_rlatt(:,jj)/dble(prime_factor(iprim))
 97        else
 98          exit
 99        end if  
100      end do
101    end do
102  end do
103  call mati3det(kptrlatt,det)
104  det=det/(common_factor(1)*common_factor(2)*common_factor(3))
105 
106  rlatt(:,:)=kptrlatt(:,:)
107  call matr3inv(rlatt,klatt)
108 !Now, klatt contains the three primitive vectors of the k lattice,
109 !in reduced coordinates. One builds all k vectors that
110 !are contained in the first Brillouin zone, with coordinates
111 !in the interval [0,1[ . First generate boundaries of a big box.
112 !In order to generate all possible vectors in the reciprocal space, 
113 !one must consider all multiples of the primitive ones, until a vector with only integers is found.
114 !The maximum bound is the scale of the corresponding kptrlatt vector, times the determinant of kptrlatt. Also consider negative vectors.
115 !On this basis, compute the bounds.
116  do jj=1,3
117 !  To accomodate the shifts, boundmin starts from -1
118 !  Well, this is not a complete solution ...
119    boundmin(jj)=-1-common_factor(jj)*abs(det)
120    boundmax(jj)=common_factor(jj)*abs(det)
121  end do
122 
123  nn=1
124  do kk=boundmin(3),boundmax(3)
125    do jj=boundmin(2),boundmax(2)
126      do ii=boundmin(1),boundmax(1)
127        do ikshft=1,nshiftk
128 
129 !        Coordinates of the trial k point with respect to the k primitive lattice
130          k1(1)=ii+shiftk(1,ikshft)
131          k1(2)=jj+shiftk(2,ikshft)
132          k1(3)=kk+shiftk(3,ikshft)
133 
134 !        Reduced coordinates of the trial k point
135          k2(:)=k1(1)*klatt(:,1)+k1(2)*klatt(:,2)+k1(3)*klatt(:,3)
136 
137 !        Eliminate the point if outside [0,1[
138          if(k2(1)<-tol10)cycle ; if(k2(1)>one-tol10)cycle
139          if(k2(2)<-tol10)cycle ; if(k2(2)>one-tol10)cycle
140          if(k2(3)<-tol10)cycle ; if(k2(3)>one-tol10)cycle
141 
142 !        Wrap the trial values in the interval ]-1/2,1/2] .
143          call wrap2_pmhalf(k2(1),k1(1),shift(1))
144          call wrap2_pmhalf(k2(2),k1(2),shift(2))
145          call wrap2_pmhalf(k2(3),k1(3),shift(3))
146          if(nn > nkpt_fullbz) then
147            write (message,'(a,i0)')' nkpt_fullbz mis-estimated, exceed nn=',nn
148            MSG_BUG(message)
149          end if
150          kpt_fullbz(:,nn)=k1(:)
151          nn=nn+1
152        end do
153      end do
154    end do
155  end do
156  nn = nn-1
157 
158  if (nn /= nkpt_fullbz) then
159    write (message,'(2(a,i0),a,a)')' nkpt_fullbz= ',nkpt_fullbz,' underestimated  nn=',nn,&
160 &   ch10, "Perhaps your k grid or shifts do not correspond to the symmetry?"
161    MSG_BUG(message)
162  end if
163 
164 end subroutine get_kpt_fullbz