File: rdmenergy.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (93 lines) | stat: -rw-r--r-- 2,092 bytes parent folder | download
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

! Copyright (C) 2005-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: rdmenergy
! !INTERFACE:
subroutine rdmenergy
! !USES:
use modmain
use modrdm
use modmpi
use modtest
! !DESCRIPTION:
!   Calculates RDMFT total energy (free energy for finite temperatures).
!
! !REVISION HISTORY:
!   Created 2008 (Sharma)
!   Updated for free energy 2009 (Baldsiefen)
!EOP
!BOC
implicit none
! local variables
integer ik,ist,is,ias
integer nr,nri,ir,i
real(8) t1
complex(8) z1
! allocatable arrays
real(8), allocatable :: rfmt(:)
complex(8), allocatable :: evecsv(:,:)
! external functions
real(8) rfmtinp
complex(8) zdotc
external rfmtinp,zdotc
allocate(rfmt(npmtmax))
! Coulomb energy from core states
engyvcl=0.d0
do ias=1,natmtot
  is=idxis(ias)
  nr=nrmt(is)
  nri=nrmti(is)
  rfmt(1:npmt(is))=0.d0
  i=1
  do ir=1,nr
    if (spincore) then
      rfmt(i)=rhocr(ir,ias,1)+rhocr(ir,ias,2)
    else
      rfmt(i)=rhocr(ir,ias,1)
    end if
    if (ir.le.nri) then
      i=i+lmmaxi
    else
      i=i+lmmaxo
    end if
  end do
  engyvcl=engyvcl+rfmtinp(nr,nri,rsp(:,is),r2sp(:,is),rfmt,vclmt(:,ias))
end do
deallocate(rfmt)
engykn=engykncr
allocate(evecsv(nstsv,nstsv))
do ik=1,nkpt
  call getevecsv(filext,ik,vkl(:,ik),evecsv)
  do ist=1,nstsv
    t1=wkpt(ik)*occsv(ist,ik)
! Coulomb energy from valence states
    engyvcl=engyvcl+t1*dble(vclmat(ist,ist,ik))
! kinetic energy from valence states
    z1=zdotc(nstsv,evecsv(:,ist),1,dkdc(:,ist,ik),1)
    engykn=engykn+t1*dble(z1)
  end do
end do
deallocate(evecsv)
! Madelung term
engymad=0.d0
do ias=1,natmtot
  is=idxis(ias)
  engymad=engymad+0.5d0*spzn(is)*(vclmt(1,ias)-vcln(1,is))*y00
end do
! exchange-correlation energy
call rdmengyxc
! total energy
engytot=0.5d0*engyvcl+engymad+engykn+engyx
if (rdmtemp.gt.0.d0) then
  call rdmentropy
  engytot=engytot-rdmtemp*rdmentrpy
end if
! write total energy to test file
call writetest(300,'RDMFT total energy',tol=1.d-6,rv=engytot)
return
end subroutine
!EOC