The Elk Code
 
Loading...
Searching...
No Matches
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:
9subroutine rdmenergy
10! !USES:
11use modmain
12use modrdm
13use 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
22implicit none
23! local variables
24integer ik,ist,is,ias
25integer nr,nri,ir,i
26real(8) wo
27complex(8) z1
28! automatic arrays
29real(8) rfmt(npmtmax)
30! allocatable arrays
31complex(8), allocatable :: evecsv(:,:)
32! external functions
33real(8), external :: rfmtinp
34complex(8), external :: zdotc
35! Coulomb energy from core states
36engyvcl=0.d0
37do 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))
65end do
67allocate(evecsv(nstsv,nstsv))
68do 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
78end do
79deallocate(evecsv)
80! Madelung term
81engymad=0.d0
82do ias=1,natmtot
83 is=idxis(ias)
84 engymad=engymad+0.5d0*spzn(is)*(vclmt(1,ias)-vcln(1,is))*y00
85end do
86! exchange-correlation energy
87call rdmengyxc
88! total energy
90if (rdmtemp > 0.d0) then
91 call rdmentropy
93end if
94! write total energy to test file
95call writetest(300,'RDMFT total energy',tol=1.d-6,rv=engytot)
96end subroutine
97!EOC
98
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00
Definition modmain.f90:1233
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
character(256) filext
Definition modmain.f90:1300
real(8) engykn
Definition modmain.f90:950
real(8) engykncr
Definition modmain.f90:952
integer lmmaxi
Definition modmain.f90:207
logical spincore
Definition modmain.f90:940
real(8) engyvcl
Definition modmain.f90:962
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
real(8) engyx
Definition modmain.f90:972
integer nstsv
Definition modmain.f90:886
real(8) engymad
Definition modmain.f90:964
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable vclmt
Definition modmain.f90:624
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
real(8) engytot
Definition modmain.f90:980
real(8), dimension(:,:), allocatable vcln
Definition modmain.f90:97
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:), allocatable rhocr
Definition modmain.f90:938
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
real(8), dimension(maxspecies) spzn
Definition modmain.f90:80
complex(8), dimension(:,:,:), allocatable vclmat
Definition modrdm.f90:13
complex(8), dimension(:,:,:), allocatable dkdc
Definition modrdm.f90:15
real(8) rdmtemp
Definition modrdm.f90:31
real(8) rdmentrpy
Definition modrdm.f90:33
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine rdmenergy
Definition rdmenergy.f90:10
subroutine rdmengyxc
Definition rdmengyxc.f90:10
subroutine rdmentropy
pure real(8) function rfmtinp(nr, nri, wr, rfmt1, rfmt2)
Definition rfmtinp.f90:10