The Elk Code
energynn.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2006 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine energynn
7 use modmain
8 implicit none
9 ! local variables
10 integer is,ias,nr,nri
11 integer iro,i0,i1
12 real(8) t1
13 ! allocatable arrays
14 complex(8), allocatable :: zvclmt(:,:),zvclir(:),zrhoir(:)
15 allocate(zvclmt(npmtmax,natmtot))
16 ! generate the nuclear monopole potentials
17 do ias=1,natmtot
18  is=idxis(ias)
19  nr=nrmt(is)
20  nri=nrmti(is)
21  iro=nri+1
22  i1=lmmaxi*(nri-1)+1
23  zvclmt(1:npmt(is),ias)=0.d0
24  zvclmt(1:i1:lmmaxi,ias)=vcln(1:nri,is)
25  i0=i1+lmmaxi
26  i1=lmmaxo*(nr-iro)+i0
27  zvclmt(i0:i1:lmmaxo,ias)=vcln(iro:nr,is)
28 end do
29 allocate(zrhoir(ngtot),zvclir(ngtot))
30 ! set the interstitial density to zero
31 zrhoir(1:ngtot)=0.d0
32 ! solve the complex Poisson's equation
34  jlgrmt,ylmg,sfacg,zrhoir,npmtmax,zvclmt,zvclir)
35 ! compute the nuclear-nuclear energy
36 engynn=0.d0
37 do ias=1,natmtot
38  is=idxis(ias)
39  t1=(dble(zvclmt(1,ias))-vcln(1,is))*y00
40  engynn=engynn+spzn(is)*t1
41 end do
42 engynn=0.5d0*engynn
43 deallocate(zvclmt,zvclir,zrhoir)
44 end subroutine
45 
complex(8), dimension(:,:), allocatable sfacg
Definition: modmain.f90:430
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer ngtot
Definition: modmain.f90:390
integer lmmaxo
Definition: modmain.f90:203
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, zrhoir, ld3, zvclmt, zvclir)
Definition: zpotcoul.f90:11
real(8), dimension(:,:), allocatable vcln
Definition: modmain.f90:97
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
real(8) engynn
Definition: modmain.f90:957
complex(8), dimension(:,:), allocatable ylmg
Definition: modmain.f90:428
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine energynn
Definition: energynn.f90:7
integer ngvec
Definition: modmain.f90:396
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(:,:,:), allocatable jlgrmt
Definition: modmain.f90:426
real(8), dimension(:), allocatable gclg
Definition: modmain.f90:424
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
integer npmtmax
Definition: modmain.f90:216
integer natmtot
Definition: modmain.f90:40
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), parameter y00
Definition: modmain.f90:1236
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150