TABLE OF CONTENTS


ABINIT/order_fs_kpts [ Functions ]

[ Top ] [ Functions ]

NAME

 order_fs_kpts

FUNCTION

 This routine re-orders the kpoints on the standard grid which belong
  to the Fermi surface: put them in increasing z, then y,  then x

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

   nkptirr = number of irreducible FS kpoints
   nkpt = input nkpt from header
   kptns = input kpt from header

OUTPUT

   FSirredtoGS = mapping of irreducible kpoints to GS set
   kptirr = irreducible FS kpoint coordinates

PARENTS

      elphon

CHILDREN

      destroy_kptrank,mkkptrank,wrap2_pmhalf

SOURCE

 35 #if defined HAVE_CONFIG_H
 36 #include "config.h"
 37 #endif
 38 
 39 #include "abi_common.h"
 40 
 41 
 42 subroutine order_fs_kpts(kptns, nkpt, kptirr,nkptirr,FSirredtoGS)
 43 
 44  use defs_basis
 45  use defs_elphon
 46  use m_profiling_abi
 47  use m_kptrank
 48 
 49  use m_numeric_tools,   only : wrap2_pmhalf
 50 
 51 !This section has been created automatically by the script Abilint (TD).
 52 !Do not modify the following lines by hand.
 53 #undef ABI_FUNC
 54 #define ABI_FUNC 'order_fs_kpts'
 55 !End of the abilint section
 56 
 57  implicit none
 58 
 59 !Arguments ------------------------------------
 60 !scalars
 61  integer,intent(in) :: nkptirr
 62  integer,intent(in) :: nkpt
 63 
 64 !arrays
 65  integer,intent(out) :: FSirredtoGS(nkptirr)
 66  real(dp),intent(in) :: kptns(3,nkpt)
 67  real(dp),intent(out) :: kptirr(3,nkptirr)
 68 
 69 !Local variables-------------------------------
 70 !scalars
 71  integer :: ikpt,jkpt,kkpt,new, ik
 72  real(dp) :: res
 73  type(kptrank_type) :: kptrank_t
 74 !arrays
 75  integer :: kptirrank(nkptirr)
 76 
 77 ! *************************************************************************
 78 
 79 !rank is used to order kpoints
 80  call mkkptrank (kptns,nkpt,kptrank_t)
 81 
 82  ik=1
 83  do ikpt=1,nkpt
 84 !  add kpt to FS kpts, in order, increasing z, then y, then x !
 85    new = 1
 86 !  look for position to insert kpt ikpt among irredkpts already found
 87    do jkpt=1,ik-1
 88      if (kptirrank(jkpt) > kptrank_t%rank(ikpt)) then
 89 !      shift all the others up
 90        do kkpt=ik-1,jkpt,-1
 91          kptirr(:,kkpt+1) = kptirr(:,kkpt)
 92          kptirrank(kkpt+1) = kptirrank(kkpt)
 93          FSirredtoGS(kkpt+1) = FSirredtoGS(kkpt)
 94        end do
 95 !      insert kpoint ikpt
 96        call wrap2_pmhalf(kptns(1,ikpt),kptirr(1,jkpt),res)
 97        call wrap2_pmhalf(kptns(2,ikpt),kptirr(2,jkpt),res)
 98        call wrap2_pmhalf(kptns(3,ikpt),kptirr(3,jkpt),res)
 99 
100        kptirrank(jkpt) = kptrank_t%rank(ikpt)
101        FSirredtoGS(jkpt) = ikpt
102        new=0
103        exit
104      end if
105    end do
106 !  ikpt not counted yet and higher rank than all previous
107    if (new == 1) then
108      call wrap2_pmhalf(kptns(1,ikpt),kptirr(1,ikpt),res)
109      call wrap2_pmhalf(kptns(2,ikpt),kptirr(2,ikpt),res)
110      call wrap2_pmhalf(kptns(3,ikpt),kptirr(3,ikpt),res)
111      kptirrank(ik) = kptrank_t%rank(ikpt)
112      FSirredtoGS(ik) = ikpt
113    end if
114    ik=ik+1
115  end do
116 
117  call destroy_kptrank (kptrank_t)
118 
119 end subroutine order_fs_kpts