v_mmp.F90 7.46 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
MODULE m_vmmp
  !     ************************************************************
  !     This subroutine calculates the potential matrix v^{s}_{m,m'}
  !     for a given atom 'n' and l-quantum number 'l'. The l,u,j's for
  !     all atoms are stored in lda_u, the density matrix is ns_mmp,
  !     and the e-e- interaction potential is u(m1,m2,m3,m4,n).
  !     For details see Eq.(16) of Shick et al. PRB 60, 10765 (1999)
  !
  !     Additionally, the total energy contribution of LDA+U (Eq.24)
  !     is calculated (e_ldau).
  !     Part of the LDA+U package                   G.B., Oct. 2000
  !     ************************************************************
CONTAINS
  SUBROUTINE v_mmp(atoms,jspins,lmaxb,ns_mmp,u,f0,f2, vs_mmp,results)

    USE m_types
    IMPLICIT NONE

    TYPE(t_results),INTENT(INOUT)   :: results
    TYPE(t_atoms),INTENT(IN)        :: atoms
    !
    ! ..  Arguments ..
    INTEGER, INTENT (IN) :: lmaxb,jspins 
    REAL,    INTENT (IN) :: u(-lmaxb:lmaxb,-lmaxb:lmaxb, -lmaxb:lmaxb,-lmaxb:lmaxb,atoms%n_u)
    REAL,    INTENT (IN) :: f0(atoms%n_u),f2(atoms%n_u)

    COMPLEX           :: ns_mmp(-lmaxb:lmaxb,-lmaxb:lmaxb,atoms%n_u,jspins)
    COMPLEX,INTENT(OUT)::vs_mmp(-lmaxb:lmaxb,-lmaxb:lmaxb,atoms%n_u,jspins)
    !
    ! ..  Local Variables ..
    INTEGER n,ispin,jspin,l ,mp,p,q,itype,m
    REAL rho_tot,u_htr,j_htr,e_ee,ns_sum,spin_deg,e_dc,e_dcc
    REAL rho_sig(jspins),v_diag(jspins),eta(0:jspins)
    !
    ! Use around-mean-field limit (true) of atomic limit (false)
    !
    !
    ! Loop over atoms
    !
    spin_deg = 1.0 / (3 - jspins)
    results%e_ldau = 0.0
    n = 0
    DO itype = 1, atoms%ntype
       l = atoms%lda_u(itype)%l
       IF (l.GE.0) THEN
          n = n + 1
          u_htr = atoms%lda_u(itype)%u / 27.21
          j_htr = atoms%lda_u(itype)%j / 27.21
          u_htr = f0(n)/27.21
          IF (l.EQ.1) THEN
             j_htr = f2(n)/(5*27.21)
          ELSEIF (l.EQ.2) THEN
             j_htr = 1.625*f2(n)/(14*27.21)
          ELSEIF (l.EQ.3) THEN
             j_htr = (286.+195*451/675+250*1001/2025)*f2(n)/(6435*27.21)
          ENDIF
          !
          ! calculate spin-density 'rho_sig' and total density 'rho_tot'
          !
          rho_tot = 0.0
          DO ispin = 1,jspins
             rho_sig(ispin) = 0.0
             DO m = -l,l
                rho_sig(ispin) = rho_sig(ispin) + REAL(ns_mmp(m,m,n,ispin))
             ENDDO
             rho_tot = rho_tot + rho_sig(ispin)
          ENDDO
          rho_sig(1) = rho_sig(1) * spin_deg  ! if jspins = 1, divide by 2
          IF (atoms%lda_u(itype)%l_amf) THEN
             eta(1) = rho_sig(1) / (2*l + 1) 
             eta(jspins) = rho_sig(jspins) / (2*l + 1) 
             eta(0) = (eta(1) + eta(jspins) ) / 2
          ELSE
             eta(0) = 1.0
             eta(1) = 1.0
             eta(jspins) = 1.0
          ENDIF
          !
          !--------------------------------------------------------------------------------------------+
          !  s       --                                        s'                    1        s   1    |
          ! V     =  >  ( <m,p|V|m',q> - <m,p|V|q,m'> d     ) n     + d    ( -U (n - -) + J (n  - -) ) |
          !  m,m'    --                                s,s'    p,q     m,m'          2            2    |
          !        p,q,s'                                                                              |
          !--------------------------------------------------------------------------------------------+     
          ! initialise vs_mmp
          !
#ifdef CPP_INVERSION
          vs_mmp(:,:,n,:) = ns_mmp(:,:,n,:)
          DO ispin = 1,jspins
             DO m = -l,l
                DO mp = -l,l
                   ns_mmp(m,mp,n,ispin) = vs_mmp(-m,-mp,n,ispin)
                ENDDO
             ENDDO
          ENDDO
#endif
          vs_mmp(:,:,n,:) = CMPLX(0.0,0.0)
          !
          ! outer spin loop - set up v_mmp
          !
          DO ispin = 1,jspins
             DO m = -l,l
                DO mp =-l,l 

                   DO jspin = 1,jspins
                      IF (ispin.EQ.jspin) THEN
                         DO p = -l,l
                            DO q = -l,l
                               vs_mmp(m,mp,n,ispin) = vs_mmp(m,mp,n,ispin) +  &
                                    ns_mmp(p, q,n,jspin) * ( u(m,p,mp,q,n) - u(m,p,q,mp,n) ) 
                            ENDDO
                         ENDDO
                      ENDIF
                      IF ((ispin.NE.jspin).OR.(jspins.EQ.1)) THEN
                         DO p = -l,l
                            DO q = -l,l
                               vs_mmp(m,mp,n,ispin) = vs_mmp(m,mp,n,ispin) +  &
                                    u(m,p,mp,q,n) * ns_mmp(p, q,n,jspin) 
                            ENDDO
                         ENDDO
                      ENDIF
                   ENDDO

                ENDDO ! m' loop
             ENDDO   ! m  loop
          ENDDO      ! outer spin loop
          !
          !  set diagonal terms and correct for non-spin-polarised case
          !
          DO ispin = 1,jspins
             v_diag(ispin) = - u_htr * ( rho_tot - 0.5*eta(0) ) + j_htr * ( rho_sig(ispin) - 0.5*eta(ispin) )
             DO m = -l,l
                DO mp = -l,l
                   vs_mmp(m,mp,n,ispin) = vs_mmp(m,mp,n,ispin) * spin_deg
                ENDDO
                vs_mmp(m,m,n,ispin) = vs_mmp(m,m,n,ispin) + v_diag(ispin)
             ENDDO
          ENDDO

          !----------------------------------------------------------------------+
          !              s                                                       !
          !  ee      1  ---   s        s                     1        s  1       !
          ! E  (n) = -  >    n      ( V     + d     ( U (n - -) - J (n - -) ))   !
          !          2  ---   m,m'     m,m'    m,m'          2           2       !
          !             m,m'                                                     !
          !----------------------------------------------------------------------+

          e_ee = 0.0
          DO ispin = 1,jspins
             DO m = -l,l
                DO mp =-l,l
                   e_ee=e_ee+REAL(vs_mmp(m,mp,n,ispin)*ns_mmp(m,mp,n,ispin))
                ENDDO
                e_ee = e_ee - v_diag(ispin) * REAL( ns_mmp(m,m,n,ispin) )
             ENDDO
          ENDDO

          !----------------------------------------------------------------------+
          !   dc       ee      U           J  --   s   s       1                 |
          !  E      = E  (n) - - n (n-1) + -  >   n  (n -1)  - - (U-J) n         |
          !   LDA+U            2           2  --               2                 |
          !                                    s                                 |
          !----------------------------------------------------------------------+

          ns_sum = 0.0
          DO ispin = 1,jspins
             ns_sum = ns_sum + rho_sig(ispin) * &
                  ( rho_sig(ispin) - eta(ispin) )
          ENDDO
          e_dc = u_htr * rho_tot * ( rho_tot - eta(0) ) - j_htr * ns_sum 
          e_dcc = (u_htr - j_htr) * rho_tot

          ns_sum = ns_sum / spin_deg
          !       e_ldau = e_ldau + (e_ee -  u_htr * rho_tot * ( rho_tot - 1. ) 
          !    +    + j_htr * ns_sum  - (u_htr - j_htr) * rho_tot) * neq(itype)
          !       write(*,*) e_ldau
          results%e_ldau = results%e_ldau + ( e_ee - e_dc - e_dcc) * atoms%neq(itype)
          !       write(*,*) e_ldau

       ENDIF
    ENDDO ! loop over atoms

    results%e_ldau = results%e_ldau / 2

  END SUBROUTINE v_mmp
END MODULE m_vmmp