TABLE OF CONTENTS


ABINIT/k_neighbors [ Functions ]

[ Top ] [ Functions ]

NAME

   k_neighbors

FUNCTION

   find 8 neighbors of given k-point on a coarse grid, and return
   them along with relative k-shift within coarse grid cell

COPYRIGHT

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

   kpt        = k-point to be interpolated to, in full BZ
   kptrlatt   = lattice vectors for coarse k-grid
   invrankkpt = rank list to find k-points

OUTPUT

   rel_kpt = k-point coordinates renormalized to coarse grid cell
   kpt_phon_indices = indices of k-points on corners of cell

PARENTS

      integrate_gamma_alt

CHILDREN

      get_rank_1kpt,interpol3d_indices,wrap2_zero_one

SOURCE

 34 #if defined HAVE_CONFIG_H
 35 #include "config.h"
 36 #endif
 37 
 38 #include "abi_common.h"
 39 
 40 
 41 subroutine k_neighbors (kpt, kptrlatt,kptrank_t, rel_kpt, kpt_phon_indices)
 42 
 43  use defs_basis
 44  use m_kptrank
 45  use m_profiling_abi
 46 
 47  use m_numeric_tools,  only : wrap2_zero_one, interpol3d_indices
 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 'k_neighbors'
 53 !End of the abilint section
 54 
 55  implicit none
 56 
 57 ! inputs
 58  real(dp), intent(in) :: kpt(3)
 59  integer, intent(in) :: kptrlatt(3,3)
 60  type(kptrank_type), intent(in) :: kptrank_t
 61 
 62 ! outputs
 63  real(dp), intent(out) :: rel_kpt(3)
 64  integer, intent(out) :: kpt_phon_indices(8)
 65 
 66 ! local vars
 67  integer :: symrankkpt
 68  integer :: ir1,ir2,ir3, pr1,pr2,pr3
 69  real(dp) :: redkpt(3), cornerkpt(3), res
 70 
 71 ! *************************************************************************
 72 
 73 !wrap fine kpt to [0,1]
 74  call wrap2_zero_one(kpt(1),redkpt(1),res)
 75  call wrap2_zero_one(kpt(2),redkpt(2),res)
 76  call wrap2_zero_one(kpt(3),redkpt(3),res)
 77 !find 8 indices of points neighboring ikpt_phon, for interpolation
 78  call interpol3d_indices (redkpt,kptrlatt(1,1),kptrlatt(2,2),kptrlatt(3,3), &
 79 & ir1,ir2,ir3, pr1,pr2,pr3)
 80 
 81 !transpose ir pr to ikpt_phon indices
 82 !order of kpt_phons:
 83 !ir1 ir2 ir3
 84  cornerkpt = (/real(ir1-1)/kptrlatt(1,1),real(ir2-1)/kptrlatt(2,2), real(ir3-1)/kptrlatt(3,3)/)
 85  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
 86  kpt_phon_indices(1) = kptrank_t%invrank(symrankkpt)
 87 !pr1 ir2 ir3 
 88  cornerkpt = (/real(pr1-1)/kptrlatt(1,1),real(ir2-1)/kptrlatt(2,2), real(ir3-1)/kptrlatt(3,3)/)
 89  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
 90  kpt_phon_indices(2) = kptrank_t%invrank(symrankkpt)
 91 !ir1 pr2 ir3 
 92  cornerkpt = (/real(ir1-1)/kptrlatt(1,1),real(pr2-1)/kptrlatt(2,2), real(ir3-1)/kptrlatt(3,3)/)
 93  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
 94  kpt_phon_indices(3) = kptrank_t%invrank(symrankkpt)
 95 !pr1 pr2 ir3 
 96  cornerkpt = (/real(pr1-1)/kptrlatt(1,1),real(pr2-1)/kptrlatt(2,2), real(ir3-1)/kptrlatt(3,3)/)
 97  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
 98  kpt_phon_indices(4) = kptrank_t%invrank(symrankkpt)
 99 !ir1 ir2 pr3 
100  cornerkpt = (/real(ir1-1)/kptrlatt(1,1),real(ir2-1)/kptrlatt(2,2), real(pr3-1)/kptrlatt(3,3)/)
101  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
102  kpt_phon_indices(5) = kptrank_t%invrank(symrankkpt)
103 !pr1 ir2 pr3 
104  cornerkpt = (/real(pr1-1)/kptrlatt(1,1),real(ir2-1)/kptrlatt(2,2), real(pr3-1)/kptrlatt(3,3)/)
105  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
106  kpt_phon_indices(6) = kptrank_t%invrank(symrankkpt)
107 !ir1 pr2 pr3 
108  cornerkpt = (/real(ir1-1)/kptrlatt(1,1),real(pr2-1)/kptrlatt(2,2), real(pr3-1)/kptrlatt(3,3)/)
109  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
110  kpt_phon_indices(7) = kptrank_t%invrank(symrankkpt)
111 !pr1 pr2 pr3 
112  cornerkpt = (/real(pr1-1)/kptrlatt(1,1),real(pr2-1)/kptrlatt(2,2), real(pr3-1)/kptrlatt(3,3)/)
113  call get_rank_1kpt (cornerkpt,symrankkpt,kptrank_t)
114  kpt_phon_indices(8) = kptrank_t%invrank(symrankkpt)
115 
116 !retrieve the gkq matrix for all q, at the neighbor k vectors
117  rel_kpt(1) = redkpt(1)*kptrlatt(1,1)-real(ir1-1)
118  rel_kpt(2) = redkpt(2)*kptrlatt(2,2)-real(ir2-1)
119  rel_kpt(3) = redkpt(3)*kptrlatt(3,3)-real(ir3-1)
120 
121 end subroutine k_neighbors