fermie.F90 7.9 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230
MODULE m_fermie
  USE m_juDFT
  !-----------------------------------------------------------------------
  !     determines the fermi energy by
  !            gaussian-integration method                          c.l.fu
  !            triangular method (or tetrahedrons)
  !            or fermi-function                                    p.kurz
  !----------------------------------------------------------------------
CONTAINS
  SUBROUTINE fermie(eig_id, mpi,kpts,obsolete,&
       input, noco,e_min,jij,cell,results)

    !---------------------------------------------------f--------------------
    !
    !     a fist (T=0) approximation to the fermi-energy is determined
    !     by:
    !           zelec = sum { spindg * we }
    !                       e=<ef
    !
    !     TREE STRUCTURE: fermie----sort
    !                           ----fergwt
    !                           ----fertri---triang
    !                                     |--dosint
    !                                     |--dosef -- trisrt
    !                                     +--doswt
    !                           ----ferhis---ef_newton
    !
    !-----------------------------------------------------------------------


    USE m_eig66_io, ONLY : read_eig
#ifdef CPP_MPI
    USE m_mpi_col_eigJ
#endif
    USE m_sort
    USE m_fertri
    USE m_ferhis
    USE m_fergwt
    USE m_types
    IMPLICIT NONE
    TYPE(t_results),INTENT(INOUT)   :: results
    TYPE(t_mpi),INTENT(IN)   :: mpi
    TYPE(t_obsolete),INTENT(IN)   :: obsolete
    TYPE(t_input),INTENT(IN)   :: input
    TYPE(t_noco),INTENT(IN)   :: noco
    TYPE(t_jij),INTENT(IN)   :: jij
    TYPE(t_cell),INTENT(IN)   :: cell
    TYPE(t_kpts),INTENT(IN)   :: kpts
    !     ..
    !     .. Scalar Arguments ..
    INTEGER, INTENT (IN) :: eig_id
    REAL,INTENT(IN)      :: e_min
    !     ..
    !     .. Array Arguments ..
    !REAL,    INTENT (OUT):: w(:,:,:) !(dimension%neigd,kpts%nkptd,dimension%jspd)
    !     ..
    !     .. Local Scalars ..
    REAL del  ,spindg,ssc ,ws,zc,tkb_1,weight
    INTEGER i,idummy,j,jsp,k,l,n,nbands,nstef,nv,nmat,nspins
    INTEGER n_help
    !     ..
    !     .. Local Arrays ..
    !
    INTEGER, ALLOCATABLE :: idxeig(:),idxjsp(:),idxkpt(:),INDEX(:)
    REAL,    ALLOCATABLE :: e(:),eig(:,:,:),we(:)
    INTEGER ne(kpts%nkptd,SIZE(results%w_iks,3))

    !--- J constants
    !--- J constants

#ifdef CPP_MPI
    INCLUDE 'mpif.h'
    INTEGER, PARAMETER :: comm = MPI_COMM_SELF
    INTEGER*4 :: nv_mpi(2),idum1d(0),idum2d(0,0)
#endif

    !     ..
    !     ..
    !***********************************************************************
    !                  ABBREVIATIONS
    !
    !     eig        : array of eigenvalues 
    !     wtkpt      : list of the weights of each k-point (from inp-file)
    !     e          : linear list of the eigenvalues 
    !     we         : list of weights of the eigenvalues in e
    !     zelec      : number of electrons 
    !     spindg     : spindegeneracy (2 in nonmagnetic calculations)
    !     seigv      : weighted sum of the occupied valence eigenvalues
    !     seigsc     : weighted sum of the semi-core eigenvalues
    !     seigscv    : sum of seigv and seigsc
    !     ts         : entropy contribution to the free energy
    !
    !***********************************************************************
    !     .. Data statements ..
    DATA del/1.0e-6/
    !     ..
    n=SIZE(results%w_iks) !size of list of all eigenvalues
    ALLOCATE (idxeig(n),idxjsp(n),idxkpt(n),INDEX(n),e(n),we(n) )
    ALLOCATE (eig(SIZE(results%w_iks,1),SIZE(results%w_iks,2),SIZE(results%w_iks,3)))

    ! initiliaze e
    e = 0
    !
    IF (jij%l_J) THEN
#ifdef CPP_MPI
       CALL mpi_col_eigJ(mpi%mpi_comm,mpi%irank,mpi%isize,kpts%nkptd,SIZE(results%w_iks,1),kpts%nkpt(1),&
            &                       jij%nkpt_l,jij%eig_l,&
            &                    kpts%bk,kpts%wtkpt,ne(1,1),eig)!keep
       IF (mpi%irank.NE.0) THEN
          DEALLOCATE( idxeig,idxjsp,idxkpt,index,e,eig,we )
          RETURN
       ENDIF
#endif
    ENDIF


    IF ( mpi%irank == 0 ) WRITE (6,FMT=8000)
