The Elk Code
k_tfvw_sp.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_tfvw_sp
8 ! !INTERFACE:
9 subroutine k_tfvw_sp(n,rhoup,rhodn,gup2,gdn2,dtdru,dtdrd,dtdgu2,dtdgd2)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! n : number of density points (in,integer)
12 ! rhoup : spin-up charge density (in,real(n))
13 ! rhodn : spin-down charge density (in,real(n))
14 ! gup2 : |grad rhoup|^2 (in,real(n))
15 ! gdn2 : |grad rhodn|^2 (in,real(n))
16 ! dtdru : dtauup/drhoup (out,real(n))
17 ! dtdrd : dtaudn/drhodn (out,real(n))
18 ! dtdgu2 : dtauup/d(|grad rhoup|^2) (out,real(n))
19 ! dtdgu2 : dtaudn/d(|grad rhodn|^2) (out,real(n))
20 ! !DESCRIPTION:
21 ! Calculates the derivatives of the spin-polarised kinetic energy density
22 ! $\partial\tau^{\uparrow}/\partial\rho^{\uparrow}$,
23 ! $\partial\tau^{\downarrow}/\partial\rho^{\downarrow}$,
24 ! $\partial\tau^{\uparrow}/\partial|\nabla\rho^{\uparrow}|^2$ and
25 ! $\partial\tau^{\downarrow}/\partial|\nabla\rho^{\downarrow}|^2$.
26 ! This is done by noting the relation for the kinetic energy functional
27 ! [G. L. Oliver and J. P. Perdew, {\it Phys. Rev. A}
28 ! {\bf 20}, 397 (1979)]
29 ! $$ T[\rho^{\uparrow},\rho^{\downarrow}]=\tfrac{1}{2}T[2\rho^{\uparrow}]
30 ! +\tfrac{1}{2}T[2\rho^{\downarrow}] $$
31 ! and taking, for example,
32 ! $$ \tau^{\uparrow}(\rho^{\uparrow},|\nabla\rho^{\uparrow}|^2)
33 ! =\tfrac{1}{2}\tau(2\rho^{\uparrow},4|\nabla\rho^{\uparrow}|^2), $$
34 ! where the gradient expansion of the unpolarised kinetic energy density is
35 ! used for $\tau$. See the routines {\tt k\_tfvw1}, {\tt ggamt\_4},
36 ! {\tt ggair\_4}, {\tt potxcmt}, and {\tt potxcir}.
37 !
38 ! !REVISION HISTORY:
39 ! Created December 2021 (JKD)
40 !EOP
41 !BOC
42 implicit none
43 ! arguments
44 integer, intent(in) :: n
45 real(8), intent(in) :: rhoup(n),rhodn(n)
46 real(8), intent(in) :: gup2(n),gdn2(n)
47 real(8), intent(out) :: dtdru(n),dtdrd(n)
48 real(8), intent(out) :: dtdgu2(n),dtdgd2(n)
49 ! local variables
50 integer i
51 do i=1,n
52  call k_tfvw1(2.d0*rhoup(i),4.d0*gup2(i),dtdru(i),dtdgu2(i))
53  dtdgu2(i)=2.d0*dtdgu2(i)
54 end do
55 do i=1,n
56  call k_tfvw1(2.d0*rhodn(i),4.d0*gdn2(i),dtdrd(i),dtdgd2(i))
57  dtdgd2(i)=2.d0*dtdgd2(i)
58 end do
59 end subroutine
60 !EOC
61 
elemental subroutine k_tfvw1(rho, grho2, dtdr, dtdgr2)
Definition: k_tfvw1.f90:10
subroutine k_tfvw_sp(n, rhoup, rhodn, gup2, gdn2, dtdru, dtdrd, dtdgu2, dtdgd2)
Definition: k_tfvw_sp.f90:10