TABLE OF CONTENTS


ABINIT/occeig [ Functions ]

[ Top ] [ Functions ]

NAME

 occeig

FUNCTION

 For each pair of active bands (m,n), generates ratios
 that depend on the difference between occupation numbers
 and eigenvalues.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 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

  doccde_k(nband_k)=derivative of occ_k wrt the energy
  doccde_kq(nband_k)=derivative of occ_kq wrt the energy
  eig0_k(nband_k)=GS eigenvalues at k
  eig0_kq(nband_k)=GS eigenvalues at k+q
  nband_k=number of bands
  occopt=option for occupancies
  occ_k(nband_k)=occupation number for each band at k
  occ_kq(nband_k)=occupation number for each band at k+q

OUTPUT

  rocceig(nband_k,nband_k)$= (occ_{k,q}(m)-occ_k(n))/(eig0_{k,q}(m)-eig0_k(n))$,
   if this ratio has been attributed to the band n, 0.0_dp otherwise

SIDE EFFECTS

NOTES

 Supposing the occupations numbers differ :
 if $abs(occ_{k,q}(m)) < abs(occ_k(n))$
  $rocceig(m,n)=(occ_{k,q}(m)-occ_k(n))/(eig0_{k,q}(m)-eig0_k(n)) $
 if $abs(occ_{k,q}(m))>abs(occ_k(n))$
  rocceig(m,n)=0.0_dp

 If the occupation numbers are close enough, then
 if the eigenvalues are also close, take the derivative
  $ rocceig(m,n)=\frac{1}{2}*docc/deig0 $
 otherwise,
  $ rocceig(m,n)=\frac{1}{2}*(occ_{k,q}(m)-occ_k(n))/(eig0_{k,q}(m)-eig0_k(n))$

PARENTS

      dfpt_nstpaw,dfpt_rhofermi,dfpt_vtorho

CHILDREN

SOURCE

 53 #if defined HAVE_CONFIG_H
 54 #include "config.h"
 55 #endif
 56 
 57 #include "abi_common.h"
 58 
 59 
 60 subroutine occeig(doccde_k,doccde_kq,eig0_k,eig0_kq,nband_k,occopt,occ_k,occ_kq,rocceig)
 61 
 62  use defs_basis
 63  use m_errors
 64  use m_profiling_abi
 65 
 66 !This section has been created automatically by the script Abilint (TD).
 67 !Do not modify the following lines by hand.
 68 #undef ABI_FUNC
 69 #define ABI_FUNC 'occeig'
 70 !End of the abilint section
 71 
 72  implicit none
 73 
 74 !Arguments ------------------------------------
 75 !scalars
 76  integer,intent(in) :: nband_k,occopt
 77 !arrays
 78  real(dp),intent(in) :: doccde_k(nband_k),doccde_kq(nband_k),eig0_k(nband_k)
 79  real(dp),intent(in) :: eig0_kq(nband_k),occ_k(nband_k),occ_kq(nband_k)
 80  real(dp),intent(out) :: rocceig(nband_k,nband_k)
 81 
 82 !Local variables-------------------------------
 83 !scalars
 84  integer :: ibandk,ibandkq
 85  real(dp) :: diffabsocc,diffeig,diffocc,ratio,sumabsocc
 86  character(len=500) :: message
 87 
 88 ! *************************************************************************
 89 
 90 !The parameter tol5 defines the treshhold for degeneracy, and the width of the step function
 91 
 92  rocceig(:,:)=0.0_dp
 93 
 94  do ibandk=1,nband_k
 95    do ibandkq=1,nband_k
 96      diffeig=eig0_kq(ibandkq)-eig0_k(ibandk)
 97      diffocc=occ_kq(ibandkq)-occ_k(ibandk)
 98 
 99      if( abs(diffeig) > tol5 ) then
100        ratio=diffocc/diffeig
101      else
102        if(occopt<3)then
103 !        In a non-metallic case, if the eigenvalues are degenerate,
104 !        the occupation numbers must also be degenerate, in which
105 !        case there is no contribution from this pair of bands
106          if( abs(diffocc) > tol5 ) then
107            write(message,'(a,a,a,a,a,a,a,2(a,i4,a,es16.6,a,es16.6,a,a),a)' ) &
108 &           'In a non-metallic case (occopt<3), for a RF calculation,',ch10,&
109 &           'if the eigenvalues are degenerate,',&
110 &           ' the occupation numbers must also be degenerate.',ch10,&
111 &           'However, the following pair of states gave :',ch10,&
112 &           'k -state, band number',ibandk,', occ=',occ_k(ibandk),&
113 &           'eigenvalue=',eig0_k(ibandk),',',ch10,&
114 &           ' kq-state, band number',ibandkq,', occ=',occ_kq(ibandkq),&
115 &           ', eigenvalue=',eig0_kq(ibandkq),'.',ch10,&
116 &           'Action : change occopt, consistently, in GS and RF calculations.'
117            MSG_ERROR(message)
118          end if
119          ratio=0.0_dp
120        else
121 !        In the metallic case, one can compute a better approximation of the
122 !        ratio by using derivatives doccde
123          ratio=0.5_dp*(doccde_kq(ibandkq)+doccde_k(ibandk))
124 !        DEBUG
125 !        write(std_out,*)' occeig : ibandkq,doccde_kq(ibandkq)',ibandkq,doccde_kq(ibandkq)
126 !        write(std_out,*)'          ibandk ,doccde_k (ibandk )',ibandk,doccde_k(ibandk)
127 !        ENDDEBUG
128        end if
129      end if
130 
131 !    Here, must pay attention to the smallness of some coefficient
132      diffabsocc=abs(occ_k(ibandk))-abs(occ_kq(ibandkq))
133      sumabsocc=abs(occ_k(ibandk))+abs(occ_kq(ibandkq))
134      if(sumabsocc>tol8)then
135        if( diffabsocc > sumabsocc*tol5 ) then
136          rocceig(ibandkq,ibandk)=ratio
137        else if ( diffabsocc >= -sumabsocc*tol5 ) then
138          rocceig(ibandkq,ibandk)=0.5_dp*ratio
139        else
140          rocceig(ibandkq,ibandk)=0.0_dp
141        end if
142      end if
143 
144    end do ! ibandkq
145  end do ! ibandk
146 
147 end subroutine occeig