TABLE OF CONTENTS


ABINIT/gaugetransfo [ Functions ]

[ Top ] [ Functions ]

NAME

 gaugetransfo

FUNCTION

 This routine allows the passage from the parallel-transport gauge
 to the diagonal gauge for the first-order wavefunctions

COPYRIGHT

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

  cg_k(2,mpw*nspinor*mband*nsppol)=planewave coefficients of wavefunctions
                                   for a particular k point.
  cwavef(2,npw1_k*nspinor)=first order wavefunction for a particular k point
                           in the parallel gauge
  eig_k(mband*nsppol)=GS eigenvalues at k (hartree)
  eig1_k(2*nsppol*mband**2)=matrix of first-order eigenvalues (hartree)
  iband=band index of the 1WF for which the transformation has to be applied
  mband=maximum number of bands
  nband_k=number of bands for this k point
  npw_k=maximum dimensioned size of npw or wfs at k
  npw1_k=number of plane waves at this k+q point
  nspinor=number of spinorial components of the wavefunctions
  nsppol=1 for unpolarized, 2 for spin-polarized
  occ_k(nband_k)=occupation number for each band (usually 2) for each k

OUTPUT

  cwavef_d(2,npw1_k*nspinor)=first order wavefunction for a particular k point
                             in the diagonal gauge

PARENTS

      dfpt_nstwf

CHILDREN

SOURCE

 44 #if defined HAVE_CONFIG_H
 45 #include "config.h"
 46 #endif
 47 
 48 #include "abi_common.h"
 49 
 50 
 51 subroutine gaugetransfo(cg_k,cwavef,cwavef_d,eig_k,eig1_k,iband,nband_k, &
 52 &                      mband,npw_k,npw1_k,nspinor,nsppol,occ_k)
 53 
 54 
 55  use defs_basis
 56  use m_profiling_abi
 57 
 58 !This section has been created automatically by the script Abilint (TD).
 59 !Do not modify the following lines by hand.
 60 #undef ABI_FUNC
 61 #define ABI_FUNC 'gaugetransfo'
 62 !End of the abilint section
 63 
 64  implicit none
 65 
 66 !Arguments ------------------------------------
 67 !scalars
 68  integer,intent(in) :: iband,mband,nband_k,npw1_k,npw_k,nspinor,nsppol
 69 !arrays
 70  real(dp),intent(in) :: cg_k(2,npw_k*nspinor*nband_k),cwavef(2,npw1_k*nspinor)
 71  real(dp),intent(in) :: eig1_k(2*nsppol*mband**2),eig_k(mband*nsppol)
 72  real(dp),intent(in) :: occ_k(nband_k)
 73  real(dp),intent(out) :: cwavef_d(2,npw1_k*nspinor)
 74 
 75 !Local variables-------------------------------
 76 !tolerance for non degenerated levels
 77 !scalars
 78  integer :: jband
 79  real(dp),parameter :: etol=1.0d-3
 80 !arrays
 81  real(dp) :: cwave0(2,npw1_k*nspinor),eig1(2)
 82 
 83 ! *********************************************************************
 84 
 85 !DEBUG
 86 !write(100,*) 'gaugetransfo: ',iband
 87 !ENDDEBUG
 88 
 89  cwavef_d(:,:) = cwavef(:,:)
 90 
 91  do jband = 1,nband_k !loop over bands
 92 
 93    if ((abs(eig_k(iband)-eig_k(jband)) > etol).and.(abs(occ_k(jband)) > tol8 )) then
 94 
 95      cwave0(:,:) = cg_k(:,1+(jband-1)*npw_k*nspinor:jband*npw_k*nspinor)
 96 
 97      eig1(1) = eig1_k(2*jband-1+(iband-1)*2*nband_k)
 98      eig1(2) = eig1_k(2*jband +(iband-1)*2*nband_k)
 99 
100      cwavef_d(1,:)=cwavef_d(1,:) &
101 &     - (eig1(1)*cwave0(1,:)-eig1(2)*cwave0(2,:))/(eig_k(jband)-eig_k(iband))
102      cwavef_d(2,:)=cwavef_d(2,:) &
103 &     - (eig1(1)*cwave0(2,:)+eig1(2)*cwave0(1,:))/(eig_k(jband)-eig_k(iband))
104 
105    end if
106 
107  end do    !loop over bands
108 
109 end subroutine gaugetransfo