TABLE OF CONTENTS
ABINIT/eigen_meandege [ Functions ]
NAME
eigen_meandege
FUNCTION
This routine takes the mean values of the responses for the eigenstates that are degenerate in energy.
COPYRIGHT
Copyright (C) 2011-2018 ABINIT group (XG) 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 .
INPUTS
eigenresp((3-option)*mband**(3-option)*nkpt*nsppol)= input eigenresp eigenrep(2*mband**2*nkpt*nsppol) for first-order derivatives of the eigenvalues eigenrep(mband*nkpt*nsppol) for Fan or Debye-Waller second-order derivatives of the eigenvalues mband= maximum number of bands natom= number of atoms in the unit cell nkpt= number of k-points nsppol= 1 for unpolarized, 2 for spin-polarized option= 1 for eigen(1), 2 for eigen(2) - Fan or Debye-Waller
OUTPUT
eigenresp_mean(mband*nkpt*nsppol)= eigenresp, averaged over degenerate states
PARENTS
dfpt_looppert,respfn
CHILDREN
SOURCE
37 #if defined HAVE_CONFIG_H 38 #include "config.h" 39 #endif 40 41 #include "abi_common.h" 42 43 44 subroutine eigen_meandege(eigen0,eigenresp,eigenresp_mean,mband,nband,nkpt,nsppol,option) 45 46 use defs_basis 47 use m_profiling_abi 48 use m_errors 49 50 !This section has been created automatically by the script Abilint (TD). 51 !Do not modify the following lines by hand. 52 #undef ABI_FUNC 53 #define ABI_FUNC 'eigen_meandege' 54 !End of the abilint section 55 56 implicit none 57 58 !Arguments ------------------------------------ 59 !scalars 60 integer,intent(in) :: mband,nkpt,nsppol,option 61 integer,intent(in) :: nband(nkpt*nsppol) 62 63 !arrays 64 real(dp),intent(in) :: eigen0(mband*nkpt*nsppol) 65 real(dp),intent(in) :: eigenresp((3-option)*mband**(3-option)*nkpt*nsppol) 66 real(dp),intent(out) :: eigenresp_mean(mband*nkpt*nsppol) 67 68 !Local variables------------------------------- 69 !scalars 70 integer :: bdtot_index,bd2tot_index,iband,ii,ikpt,isppol,nband_k 71 real(dp) :: eig0,mean 72 character(len=500) :: message 73 !arrays 74 75 ! ********************************************************************* 76 77 !DEBUG 78 !write(std_out,*)' eigen_meandege : enter ' 79 !ENDDEBUG 80 81 if(option/=1 .and. option/=2)then 82 write(message, '(a,i0)' )' The argument option should be 1 or 2, while it is found that option=',option 83 MSG_BUG(message) 84 end if 85 86 bdtot_index=0 ; bd2tot_index=0 87 do isppol=1,nsppol 88 do ikpt=1,nkpt 89 nband_k=nband(ikpt+(isppol-1)*nkpt) 90 if(option==1)then 91 do iband=1,nband_k 92 eigenresp_mean(iband+bdtot_index)=& 93 & eigenresp(2*iband-1 + (iband-1)*2*nband_k + bd2tot_index) 94 end do 95 else if(option==2)then 96 do iband=1,nband_k 97 eigenresp_mean(iband+bdtot_index)=eigenresp(iband+bdtot_index) 98 end do 99 end if 100 101 ! Treat the case of degeneracies : take the mean of degenerate states 102 if(nband_k>1)then 103 eig0=eigen0(1+bdtot_index) 104 ii=1 105 do iband=2,nband_k 106 if(eigen0(iband+bdtot_index)-eig0<tol8)then 107 ii=ii+1 108 else 109 mean=sum(eigenresp_mean(iband-ii+bdtot_index:iband-1+bdtot_index))/ii 110 eigenresp_mean(iband-ii+bdtot_index:iband-1+bdtot_index)=mean 111 ii=1 112 end if 113 eig0=eigen0(iband+bdtot_index) 114 if(iband==nband_k)then 115 mean=sum(eigenresp_mean(iband-ii+1+bdtot_index:iband+bdtot_index))/ii 116 eigenresp_mean(iband-ii+1+bdtot_index:iband+bdtot_index)=mean 117 end if 118 end do 119 end if 120 121 bdtot_index=bdtot_index+nband_k 122 bd2tot_index=bd2tot_index+2*nband_k**2 123 end do 124 end do 125 126 !DEBUG 127 !write(std_out,*)' eigen_meandege : exit' 128 !ENDDEBUG 129 130 end subroutine eigen_meandege