TABLE OF CONTENTS


ABINIT/xcxalp [ Functions ]

[ Top ] [ Functions ]

NAME

 xcxalp

FUNCTION

 Returns exc, vxc, and eventually d(vxc)/d($\rho$) from input $\rho$.
 "X$\alpha$" method is used in this subroutine:
 a single fixed value is chosen for "alpha", set below.
 Expression is exc=-alpha*efac/rs (hartree), efac below.
 rs = $(3/(4\pi))^{1/3}* \rho (r)^{-1/3}$.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  npt=number of real space points on which density is provided
  order=gives the maximal derivative of Exc computed.
  rspts(npt)=Wigner-Seitz radii, at each point

OUTPUT

  exc(npt)=exchange-correlation energy density (hartree)
  vxc(npt)=xc potential (d($\rho$*exc)/d($\rho$)) (hartree)
  if(order>1) dvxc(npt)=derivative d(vxc)/d($\rho$) (hartree*bohr^3)

PARENTS

      drivexc

CHILDREN

SOURCE

 37 #if defined HAVE_CONFIG_H
 38 #include "config.h"
 39 #endif
 40 
 41 #include "abi_common.h"
 42 
 43 
 44 subroutine xcxalp(exc,npt,order,rspts,vxc, dvxc)  ! dvxc is optional
 45 
 46  use defs_basis
 47  use m_errors
 48  use m_profiling_abi
 49 
 50 !This section has been created automatically by the script Abilint (TD).
 51 !Do not modify the following lines by hand.
 52 #undef ABI_FUNC
 53 #define ABI_FUNC 'xcxalp'
 54 !End of the abilint section
 55 
 56  implicit none
 57 
 58 !Arguments ------------------------------------
 59 !scalars
 60  integer,intent(in) :: npt,order
 61 !arrays
 62  real(dp),intent(in) :: rspts(npt)
 63  real(dp),intent(out) :: exc(npt),vxc(npt)
 64  real(dp),intent(out),optional :: dvxc(npt)
 65 
 66 !Local variables-------------------------------
 67 !Set value of alpha in "X-alpha" method
 68 !scalars
 69  integer :: ipt
 70  real(dp),parameter :: alpha=1.0_dp
 71  real(dp) :: dfac,efac,rs,rsm1,vfac
 72  character(len=500) :: message
 73 
 74 ! *************************************************************************
 75 
 76 !Checks the values of order
 77  if(order<0 .or. order>2)then
 78    write(message, '(a,a,a,i3)' )&
 79 &   'With X-alpha xc functional, the only',ch10,&
 80 &   'allowed values for order are 0, 1 or 2, while it is found to be',order
 81    MSG_BUG(message)
 82  end if
 83 
 84 !Compute vfac=(3/(2*Pi))^(2/3)
 85  vfac=(1.5_dp/pi)**(2.0_dp/3.0_dp)
 86 !Compute efac=(3/4)*vfac
 87  efac=0.75_dp*vfac
 88 !Compute dfac=(4*Pi/9)*vfac
 89  dfac=(4.0_dp*pi/9.0_dp)*vfac
 90 
 91 !separate cases with respect to order
 92  if(order==2) then
 93 !  Loop over grid points
 94    do ipt=1,npt
 95      rs=rspts(ipt)
 96      rsm1=1.0_dp/rs
 97 !    compute energy density (hartree)
 98      exc(ipt)=-alpha*efac*rsm1
 99 !    compute potential (hartree)
100      vxc(ipt)=-alpha*vfac*rsm1
101 !    compute d(vxc)/d(rho) (hartree*bohr^3)
102      dvxc(ipt)=-alpha*dfac*rs**2
103    end do
104  else
105 !  Loop over grid points
106    do ipt=1,npt
107      rs=rspts(ipt)
108      rsm1=1.0_dp/rs
109 !    compute energy density (hartree)
110      exc(ipt)=-alpha*efac*rsm1
111 !    compute potential (hartree)
112      vxc(ipt)=-alpha*vfac*rsm1
113    end do
114  end if
115 !
116 end subroutine xcxalp