The Elk Code
rdmenergy.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2005-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: rdmenergy
8 ! !INTERFACE:
9 subroutine rdmenergy
10 ! !USES:
11 use modmain
12 use modrdm
13 use modtest
14 ! !DESCRIPTION:
15 ! Calculates RDMFT total energy (free energy for finite temperatures).
16 !
17 ! !REVISION HISTORY:
18 ! Created 2008 (Sharma)
19 ! Updated for free energy 2009 (Baldsiefen)
20 !EOP
21 !BOC
22 implicit none
23 ! local variables
24 integer ik,ist,is,ias
25 integer nr,nri,ir,i
26 real(8) wo
27 complex(8) z1
28 ! automatic arrays
29 real(8) rfmt(npmtmax)
30 ! allocatable arrays
31 complex(8), allocatable :: evecsv(:,:)
32 ! external functions
33 real(8), external :: rfmtinp
34 complex(8), external :: zdotc
35 ! Coulomb energy from core states
36 engyvcl=0.d0
37 do ias=1,natmtot
38  is=idxis(ias)
39  nr=nrmt(is)
40  nri=nrmti(is)
41  rfmt(1:npmt(is))=0.d0
42  i=1
43  if (spincore) then
44 ! spin-polarised core
45  do ir=1,nri
46  rfmt(i)=rhocr(ir,ias,1)+rhocr(ir,ias,2)
47  i=i+lmmaxi
48  end do
49  do ir=nri+1,nr
50  rfmt(i)=rhocr(ir,ias,1)+rhocr(ir,ias,2)
51  i=i+lmmaxo
52  end do
53  else
54 ! spin-unpolarised core
55  do ir=1,nri
56  rfmt(i)=rhocr(ir,ias,1)
57  i=i+lmmaxi
58  end do
59  do ir=nri+1,nr
60  rfmt(i)=rhocr(ir,ias,1)
61  i=i+lmmaxo
62  end do
63  end if
64  engyvcl=engyvcl+rfmtinp(nr,nri,wr2mt(:,is),rfmt,vclmt(:,ias))
65 end do
67 allocate(evecsv(nstsv,nstsv))
68 do ik=1,nkpt
69  call getevecsv(filext,ik,vkl(:,ik),evecsv)
70  do ist=1,nstsv
71  wo=wkpt(ik)*occsv(ist,ik)
72 ! Coulomb energy from valence states
73  engyvcl=engyvcl+wo*dble(vclmat(ist,ist,ik))
74 ! kinetic energy from valence states
75  z1=zdotc(nstsv,evecsv(:,ist),1,dkdc(:,ist,ik),1)
76  engykn=engykn+wo*z1%re
77  end do
78 end do
79 deallocate(evecsv)
80 ! Madelung term
81 engymad=0.d0
82 do ias=1,natmtot
83  is=idxis(ias)
84  engymad=engymad+0.5d0*spzn(is)*(vclmt(1,ias)-vcln(1,is))*y00
85 end do
86 ! exchange-correlation energy
87 call rdmengyxc
88 ! total energy
90 if (rdmtemp > 0.d0) then
91  call rdmentropy
93 end if
94 ! write total energy to test file
95 call writetest(300,'RDMFT total energy',tol=1.d-6,rv=engytot)
96 end subroutine
97 !EOC
98 
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
real(8) engyx
Definition: modmain.f90:975
integer lmmaxo
Definition: modmain.f90:203
integer nkpt
Definition: modmain.f90:461
real(8), dimension(:,:), allocatable vcln
Definition: modmain.f90:97
real(8) engykn
Definition: modmain.f90:953
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
real(8), dimension(:,:,:), allocatable rhocr
Definition: modmain.f90:941
subroutine rdmengyxc
Definition: rdmengyxc.f90:10
complex(8), dimension(:,:,:), allocatable dkdc
Definition: modrdm.f90:15
real(8) rdmentrpy
Definition: modrdm.f90:33
real(8) rdmtemp
Definition: modrdm.f90:31
Definition: modrdm.f90:6
integer nstsv
Definition: modmain.f90:889
real(8), dimension(:), allocatable wkpt
Definition: modmain.f90:475
real(8) engytot
Definition: modmain.f90:983
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
real(8) engyvcl
Definition: modmain.f90:965
subroutine rdmentropy
Definition: rdmentropy.f90:10
real(8) engykncr
Definition: modmain.f90:955
real(8) engymad
Definition: modmain.f90:967
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
integer natmtot
Definition: modmain.f90:40
complex(8), dimension(:,:,:), allocatable vclmat
Definition: modrdm.f90:13
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
subroutine rdmenergy
Definition: rdmenergy.f90:10
real(8), parameter y00
Definition: modmain.f90:1236
logical spincore
Definition: modmain.f90:943
real(8), dimension(:,:), allocatable wr2mt
Definition: modmain.f90:183
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150