8000 FORMAT (/,/,1x,'fermi energy and band-weighting factors:')
    !
    !---> READ IN EIGENVALUES
    !
    spindg = 2.0/REAL(input%jspins)
    n = 0
    results%seigsc = 0.0
    ssc = 0.0
    n_help = 0
    !
    !---> pk non-collinear
    IF (noco%l_noco) THEN
       nspins = 1
    ELSE
       nspins = input%jspins
    ENDIF
    !---> pk non-collinear
    !
    DO  jsp = 1,nspins
       DO  k = 1,kpts%nkpt
          CALL read_eig(eig_id,k,jsp,neig=ne(k,jsp),eig=eig(:,k,jsp))
          IF ( mpi%irank == 0 ) THEN
             WRITE (6,'(a2,3f10.5,f12.6)') 'at',kpts%bk(:,k),kpts%wtkpt(k)
             WRITE (6,'(i5,a14)') ne(k,jsp),' eigenvalues :' 
             WRITE (6,'(8f12.6)') (eig(i,k,jsp),i=1,ne(k,jsp))
          END IF
          nv= -1 
          !
          !--->          STORE EIGENVALUES AND WEIGHTS IN A LINEAR LIST. AND MEMORIZE 
          !--->          CONECTION TO THE ORIGINAL ARRAYS
          !
          DO  j = 1,ne(k,jsp)
             e(n+j) = eig(j,k,jsp)
             we(n+j) = kpts%wtkpt(k)
             idxeig(n+j) = j+n_help
             idxkpt(n+j) = k
             idxjsp(n+j) = jsp
 	  ENDDO
          !--->          COUNT THE NUMBER OF EIGENVALUES
          n = n + ne(k,jsp)
       ENDDO
    ENDDO
    CALL sort(n,e,index)

    !     Check if no deep eigenvalue is found
    IF (e_min-MINVAL(e(1:n))>1.0) THEN
       WRITE(6,*) 'WARNING: Too low eigenvalue detected:'
       WRITE(6,*) 'min E=', MINVAL(e(1:n)),' min(enpara)=',&
            &             e_min
       CALL juDFT_warn("Too low eigenvalue detected",calledby="fermi" &
            &     ,hint ="If the lowest eigenvalue is more than 1Htr below "//&
            &     "the lowest energy parameter, you probably have picked up"//&
            &     " a ghoststate")
    END IF
    !
    !---> DETERMINE EF BY SUMMING WEIGHTS
    !
    weight = input%zelec/spindg
    results%seigv = 0.0e0
    ws = 0.0e0
    l = 0
    DO WHILE ((ws+del).LT.weight)
       l = l + 1
       IF (l.GT.n) THEN
          IF ( mpi%irank == 0 ) THEN
             WRITE (16,FMT=8010) n,ws,weight
             WRITE (6,FMT=8010) n,ws,weight
          END IF
          CALL juDFT_error("fermi",calledby="fermie")
8010      FORMAT (/,10x,'error: not enough wavefunctions.',i10,&
               &             2d20.10)
       END IF
       ws = ws + we(INDEX(l))
       results%seigv = results%seigv + e(INDEX(l))*we(INDEX(l))*spindg
       !         WRITE (6,FMT='(2f10.7)') e(index(l)),we(index(l))
    END DO
    results%ef = e(INDEX(l))
    nstef = l
    zc = input%zelec
    IF ( mpi%irank == 0 ) WRITE (6,FMT=8020) results%ef,nstef,results%seigv,ws,results%seigsc,ssc

    !+po
    results%ts = 0.0
    !-po
    results%w_iks = 0.0
    IF (input%gauss) THEN
       CALL fergwt(kpts,input,mpi,ne, eig,results)
    ELSE IF (input%tria) THEN
       CALL fertri(mpi%irank, ne,kpts%nkpt,nspins,zc,eig,kpts%bk,spindg,&
            results%ef,results%seigv,results%w_iks)
    ELSE
       nspins = input%jspins
       IF (noco%l_noco) nspins = 1
       tkb_1 = input%tkb
       CALL ferhis(input,kpts,mpi,results,index,idxeig,idxkpt,idxjsp, n,&
            nstef,ws,spindg,weight,e,ne,we, noco,jij,cell)
    END IF
    !     7.12.95 r.pentcheva seigscv must be calculated outside if (gauss)
    results%seigscv = results%seigsc + results%seigv
    !
    DEALLOCATE ( idxeig,idxjsp,idxkpt,index,e,eig,we )
    !
    RETURN
8020 FORMAT (/,'FERMIE:',/,&
         &       10x,'first approx. to results%ef    (T=0)  :',f10.6,' htr',&
         &       '   (energy of the highest occ. eigenvalue)',/,&
         &       10x,'number of occ. states  (T=0)  :',i10,/,&
         &       10x,'first approx. to results%seigv (T=0)  :',f10.6,' htr',/,&
         &       10x,'sum of weights of occ. states :',f10.6,/,&
         &       10x,'sum of semicore eigenvalues   :',f10.6,' htr',/,&
         &       10x,'sum of semicore charge        :',f10.6,' e',/)
  END SUBROUTINE fermie
     END MODULE m_fermie