TABLE OF CONTENTS


ABINIT/overlap_g [ Functions ]

[ Top ] [ Functions ]

NAME

 overlap_g

FUNCTION

 Compute the scalar product between WF at two different k-points
 < u_{n,k1} | u_{n,k2} >

COPYRIGHT

 Copyright (C) 1999-2017 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 .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

 mpw = maximum dimensioned size of npw
 npw_k1 = number of plane waves at k1
 npw_k2 = number of plane waves at k2
 nspinor = 1 for scalar, 2 for spinor wavefunctions
 pwind_k = array required to compute the scalar product (see initberry.f)
 vect1 = wavefunction at k1: | u_{n,k1} >
 vect2 = wavefunction at k1: | u_{n,k2} >

OUTPUT

 doti = imaginary part of the scalarproduct
 dotr = real part of the scalarproduct

NOTES

 In case a G-vector of the basis sphere of plane waves at k1
 does not belong to the basis sphere of plane waves at k2,
 pwind = 0. Therefore, the dimensions of vect1 &
 vect2 are (1:2,0:mpw) and the element (1:2,0) MUST be
 set to zero.

 The current implementation if not compatible with TR-symmetry (i.e. istwfk/=1) !

PARENTS

      dfptff_bec,dfptff_die,dfptff_ebp,dfptff_edie,dfptff_gbefd
      dfptff_gradberry,qmatrix,smatrix

CHILDREN

SOURCE

 47 #if defined HAVE_CONFIG_H
 48 #include "config.h"
 49 #endif
 50 
 51 #include "abi_common.h"
 52 
 53 subroutine overlap_g(doti,dotr,mpw,npw_k1,npw_k2,nspinor,pwind_k,vect1,vect2)
 54 
 55  use defs_basis
 56  use m_profiling_abi
 57  use m_errors
 58 
 59 !This section has been created automatically by the script Abilint (TD).
 60 !Do not modify the following lines by hand.
 61 #undef ABI_FUNC
 62 #define ABI_FUNC 'overlap_g'
 63 !End of the abilint section
 64 
 65  implicit none
 66 
 67 !Arguments ------------------------------------
 68 !scalars
 69  integer,intent(in) :: mpw,npw_k1,npw_k2,nspinor
 70  real(dp),intent(out) :: doti,dotr
 71 !arrays
 72  integer,intent(in) :: pwind_k(mpw)
 73  real(dp),intent(in) :: vect1(1:2,0:mpw*nspinor),vect2(1:2,0:mpw*nspinor)
 74 
 75 !Local variables-------------------------------
 76 !scalars
 77  integer :: ipw,ispinor,jpw,spnshft1,spnshft2
 78  character(len=500) :: message
 79 
 80 ! *************************************************************************
 81 
 82 !Check if vect1(:,0) = 0 and vect2(:,0) = 0
 83  if ((abs(vect1(1,0)) > tol12).or.(abs(vect1(2,0)) > tol12).or. &
 84 & (abs(vect2(1,0)) > tol12).or.(abs(vect2(2,0)) > tol12)) then
 85    message = ' vect1(:,0) and/or vect2(:,0) are not equal to zero'
 86    MSG_BUG(message)
 87  end if
 88 
 89 !Compute the scalar product
 90  dotr = zero; doti = zero
 91  do ispinor = 1, nspinor
 92    spnshft1 = (ispinor-1)*npw_k1
 93    spnshft2 = (ispinor-1)*npw_k2
 94 !$OMP PARALLEL DO PRIVATE(jpw) REDUCTION(+:doti,dotr) 
 95    do ipw = 1, npw_k1
 96      jpw = pwind_k(ipw)
 97      dotr = dotr + vect1(1,spnshft1+ipw)*vect2(1,spnshft2+jpw) + vect1(2,spnshft1+ipw)*vect2(2,spnshft2+jpw)
 98      doti = doti + vect1(1,spnshft1+ipw)*vect2(2,spnshft2+jpw) - vect1(2,spnshft1+ipw)*vect2(1,spnshft2+jpw)
 99    end do
100  end do
101 
102 end subroutine overlap_g