~siesta-pseudos-bases/siesta/trunk-psml

« back to all changes in this revision

Viewing changes to Src/m_energies.F90

  • Committer: Alberto Garcia
  • Date: 2019-09-02 14:09:43 UTC
  • mfrom: (427.6.323 trunk)
  • Revision ID: albertog@icmab.es-20190902140943-mzmbe1jacgefpgxw
Sync to trunk-776 (notably nc/soc wavefunction support)


Show diffs side-by-side

added added

removed removed

Lines of Context:
47
47
  real(dp):: Eldau      
48
48
  real(dp):: DEldau
49
49
 
50
 
  real(dp) :: DE_NEGF  ! NEGF total energy contribution = - e * \sum_i N_i \mu_i
 
50
  real(dp) :: NEGF_DE  ! NEGF total energy contribution = - e * \sum_i N_i \mu_i
 
51
  ! Generally we should only calculate energies in the regions where we are updating elements
 
52
  ! As such we need a specific energy for the NEGF part
 
53
  ! Their meaning is directly transferable to the above listed energies (in the updating regions)
 
54
  real(dp) :: NEGF_Ebs
 
55
  real(dp) :: NEGF_Ekin
 
56
  real(dp) :: NEGF_Enl
 
57
  real(dp) :: NEGF_DEharr
 
58
  real(dp) :: NEGF_Eharrs
 
59
  real(dp) :: NEGF_Etot
 
60
  real(dp) :: NEGF_FreeE
51
61
 
52
62
contains
53
63
 
85
95
    Eso = 0._dp
86
96
    Eldau = 0._dp      
87
97
    DEldau = 0._dp
88
 
    DE_NEGF = 0._dp
 
98
 
 
99
    ! NEGF part
 
100
    NEGF_DE = 0._dp
 
101
    NEGF_Ebs = 0._dp
 
102
    NEGF_Ekin = 0._dp
 
103
    NEGF_Enl = 0._dp
 
104
    NEGF_DEharr = 0._dp
 
105
    NEGF_Eharrs = 0._dp
 
106
    NEGF_Etot = 0._dp
 
107
    NEGF_FreeE = 0._dp
89
108
 
90
109
  end subroutine init_Energies
91
110
 
102
121
  subroutine update_E0()
103
122
 
104
123
    E0 = Ena + Ekin + Enl + Eso - Eions
105
 
    
 
124
 
106
125
  end subroutine update_E0
107
126
  
108
127
  subroutine update_Etot()
 
128
    use m_ts_global_vars, only: TSrun
109
129
    
110
130
    ! DUext (external electric field) -- should it be in or out?
111
131
    Etot = Ena + Ekin + Enl + Eso - Eions + &
112
 
         DEna + DUscf + DUext + Exc + &
113
 
         Ecorrec + Emad + Emm + Emeta + Eldau
114
 
    ! Commented out the NEGF contribution to the total energy
115
 
    ! We know it is wrong, but we estimate it. See output
116
 
    !    Etot = Etot + DE_NEGF
 
132
        DEna + DUscf + DUext + Exc + &
 
133
        Ecorrec + Emad + Emm + Emeta + Eldau
117
134
 
 
135
    if ( TSrun ) then
 
136
      NEGF_Etot = Ena + NEGF_Ekin + NEGF_Enl - Eions + &
 
137
          DEna + DUscf + DUext + Exc + Ecorrec + Emad + Emm + Emeta + &
 
138
          Eldau + NEGF_DE
 
139
    end if
 
140
    
118
141
  end subroutine update_Etot
119
142
 
120
143
  !> @param kBT the temperature in energy
121
144
  subroutine update_FreeE( kBT )
 
145
    use m_ts_global_vars, only: TSrun
122
146
    real(dp), intent(in) :: kBT
123
147
 
124
148
    FreeE = Etot - kBT * Entropy
125
149
 
 
150
    if ( TSrun ) then
 
151
      NEGF_FreeE = NEGF_Etot - kBT * Entropy
 
152
    end if
 
153
 
126
154
  end subroutine update_FreeE
127
155
 
128
156
  !> @param kBT the temperature in energy