The Elk Code
energytd.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
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 energytd
7 use modmain
8 use modtddft
9 implicit none
10 ! local variables
11 integer is,ias
12 real(8) ca,engya,t1,t2
13 ! external functions
14 real(8), external :: rfinp
15 ! Coulomb potential energy
17 ! Madelung term
18 engymad=0.d0
19 do ias=1,natmtot
20  is=idxis(ias)
21  engymad=engymad+0.5d0*spzn(is)*(vclmt(1,ias)-vcln(1,is))*y00
22 end do
23 ! exchange and correlation energy
24 engyx=rfinp(rhomt,rhoir,exmt,exir)
25 engyc=rfinp(rhomt,rhoir,ecmt,ecir)
26 ! vector potential contributions to energy
27 ca=-1.d0/solsc
28 ! coupling term -1/c A(t)⋅J(t)
29 engya=ca*dot_product(afieldt(1:3,itimes),jtot(1:3))
30 ! constant term 1/2c² A²(t)Q where Q is the 'dynamic' charge
31 ca=0.5d0/solsc**2
32 t1=sum(afieldt(1:3,itimes)**2)
33 t2=sum(chgstot(1:3))/3.d0
34 engya=engya+ca*t1*(chgtot-t2)
35 ! total energy
37 end subroutine
38 
subroutine energytd
Definition: energytd.f90:7
real(8) engyx
Definition: modmain.f90:975
real(8), dimension(3) jtot
Definition: modmain.f90:748
real(8), dimension(:), allocatable ecir
Definition: modmain.f90:632
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
real(8), dimension(:,:), allocatable vcln
Definition: modmain.f90:97
real(8) engykn
Definition: modmain.f90:953
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
real(8), dimension(:,:), allocatable ecmt
Definition: modmain.f90:632
real(8), dimension(:,:), allocatable exmt
Definition: modmain.f90:630
real(8) engytot
Definition: modmain.f90:983
real(8) engyvcl
Definition: modmain.f90:965
real(8), dimension(:), allocatable vclir
Definition: modmain.f90:624
real(8) engyc
Definition: modmain.f90:977
real(8) solsc
Definition: modmain.f90:1253
real(8) engymad
Definition: modmain.f90:967
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
integer itimes
Definition: modtddft.f90:46
real(8) chgtot
Definition: modmain.f90:726
real(8), dimension(3) chgstot
Definition: modtddft.f90:84
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
integer natmtot
Definition: modmain.f90:40
real(8), parameter y00
Definition: modmain.f90:1236
real(8), dimension(:), allocatable exir
Definition: modmain.f90:630