TABLE OF CONTENTS


ABINIT/kptfine_av [ Functions ]

[ Top ] [ Functions ]

NAME

 kptfine_av

FUNCTION

 This routine compute the k-points of a fine grid that are around a 
 define k-point of a coarse mesh.  

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (SP)
 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 .

INPUTS

  center(3) = the point of the coarse mesh around which you want know which
              k-points of the fine mesh belong to.
  qptrlatt(3,3) = qptrlatt of the considered calculation (this is obtained
              from the input variable ngqpt and shiftq.
  kpt_fine(3,nkpt_fine) = this table contain all the k-points of the fine grid
              in the full BZ (no sym op. allowed) and is read from the header
              of the dense WF file.
  nkpt_fine = number of k-points of the fine grid read from the header of the
              dense WF file. 

OUTPUT

  kpt_fine_sub(nkpt_sub) = k-points of the fine grid that are around center(3)
  nkpt_sub = number of k-points of the fine grid that are around center(3)
  wgt_sub(nkpt_sub) = weight of the k-points of the fine grid that are around
              center(3).  

PARENTS

      eig2stern,eig2tot

CHILDREN

SOURCE

 40 #if defined HAVE_CONFIG_H
 41 #include "config.h"
 42 #endif
 43 
 44 #include "abi_common.h"
 45 
 46 subroutine kptfine_av(center,qptrlatt,kpt_fine,nkpt_fine,kpt_fine_sub,&
 47 &                nkpt_sub,wgt_sub)
 48 
 49  use defs_basis
 50  use m_profiling_abi
 51  use m_xmpi
 52  use m_errors
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'kptfine_av'
 58 !End of the abilint section
 59 
 60  implicit none
 61 
 62 !Arguments ------------------------------------
 63 !scalars
 64  integer,intent(in)   :: nkpt_fine
 65  integer,intent(out)  :: nkpt_sub
 66 !arrays
 67  integer,intent(in)   :: qptrlatt(3,3)
 68  real(dp),intent(in)  :: kpt_fine(3,nkpt_fine)
 69  real(dp),intent(in)  :: center(3)
 70  integer,pointer      :: kpt_fine_sub(:)
 71  real(dp),pointer     :: wgt_sub(:)
 72 
 73 !Local variables-------------------------------
 74 !scalars
 75  integer :: ikpt,aa,bb,cc
 76  integer :: ii,jj
 77 !arrays
 78  real(dp) :: center_ref(3)
 79  real(dp) :: kpt_fine_ref(3)
 80  real(dp) :: kpt_tmp(3),kpt_tmp2(3)
 81  integer,allocatable  :: kpt_fine_sub_tmp(:)
 82  real(dp),allocatable :: wgt_sub_tmp(:)
 83  logical :: found(3)
 84 
 85 ! *************************************************************************
 86 
 87  ABI_ALLOCATE(kpt_fine_sub_tmp,(nkpt_fine))
 88  ABI_ALLOCATE(wgt_sub_tmp,(nkpt_fine))
 89 
 90 !It is easier to work in real space using the qptrlatt matrices because in this
 91 !referential any k-points sampling will be cast into an orthorhombic shape.
 92 !In that space we can simply take all k-points of the fine grid that between
 93 !center_ref-0.5 and center_ref+0.5
 94 
 95  center_ref = MATMUL(qptrlatt,center)
 96 
 97 !When considering points center(3) that lying close or on a BZ edge we need to
 98 !take the k-points of the fine grid taking into account unklamp vectors. This
 99 !is done with the aa, bb and cc loops. 
100 
101  ii = 1
102  do ikpt=1,nkpt_fine
103    kpt_tmp = kpt_fine(:,ikpt)
104    do aa=-1,1
105      kpt_tmp2(1) = kpt_tmp(1)+aa
106      do bb=-1,1
107        kpt_tmp2(2) = kpt_tmp(2)+bb
108        do cc=-1,1
109          kpt_tmp2(3) = kpt_tmp(3)+cc
110          kpt_fine_ref = MATMUL(qptrlatt,kpt_tmp2)
111          if((kpt_fine_ref(1)>=center_ref(1)-0.5-tol8).and.&
112 &         (kpt_fine_ref(1)<=center_ref(1)+0.5+tol8)) then
113            if((kpt_fine_ref(2)>=center_ref(2)-0.5-tol8).and.&
114 &           (kpt_fine_ref(2)<=center_ref(2)+0.5+tol8)) then
115              if((kpt_fine_ref(3)>=center_ref(3)-0.5-tol8).and.&
116 &             (kpt_fine_ref(3)<=center_ref(3)+0.5+tol8)) then
117                kpt_fine_sub_tmp(ii) = ikpt
118                ii = ii +1
119              end if
120            end if
121          end if
122        end do
123      end do
124    end do
125  end do
126 
127  nkpt_sub = ii-1
128  ABI_ALLOCATE(kpt_fine_sub,(nkpt_sub))
129  ABI_ALLOCATE(wgt_sub,(nkpt_sub))
130 
131  do jj=1,nkpt_sub
132    kpt_fine_sub(jj) = kpt_fine_sub_tmp(jj)
133  end do
134 
135 !We then compute a weight function. This weight function is simply a
136 !rectangular weight function that take the value 1 for k-points of the fine
137 !grid inside the cube, 0.5 for k-points that are lying on one face of the cube,
138 !0.25 for k-points that are lying on an edge of the cube and 0.125 for k-points
139 !that are lying on a peak of the cube.  
140 
141  wgt_sub(:) = 1.0
142 
143  do ikpt=1,nkpt_sub
144    found(:) = .True.
145    kpt_tmp = kpt_fine(:,kpt_fine_sub(ikpt))
146    do aa=-1,1
147      kpt_tmp2(1) = kpt_tmp(1)+aa
148      do bb=-1,1
149        kpt_tmp2(2) = kpt_tmp(2)+bb
150        do cc=-1,1
151          kpt_tmp2(3) = kpt_tmp(3)+cc
152          kpt_fine_ref = MATMUL(qptrlatt,kpt_tmp2)
153          if((ABS(kpt_fine_ref(1)-center_ref(1)-0.5)< tol8) .or.&
154 &         (ABS(kpt_fine_ref(1)-center_ref(1)+0.5) < tol8)) then
155            if(found(1)) then
156              wgt_sub(ikpt) = wgt_sub(ikpt)*0.5
157              found(1) = .False.
158            end if
159          end if
160          if((ABS(kpt_fine_ref(2)-center_ref(2)-0.5) < tol8) .or.&
161 &         (ABS(kpt_fine_ref(2)-center_ref(2)+0.5) < tol8)) then
162            if(found(2)) then
163              wgt_sub(ikpt) = wgt_sub(ikpt)*0.5
164              found(2) = .False.
165            end if
166          end if
167          if((ABS(kpt_fine_ref(3)-center_ref(3)-0.5)< tol8) .or.&
168 &         (ABS(kpt_fine_ref(3)-center_ref(3)+0.5) < tol8)) then
169            if(found(3)) then
170              wgt_sub(ikpt) = wgt_sub(ikpt)*0.5
171              found(3) = .False.
172            end if
173          end if
174        end do
175      end do
176    end do
177  end do
178 
179  ABI_DEALLOCATE(kpt_fine_sub_tmp)
180  ABI_DEALLOCATE(wgt_sub_tmp)
181 
182 end subroutine kptfine_av