The Elk Code
Loading...
Searching...
No Matches
dengyeph.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2020 Chung-Yu Wang, J. K. Dewhurst, S. Sharma and
3
! E. K. U. Gross. This file is distributed under the terms of the GNU General
4
! Public License. See the file COPYING for license details.
5
6
subroutine
dengyeph
7
use
modmain
8
use
modphonon
9
use
modbog
10
implicit none
11
! local variables
12
integer
ik,iq,i
13
real
(8) w
14
! change in electron energy per unit cell
15
dengye
=0.d0
16
do
ik=1,
nkpt
17
w=
wkpt
(ik)
18
do
i=1,
nstsv
19
dengye
=
dengye
+w*abs((
evalsv
(i,ik)-
efermi
)*dble(
dvv
(i,i,ik)))
20
end do
21
end do
22
dengye
=abs(
occmax
*
dengye
)
23
! change in phonon energy per unit cell
24
dengyph
=0.d0
25
do
iq=1,
nqpt
26
w=
wqpt
(iq)
27
do
i=1,
nbph
28
dengyph
=
dengyph
+w*abs(
wphq
(i,iq)*dble(
dxx
(i,i,iq)))
29
end do
30
end do
31
! sum of both
32
dengy
=
dengye
+
dengyph
33
end subroutine
34
dengyeph
subroutine dengyeph
Definition
dengyeph.f90:7
modbog
Definition
modbog.f90:6
modbog::dvv
complex(8), dimension(:,:,:), pointer, contiguous dvv
Definition
modbog.f90:25
modbog::dxx
complex(8), dimension(:,:,:), pointer, contiguous dxx
Definition
modbog.f90:43
modmain
Definition
modmain.f90:6
modmain::wkpt
real(8), dimension(:), allocatable wkpt
Definition
modmain.f90:475
modmain::efermi
real(8) efermi
Definition
modmain.f90:904
modmain::nqpt
integer nqpt
Definition
modmain.f90:525
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::wqpt
real(8), dimension(:), allocatable wqpt
Definition
modmain.f90:549
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition
modmain.f90:918
modphonon
Definition
modphonon.f90:6
modphonon::nbph
integer nbph
Definition
modphonon.f90:13
modphonon::dengye
real(8) dengye
Definition
modphonon.f90:140
modphonon::dengy
real(8) dengy
Definition
modphonon.f90:144
modphonon::dengyph
real(8) dengyph
Definition
modphonon.f90:142
modphonon::wphq
real(8), dimension(:,:), allocatable wphq
Definition
modphonon.f90:31
dengyeph.f90
Generated by
1.9.8