The Elk Code
k_tfvw1.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2021 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: k_tfvw1
8 ! !INTERFACE:
9 elemental subroutine k_tfvw1(rho,grho2,dtdr,dtdgr2)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! rho : spin-unpolarised charge density (in,real)
12 ! grho2 : |grad rho|^2 (in,real)
13 ! dtdr : dtau/drho (out,real)
14 ! dtdgr2 : dtau/d(|grad rho|^2) (out,real)
15 ! !DESCRIPTION:
16 ! Calculates the derivatives $\partial\tau/\partial\rho$ and
17 ! $\partial\tau/\partial|\nabla\rho|^2$ of the gradient expansion of
18 ! the kinetic energy density $\tau$. This includes the Thomas-Fermi and
19 ! von Weizsacker terms:
20 ! $$ \tau=\frac{3}{10}(3\pi^2)^{2/3}\rho^{5/3}
21 ! +\frac{1}{72}\frac{|\nabla\rho|^2}{\rho}. $$
22 !
23 ! !REVISION HISTORY:
24 ! Created December 2021 (JKD)
25 !EOP
26 !BOC
27 implicit none
28 ! arguments
29 real(8), intent(in) :: rho,grho2
30 real(8), intent(out) :: dtdr,dtdgr2
31 ! local variables
32 real(8), parameter :: pi=3.1415926535897932385d0
33 ! Thomas-Fermi coefficient
34 real(8), parameter :: ctf=(3.d0/10.d0)*(3.d0*pi**2)**(2.d0/3.d0)
35 ! von Weizsacker coefficient
36 real(8), parameter :: cvw=1.d0/72.d0
37 real(8) ri,t1,t2
38 if ((rho < 1.d-20).or.(grho2 < 0.d0)) then
39  dtdr=0.d0
40  dtdgr2=0.d0
41  return
42 end if
43 ri=1.d0/rho
44 t1=ctf*(5.d0/3.d0)*rho**(2.d0/3.d0)
45 t2=cvw*ri
46 dtdr=t1-t2*grho2*ri
47 dtdgr2=t2
48 end subroutine
49 !EOC
50 
elemental subroutine k_tfvw1(rho, grho2, dtdr, dtdgr2)
Definition: k_tfvw1.f90:10