TABLE OF CONTENTS


ABINIT/expibi [ Functions ]

[ Top ] [ Functions ]

NAME

 expibi

FUNCTION

 Routine that computes exp(i (-b_ket).R) at each site.

COPYRIGHT

 Copyright (C) 2005-2018 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 .

INPUTS

  dkvecs(3) :: $\Delta k$ increment
  natom :: number of atoms in unit cell
  xred(natom,3) :: reduced coordinates of atoms in unit cell

OUTPUT

  calc_expibi(2,natom) :: phase factors at each atom for vector shift

SIDE EFFECTS

NOTES

PARENTS

      initberry,overlap_k1k2_paw

CHILDREN

SOURCE

34 #if defined HAVE_CONFIG_H
35 #include "config.h"
36 #endif
37 
38 #include "abi_common.h"
39 
40  subroutine expibi(calc_expibi,dkvecs,natom,xred)
41 
42  use m_profiling_abi
43  use m_errors
44  use defs_basis
45 
46 !This section has been created automatically by the script Abilint (TD).
47 !Do not modify the following lines by hand.
48 #undef ABI_FUNC
49 #define ABI_FUNC 'expibi'
50 !End of the abilint section
51 
52  implicit none
53 
54 !Arguments---------------------------
55 !scalars
56  integer,intent(in) :: natom
57  real(dp),intent(out) :: calc_expibi(2,natom)
58 !arrays
59  real(dp),intent(in) :: dkvecs(3),xred(3,natom)
60 
61 !Local variables---------------------------
62 !scalars
63  integer :: iatom
64  real(dp) :: bdotr
65 
66 ! *************************************************************************
67 
68  calc_expibi(:,:) = zero
69 
70  !calc_expibi(2,natom)
71  !used for PAW field calculations (distributed over atomic sites)
72  !stores the on-site phase factors arising from
73  !$\langle\phi_{i,k}|\phi_{j,k+\sigma_k k_k}\rangle$
74  !where $\sigma = \pm 1$. These overlaps arise in various Berry
75  !phase calculations of electric and magnetic polarization. The on-site
76  !phase factor is $\exp[-i\sigma_k k_k)\cdot I]$ where
77  !$I$ is the nuclear position. 
78 
79  do iatom = 1, natom
80 
81     !    note the definition used for the k-dependence of the PAW basis functions:
82     !$|\phi_{i,k}\rangle = exp(-i k\cdot r)|\phi_i\rangle
83     !    see Umari, Gonze, and Pasquarello, PRB 69,235102 Eq. 23. 
84    bdotr = DOT_PRODUCT(xred(1:3,iatom),-dkvecs(1:3))
85     !    here is exp(i b.R) for the given site
86    calc_expibi(1,iatom) = cos(two_pi*bdotr)
87    calc_expibi(2,iatom) = sin(two_pi*bdotr)
88 
89  end do ! end loop on natom
90 
91 end subroutine expibi