The Elk Code
 
Loading...
Searching...
No Matches
energykncr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6subroutine energykncr
7use modmain
8implicit none
9integer ist,is,ias
10integer nr,nri,ir,i
11! automatic arrays
12real(8) rfmt(npmtmax)
13! external functions
14real(8), external :: rfmtinp
15! calculate the kinetic energy for core states
16engykncr=0.d0
17do ias=1,natmtot
18 is=idxis(ias)
19 nr=nrmt(is)
20 nri=nrmti(is)
21! sum of core eigenvalues
22 do ist=1,nstsp(is)
23 if (spcore(ist,is)) engykncr=engykncr+occcr(ist,ias)*evalcr(ist,ias)
24 end do
25! core density
26 rfmt(1:npmt(is))=0.d0
27 i=1
28 if (spincore) then
29! spin-polarised core
30 do ir=1,nri
31 rfmt(i)=rhocr(ir,ias,1)+rhocr(ir,ias,2)
32 i=i+lmmaxi
33 end do
34 do ir=nri+1,nr
35 rfmt(i)=rhocr(ir,ias,1)+rhocr(ir,ias,2)
36 i=i+lmmaxo
37 end do
38 else
39! spin-unpolarised core
40 do ir=1,nri
41 rfmt(i)=rhocr(ir,ias,1)
42 i=i+lmmaxi
43 end do
44 do ir=nri+1,nr
45 rfmt(i)=rhocr(ir,ias,1)
46 i=i+lmmaxo
47 end do
48 end if
49 engykncr=engykncr-rfmtinp(nr,nri,wr2mt(:,is),rfmt,vsmt(:,ias))
50end do
51end subroutine
52
subroutine energykncr
Definition energykncr.f90:7
real(8), dimension(:,:), allocatable evalcr
Definition modmain.f90:934
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
real(8) engykncr
Definition modmain.f90:952
integer lmmaxi
Definition modmain.f90:207
logical spincore
Definition modmain.f90:940
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:), pointer, contiguous vsmt
Definition modmain.f90:649
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:), allocatable rhocr
Definition modmain.f90:938
real(8), dimension(:,:), allocatable occcr
Definition modmain.f90:932
pure real(8) function rfmtinp(nr, nri, wr, rfmt1, rfmt2)
Definition rfmtinp.f90:10