TABLE OF CONTENTS


ABINIT/m_xcpbe [ Modules ]

[ Top ] [ Modules ]

NAME

  m_xcpbe

FUNCTION

 Treat XC functionals closely linked with the Perdew-Wang 92 LSD and the PBE GGA.

COPYRIGHT

  Copyright (C) 1998-2024 ABINIT group (XG,MF,LG,CE)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_xcpbe
23 
24  use defs_basis
25  use m_abicore
26  use m_errors
27 
28  use m_numeric_tools,      only : invcb
29 
30  implicit none
31 
32  private

ABINIT/xcpbe [ Functions ]

[ Top ] [ Functions ]

NAME

 xcpbe

FUNCTION

 Treat XC functionals closely linked with the Perdew-Wang 92 LSD and the PBE GGA.

 For a series of values of the density and, if GGA, the square of the
 gradient of the density, return the associated Exc energy,
 potential, and, in case of response-function, functions needed
 to build the XC kernel.

 If option==2, Exchange-correlation functional from Perdew-Burke-Ernzerhof,
 Phys.Rev.Lett. 77, 3866 (1996) [[cite:Perdew1996]].
 If option==1, Reduces to Perdew-Wang LSD , PRB45,13244 (1992) [[cite:Perdew1992]].
 If option==-1 or -2, take only exchange part of PW (-1) or PBE (-2) functionals.
 If option==-4, C09x exchange functional of V. R. Cooper, PRB 81, 161104(R) (2010) [[cite:Cooper2010]].
 If option==3 or 4, take exchange plus RPA correlation
   part of LSD PW (3) or GGA PBE (4) functionals.
 If option==5, revPBE functional of Zhang and Yang, PRL 80, 890 (1998) [[cite:Zhang1998]]
 If option==6, RPBE functional of Hammer, Hansen and Norskov, PRB 59, 7413 (1999) [[cite:Hammer1999]]
 If option==7, WC functional of Wu and Cohen, PRB 73, 235116 (2006) [[cite:Wu2006]]

INPUTS

  exexch= choice of local exact exchange. Active if exexch=1
  npts= number of points to be computed
  nspden=1 for unpolarized, 2 for spin-polarized
  grho2_updn(npts,2*nspden-1)=square of the gradient of the spin-up,
     and, if nspden==2, spin-down, and total density (Hartree/Bohr**2),
     only used if gradient corrected functional (option=2,-2,-4 and 4 or beyond)
  option= see above
  order=its absolute value gives the maximal derivative of Exc to be computed.
  rho_updn(npts,nspden)=spin-up and spin-down density (Hartree/bohr**3)
  ndvxci= size of dvxci(npts,ndvxci)
  nd2vxci=size of d2vxci(npts,nd2vxci)

OUTPUT

  d2vxci=third derivative of the xc energy with respect to the density, only
    only if local-density approximation
   calculated if order==3
   In case of local energy functional (option=1,-1 or 3):
    d2vxci(npts,nd2vxc)=              (Hartree*bohr^3)
     if(nspden=1): d2vxci(:,1)=-(2/3)*dvxci/d$\rho$
                                  (dvxci is the second derivative, see below)
     if(nspden=2): d2vxci(:,1)=3rd order derivative of XC energy with respect to rhouprhouprhoup,
                   d2vxci(:,2)=3rd order derivative of XC energy with respect to rhouprhouprhodn
                   d2vxci(:,3)=3rd order derivative of XC energy with respect to rhodnrhouprhodn
                   d2vxci(:,4)=3rd order derivative of XC energy with respect to rhodnrhodnrhodn
  dvxcdgr(npts,3)=partial derivative of the exchange-correlation
    energy (exci*$\rho$) with respect to the spin-up (dvxcdgr(:,1)),
    spin-down (dvxcdgr(:,2)), or total spin (dvxcdgr(:,3)) gradients of the density
    divided by the norm of the gradient (the definition changed in v3.3)

  dvxci=partial second derivatives of the xc energy, only if abs(order)>1
   In case of local energy functional (option=1,-1 or 3):
    dvxci(npts,1+nspden)=              (Hartree*bohr^3)
     if(nspden=1 .and. order==2): dvxci(:,1)=dvxc/d$\rho$ , dvxc(:,2) empty
     if(nspden=1 .and. order==-2): also compute dvxci(:,2)=dvxc($\uparrow$)/d$\rho(\downarrow)$
     if(nspden=2): dvxci(:,1)=dvxc($\uparrow$)/d$\rho(\uparrow)$,
                   dvxci(:,2)=dvxc($\uparrow$)/d$\rho(\downarrow)$,
                   dvxci(:,3)=dvxc($\downarrow$)/d$\rho(\downarrow)$
   In case of gradient corrected functional (option=2,-2, 4, 5, 6, -4):
    dvxci(npts,15)=
     dvxci(:,1)= d2Ex/drho_up drho_up
     dvxci(:,2)= d2Ex/drho_dn drho_dn
     dvxci(:,3)= dEx/d(abs(grad(rho_up))) / abs(grad(rho_up))
     dvxci(:,4)= dEx/d(abs(grad(rho_dn))) / abs(grad(rho_dn))
     dvxci(:,5)= d2Ex/d(abs(grad(rho_up))) drho_up / abs(grad(rho_up))
     dvxci(:,6)= d2Ex/d(abs(grad(rho_dn))) drho_dn / abs(grad(rho_dn))
     dvxci(:,7)= 1/abs(grad(rho_up)) * d/drho_up (dEx/d(abs(grad(rho_up))) /abs(grad(rho_up)))
     dvxci(:,8)= 1/abs(grad(rho_dn)) * d/drho_dn (dEx/d(abs(grad(rho_dn))) /abs(grad(rho_dn)))
     dvxci(:,9)= d2Ec/drho_up drho_up
     dvxci(:,10)=d2Ec/drho_up drho_dn
     dvxci(:,11)=d2Ec/drho_dn drho_dn
     dvxci(:,12)=dEc/d(abs(grad(rho))) / abs(grad(rho))
     dvxci(:,13)=d2Ec/d(abs(grad(rho))) drho_up / abs(grad(rho))
     dvxci(:,14)=d2Ec/d(abs(grad(rho))) drho_dn / abs(grad(rho))
     dvxci(:,15)=1/abs(grad(rho)) * d/drho (dEc/d(abs(grad(rho))) /abs(grad(rho)))

  exci(npts)=exchange-correlation energy density (hartree)
  vxci(npts,nspden)=partial derivative of the exchange-correlation energy (exci*$\rho$)
    with respect to the spin-down (vxci(:,1)) and spin-up (vxci(:,2) densities
 Normalization: Exc=$\int (exc(r)*\rho (r) d^3 r)$ for $\rho$(r)=electron density.

TODO

  WARNING: option=4 not yet implemented.

SOURCE

 130 subroutine xcpbe(exci,npts,nspden,option,order,rho_updn,vxci,ndvxci,nd2vxci, & !Mandatory Arguments
 131 &                d2vxci,dvxcdgr,dvxci,exexch,grho2_updn)                       !Optional Arguments
 132 
 133 !Arguments ------------------------------------
 134 !scalars
 135  integer,intent(in) :: ndvxci,nd2vxci,npts,nspden,option,order
 136  integer,intent(in),optional :: exexch
 137 !arrays
 138  real(dp),intent(in) :: rho_updn(npts,nspden)
 139  real(dp),intent(in),optional :: grho2_updn(npts,2*nspden-1)
 140  real(dp),intent(out) :: exci(npts),vxci(npts,nspden)
 141  real(dp),intent(out),optional :: d2vxci(npts,nd2vxci),dvxcdgr(npts,3)
 142  real(dp),intent(out),optional :: dvxci(npts,ndvxci)
 143 
 144 !Local variables-------------------------------
 145 ! The "accurate" value of mu is taken from the PBE code.
 146 ! The value of mu_c09 is taken from the paper (see above) in order to recover
 147 ! the GEA behaviour of the enhancement factor for small values of s rather than
 148 ! the LSD linear response limit used in revPBE.
 149 !scalars
 150  integer,save :: initialized=0
 151  integer :: exexch_,ipts,ispden
 152  real(dp),parameter :: alpha_c09=0.0483_dp
 153 !GMR (to match Libxc one should comment the next line and uncomment the following )
 154  real(dp),parameter :: alpha_zeta2=1.0_dp-1.0e-6_dp,alpha_zeta=1.0_dp-1.0e-6_dp
 155 !real(dp),parameter :: alpha_zeta2=1.0_dp,alpha_zeta=1.0_dp
 156 !GMR
 157 !GMR (to match Libxc one should comment the next line and uncomment the following )
 158  real(dp),parameter :: b_wc=0.123456790123_dp,beta=0.066725_dp
 159 !real(dp),parameter :: b_wc=0.123456790123_dp,beta=0.06672455060314922_dp
 160 !GMR
 161  real(dp),parameter :: beta_inv=1.0_dp/beta,c_wc=0.00793746933516_dp
 162 !GMR (to match Libxc one should comment the next line and uncomment the following two)
 163  real(dp),parameter :: fsec_inv=1.0_dp/1.709921_dp,kappa_pbe=0.804_dp
 164 !real(dp),parameter :: fsec_inv=1.0_dp/1.709920934161365617563962776245_dp
 165 !real(dp),parameter :: kappa_pbe=0.804_dp
 166 !GMR
 167  real(dp),parameter :: kappa_revpbe=1.245_dp,mu=0.2195149727645171_dp
 168  real(dp),parameter :: kappa_c09=1.245_dp, mu_c09=0.0617_dp
 169  real(dp),parameter :: mu_divkappa_pbe=mu/kappa_pbe
 170  real(dp),parameter :: mu_divkappa_revpbe=mu/kappa_revpbe
 171  real(dp),parameter :: rsfac=0.6203504908994000_dp,tolgrad=tol10
 172  real(dp),save :: beta_gamma,coeff_tt,factf_zeta,factfp_zeta,gamma,gamma_inv
 173  real(dp),save :: sixpi2_1_3,sixpi2m1_3,sq_rsfac,sq_rsfac_inv,threefourth_divpi,twom1_3
 174  real(dp) :: aa,arg_rr,bb,cc,coeff_aa,alphs2,alphmu,coeff_qq,coeffss,d2aa_drs2,d2aa_drsdzeta
 175  real(dp) :: d2aa_dzeta2,d2bb_drs2,d2bb_drsdzeta,d2bb_dzeta2,d2cc_dbb2
 176  real(dp) :: d2cc_drs2,d2cc_drsdzeta,d2cc_dzeta2,d2ecrs0_drs2,d2ecrs1_drs2
 177  real(dp) :: d2ecrs_drdn2,d2ecrs_drdndrup,d2ecrs_drho2,d2ecrs_drs2
 178  real(dp) :: d2ecrs_drsdzeta,d2ecrs_drup2,d2ecrs_dzeta2,d2fxdg2,d2fxdn2
 179  real(dp) :: d2fxdndg,d2fxdss2,d2fzeta4_dzeta2,d2gcrs_drs2,d2hh_drs2
 180  real(dp) :: d2hh_drsdtt,d2hh_drsdzeta,d2hh_dtt2,d2hh_dttdzeta,d2hh_dzeta2
 181  real(dp) :: d2macrs_drs2,d2pade_drs2,d2pade_drsdtt,d2pade_drsdzeta,d2pade_dtt2
 182  real(dp) :: d2pade_dttdzeta,d2pade_dxx2,d2pade_dzeta2,d2qq_drs2,d2qq_drsdtt
 183  real(dp) :: d2qq_drsdzeta,d2qq_dtt2,d2qq_dttdzeta,d2qq_dzeta2,d2rhohh_drho2
 184  real(dp) :: d2rhohh_drhodg,d2rr_dqq2,d2rr_drs2,d2rr_drsdtt,d2rr_drsdzeta
 185  real(dp) :: d2rr_dtt2,d2rr_dttdzeta,d2rr_dzeta2,d2rs_dn2,d2ssdn2,d2ssdndg
 186  real(dp) :: d2vcrs_drs2,d2xx_drs2,d2xx_drsdtt,d2xx_drsdzeta,d2xx_dttdzeta
 187  real(dp) :: d2xx_dzeta2,d3ecrs0_drs3,d3ecrs_drup3,d3ecrs_drup2drdn
 188  real(dp) :: d3ecrs_drupdrdn2,d3ecrs_drdn3,d3ecrs_dzeta3
 189  real(dp) :: d3ecrs_drs2dzeta,d3ecrs_dzeta2drs,d3ecrs1_drs3,d3gcrs_drs3
 190  real(dp) :: d3ecrs_drs3,d3macrs_drs3
 191  real(dp) :: d_wc,daa_drs,daa_dzeta,dbb_drs,dbb_dzeta
 192  real(dp) :: dcc_dbb,dcc_drs,dcc_dzeta,decrs0_drs,decrs1_drs,decrs_drs
 193  real(dp) :: decrs_dzeta,dfxdg,dfxdn,dfxdss,dfzeta4_dzeta,dgcrs_drs
 194  real(dp) :: dhh_drs,dhh_dtt,dhh_dzeta,div_rr,divss,dmacrs_drs,dpade_drs
 195  real(dp) :: dpade_dtt,dpade_dxx,dpade_dzeta,dqq_drs,dqq_dtt,dqq_dzeta
 196  real(dp) :: drhohh_drho,drr_dqq,drr_drs,drr_dtt,drr_dzeta,drs_dn,dssdg,dssdn
 197  real(dp) :: dtt_dg,dvcrs_drs,dxx_drs,dxx_dtt,dxx_dzeta,ec0_a1,ec0_aa,ec0_b1
 198  real(dp) :: ec0_b2,ec0_b3,ec0_b4,ec0_den,ec0_f1,ec0_f2,ec0_log,ec0_q0,ec0_q1
 199  real(dp) :: ec0_q1p,ec0_q1pp,ec0_q1ppp,ec1_a1,ec1_aa,ec1_b1,ec1_b2,ec1_b3
 200  real(dp) :: ec1_b4,ec1_den,ec1_f1,ec1_f2,ec1_log,ec1_q0,ec1_q1,ec1_q1p,ec1_q1pp
 201  real(dp) :: ec1_q1ppp,ecrs,ecrs0,factfppp_zeta
 202  real(dp) :: ecrs1,ex_gga,ex_lsd,exc,exp_pbe,expss,f_zeta,factfpp_zeta
 203  real(dp) :: fp_zeta,fpp_zeta,fppp_zeta,fx,gamphi3inv,gcrs,grrho2,hh,kappa
 204  real(dp) :: mac_a1,mac_aa,mac_b1,mac_b2,mac_b3,mac_b4,mac_den,mac_f1,mac_f2,mac_log,mac_q0
 205  real(dp) :: mac_q1,mac_q1ppp
 206  real(dp) :: mac_q1p,mac_q1pp,macrs,mu_divkappa,p1_wc,p2_wc,pade,pade_den
 207  real(dp) :: phi3_zeta,phi_logder,phi_zeta,phi_zeta_inv,phip_zeta,phipp_zeta,qq
 208  real(dp) :: rho,rho_inv,rhomot
 209  real(dp) :: rhotmo6,rhotmot,rhoto6,rhotot,rhotot_inv,rr,rs,rsm1_2,sqr_rs
 210  real(dp) :: sqr_sqr_rs,ss,tt,vxcadd,xx,zeta,zeta4,zetm_1_3,zetp_1_3
 211  real(dp) :: a1fa,a2fa,b1fa,b2fa,c1fa,c2fa,e1fa,e2fa,f1fa,f2fa,g1fa,g2fa,h1fa,h2fa
 212  real(dp) :: i1fa,i2fa,m1fa,m2fa,n1fa,n2fa
 213  real(dp) :: sp1_up3,sp1_up2dn,sp1_updn2,sp1_dn3
 214  real(dp) :: sp2_up3,sp2_up2dn,sp2_updn2,sp2_dn3
 215  real(dp) :: sp3_up3,sp3_up2dn,sp3_updn2,sp3_dn3
 216  real(dp) :: d3ecrs_sp0,d3ecrs_sp1,d3ecrs_sp2,d3ecrs_sp3
 217  character(len=500) :: message
 218 !arrays
 219  real(dp),allocatable :: rho_updnm1_3(:,:),rhoarr(:),rhom1_3(:),zetm(:)
 220  real(dp),allocatable :: zetmm1_3(:),zetp(:),zetpm1_3(:)
 221 !no_abirules
 222 !integer :: debug
 223 !real(dp) :: delta,factor,grr,rho_dn,rho_dnm,rho_dnp,rho_up,rho_upm,rho_upp,zeta_mean
 224 !real(dp), allocatable :: wecrsz(:,:),d1wecrsz(:,:),d2wecrsz(:,:),d3wecrsz(:,:)
 225 !real(dp) :: d3ecrs_drho3,d3ecrs_drhodndrho2,d3ecrs_drhoupdrho2
 226 !real(dp) :: ec1_q0p,mac_q0p,sigma1,sigma2,sigma3
 227 
 228 ! *************************************************************************
 229 
 230 !DEBUG
 231 !write(std_out,*)' xcpbe : enter'
 232 !ENDDEBUG
 233 
 234  d_wc=mu-b_wc
 235  exexch_=0;if(present(exexch)) exexch_=exexch
 236 
 237 !DEBUG
 238 !allocate(wecrsz(npts,8),d1wecrsz(npts,8),d2wecrsz(npts,8),d3wecrsz(npts,8))
 239 !ENDDEBUG
 240 
 241  if (option<=-4 .or. option==0 .or. option==4 .or. option>=8 ) then
 242    write(message, '(a,a,a,a,i12,a)' ) ch10,&
 243 &   ' xcpbe : BUG -',ch10,&
 244 &   '  Option must be 1, 2, 3, 5, 6, 7, 8, -1 or -2 ; argument was ',option,'.'
 245 !  ABI_BUG(message)
 246  end if
 247 
 248 !Checks the compatibility between the presence of dvxci and ndvxci
 249  if(ndvxci /=0 .neqv. present(dvxci))then
 250    message = ' If ndvxci/=0 there must the optional argument dvxci'
 251    ABI_BUG(message)
 252  end if
 253 
 254 !Checks the compatibility between the inputs and the presence of the optional arguments
 255  if(ndvxci /= 0 .and. abs(order) <= 1)then
 256    write(message, '(3a,i8,a)' )&
 257 &   'The order does not require the presence of dvxci',ch10,&
 258 &   'that is allowed when |order|>1, while we have',order,'.'
 259    ABI_BUG(message)
 260  end if
 261 
 262  if(ndvxci /= 0 .and. (&
 263 & ((option == 1 .or. option == -1 .or. option == 3) .and. ndvxci /= nspden + 1)&
 264 & .or. (option == -2 .and. ndvxci /= 8)&
 265 & .or. ((option == 2 .or. option == 5 .or. option == 6 .or. option == 7 .or. option == 8) .and. ndvxci /= 15)&
 266 & ))then
 267    write(message, '(12a,4(a,i5))' )&
 268 &   '  The option is not consistent with the value of ndvxci',ch10,&
 269 &   '  Allowed values are:',ch10,&
 270 &   '  ndvxci     option',ch10,&
 271 &   ' nspden+1    1,-1,3',ch10,&
 272 &   '    8          -2',ch10,&
 273 &   '    15       2, 5,6,7',ch10,&
 274 &   '  While we have: ndvxc=',ndvxci,', option=',option,', nspden=',nspden,', order=',order
 275    ABI_BUG(message)
 276  end if
 277 
 278  if ((option == 1 .or. option == -1 .or. option ==3) .and.  (present(grho2_updn) .or. present(dvxcdgr))) then
 279    write(message, '(a,a,a,i6,a)' )&
 280 &   'The option chosen does not need the presence',ch10,&
 281 &   'of the gradient, or of the array dvxcdgr in the input, needed if option/=1,-1,3 , while we have',option,'.'
 282    ABI_BUG(message)
 283  end if
 284 
 285  if (order /= 3 .and. present(d2vxci)) then
 286    write(message, '(a,a,a,i6,a)' )&
 287 &   'The order chosen does not need the presence',ch10,&
 288 &   'of the array d2vxci, needed if order=3 , while we have',order,'.'
 289    ABI_BUG(message)
 290  end if
 291 
 292  if(initialized==0)then
 293    twom1_3=two**(-third)
 294    sixpi2_1_3=(six*pi**2)**third
 295    sixpi2m1_3=one/sixpi2_1_3
 296    threefourth_divpi=three_quarters*piinv
 297    gamma=(one-log(two))*piinv**2
 298    gamma_inv=one/gamma
 299    beta_gamma=beta*gamma_inv
 300    factf_zeta= one / ( two**(four/three)-two )
 301    factfp_zeta= four_thirds * factf_zeta * alpha_zeta2
 302    coeff_tt= one/ (four*four*piinv*(three*pi**2)**third)
 303 !  coeff_tt= two * sqrt(four*piinv*(three*pi**2)**third)
 304    sq_rsfac=sqrt(rsfac)
 305    sq_rsfac_inv=one/sq_rsfac
 306    initialized=1
 307  end if
 308 
 309 !Parameters for the Perdew-Wang 92 LSD as well as LSD-RPA,
 310 !see Table I of Phys.Rev.B 45,13244 (1992) [[cite:Perdew1992]]
 311 !GMR (to match Libxc one should comment the next line and uncomment the following two)
 312  ec0_aa=0.031091_dp  ; ec1_aa=0.015545_dp ; mac_aa=0.016887_dp
 313 !ec0_aa=0.0310907_dp  ; ec1_aa=0.01554535_dp ; mac_aa=0.0168869_dp
 314 !GMR
 315  if(option/=3 .and. option/=4)then
 316    ec0_a1=0.21370_dp  ; ec1_a1=0.20548_dp  ; mac_a1=0.11125_dp
 317    ec0_b1=7.5957_dp   ; ec1_b1=14.1189_dp  ; mac_b1=10.357_dp
 318    ec0_b2=3.5876_dp   ; ec1_b2=6.1977_dp   ; mac_b2=3.6231_dp
 319    ec0_b3=1.6382_dp   ; ec1_b3=3.3662_dp   ; mac_b3=0.88026_dp
 320    ec0_b4=0.49294_dp  ; ec1_b4=0.62517_dp  ; mac_b4=0.49671_dp
 321  else  ! RPA values
 322    ec0_a1=0.082477_dp ; ec1_a1=0.035374_dp ; mac_a1=0.028829_dp
 323    ec0_b1=5.1486_dp   ; ec1_b1=6.4869_dp   ; mac_b1=10.357_dp
 324    ec0_b2=1.6483_dp   ; ec1_b2=1.3083_dp   ; mac_b2=3.6231_dp
 325    ec0_b3=0.23647_dp  ; ec1_b3=0.11518_dp  ; mac_b3=0.479_dp
 326    ec0_b4=0.20614_dp  ; ec1_b4=0.082349_dp ; mac_b4=0.112279_dp
 327  end if
 328 
 329  if(option/=5 .and. option/=-4)then
 330    kappa=kappa_pbe
 331    mu_divkappa=mu_divkappa_pbe
 332  end if
 333  if(option==5)then
 334    kappa=kappa_revpbe
 335    mu_divkappa=mu_divkappa_revpbe
 336  end if
 337  if(option==-4)then
 338    kappa=kappa_c09
 339  end if
 340 !DEBUG
 341 !Finite-difference debugging, do not take away
 342 !Note : here work with collinear gradients. Might be generalized ...
 343 !debug=2  ! Choose 1 (rho grads) or 2 (grho grads)
 344 !if(order==3)debug=1
 345 !factor=1.0_dp
 346 !zeta_mean=0.98_dp
 347 !!zeta_mean=zero
 348 !delta=0.000025*factor
 349 !delta=0.0000125*factor
 350 !if(debug/=0)then
 351 !do ipts=1,npts-4,5
 352 !rho=ipts*0.01_dp*factor
 353 !rho_up=rho*(1.0_dp+zeta_mean)*0.5_dp
 354 !rho_dn=rho*(1.0_dp-zeta_mean)*0.5_dp
 355 !rho_upp=rho_up+delta
 356 !rho_upm=rho_up-delta
 357 !rho_dnp=rho_dn+delta
 358 !rho_dnm=rho_dn-delta
 359 !! Here, vary rho
 360 !if(debug==1)then
 361 !rho_updn(ipts  ,1)=rho_up ; rho_updn(ipts  ,2)=rho_dn
 362 !rho_updn(ipts+1,1)=rho_upp; rho_updn(ipts+1,2)=rho_dn
 363 !rho_updn(ipts+2,1)=rho_upm; rho_updn(ipts+2,2)=rho_dn
 364 !rho_updn(ipts+3,1)=rho_up ; rho_updn(ipts+3,2)=rho_dnp
 365 !rho_updn(ipts+4,1)=rho_up ; rho_updn(ipts+4,2)=rho_dnm
 366 !grho2_updn(ipts:ipts+4,1)=(0.2_dp*factor)**2     ! grad2 of spin up density
 367 !grho2_updn(ipts:ipts+4,2)=(0.2_dp*factor)**2     ! grad2 of spin down density
 368 !grho2_updn(ipts:ipts+4,3)=(0.3_dp*factor)**2     ! grad2 of total density
 369 !else
 370 !!  Here, vary grho (interchange rho and grho)
 371 !grho2_updn(ipts  ,1)=rho_up**2 ; grho2_updn(ipts  ,2)=rho_dn**2
 372 !grho2_updn(ipts+1,1)=rho_upp**2; grho2_updn(ipts+1,2)=rho_dn**2
 373 !grho2_updn(ipts+2,1)=rho_upm**2; grho2_updn(ipts+2,2)=rho_dn**2
 374 !grho2_updn(ipts+3,1)=rho_up**2 ; grho2_updn(ipts+3,2)=rho_dnp**2
 375 !grho2_updn(ipts+4,1)=rho_up**2 ; grho2_updn(ipts+4,2)=rho_dnm**2
 376 !grho2_updn(ipts  ,3)=(ipts*0.01_dp*factor)**2
 377 !grho2_updn(ipts+1,3)=(ipts*0.01_dp*factor+delta)**2
 378 !grho2_updn(ipts+2,3)=(ipts*0.01_dp*factor-delta)**2
 379 !grho2_updn(ipts+3,3)=(ipts*0.01_dp*factor+delta)**2   ! identical to ipts+1
 380 !grho2_updn(ipts+4,3)=(ipts*0.01_dp*factor-delta)**2   ! identical to ipts+2
 381 !rho_updn(ipts:ipts+4,1)=0.2_dp*factor*(1.0_dp+zeta_mean)*0.5_dp    ! spin up density
 382 !rho_updn(ipts:ipts+4,2)=0.2_dp*factor*(1.0_dp-zeta_mean)*0.5_dp    ! spin down density
 383 !end if
 384 !end do
 385 !end if
 386 !Usual option :
 387 !nspden=2 ; order=2
 388 !GGA
 389 !nspden=2 ; order=1
 390 !Might take also, although finite difference later is meaningless
 391 !nspden=1 ; order=-2
 392 !Here, alternative specification, in terms of defined rs and zeta
 393 !do ipts=1,5
 394 !if(ipts==1)then ;rs=0.01_dp ; zeta=0.98_dp ; endif
 395 !if(ipts==2)then ;rs=0.01_dp+delta ; zeta=0.98_dp ; endif
 396 !if(ipts==3)then ;rs=0.01_dp-delta ; zeta=0.98_dp ; endif
 397 !if(ipts==4)then ;rs=0.01_dp ; zeta=0.98_dp+delta ; endif
 398 !if(ipts==5)then ;rs=0.01_dp ; zeta=0.98_dp-delta ; endif
 399 !rho=(rsfac/rs)**3
 400 !rho_up=rho*(1.0_dp+zeta)*0.5_dp
 401 !rho_dn=rho*(1.0_dp-zeta)*0.5_dp
 402 !rho_updn(ipts  ,1)=rho_up ; rho_updn(ipts  ,2)=rho_dn
 403 !enddo
 404 !ENDDEBUG
 405 
 406  if(order**2 >1)then
 407    factfpp_zeta= third * factfp_zeta * alpha_zeta2
 408  end if
 409 
 410 
 411  ABI_MALLOC(rhoarr,(npts))
 412  ABI_MALLOC(rhom1_3,(npts))
 413  ABI_MALLOC(rho_updnm1_3,(npts,2))
 414  ABI_MALLOC(zetm,(npts))
 415  ABI_MALLOC(zetmm1_3,(npts))
 416  ABI_MALLOC(zetp,(npts))
 417  ABI_MALLOC(zetpm1_3,(npts))
 418 
 419  do ispden=1,nspden
 420    call invcb(rho_updn(:,ispden),rho_updnm1_3(:,ispden),npts)
 421  end do
 422 
 423 
 424  if(nspden==1)then
 425    rhoarr(:)=two*rho_updn(:,1)
 426    rhom1_3(:)=twom1_3*rho_updnm1_3(:,1)
 427    rho_updnm1_3(:,2)=rho_updnm1_3(:,1)
 428  else
 429    rhoarr(:)=rho_updn(:,1)+rho_updn(:,2)
 430    call invcb(rhoarr,rhom1_3,npts)
 431    do ipts=1,npts
 432      rhotmot=rhom1_3(ipts)
 433      rhotot_inv=rhotmot*rhotmot*rhotmot
 434      zeta=(rho_updn(ipts,1)-rho_updn(ipts,2))*rhotot_inv
 435      zetp(ipts)=1.0_dp+zeta*alpha_zeta
 436      zetm(ipts)=1.0_dp-zeta*alpha_zeta
 437    end do
 438    call invcb(zetp,zetpm1_3,npts)
 439    call invcb(zetm,zetmm1_3,npts)
 440  end if
 441 
 442 
 443 !fab: eliminate the following restriction
 444 
 445 !if (order==3 .and. nspden == 1) d2vxci(:,:)=0._dp
 446 
 447 
 448  if (order==3) d2vxci(:,:)=0._dp
 449 
 450 !!!Loop unrolling summary
 451 !Completely unrolled for spin non-polarized case
 452 !To be optimized for spin-polarized cases
 453 !The loops are unrolled as follows:
 454 !nspden=1     (line 433)
 455 !order^2<=1  (line 460)
 456 !option=2,5 (line 462)
 457 !option=6,7 (line 630)
 458 !option=-1  (line 825)
 459 !option=-2  (line 853)
 460 !option=1   (line 904)
 461 !option=3   (line 963)
 462 !order=3     (line 1024)
 463 !option=2,5
 464 !option=6,7
 465 !option=-1
 466 !option=-2
 467 !option=1
 468 !option=3
 469 !order=-2    (line 1983)
 470 !option=2,5
 471 !option=6,7
 472 !option=-1
 473 !option=-2
 474 !option=1
 475 !option=3
 476 !order^2>1   (line 2875)
 477 !option=2,5
 478 !option=6,7
 479 !option=-1
 480 !option=-2
 481 !option=1
 482 !option=3
 483 !nspden=2     (line 3750)
 484 !order^2<=1  (line 3754)
 485 !order^2>1 (with if statements inside distinguishing between order=3 or -2)   (line 4000)
 486 !!!End loop unrolling summary
 487 
 488 !we separate different cases, depending on nspden
 489  if (nspden==1) then
 490 !  we separate different cases, depending on order
 491    if (order**2<=1) then
 492 !    we separate different cases, depending on option
 493      if(option==2 .or. option==5)then
 494 
 495        do ipts=1,npts
 496 
 497          rhotot=rhoarr(ipts)
 498          rhotmot=rhom1_3(ipts)
 499          rhotot_inv=rhotmot*rhotmot*rhotmot
 500          rhotmo6=sqrt(rhotmot)
 501          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
 502 !        -----------------------------------------------------------------------
 503 !        First take care of the exchange part of the functional
 504 
 505          exc=zero
 506          dvxcdgr(ipts,3)=zero
 507 !        loop over the spin
 508          ispden=1
 509          rho   =rho_updn(ipts,ispden)
 510          rhomot=rho_updnm1_3(ipts,ispden)
 511          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
 512 !        Perdew-Burke-Ernzerhof GGA, exchange part
 513          rho_inv=rhomot*rhomot*rhomot
 514          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
 515          ss=grho2_updn(ipts,ispden)*coeffss
 516          divss=one/(one+mu_divkappa*ss)
 517          dfxdss= mu*divss*divss
 518          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
 519          fx    = one+kappa*(one-divss)
 520          ex_gga= ex_lsd*fx
 521          dssdn=-eight*third*ss*rho_inv
 522          dfxdn  = dfxdss*dssdn
 523          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
 524 !        The new definition (v3.3) includes the division by the norm of the gradient
 525          dssdg =two*coeffss
 526          dfxdg=dfxdss*dssdg
 527          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
 528          exc=exc+ex_gga*rho
 529 
 530 !        end of loop over the spin
 531 !        If non spin-polarized, treat spin down contribution now, similar to spin up
 532          exc=exc*2
 533          exci(ipts)=exc*rhotot_inv
 534          if(exexch_==1) cycle
 535 !        -----------------------------------------------------------------------------
 536 !        Then takes care of the LSD correlation part of the functional
 537 
 538 
 539          rs=rsfac*rhotmot
 540          sqr_rs=sq_rsfac*rhotmo6
 541          rsm1_2=sq_rsfac_inv*rhoto6
 542 
 543 !        Formulas A6-A8 of PW92LSD
 544          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
 545          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
 546          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
 547          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
 548 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
 549          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
 550          ecrs0=ec0_q0*ec0_log
 551          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
 552 
 553          ecrs=ecrs0
 554          decrs_drs=decrs0_drs
 555          decrs_dzeta=0.0_dp
 556          zeta=0.0_dp
 557 
 558 !        Add LSD correlation functional to GGA exchange functional
 559          exci(ipts)=exci(ipts)+ecrs
 560          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
 561 
 562 
 563 !        -----------------------------------------------------------------------------
 564 !        Eventually add the GGA correlation part of the PBE functional
 565 !        Note : the computation of the potential in the spin-unpolarized
 566 !        case could be optimized much further. Other optimizations are left to do.
 567 
 568          phi_zeta=1.0_dp
 569          phip_zeta=0.0_dp
 570          phi_zeta_inv=1.0_dp
 571          phi_logder=0.0_dp
 572          phi3_zeta=1.0_dp
 573          gamphi3inv=gamma_inv
 574          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
 575 
 576 !        From ec to bb
 577          bb=ecrs*gamphi3inv
 578          dbb_drs=decrs_drs*gamphi3inv
 579          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
 580 
 581 !        From bb to cc
 582          exp_pbe=exp(-bb)
 583          cc=one/(exp_pbe-one)
 584          dcc_dbb=cc*cc*exp_pbe
 585          dcc_drs=dcc_dbb*dbb_drs
 586          dcc_dzeta=dcc_dbb*dbb_dzeta
 587 
 588 !        From cc to aa
 589          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
 590          aa=coeff_aa*cc
 591          daa_drs=coeff_aa*dcc_drs
 592          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
 593 
 594 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
 595          grrho2=four*grho2_updn(ipts,1)
 596          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
 597 !        Note that tt is (the t variable of PBE divided by phi) squared
 598          tt=half*grrho2*dtt_dg
 599 
 600 !        Get xx from aa and tt
 601          xx=aa*tt
 602          dxx_drs=daa_drs*tt
 603          dxx_dzeta=daa_dzeta*tt
 604          dxx_dtt=aa
 605 
 606 !        From xx to pade
 607          pade_den=one/(one+xx*(one+xx))
 608          pade=(one+xx)*pade_den
 609          dpade_dxx=-xx*(two+xx)*pade_den**2
 610          dpade_drs=dpade_dxx*dxx_drs
 611          dpade_dtt=dpade_dxx*dxx_dtt
 612          dpade_dzeta=dpade_dxx*dxx_dzeta
 613 
 614 !        From pade to qq
 615          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
 616          qq=coeff_qq*pade
 617          dqq_drs=coeff_qq*dpade_drs
 618          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
 619          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
 620 
 621 !        From qq to rr
 622          arg_rr=one+beta*gamma_inv*qq
 623          div_rr=one/arg_rr
 624          rr=gamma*log(arg_rr)
 625          drr_dqq=beta*div_rr
 626          drr_drs=drr_dqq*dqq_drs
 627          drr_dtt=drr_dqq*dqq_dtt
 628          drr_dzeta=drr_dqq*dqq_dzeta
 629 
 630 !        From rr to hh
 631          hh=phi3_zeta*rr
 632          dhh_drs=phi3_zeta*drr_drs
 633          dhh_dtt=phi3_zeta*drr_dtt
 634          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
 635 
 636 !        The GGA correlation energy is added
 637          exci(ipts)=exci(ipts)+hh
 638 
 639 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
 640 
 641 !        From hh to the derivative of the energy wrt the density
 642          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
 643          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
 644 
 645 !        From hh to the derivative of the energy wrt to the gradient of the
 646 !        density, divided by the gradient of the density
 647 !        (The v3.3 definition includes the division by the norm of the gradient)
 648          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
 649 
 650 !        End condition of GGA
 651 
 652 !        Correlation has been added
 653 !        -----------------------------------------------------------------------------
 654 
 655 !        vxci(ipts,2)=vxci(ipts,1)
 656          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
 657 
 658        end do
 659      else if((option==6) .or. (option==7)) then
 660 
 661        do ipts=1,npts
 662 
 663          rhotot=rhoarr(ipts)
 664          rhotmot=rhom1_3(ipts)
 665          rhotot_inv=rhotmot*rhotmot*rhotmot
 666          rhotmo6=sqrt(rhotmot)
 667          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
 668 !        -----------------------------------------------------------------------
 669 !        First take care of the exchange part of the functional
 670 
 671          exc=zero
 672          dvxcdgr(ipts,3)=zero
 673 !        loop over the spin
 674          ispden=1
 675          rho   =rho_updn(ipts,ispden)
 676          rhomot=rho_updnm1_3(ipts,ispden)
 677          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
 678 !        Perdew-Burke-Ernzerhof GGA, exchange part
 679          rho_inv=rhomot*rhomot*rhomot
 680          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
 681          ss=grho2_updn(ipts,ispden)*coeffss
 682 
 683 !        This is RPBE modification
 684          if (option==6) then
 685            divss=exp(-mu_divkappa*ss)
 686            dfxdss= mu*divss
 687            d2fxdss2=-mu*mu_divkappa*divss
 688 
 689            fx    = one+kappa*(one-divss)
 690            ex_gga= ex_lsd*fx
 691            dssdn=-eight*third*ss*rho_inv
 692            dfxdn  = dfxdss*dssdn
 693            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
 694 !          The new definition (v3.3) includes the division by the norm of the gradient
 695            dssdg =two*coeffss
 696            dfxdg=dfxdss*dssdg
 697            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
 698            exc=exc+ex_gga*rho
 699 !          This is the Wu and Cohen modification
 700          else
 701            expss=exp(-ss)
 702            p1_wc=b_wc+(mu-b_wc)*(one-ss)*expss+two*c_wc*ss/(one+c_wc*ss*ss)
 703            p2_wc=d_wc*(ss-two)*expss+two*c_wc/(one+c_wc*ss*ss)-&
 704 &           four*c_wc*c_wc*ss*ss/((one+c_wc*ss*ss)*(one+c_wc*ss*ss))
 705            divss=one/(one+(b_wc*ss+d_wc*ss*expss+log(one+c_wc*ss*ss))/kappa)
 706            dfxdss=p1_wc*divss*divss
 707            d2fxdss2=p2_wc*divss*divss-two*divss*divss*divss*p1_wc*p1_wc/kappa
 708 
 709            fx    = one+kappa*(one-divss)
 710            ex_gga= ex_lsd*fx
 711            dssdn=-eight*third*ss*rho_inv
 712            dfxdn  = dfxdss*dssdn
 713            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
 714 !          The new definition (v3.3) includes the division by the norm of the gradient
 715            dssdg =two*coeffss
 716            dfxdg=dfxdss*dssdg
 717            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
 718            exc=exc+ex_gga*rho
 719          end if
 720 
 721 !        end of loop over the spin
 722 !        If non spin-polarized, treat spin down contribution now, similar to spin up
 723          exc=exc*2
 724          exci(ipts)=exc*rhotot_inv
 725          if(exexch_==1) cycle
 726 !        -----------------------------------------------------------------------------
 727 !        Then takes care of the LSD correlation part of the functional
 728 
 729 
 730          rs=rsfac*rhotmot
 731          sqr_rs=sq_rsfac*rhotmo6
 732          rsm1_2=sq_rsfac_inv*rhoto6
 733 
 734 !        Formulas A6-A8 of PW92LSD
 735          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
 736          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
 737          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
 738          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
 739 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
 740          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
 741          ecrs0=ec0_q0*ec0_log
 742          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
 743 
 744          ecrs=ecrs0
 745          decrs_drs=decrs0_drs
 746          decrs_dzeta=0.0_dp
 747          zeta=0.0_dp
 748 
 749 !        Add LSD correlation functional to GGA exchange functional
 750          exci(ipts)=exci(ipts)+ecrs
 751          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
 752 
 753 
 754 !        -----------------------------------------------------------------------------
 755 !        Eventually add the GGA correlation part of the PBE functional
 756 !        Note : the computation of the potential in the spin-unpolarized
 757 !        case could be optimized much further. Other optimizations are left to do.
 758 
 759          phi_zeta=1.0_dp
 760          phip_zeta=0.0_dp
 761          phi_zeta_inv=1.0_dp
 762          phi_logder=0.0_dp
 763          phi3_zeta=1.0_dp
 764          gamphi3inv=gamma_inv
 765          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
 766 
 767 !        From ec to bb
 768          bb=ecrs*gamphi3inv
 769          dbb_drs=decrs_drs*gamphi3inv
 770          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
 771 
 772 !        From bb to cc
 773          exp_pbe=exp(-bb)
 774          cc=one/(exp_pbe-one)
 775          dcc_dbb=cc*cc*exp_pbe
 776          dcc_drs=dcc_dbb*dbb_drs
 777          dcc_dzeta=dcc_dbb*dbb_dzeta
 778 
 779 !        From cc to aa
 780          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
 781          aa=coeff_aa*cc
 782          daa_drs=coeff_aa*dcc_drs
 783          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
 784 
 785 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
 786          grrho2=four*grho2_updn(ipts,1)
 787          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
 788 !        Note that tt is (the t variable of PBE divided by phi) squared
 789          tt=half*grrho2*dtt_dg
 790 
 791 !        Get xx from aa and tt
 792          xx=aa*tt
 793          dxx_drs=daa_drs*tt
 794          dxx_dzeta=daa_dzeta*tt
 795          dxx_dtt=aa
 796 
 797 !        From xx to pade
 798          pade_den=one/(one+xx*(one+xx))
 799          pade=(one+xx)*pade_den
 800          dpade_dxx=-xx*(two+xx)*pade_den**2
 801          dpade_drs=dpade_dxx*dxx_drs
 802          dpade_dtt=dpade_dxx*dxx_dtt
 803          dpade_dzeta=dpade_dxx*dxx_dzeta
 804 
 805 !        From pade to qq
 806          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
 807          qq=coeff_qq*pade
 808          dqq_drs=coeff_qq*dpade_drs
 809          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
 810          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
 811 
 812 !        From qq to rr
 813          arg_rr=one+beta*gamma_inv*qq
 814          div_rr=one/arg_rr
 815          rr=gamma*log(arg_rr)
 816          drr_dqq=beta*div_rr
 817          drr_drs=drr_dqq*dqq_drs
 818          drr_dtt=drr_dqq*dqq_dtt
 819          drr_dzeta=drr_dqq*dqq_dzeta
 820 
 821 !        From rr to hh
 822          hh=phi3_zeta*rr
 823          dhh_drs=phi3_zeta*drr_drs
 824          dhh_dtt=phi3_zeta*drr_dtt
 825          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
 826 
 827 !        The GGA correlation energy is added
 828          exci(ipts)=exci(ipts)+hh
 829 
 830 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
 831 
 832 !        From hh to the derivative of the energy wrt the density
 833          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
 834          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
 835 
 836 !        From hh to the derivative of the energy wrt to the gradient of the
 837 !        density, divided by the gradient of the density
 838 !        (The v3.3 definition includes the division by the norm of the gradient)
 839          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
 840 
 841 !        End condition of GGA
 842 
 843 !        Correlation has been added
 844 !        -----------------------------------------------------------------------------
 845 
 846 !        vxci(ipts,2)=vxci(ipts,1)
 847          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
 848 
 849        end do
 850 
 851 
 852      else if (option==-1) then
 853 
 854        do ipts=1,npts
 855 
 856          rhotot=rhoarr(ipts)
 857          rhotmot=rhom1_3(ipts)
 858          rhotot_inv=rhotmot*rhotmot*rhotmot
 859          rhotmo6=sqrt(rhotmot)
 860          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
 861 !        -----------------------------------------------------------------------
 862 !        First take care of the exchange part of the functional
 863 
 864          exc=zero
 865 !        loop over the spin
 866          ispden=1
 867          rho   =rho_updn(ipts,ispden)
 868          rhomot=rho_updnm1_3(ipts,ispden)
 869          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
 870 !        Perdew_Wang 91 LSD
 871          vxci(ipts,ispden)=four_thirds*ex_lsd
 872          exc=exc+ex_lsd*rho
 873 
 874 !        end of loop over the spin
 875 !        If non spin-polarized, treat spin down contribution now, similar to spin up
 876          exc=exc*2
 877          exci(ipts)=exc*rhotot_inv
 878        end do
 879 
 880      else if(option==-2) then
 881 
 882 
 883        do ipts=1,npts
 884 
 885          rhotot=rhoarr(ipts)
 886          rhotmot=rhom1_3(ipts)
 887          rhotot_inv=rhotmot*rhotmot*rhotmot
 888          rhotmo6=sqrt(rhotmot)
 889          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
 890 !        -----------------------------------------------------------------------
 891 !        First take care of the exchange part of the functional
 892 
 893          exc=zero
 894          dvxcdgr(ipts,3)=zero
 895 !        loop over the spin
 896          ispden=1
 897          rho   =rho_updn(ipts,ispden)
 898          rhomot=rho_updnm1_3(ipts,ispden)
 899          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
 900 !        Perdew-Burke-Ernzerhof GGA, exchange part
 901          rho_inv=rhomot*rhomot*rhomot
 902          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
 903          ss=grho2_updn(ipts,ispden)*coeffss
 904          divss=one/(one+mu_divkappa*ss)
 905          dfxdss= mu*divss*divss
 906          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
 907          fx    = one+kappa*(one-divss)
 908          ex_gga= ex_lsd*fx
 909          dssdn=-eight*third*ss*rho_inv
 910          dfxdn  = dfxdss*dssdn
 911          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
 912 !        The new definition (v3.3) includes the division by the norm of the gradient
 913          dssdg =two*coeffss
 914          dfxdg=dfxdss*dssdg
 915          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
 916          exc=exc+ex_gga*rho
 917 
 918 !        end of loop over the spin
 919 !        If non spin-polarized, treat spin down contribution now, similar to spin up
 920          exc=exc*2
 921          exci(ipts)=exc*rhotot_inv
 922 
 923 !        Correlation has been added
 924 !        -----------------------------------------------------------------------------
 925 
 926 !        vxci(ipts,2)=vxci(ipts,1)
 927          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
 928 
 929        end do
 930 
 931 
 932      else if(option==-4) then
 933 
 934 
 935        do ipts=1,npts
 936 
 937          rhotot=rhoarr(ipts)
 938          rhotmot=rhom1_3(ipts)
 939          rhotot_inv=rhotmot*rhotmot*rhotmot
 940          rhotmo6=sqrt(rhotmot)
 941          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
 942 !        -----------------------------------------------------------------------
 943 !        First take care of the exchange part of the functional
 944 
 945          exc=zero
 946          dvxcdgr(ipts,3)=zero
 947 !        loop over the spin
 948          ispden=1
 949          rho   =rho_updn(ipts,ispden)
 950          rhomot=rho_updnm1_3(ipts,ispden)
 951          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
 952 !        VALENTINO R. COOPER C09x GGA, This is an exchange term proposed
 953 !        to use together with vdw-DF (see above).
 954          rho_inv=rhomot*rhomot*rhomot
 955          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
 956 !        the quarter that is lacking is compensated by the grho2_updn in the
 957 !        next line.
 958          ss=grho2_updn(ipts,ispden)*coeffss
 959          alphs2=alpha_c09*ss
 960 ! jmb : overflow if alphs2 > 100
 961          if (alphs2 > 100.0 ) alphs2=100.0
 962          alphmu=alpha_c09*mu_c09
 963          dfxdss= mu_c09*exp(-alphs2)*(one-alphs2)+&
 964 &         kappa*alpha_c09*exp(-alphs2/two)/two
 965          d2fxdss2=-alphmu*exp(-alphs2)*(two-alphs2)-&
 966 &         kappa*(alpha_c09**two)*exp(alphs2/two)/four
 967          fx    = one+mu_c09*ss*exp(-alphs2)+kappa*(one-exp(-alphs2/two))
 968          ex_gga= ex_lsd*fx
 969          dssdn=-eight*third*ss*rho_inv
 970          dfxdn  = dfxdss*dssdn
 971          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
 972 !
 973 !        The new definition (v3.3) includes the division by the norm of the gradient
 974          dssdg =two*coeffss
 975          dfxdg=dfxdss*dssdg
 976          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg !here also
 977          exc=exc+ex_gga*rho
 978 
 979 !        end of loop over the spin
 980 !        If non spin-polarized, treat spin down contribution now, similar to spin up
 981          exc=exc*2
 982          exci(ipts)=exc*rhotot_inv
 983 
 984 !        Correlation has been added
 985 !        -----------------------------------------------------------------------------
 986 
 987 !        vxci(ipts,2)=vxci(ipts,1)
 988          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
 989 
 990        end do
 991 
 992 
 993      else if(option==1)then
 994 
 995 
 996        do ipts=1,npts
 997 
 998          rhotot=rhoarr(ipts)
 999          rhotmot=rhom1_3(ipts)
1000          rhotot_inv=rhotmot*rhotmot*rhotmot
1001          rhotmo6=sqrt(rhotmot)
1002          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1003 !        -----------------------------------------------------------------------
1004 !        First take care of the exchange part of the functional
1005 
1006          exc=zero
1007 !        loop over the spin
1008          ispden=1
1009          rho   =rho_updn(ipts,ispden)
1010          rhomot=rho_updnm1_3(ipts,ispden)
1011          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1012 !        Perdew_Wang 91 LSD
1013          vxci(ipts,ispden)=four_thirds*ex_lsd
1014          exc=exc+ex_lsd*rho
1015 
1016 !        end of loop over the spin
1017 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1018          exc=exc*2
1019          exci(ipts)=exc*rhotot_inv
1020 !        -----------------------------------------------------------------------------
1021 !        Then takes care of the LSD correlation part of the functional
1022 
1023 
1024          rs=rsfac*rhotmot
1025          sqr_rs=sq_rsfac*rhotmo6
1026          rsm1_2=sq_rsfac_inv*rhoto6
1027 
1028 !        Formulas A6-A8 of PW92LSD
1029          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
1030          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
1031          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
1032          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
1033 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
1034          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
1035          ecrs0=ec0_q0*ec0_log
1036          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
1037 
1038          ecrs=ecrs0
1039          decrs_drs=decrs0_drs
1040          decrs_dzeta=0.0_dp
1041          zeta=0.0_dp
1042 
1043 !        Add LSD correlation functional to GGA exchange functional
1044          exci(ipts)=exci(ipts)+ecrs
1045          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
1046 
1047 !        Correlation has been added
1048 !        -----------------------------------------------------------------------------
1049 
1050        end do
1051 
1052      else if (option==3) then
1053        do ipts=1,npts
1054 
1055          rhotot=rhoarr(ipts)
1056          rhotmot=rhom1_3(ipts)
1057          rhotot_inv=rhotmot*rhotmot*rhotmot
1058          rhotmo6=sqrt(rhotmot)
1059          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1060 !        -----------------------------------------------------------------------
1061 !        First take care of the exchange part of the functional
1062 
1063          exc=zero
1064 !        loop over the spin
1065          ispden=1
1066          rho   =rho_updn(ipts,ispden)
1067          rhomot=rho_updnm1_3(ipts,ispden)
1068          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1069 !        Perdew_Wang 91 LSD
1070          vxci(ipts,ispden)=four_thirds*ex_lsd
1071          exc=exc+ex_lsd*rho
1072 
1073 !        end of loop over the spin
1074 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1075          exc=exc*2
1076          exci(ipts)=exc*rhotot_inv
1077 !        -----------------------------------------------------------------------------
1078 !        Then takes care of the LSD correlation part of the functional
1079 
1080 
1081          rs=rsfac*rhotmot
1082          sqr_rs=sq_rsfac*rhotmo6
1083          rsm1_2=sq_rsfac_inv*rhoto6
1084 
1085 !        Formulas A6-A8 of PW92LSD
1086          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
1087          sqr_sqr_rs=max(1.e-15_dp,sqrt(sqr_rs))
1088          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs/sqr_sqr_rs)
1089          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+3.5_dp*ec0_b4*rs/sqr_sqr_rs)
1090          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
1091 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
1092          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
1093          ecrs0=ec0_q0*ec0_log
1094          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
1095 
1096          ecrs=ecrs0
1097          decrs_drs=decrs0_drs
1098          decrs_dzeta=0.0_dp
1099          zeta=0.0_dp
1100 
1101 !        Add LSD correlation functional to GGA exchange functional
1102          exci(ipts)=exci(ipts)+ecrs
1103          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
1104 
1105 !        Correlation has been added
1106 !        -----------------------------------------------------------------------------
1107 
1108        end do
1109 
1110 
1111      end if
1112 
1113    else if (order==3) then
1114 !    separate cases with respect to option
1115      if(option==2 .or. option==5) then
1116 
1117        do ipts=1,npts
1118 
1119          rhotot=rhoarr(ipts)
1120          rhotmot=rhom1_3(ipts)
1121          rhotot_inv=rhotmot*rhotmot*rhotmot
1122          rhotmo6=sqrt(rhotmot)
1123          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1124 !        -----------------------------------------------------------------------
1125 !        First take care of the exchange part of the functional
1126 
1127          exc=zero
1128          dvxcdgr(ipts,3)=zero
1129 !        loop over the spin
1130          ispden=1
1131          rho   =rho_updn(ipts,ispden)
1132          rhomot=rho_updnm1_3(ipts,ispden)
1133          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1134 !        Perdew-Burke-Ernzerhof GGA, exchange part
1135          rho_inv=rhomot*rhomot*rhomot
1136          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
1137          ss=grho2_updn(ipts,ispden)*coeffss
1138          divss=one/(one+mu_divkappa*ss)
1139          dfxdss= mu*divss*divss
1140          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
1141          fx    = one+kappa*(one-divss)
1142          ex_gga= ex_lsd*fx
1143          dssdn=-eight*third*ss*rho_inv
1144          dfxdn  = dfxdss*dssdn
1145          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
1146 !        The new definition (v3.3) includes the division by the norm of the gradient
1147          dssdg =two*coeffss
1148          dfxdg=dfxdss*dssdg
1149          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
1150          exc=exc+ex_gga*rho
1151 
1152 !        Perdew-Burke-Ernzerhof GGA, exchange part
1153 !        Components 3 or 4
1154          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
1155 !        Components 1 or 2
1156          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
1157          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
1158          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
1159 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
1160 !        Components 5 or 6
1161          d2ssdndg=-eight*third*dssdg*rho_inv
1162          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
1163          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
1164 !        Components 7 or 8
1165          d2fxdg2=d2fxdss2*dssdg**2
1166          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
1167 !        For the time being, treat non-spin-polarized like spin-polarized
1168          dvxci(ipts,2)=dvxci(ipts,1)
1169          dvxci(ipts,4)=dvxci(ipts,3)
1170          dvxci(ipts,6)=dvxci(ipts,5)
1171          dvxci(ipts,8)=dvxci(ipts,7)
1172 
1173 !        end of loop over the spin
1174 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1175          exc=exc*2
1176          exci(ipts)=exc*rhotot_inv
1177 !        -----------------------------------------------------------------------------
1178 !        Then takes care of the LSD correlation part of the functional
1179 
1180 
1181          rs=rsfac*rhotmot
1182          sqr_rs=sq_rsfac*rhotmo6
1183          rsm1_2=sq_rsfac_inv*rhoto6
1184 
1185 !        Formulas A6-A8 of PW92LSD
1186          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
1187          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
1188          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
1189          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
1190 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
1191          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
1192          ecrs0=ec0_q0*ec0_log
1193          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
1194          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
1195          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
1196 &         -ec0_q0*ec0_q1pp*ec0_den                        &
1197 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
1198          ec0_q1ppp = 0.75_dp*ec0_aa*(rsm1_2**5)*(ec0_b1-ec0_b3*rs)
1199          ec0_f1 = 1._dp/(ec0_q1*ec0_q1*(1._dp + ec0_q1))
1200          ec0_f2 = 1._dp/(ec0_q1*(1+ec0_q1))
1201          d3ecrs0_drs3 = 6._dp*ec0_q1p*ec0_f1*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + &
1202 &         ec0_q0*ec0_q1pp) - &
1203 &         ec0_f2*(-6._dp*ec0_aa*ec0_a1*ec0_q1pp + ec0_q0*ec0_q1ppp + &
1204 &         ec0_f2*(3._dp*ec0_q1p*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + ec0_q0*ec0_q1pp) + &
1205 &         ec0_f2*2._dp*ec0_q0*(ec0_q1p**3)*(1._dp + 3._dp*ec0_q1*(1._dp + ec0_q1))))
1206 
1207          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
1208          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
1209          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
1210          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
1211          mac_log=-log( mac_q1*mac_q1*mac_den )
1212          macrs=mac_q0*mac_log
1213          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
1214 
1215          ecrs=ecrs0
1216          decrs_drs=decrs0_drs
1217          decrs_dzeta=0.0_dp
1218          d2ecrs_drs2=d2ecrs0_drs2
1219          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
1220          d2ecrs_drsdzeta=zero
1221          zeta=0.0_dp
1222 
1223 
1224 !        Add LSD correlation functional to GGA exchange functional
1225          exci(ipts)=exci(ipts)+ecrs
1226          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
1227 
1228          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
1229 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
1230          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
1231          dvxci(ipts,9)=d2ecrs_drho2
1232          dvxci(ipts,10)=d2ecrs_drho2
1233          dvxci(ipts,11)=d2ecrs_drho2
1234 
1235 !        -----------------------------------------------------------------------------
1236 !        Eventually add the GGA correlation part of the PBE functional
1237 !        Note : the computation of the potential in the spin-unpolarized
1238 !        case could be optimized much further. Other optimizations are left to do.
1239 
1240          phi_zeta=1.0_dp
1241          phip_zeta=0.0_dp
1242          phi_zeta_inv=1.0_dp
1243          phi_logder=0.0_dp
1244          phi3_zeta=1.0_dp
1245          gamphi3inv=gamma_inv
1246          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
1247 
1248 !        From ec to bb
1249          bb=ecrs*gamphi3inv
1250          dbb_drs=decrs_drs*gamphi3inv
1251          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
1252          d2bb_drs2=d2ecrs_drs2*gamphi3inv
1253          d2bb_drsdzeta=gamphi3inv*(d2ecrs_drsdzeta-three*decrs_drs*phi_logder)
1254          d2bb_dzeta2=gamphi3inv*(d2ecrs_dzeta2-six*decrs_dzeta*phi_logder+&
1255 &         12.0_dp*ecrs*phi_logder*phi_logder-three*ecrs*phi_zeta_inv*phipp_zeta)
1256 
1257 !        From bb to cc
1258          exp_pbe=exp(-bb)
1259          cc=one/(exp_pbe-one)
1260          dcc_dbb=cc*cc*exp_pbe
1261          dcc_drs=dcc_dbb*dbb_drs
1262          dcc_dzeta=dcc_dbb*dbb_dzeta
1263          d2cc_dbb2=cc*cc*exp_pbe*(two*cc*exp_pbe-one)
1264          d2cc_drs2=d2cc_dbb2*dbb_drs*dbb_drs+dcc_dbb*d2bb_drs2
1265          d2cc_drsdzeta=d2cc_dbb2*dbb_drs*dbb_dzeta+dcc_dbb*d2bb_drsdzeta
1266          d2cc_dzeta2=d2cc_dbb2*dbb_dzeta*dbb_dzeta+dcc_dbb*d2bb_dzeta2
1267 
1268 !        From cc to aa
1269          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
1270          aa=coeff_aa*cc
1271          daa_drs=coeff_aa*dcc_drs
1272          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
1273          d2aa_drs2=coeff_aa*d2cc_drs2
1274          d2aa_drsdzeta=-two*daa_drs*phi_logder+coeff_aa*d2cc_drsdzeta
1275          d2aa_dzeta2=aa*(-two*phi_zeta_inv*phipp_zeta+six*phi_logder*phi_logder)+&
1276 &         coeff_aa*(-four*dcc_dzeta*phi_logder+d2cc_dzeta2)
1277 
1278 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
1279          grrho2=four*grho2_updn(ipts,1)
1280          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
1281 !        Note that tt is (the t variable of PBE divided by phi) squared
1282          tt=half*grrho2*dtt_dg
1283 
1284 !        Get xx from aa and tt
1285          xx=aa*tt
1286          dxx_drs=daa_drs*tt
1287          dxx_dzeta=daa_dzeta*tt
1288          dxx_dtt=aa
1289          d2xx_drs2=d2aa_drs2*tt
1290          d2xx_drsdzeta=d2aa_drsdzeta*tt
1291          d2xx_drsdtt=daa_drs
1292          d2xx_dttdzeta=daa_dzeta
1293          d2xx_dzeta2=d2aa_dzeta2*tt
1294 
1295 !        From xx to pade
1296          pade_den=one/(one+xx*(one+xx))
1297          pade=(one+xx)*pade_den
1298          dpade_dxx=-xx*(two+xx)*pade_den**2
1299          dpade_drs=dpade_dxx*dxx_drs
1300          dpade_dtt=dpade_dxx*dxx_dtt
1301          dpade_dzeta=dpade_dxx*dxx_dzeta
1302          d2pade_dxx2=two*(-one+xx*xx*(three+xx))*pade_den*pade_den*pade_den
1303          d2pade_drs2=d2pade_dxx2*dxx_drs*dxx_drs+dpade_dxx*d2xx_drs2
1304          d2pade_drsdtt=d2pade_dxx2*dxx_drs*dxx_dtt+dpade_dxx*d2xx_drsdtt
1305          d2pade_drsdzeta=d2pade_dxx2*dxx_drs*dxx_dzeta+dpade_dxx*d2xx_drsdzeta
1306          d2pade_dtt2=d2pade_dxx2*dxx_dtt*dxx_dtt
1307          d2pade_dttdzeta=d2pade_dxx2*dxx_dtt*dxx_dzeta+dpade_dxx*d2xx_dttdzeta
1308          d2pade_dzeta2=d2pade_dxx2*dxx_dzeta*dxx_dzeta+dpade_dxx*d2xx_dzeta2
1309 
1310 !        From pade to qq
1311          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
1312          qq=coeff_qq*pade
1313          dqq_drs=coeff_qq*dpade_drs
1314          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
1315          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
1316          d2qq_drs2=coeff_qq*d2pade_drs2
1317          d2qq_drsdtt=phi_zeta_inv*phi_zeta_inv*(dpade_drs+tt*d2pade_drsdtt)
1318          d2qq_drsdzeta=coeff_qq*(d2pade_drsdzeta-two*dpade_drs*phi_logder)
1319          d2qq_dtt2=phi_zeta_inv*phi_zeta_inv*(two*dpade_dtt+tt*d2pade_dtt2)
1320          d2qq_dttdzeta=phi_zeta_inv*phi_zeta_inv*(dpade_dzeta-two*pade*phi_logder)+&
1321 &         coeff_qq*(d2pade_dttdzeta-two*dpade_dtt*phi_logder)
1322          d2qq_dzeta2=coeff_qq*( d2pade_dzeta2-four*dpade_dzeta*phi_logder &
1323 &         +six*pade*phi_logder*phi_logder            &
1324 &         -two*pade*phi_zeta_inv*phipp_zeta)
1325 
1326 !        From qq to rr
1327          arg_rr=one+beta*gamma_inv*qq
1328          div_rr=one/arg_rr
1329          rr=gamma*log(arg_rr)
1330          drr_dqq=beta*div_rr
1331          drr_drs=drr_dqq*dqq_drs
1332          drr_dtt=drr_dqq*dqq_dtt
1333          drr_dzeta=drr_dqq*dqq_dzeta
1334          d2rr_dqq2=-div_rr**2*beta*beta*gamma_inv
1335          d2rr_drs2=d2rr_dqq2*dqq_drs*dqq_drs+drr_dqq*d2qq_drs2
1336          d2rr_drsdtt=d2rr_dqq2*dqq_drs*dqq_dtt+drr_dqq*d2qq_drsdtt
1337          d2rr_drsdzeta=d2rr_dqq2*dqq_drs*dqq_dzeta+drr_dqq*d2qq_drsdzeta
1338          d2rr_dtt2=d2rr_dqq2*dqq_dtt*dqq_dtt+drr_dqq*d2qq_dtt2
1339          d2rr_dttdzeta=d2rr_dqq2*dqq_dtt*dqq_dzeta+drr_dqq*d2qq_dttdzeta
1340          d2rr_dzeta2=d2rr_dqq2*dqq_dzeta*dqq_dzeta+drr_dqq*d2qq_dzeta2
1341 
1342 !        From rr to hh
1343          hh=phi3_zeta*rr
1344          dhh_drs=phi3_zeta*drr_drs
1345          dhh_dtt=phi3_zeta*drr_dtt
1346          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
1347          d2hh_drs2=phi3_zeta*d2rr_drs2
1348          d2hh_drsdtt=phi3_zeta*d2rr_drsdtt
1349          d2hh_drsdzeta=phi3_zeta*(d2rr_drsdzeta+three*drr_drs*phi_logder)
1350          d2hh_dtt2=phi3_zeta*d2rr_dtt2
1351          d2hh_dttdzeta=phi3_zeta*(d2rr_dttdzeta+three*drr_dtt*phi_logder)
1352          d2hh_dzeta2=phi3_zeta*(six*rr*phi_logder*phi_logder+&
1353 &         six*phi_logder*drr_dzeta+d2rr_dzeta2)  &
1354 &         +three*phi_zeta*phi_zeta*rr*phipp_zeta
1355 
1356 
1357 !        The GGA correlation energy is added
1358          exci(ipts)=exci(ipts)+hh
1359 
1360 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
1361 
1362 
1363 
1364 !        From hh to the derivative of the energy wrt the density
1365          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
1366          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
1367 
1368 !        From hh to the derivative of the energy wrt to the gradient of the
1369 !        density, divided by the gradient of the density
1370 !        (The v3.3 definition includes the division by the norm of the gradient)
1371          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
1372 
1373          d2rhohh_drho2=rhotot_inv*&
1374 &         (-two*ninth*rs*dhh_drs +seven*four*ninth*tt*dhh_dtt &
1375 &         +ninth*rs*rs*d2hh_drs2+zeta*zeta*d2hh_dzeta2+(seven*third*tt)**2*d2hh_dtt2 &
1376 &         +two*third*rs*zeta*d2hh_drsdzeta+two*seven*ninth*rs*tt*d2hh_drsdtt &
1377 &         +two*seven*third*tt*zeta*d2hh_dttdzeta)
1378          d2rhohh_drhodg=dtt_dg*(-four*third*dhh_dtt-third*rs*d2hh_drsdtt &
1379 &         -zeta*d2hh_dttdzeta-seven*third*tt*d2hh_dtt2)
1380 
1381 !        Component 12 : first derivative with respect to the gradient
1382 !        of the density, div by the grad of the density
1383          dvxci(ipts,12)=dvxcdgr(ipts,3)
1384 !        Components 9, 10 and 11 : second derivatives with respect to the spin-density
1385 !        Note that there is already a contribution from LSDA
1386          dvxci(ipts,9)=dvxci(ipts,9)+d2rhohh_drho2+rhotot_inv*           &
1387 &         ( d2hh_dzeta2*(one-two*zeta) &
1388 &         -two*third*rs*d2hh_drsdzeta-14.0_dp*third*tt*d2hh_dttdzeta)
1389          dvxci(ipts,10)=dvxci(ipts,10)+d2rhohh_drho2-rhotot_inv*d2hh_dzeta2
1390          dvxci(ipts,11)=dvxci(ipts,11)+d2rhohh_drho2+rhotot_inv*           &
1391 &         ( d2hh_dzeta2*(one+two*zeta) &
1392 &         +two*third*rs*d2hh_drsdzeta+14.0_dp*third*tt*d2hh_dttdzeta)
1393 !        Components 13 and 14 : second derivatives with respect to spin density
1394 !        and gradient, divided by the gradient
1395          dvxci(ipts,13)=d2rhohh_drhodg+dtt_dg*d2hh_dttdzeta
1396          dvxci(ipts,14)=d2rhohh_drhodg-dtt_dg*d2hh_dttdzeta
1397 !        Component 15 : derivative of the (derivative wrt the gradient div by the grad),
1398 !        divided by the grad
1399          dvxci(ipts,15)=rhotot*d2hh_dtt2*dtt_dg*dtt_dg
1400 
1401 
1402 !        End condition of GGA
1403 
1404 !        Correlation has been added
1405 !        -----------------------------------------------------------------------------
1406 
1407 !        vxci(ipts,2)=vxci(ipts,1)
1408          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
1409 
1410        end do
1411      else if ((option==6) .or. (option==7)) then
1412 
1413        do ipts=1,npts
1414 
1415          rhotot=rhoarr(ipts)
1416          rhotmot=rhom1_3(ipts)
1417          rhotot_inv=rhotmot*rhotmot*rhotmot
1418          rhotmo6=sqrt(rhotmot)
1419          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1420 !        -----------------------------------------------------------------------
1421 !        First take care of the exchange part of the functional
1422 
1423          exc=zero
1424          dvxcdgr(ipts,3)=zero
1425 !        loop over the spin
1426          ispden=1
1427          rho   =rho_updn(ipts,ispden)
1428          rhomot=rho_updnm1_3(ipts,ispden)
1429          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1430 !        Perdew-Burke-Ernzerhof GGA, exchange part
1431          rho_inv=rhomot*rhomot*rhomot
1432          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
1433          ss=grho2_updn(ipts,ispden)*coeffss
1434 
1435          if (option==6) then
1436            divss=exp(-mu_divkappa*ss)
1437            dfxdss= mu*divss
1438            d2fxdss2=-mu*mu_divkappa*divss
1439 
1440            fx    = one+kappa*(one-divss)
1441            ex_gga= ex_lsd*fx
1442            dssdn=-eight*third*ss*rho_inv
1443            dfxdn  = dfxdss*dssdn
1444            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
1445 !          The new definition (v3.3) includes the division by the norm of the gradient
1446            dssdg =two*coeffss
1447            dfxdg=dfxdss*dssdg
1448            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
1449            exc=exc+ex_gga*rho
1450 !          This is the Wu and Cohen modification
1451          else
1452            expss=exp(-ss)
1453            p1_wc=b_wc+(mu-b_wc)*(one-ss)*expss+two*c_wc*ss/(one+c_wc*ss*ss)
1454            p2_wc=d_wc*(ss-two)*expss+two*c_wc/(one+c_wc*ss*ss)-&
1455 &           four*c_wc*c_wc*ss*ss/((one+c_wc*ss*ss)*(one+c_wc*ss*ss))
1456            divss=one/(one+(b_wc*ss+d_wc*ss*expss+log(one+c_wc*ss*ss))/kappa)
1457            dfxdss=p1_wc*divss*divss
1458            d2fxdss2=p2_wc*divss*divss-two*divss*divss*divss*p1_wc*p1_wc/kappa
1459 
1460            fx    = one+kappa*(one-divss)
1461            ex_gga= ex_lsd*fx
1462            dssdn=-eight*third*ss*rho_inv
1463            dfxdn  = dfxdss*dssdn
1464            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
1465 !          The new definition (v3.3) includes the division by the norm of the gradient
1466            dssdg =two*coeffss
1467            dfxdg=dfxdss*dssdg
1468            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
1469            exc=exc+ex_gga*rho
1470          end if
1471 
1472 !        Perdew-Burke-Ernzerhof GGA, exchange part
1473 !        Components 3 or 4
1474          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
1475 !        Components 1 or 2
1476          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
1477          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
1478          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
1479 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
1480 !        Components 5 or 6
1481          d2ssdndg=-eight*third*dssdg*rho_inv
1482          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
1483          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
1484 !        Components 7 or 8
1485          d2fxdg2=d2fxdss2*dssdg**2
1486          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
1487 !        For the time being, treat non-spin-polarized like spin-polarized
1488          dvxci(ipts,2)=dvxci(ipts,1)
1489          dvxci(ipts,4)=dvxci(ipts,3)
1490          dvxci(ipts,6)=dvxci(ipts,5)
1491          dvxci(ipts,8)=dvxci(ipts,7)
1492 
1493 !        end of loop over the spin
1494 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1495          exc=exc*2
1496          exci(ipts)=exc*rhotot_inv
1497 !        -----------------------------------------------------------------------------
1498 !        Then takes care of the LSD correlation part of the functional
1499 
1500 
1501          rs=rsfac*rhotmot
1502          sqr_rs=sq_rsfac*rhotmo6
1503          rsm1_2=sq_rsfac_inv*rhoto6
1504 
1505 !        Formulas A6-A8 of PW92LSD
1506          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
1507          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
1508          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
1509          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
1510 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
1511          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
1512          ecrs0=ec0_q0*ec0_log
1513          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
1514          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
1515          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
1516 &         -ec0_q0*ec0_q1pp*ec0_den                        &
1517 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
1518          ec0_q1ppp = 0.75_dp*ec0_aa*(rsm1_2**5)*(ec0_b1-ec0_b3*rs)
1519          ec0_f1 = 1._dp/(ec0_q1*ec0_q1*(1._dp + ec0_q1))
1520          ec0_f2 = 1._dp/(ec0_q1*(1+ec0_q1))
1521          d3ecrs0_drs3 = 6._dp*ec0_q1p*ec0_f1*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + &
1522 &         ec0_q0*ec0_q1pp) - &
1523 &         ec0_f2*(-6._dp*ec0_aa*ec0_a1*ec0_q1pp + ec0_q0*ec0_q1ppp + &
1524 &         ec0_f2*(3._dp*ec0_q1p*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + ec0_q0*ec0_q1pp) + &
1525 &         ec0_f2*2._dp*ec0_q0*(ec0_q1p**3)*(1._dp + 3._dp*ec0_q1*(1._dp + ec0_q1))))
1526 
1527          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
1528          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
1529          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
1530          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
1531          mac_log=-log( mac_q1*mac_q1*mac_den )
1532          macrs=mac_q0*mac_log
1533          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
1534 
1535          ecrs=ecrs0
1536          decrs_drs=decrs0_drs
1537          decrs_dzeta=0.0_dp
1538          d2ecrs_drs2=d2ecrs0_drs2
1539          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
1540          d2ecrs_drsdzeta=zero
1541          zeta=0.0_dp
1542 
1543 
1544 !        Add LSD correlation functional to GGA exchange functional
1545          exci(ipts)=exci(ipts)+ecrs
1546          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
1547 
1548          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
1549 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
1550          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
1551          dvxci(ipts,9)=d2ecrs_drho2
1552          dvxci(ipts,10)=d2ecrs_drho2
1553          dvxci(ipts,11)=d2ecrs_drho2
1554 
1555 !        -----------------------------------------------------------------------------
1556 !        Eventually add the GGA correlation part of the PBE functional
1557 !        Note : the computation of the potential in the spin-unpolarized
1558 !        case could be optimized much further. Other optimizations are left to do.
1559 
1560          phi_zeta=1.0_dp
1561          phip_zeta=0.0_dp
1562          phi_zeta_inv=1.0_dp
1563          phi_logder=0.0_dp
1564          phi3_zeta=1.0_dp
1565          gamphi3inv=gamma_inv
1566          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
1567 
1568 !        From ec to bb
1569          bb=ecrs*gamphi3inv
1570          dbb_drs=decrs_drs*gamphi3inv
1571          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
1572          d2bb_drs2=d2ecrs_drs2*gamphi3inv
1573          d2bb_drsdzeta=gamphi3inv*(d2ecrs_drsdzeta-three*decrs_drs*phi_logder)
1574          d2bb_dzeta2=gamphi3inv*(d2ecrs_dzeta2-six*decrs_dzeta*phi_logder+&
1575 &         12.0_dp*ecrs*phi_logder*phi_logder-three*ecrs*phi_zeta_inv*phipp_zeta)
1576 
1577 !        From bb to cc
1578          exp_pbe=exp(-bb)
1579          cc=one/(exp_pbe-one)
1580          dcc_dbb=cc*cc*exp_pbe
1581          dcc_drs=dcc_dbb*dbb_drs
1582          dcc_dzeta=dcc_dbb*dbb_dzeta
1583          d2cc_dbb2=cc*cc*exp_pbe*(two*cc*exp_pbe-one)
1584          d2cc_drs2=d2cc_dbb2*dbb_drs*dbb_drs+dcc_dbb*d2bb_drs2
1585          d2cc_drsdzeta=d2cc_dbb2*dbb_drs*dbb_dzeta+dcc_dbb*d2bb_drsdzeta
1586          d2cc_dzeta2=d2cc_dbb2*dbb_dzeta*dbb_dzeta+dcc_dbb*d2bb_dzeta2
1587 
1588 !        From cc to aa
1589          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
1590          aa=coeff_aa*cc
1591          daa_drs=coeff_aa*dcc_drs
1592          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
1593          d2aa_drs2=coeff_aa*d2cc_drs2
1594          d2aa_drsdzeta=-two*daa_drs*phi_logder+coeff_aa*d2cc_drsdzeta
1595          d2aa_dzeta2=aa*(-two*phi_zeta_inv*phipp_zeta+six*phi_logder*phi_logder)+&
1596 &         coeff_aa*(-four*dcc_dzeta*phi_logder+d2cc_dzeta2)
1597 
1598 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
1599          grrho2=four*grho2_updn(ipts,1)
1600          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
1601 !        Note that tt is (the t variable of PBE divided by phi) squared
1602          tt=half*grrho2*dtt_dg
1603 
1604 !        Get xx from aa and tt
1605          xx=aa*tt
1606          dxx_drs=daa_drs*tt
1607          dxx_dzeta=daa_dzeta*tt
1608          dxx_dtt=aa
1609          d2xx_drs2=d2aa_drs2*tt
1610          d2xx_drsdzeta=d2aa_drsdzeta*tt
1611          d2xx_drsdtt=daa_drs
1612          d2xx_dttdzeta=daa_dzeta
1613          d2xx_dzeta2=d2aa_dzeta2*tt
1614 
1615 !        From xx to pade
1616          pade_den=one/(one+xx*(one+xx))
1617          pade=(one+xx)*pade_den
1618          dpade_dxx=-xx*(two+xx)*pade_den**2
1619          dpade_drs=dpade_dxx*dxx_drs
1620          dpade_dtt=dpade_dxx*dxx_dtt
1621          dpade_dzeta=dpade_dxx*dxx_dzeta
1622          d2pade_dxx2=two*(-one+xx*xx*(three+xx))*pade_den*pade_den*pade_den
1623          d2pade_drs2=d2pade_dxx2*dxx_drs*dxx_drs+dpade_dxx*d2xx_drs2
1624          d2pade_drsdtt=d2pade_dxx2*dxx_drs*dxx_dtt+dpade_dxx*d2xx_drsdtt
1625          d2pade_drsdzeta=d2pade_dxx2*dxx_drs*dxx_dzeta+dpade_dxx*d2xx_drsdzeta
1626          d2pade_dtt2=d2pade_dxx2*dxx_dtt*dxx_dtt
1627          d2pade_dttdzeta=d2pade_dxx2*dxx_dtt*dxx_dzeta+dpade_dxx*d2xx_dttdzeta
1628          d2pade_dzeta2=d2pade_dxx2*dxx_dzeta*dxx_dzeta+dpade_dxx*d2xx_dzeta2
1629 
1630 !        From pade to qq
1631          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
1632          qq=coeff_qq*pade
1633          dqq_drs=coeff_qq*dpade_drs
1634          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
1635          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
1636          d2qq_drs2=coeff_qq*d2pade_drs2
1637          d2qq_drsdtt=phi_zeta_inv*phi_zeta_inv*(dpade_drs+tt*d2pade_drsdtt)
1638          d2qq_drsdzeta=coeff_qq*(d2pade_drsdzeta-two*dpade_drs*phi_logder)
1639          d2qq_dtt2=phi_zeta_inv*phi_zeta_inv*(two*dpade_dtt+tt*d2pade_dtt2)
1640          d2qq_dttdzeta=phi_zeta_inv*phi_zeta_inv*(dpade_dzeta-two*pade*phi_logder)+&
1641 &         coeff_qq*(d2pade_dttdzeta-two*dpade_dtt*phi_logder)
1642          d2qq_dzeta2=coeff_qq*( d2pade_dzeta2-four*dpade_dzeta*phi_logder &
1643 &         +six*pade*phi_logder*phi_logder            &
1644 &         -two*pade*phi_zeta_inv*phipp_zeta)
1645 
1646 !        From qq to rr
1647          arg_rr=one+beta*gamma_inv*qq
1648          div_rr=one/arg_rr
1649          rr=gamma*log(arg_rr)
1650          drr_dqq=beta*div_rr
1651          drr_drs=drr_dqq*dqq_drs
1652          drr_dtt=drr_dqq*dqq_dtt
1653          drr_dzeta=drr_dqq*dqq_dzeta
1654          d2rr_dqq2=-div_rr**2*beta*beta*gamma_inv
1655          d2rr_drs2=d2rr_dqq2*dqq_drs*dqq_drs+drr_dqq*d2qq_drs2
1656          d2rr_drsdtt=d2rr_dqq2*dqq_drs*dqq_dtt+drr_dqq*d2qq_drsdtt
1657          d2rr_drsdzeta=d2rr_dqq2*dqq_drs*dqq_dzeta+drr_dqq*d2qq_drsdzeta
1658          d2rr_dtt2=d2rr_dqq2*dqq_dtt*dqq_dtt+drr_dqq*d2qq_dtt2
1659          d2rr_dttdzeta=d2rr_dqq2*dqq_dtt*dqq_dzeta+drr_dqq*d2qq_dttdzeta
1660          d2rr_dzeta2=d2rr_dqq2*dqq_dzeta*dqq_dzeta+drr_dqq*d2qq_dzeta2
1661 
1662 !        From rr to hh
1663          hh=phi3_zeta*rr
1664          dhh_drs=phi3_zeta*drr_drs
1665          dhh_dtt=phi3_zeta*drr_dtt
1666          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
1667          d2hh_drs2=phi3_zeta*d2rr_drs2
1668          d2hh_drsdtt=phi3_zeta*d2rr_drsdtt
1669          d2hh_drsdzeta=phi3_zeta*(d2rr_drsdzeta+three*drr_drs*phi_logder)
1670          d2hh_dtt2=phi3_zeta*d2rr_dtt2
1671          d2hh_dttdzeta=phi3_zeta*(d2rr_dttdzeta+three*drr_dtt*phi_logder)
1672          d2hh_dzeta2=phi3_zeta*(six*rr*phi_logder*phi_logder+&
1673 &         six*phi_logder*drr_dzeta+d2rr_dzeta2)  &
1674 &         +three*phi_zeta*phi_zeta*rr*phipp_zeta
1675 
1676 
1677 !        The GGA correlation energy is added
1678          exci(ipts)=exci(ipts)+hh
1679 
1680 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
1681 
1682 !        From hh to the derivative of the energy wrt the density
1683          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
1684          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
1685 
1686 !        From hh to the derivative of the energy wrt to the gradient of the
1687 !        density, divided by the gradient of the density
1688 !        (The v3.3 definition includes the division by the norm of the gradient)
1689          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
1690 
1691          d2rhohh_drho2=rhotot_inv*&
1692 &         (-two*ninth*rs*dhh_drs +seven*four*ninth*tt*dhh_dtt &
1693 &         +ninth*rs*rs*d2hh_drs2+zeta*zeta*d2hh_dzeta2+(seven*third*tt)**2*d2hh_dtt2 &
1694 &         +two*third*rs*zeta*d2hh_drsdzeta+two*seven*ninth*rs*tt*d2hh_drsdtt &
1695 &         +two*seven*third*tt*zeta*d2hh_dttdzeta)
1696          d2rhohh_drhodg=dtt_dg*(-four*third*dhh_dtt-third*rs*d2hh_drsdtt &
1697 &         -zeta*d2hh_dttdzeta-seven*third*tt*d2hh_dtt2)
1698 
1699 !        Component 12 : first derivative with respect to the gradient
1700 !        of the density, div by the grad of the density
1701          dvxci(ipts,12)=dvxcdgr(ipts,3)
1702 !        Components 9, 10 and 11 : second derivatives with respect to the spin-density
1703 !        Note that there is already a contribution from LSDA
1704          dvxci(ipts,9)=dvxci(ipts,9)+d2rhohh_drho2+rhotot_inv*           &
1705 &         ( d2hh_dzeta2*(one-two*zeta) &
1706 &         -two*third*rs*d2hh_drsdzeta-14.0_dp*third*tt*d2hh_dttdzeta)
1707          dvxci(ipts,10)=dvxci(ipts,10)+d2rhohh_drho2-rhotot_inv*d2hh_dzeta2
1708          dvxci(ipts,11)=dvxci(ipts,11)+d2rhohh_drho2+rhotot_inv*           &
1709 &         ( d2hh_dzeta2*(one+two*zeta) &
1710 &         +two*third*rs*d2hh_drsdzeta+14.0_dp*third*tt*d2hh_dttdzeta)
1711 !        Components 13 and 14 : second derivatives with respect to spin density
1712 !        and gradient, divided by the gradient
1713          dvxci(ipts,13)=d2rhohh_drhodg+dtt_dg*d2hh_dttdzeta
1714          dvxci(ipts,14)=d2rhohh_drhodg-dtt_dg*d2hh_dttdzeta
1715 !        Component 15 : derivative of the (derivative wrt the gradient div by the grad),
1716 !        divided by the grad
1717          dvxci(ipts,15)=rhotot*d2hh_dtt2*dtt_dg*dtt_dg
1718 
1719 
1720 !        End condition of GGA
1721 
1722 
1723 !        Correlation has been added
1724 !        -----------------------------------------------------------------------------
1725 
1726 !        vxci(ipts,2)=vxci(ipts,1)
1727          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
1728        end do
1729 
1730      else if (option==-1) then
1731 
1732        do ipts=1,npts
1733 
1734          rhotot=rhoarr(ipts)
1735          rhotmot=rhom1_3(ipts)
1736          rhotot_inv=rhotmot*rhotmot*rhotmot
1737          rhotmo6=sqrt(rhotmot)
1738          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1739 !        -----------------------------------------------------------------------
1740 !        First take care of the exchange part of the functional
1741 
1742          exc=zero
1743 !        loop over the spin
1744          ispden=1
1745          rho   =rho_updn(ipts,ispden)
1746          rhomot=rho_updnm1_3(ipts,ispden)
1747          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1748 !        Perdew_Wang 91 LSD
1749          vxci(ipts,ispden)=four_thirds*ex_lsd
1750          exc=exc+ex_lsd*rho
1751 !        Perdew_Wang 91 LSD
1752          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
1753 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
1754          dvxci(ipts,2)=zero
1755 !        If non-spin-polarized, first component of dvxci is second
1756 !        derivative with respect to TOTAL density.
1757          dvxci(ipts,1)=dvxci(ipts,1)*half
1758          if(order==3)then
1759 !          Compute the second derivative of vx
1760 !          vx^(2) = -2*vx^(1)/(3*rhotot)
1761            d2vxci(ipts,1) = -2._dp*dvxci(ipts,1)/(3._dp*rhotot)
1762          end if
1763 !        end of loop over the spin
1764 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1765          exc=exc*2
1766          exci(ipts)=exc*rhotot_inv
1767 !        -----------------------------------------------------------------------------
1768 !        Then takes care of the LSD correlation part of the functional
1769 
1770 !        Correlation has been added
1771 !        -----------------------------------------------------------------------------
1772 
1773        end do
1774      else if (option==-2) then
1775 
1776        do ipts=1,npts
1777 
1778          rhotot=rhoarr(ipts)
1779          rhotmot=rhom1_3(ipts)
1780          rhotot_inv=rhotmot*rhotmot*rhotmot
1781          rhotmo6=sqrt(rhotmot)
1782          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1783 !        -----------------------------------------------------------------------
1784 !        First take care of the exchange part of the functional
1785 
1786          exc=zero
1787          dvxcdgr(ipts,3)=zero
1788 !        loop over the spin
1789          ispden=1
1790          rho   =rho_updn(ipts,ispden)
1791          rhomot=rho_updnm1_3(ipts,ispden)
1792          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1793 !        Perdew-Burke-Ernzerhof GGA, exchange part
1794          rho_inv=rhomot*rhomot*rhomot
1795          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
1796          ss=grho2_updn(ipts,ispden)*coeffss
1797          divss=one/(one+mu_divkappa*ss)
1798          dfxdss= mu*divss*divss
1799          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
1800          fx    = one+kappa*(one-divss)
1801          ex_gga= ex_lsd*fx
1802          dssdn=-eight*third*ss*rho_inv
1803          dfxdn  = dfxdss*dssdn
1804          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
1805 !        The new definition (v3.3) includes the division by the norm of the gradient
1806          dssdg =two*coeffss
1807          dfxdg=dfxdss*dssdg
1808          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
1809          exc=exc+ex_gga*rho
1810 
1811 !        Perdew-Burke-Ernzerhof GGA, exchange part
1812 !        Components 3 or 4
1813          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
1814 !        Components 1 or 2
1815          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
1816          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
1817          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
1818 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
1819 !        Components 5 or 6
1820          d2ssdndg=-eight*third*dssdg*rho_inv
1821          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
1822          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
1823 !        Components 7 or 8
1824          d2fxdg2=d2fxdss2*dssdg**2
1825          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
1826 !        For the time being, treat non-spin-polarized like spin-polarized
1827          dvxci(ipts,2)=dvxci(ipts,1)
1828          dvxci(ipts,4)=dvxci(ipts,3)
1829          dvxci(ipts,6)=dvxci(ipts,5)
1830          dvxci(ipts,8)=dvxci(ipts,7)
1831 
1832 !        end of loop over the spin
1833 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1834          exc=exc*2
1835          exci(ipts)=exc*rhotot_inv
1836 !        -----------------------------------------------------------------------------
1837 !        Then takes care of the LSD correlation part of the functional
1838 
1839 !        Correlation has been added
1840 !        -----------------------------------------------------------------------------
1841 
1842 !        vxci(ipts,2)=vxci(ipts,1)
1843          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
1844 
1845        end do
1846 
1847 
1848      else if(option==-4) then
1849 
1850 
1851        do ipts=1,npts
1852 
1853          rhotot=rhoarr(ipts)
1854          rhotmot=rhom1_3(ipts)
1855          rhotot_inv=rhotmot*rhotmot*rhotmot
1856          rhotmo6=sqrt(rhotmot)
1857          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1858 !        -----------------------------------------------------------------------
1859 !        First take care of the exchange part of the functional
1860 
1861          exc=zero
1862          dvxcdgr(ipts,3)=zero
1863 !        loop over the spin
1864          ispden=1
1865          rho   =rho_updn(ipts,ispden)
1866          rhomot=rho_updnm1_3(ipts,ispden)
1867          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1868 !        VALENTINO R. COOPER C09x GGA, This is an exchange term proposed
1869 !        to use together with vdw-DF (see above).
1870          rho_inv=rhomot*rhomot*rhomot
1871          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
1872 !        the quarter that is lacking is compensated by the grho2_updn in the
1873 !        next line.
1874          ss=grho2_updn(ipts,ispden)*coeffss
1875          alphs2=alpha_c09*ss
1876          alphmu=alpha_c09*mu_c09
1877          dfxdss= mu_c09*exp(-alphs2)*(one-alphs2)+&
1878 &         kappa*alpha_c09*exp(-alphs2/two)/two
1879          d2fxdss2=-alphmu*exp(-alphs2)*(two-alphs2)-&
1880 &         kappa*(alpha_c09**two)*exp(alphs2/two)/four
1881          fx    = one+mu_c09*ss*exp(-alphs2)+kappa*(one-exp(-alphs2/two))
1882          ex_gga= ex_lsd*fx
1883          dssdn=-eight*third*ss*rho_inv
1884          dfxdn  = dfxdss*dssdn
1885          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
1886 !        The new definition (v3.3) includes the division by the norm of the gradient
1887          dssdg =two*coeffss
1888          dfxdg=dfxdss*dssdg
1889          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
1890          exc=exc+ex_gga*rho
1891 
1892 !        Cooper C09x GGA exchange
1893 !        Components 3 or 4
1894          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
1895 !        Components 1 or 2
1896          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
1897          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
1898          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
1899 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
1900 !        Components 5 or 6
1901          d2ssdndg=-eight*third*dssdg*rho_inv
1902          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
1903          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
1904 !        Components 7 or 8
1905          d2fxdg2=d2fxdss2*dssdg**2
1906          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
1907 !        For the time being, treat non-spin-polarized like spin-polarized
1908          dvxci(ipts,2)=dvxci(ipts,1)
1909          dvxci(ipts,4)=dvxci(ipts,3)
1910          dvxci(ipts,6)=dvxci(ipts,5)
1911          dvxci(ipts,8)=dvxci(ipts,7)
1912 
1913 !        end of loop over the spin
1914 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1915          exc=exc*2
1916          exci(ipts)=exc*rhotot_inv
1917 !        -----------------------------------------------------------------------------
1918 !        Then takes care of the LSD correlation part of the functional
1919 
1920 !        Correlation has been added
1921 !        -----------------------------------------------------------------------------
1922 
1923 !        vxci(ipts,2)=vxci(ipts,1)
1924          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
1925 
1926        end do
1927 
1928 
1929      else if(option==1) then
1930 
1931        do ipts=1,npts
1932 
1933 
1934          rhotot=rhoarr(ipts)
1935          rhotmot=rhom1_3(ipts)
1936          rhotot_inv=rhotmot*rhotmot*rhotmot
1937          rhotmo6=sqrt(rhotmot)
1938          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
1939 !        -----------------------------------------------------------------------
1940 !        First take care of the exchange part of the functional
1941 
1942          exc=zero
1943 !        loop over the spin
1944          ispden=1
1945          rho   =rho_updn(ipts,ispden)
1946          rhomot=rho_updnm1_3(ipts,ispden)
1947          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
1948 !        Perdew_Wang 91 LSD
1949          vxci(ipts,ispden)=four_thirds*ex_lsd
1950          exc=exc+ex_lsd*rho
1951 
1952 !        Perdew_Wang 91 LSD
1953          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
1954 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
1955          dvxci(ipts,2)=zero
1956 !        If non-spin-polarized, first component of dvxci is second
1957 !        derivative with respect to TOTAL density.
1958          dvxci(ipts,1)=dvxci(ipts,1)*half
1959          if(order==3)then
1960 !          Compute the second derivative of vx
1961 !          vx^(2) = -2*vx^(1)/(3*rhotot)
1962            d2vxci(ipts,1) = -2._dp*dvxci(ipts,1)/(3._dp*rhotot)
1963          end if
1964 !        end of loop over the spin
1965 !        If non spin-polarized, treat spin down contribution now, similar to spin up
1966          exc=exc*2
1967          exci(ipts)=exc*rhotot_inv
1968 !        -----------------------------------------------------------------------------
1969 !        Then takes care of the LSD correlation part of the functional
1970 
1971          rs=rsfac*rhotmot
1972          sqr_rs=sq_rsfac*rhotmo6
1973          rsm1_2=sq_rsfac_inv*rhoto6
1974 
1975 !        Formulas A6-A8 of PW92LSD
1976          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
1977          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
1978          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
1979          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
1980 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
1981          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
1982          ecrs0=ec0_q0*ec0_log
1983          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
1984          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
1985          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
1986 &         -ec0_q0*ec0_q1pp*ec0_den                        &
1987 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
1988          ec0_q1ppp = 0.75_dp*ec0_aa*(rsm1_2**5)*(ec0_b1-ec0_b3*rs)
1989          ec0_f1 = 1._dp/(ec0_q1*ec0_q1*(1._dp + ec0_q1))
1990          ec0_f2 = 1._dp/(ec0_q1*(1+ec0_q1))
1991          d3ecrs0_drs3 = 6._dp*ec0_q1p*ec0_f1*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + &
1992 &         ec0_q0*ec0_q1pp) - &
1993 &         ec0_f2*(-6._dp*ec0_aa*ec0_a1*ec0_q1pp + ec0_q0*ec0_q1ppp + &
1994 &         ec0_f2*(3._dp*ec0_q1p*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + ec0_q0*ec0_q1pp) + &
1995 &         ec0_f2*2._dp*ec0_q0*(ec0_q1p**3)*(1._dp + 3._dp*ec0_q1*(1._dp + ec0_q1))))
1996 
1997          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
1998          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
1999          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
2000          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
2001          mac_log=-log( mac_q1*mac_q1*mac_den )
2002          macrs=mac_q0*mac_log
2003          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
2004 
2005          ecrs=ecrs0
2006          decrs_drs=decrs0_drs
2007          decrs_dzeta=0.0_dp
2008          d2ecrs_drs2=d2ecrs0_drs2
2009          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
2010          d2ecrs_drsdzeta=zero
2011          zeta=0.0_dp
2012 
2013 
2014 !        Add LSD correlation functional to GGA exchange functional
2015          exci(ipts)=exci(ipts)+ecrs
2016          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
2017 
2018          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
2019 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
2020          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
2021          dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drho2
2022          d2vcrs_drs2 = third*(d2ecrs_drs2 - rs*d3ecrs0_drs3)
2023          drs_dn = -1._dp*four_pi*ninth*rs**4
2024          d2rs_dn2 = 64._dp*pi*pi*(rs**7)/81._dp
2025          if(order==3)then
2026            d2vxci(ipts,1) = d2vxci(ipts,1) + d2vcrs_drs2*drs_dn*drs_dn + &
2027 &           dvcrs_drs*d2rs_dn2
2028 
2029 
2030 
2031 
2032 
2033          end if
2034 
2035 !        -----------------------------------------------------------------------------
2036 !        Eventually add the GGA correlation part of the PBE functional
2037 !        Note : the computation of the potential in the spin-unpolarized
2038 !        case could be optimized much further. Other optimizations are left to do.
2039 
2040 
2041 !        Correlation has been added
2042 !        -----------------------------------------------------------------------------
2043        end do
2044      else if (option==3) then
2045 
2046        do ipts=1,npts
2047 
2048          rhotot=rhoarr(ipts)
2049          rhotmot=rhom1_3(ipts)
2050          rhotot_inv=rhotmot*rhotmot*rhotmot
2051          rhotmo6=sqrt(rhotmot)
2052          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
2053 !        -----------------------------------------------------------------------
2054 !        First take care of the exchange part of the functional
2055 
2056          exc=zero
2057 !        loop over the spin
2058          ispden=1
2059          rho   =rho_updn(ipts,ispden)
2060          rhomot=rho_updnm1_3(ipts,ispden)
2061          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
2062 !        Perdew_Wang 91 LSD
2063          vxci(ipts,ispden)=four_thirds*ex_lsd
2064          exc=exc+ex_lsd*rho
2065 
2066 
2067 !        Perdew_Wang 91 LSD
2068          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
2069 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
2070          dvxci(ipts,2)=zero
2071 !        If non-spin-polarized, first component of dvxci is second
2072 !        derivative with respect to TOTAL density.
2073          dvxci(ipts,1)=dvxci(ipts,1)*half
2074          if(order==3)then
2075 !          Compute the second derivative of vx
2076 !          vx^(2) = -2*vx^(1)/(3*rhotot)
2077            d2vxci(ipts,1) = -2._dp*dvxci(ipts,1)/(3._dp*rhotot)
2078          end if
2079 !        end of loop over the spin
2080 !        If non spin-polarized, treat spin down contribution now, similar to spin up
2081          exc=exc*2
2082          exci(ipts)=exc*rhotot_inv
2083 !        -----------------------------------------------------------------------------
2084 !        Then takes care of the LSD correlation part of the functional
2085 
2086          rs=rsfac*rhotmot
2087          sqr_rs=sq_rsfac*rhotmo6
2088          rsm1_2=sq_rsfac_inv*rhoto6
2089 
2090 !        Formulas A6-A8 of PW92LSD
2091          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
2092          sqr_sqr_rs=max(1.e-15_dp,sqrt(sqr_rs))
2093          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs/sqr_sqr_rs)
2094          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+3.5_dp*ec0_b4*rs/sqr_sqr_rs)
2095          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
2096 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
2097          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
2098          ecrs0=ec0_q0*ec0_log
2099          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
2100          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
2101          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
2102 &         -ec0_q0*ec0_q1pp*ec0_den                        &
2103 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
2104          ec0_q1ppp = 0.75_dp*ec0_aa*(rsm1_2**5)*(ec0_b1-ec0_b3*rs)
2105          ec0_f1 = 1._dp/(ec0_q1*ec0_q1*(1._dp + ec0_q1))
2106          ec0_f2 = 1._dp/(ec0_q1*(1+ec0_q1))
2107          d3ecrs0_drs3 = 6._dp*ec0_q1p*ec0_f1*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + &
2108 &         ec0_q0*ec0_q1pp) - &
2109 &         ec0_f2*(-6._dp*ec0_aa*ec0_a1*ec0_q1pp + ec0_q0*ec0_q1ppp + &
2110 &         ec0_f2*(3._dp*ec0_q1p*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + ec0_q0*ec0_q1pp) + &
2111 &         ec0_f2*2._dp*ec0_q0*(ec0_q1p**3)*(1._dp + 3._dp*ec0_q1*(1._dp + ec0_q1))))
2112 
2113          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
2114          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
2115          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
2116          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
2117          mac_log=-log( mac_q1*mac_q1*mac_den )
2118          macrs=mac_q0*mac_log
2119          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
2120 
2121          ecrs=ecrs0
2122          decrs_drs=decrs0_drs
2123          decrs_dzeta=0.0_dp
2124          d2ecrs_drs2=d2ecrs0_drs2
2125          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
2126          d2ecrs_drsdzeta=zero
2127          zeta=0.0_dp
2128 
2129 
2130 !        Add LSD correlation functional to GGA exchange functional
2131          exci(ipts)=exci(ipts)+ecrs
2132          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
2133 
2134          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
2135 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
2136          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
2137 
2138          dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drho2
2139          d2vcrs_drs2 = third*(d2ecrs_drs2 - rs*d3ecrs0_drs3)
2140          drs_dn = -1._dp*four_pi*ninth*rs**4
2141          d2rs_dn2 = 64._dp*pi*pi*(rs**7)/81._dp
2142          if(order==3)then
2143            d2vxci(ipts,1) = d2vxci(ipts,1) + d2vcrs_drs2*drs_dn*drs_dn + &
2144 &           dvcrs_drs*d2rs_dn2
2145 
2146 
2147 
2148          end if
2149 
2150 
2151        end do
2152 
2153      end if
2154 
2155    else if(order==-2) then
2156      if(option==2 .or. option==5) then
2157        do ipts=1,npts
2158 
2159          rhotot=rhoarr(ipts)
2160          rhotmot=rhom1_3(ipts)
2161          rhotot_inv=rhotmot*rhotmot*rhotmot
2162          rhotmo6=sqrt(rhotmot)
2163          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
2164 !        -----------------------------------------------------------------------
2165 !        First take care of the exchange part of the functional
2166 
2167          exc=zero
2168          dvxcdgr(ipts,3)=zero
2169 !        loop over the spin
2170          ispden=1
2171          rho   =rho_updn(ipts,ispden)
2172          rhomot=rho_updnm1_3(ipts,ispden)
2173          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
2174 !        Perdew-Burke-Ernzerhof GGA, exchange part
2175          rho_inv=rhomot*rhomot*rhomot
2176          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
2177          ss=grho2_updn(ipts,ispden)*coeffss
2178          divss=one/(one+mu_divkappa*ss)
2179          dfxdss= mu*divss*divss
2180          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
2181          fx    = one+kappa*(one-divss)
2182          ex_gga= ex_lsd*fx
2183          dssdn=-eight*third*ss*rho_inv
2184          dfxdn  = dfxdss*dssdn
2185          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
2186 !        The new definition (v3.3) includes the division by the norm of the gradient
2187          dssdg =two*coeffss
2188          dfxdg=dfxdss*dssdg
2189          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
2190          exc=exc+ex_gga*rho
2191 
2192 !        Perdew-Burke-Ernzerhof GGA, exchange part
2193 !        Components 3 or 4
2194          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
2195 !        Components 1 or 2
2196          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
2197          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
2198          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
2199 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
2200 !        Components 5 or 6
2201          d2ssdndg=-eight*third*dssdg*rho_inv
2202          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
2203          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
2204 !        Components 7 or 8
2205          d2fxdg2=d2fxdss2*dssdg**2
2206          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
2207 !        For the time being, treat non-spin-polarized like spin-polarized
2208          dvxci(ipts,2)=dvxci(ipts,1)
2209          dvxci(ipts,4)=dvxci(ipts,3)
2210          dvxci(ipts,6)=dvxci(ipts,5)
2211          dvxci(ipts,8)=dvxci(ipts,7)
2212 !        end of loop over the spin
2213 !        If non spin-polarized, treat spin down contribution now, similar to spin up
2214          exc=exc*2
2215          exci(ipts)=exc*rhotot_inv
2216 !        -----------------------------------------------------------------------------
2217 !        Then takes care of the LSD correlation part of the functional
2218 
2219 
2220          rs=rsfac*rhotmot
2221          sqr_rs=sq_rsfac*rhotmo6
2222          rsm1_2=sq_rsfac_inv*rhoto6
2223 
2224 !        Formulas A6-A8 of PW92LSD
2225          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
2226          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
2227          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
2228          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
2229 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
2230          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
2231          ecrs0=ec0_q0*ec0_log
2232          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
2233          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
2234          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
2235 &         -ec0_q0*ec0_q1pp*ec0_den                        &
2236 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
2237 
2238          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
2239          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
2240          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
2241          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
2242          mac_log=-log( mac_q1*mac_q1*mac_den )
2243          macrs=mac_q0*mac_log
2244          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
2245 
2246          ecrs=ecrs0
2247          decrs_drs=decrs0_drs
2248          decrs_dzeta=0.0_dp
2249          d2ecrs_drs2=d2ecrs0_drs2
2250          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
2251          d2ecrs_drsdzeta=zero
2252          zeta=0.0_dp
2253 
2254 
2255 !        Add LSD correlation functional to GGA exchange functional
2256          exci(ipts)=exci(ipts)+ecrs
2257          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
2258 
2259          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
2260 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
2261          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
2262          dvxci(ipts,9)=d2ecrs_drho2
2263          dvxci(ipts,10)=d2ecrs_drho2
2264          dvxci(ipts,11)=d2ecrs_drho2
2265 
2266 !        -----------------------------------------------------------------------------
2267 !        Eventually add the GGA correlation part of the PBE functional
2268 !        Note : the computation of the potential in the spin-unpolarized
2269 !        case could be optimized much further. Other optimizations are left to do.
2270 
2271          phi_zeta=1.0_dp
2272          phip_zeta=0.0_dp
2273          phi_zeta_inv=1.0_dp
2274          phi_logder=0.0_dp
2275          phi3_zeta=1.0_dp
2276          gamphi3inv=gamma_inv
2277          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
2278 
2279 !        From ec to bb
2280          bb=ecrs*gamphi3inv
2281          dbb_drs=decrs_drs*gamphi3inv
2282          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
2283          d2bb_drs2=d2ecrs_drs2*gamphi3inv
2284          d2bb_drsdzeta=gamphi3inv*(d2ecrs_drsdzeta-three*decrs_drs*phi_logder)
2285          d2bb_dzeta2=gamphi3inv*(d2ecrs_dzeta2-six*decrs_dzeta*phi_logder+&
2286 &         12.0_dp*ecrs*phi_logder*phi_logder-three*ecrs*phi_zeta_inv*phipp_zeta)
2287 
2288 !        From bb to cc
2289          exp_pbe=exp(-bb)
2290          cc=one/(exp_pbe-one)
2291          dcc_dbb=cc*cc*exp_pbe
2292          dcc_drs=dcc_dbb*dbb_drs
2293          dcc_dzeta=dcc_dbb*dbb_dzeta
2294          d2cc_dbb2=cc*cc*exp_pbe*(two*cc*exp_pbe-one)
2295          d2cc_drs2=d2cc_dbb2*dbb_drs*dbb_drs+dcc_dbb*d2bb_drs2
2296          d2cc_drsdzeta=d2cc_dbb2*dbb_drs*dbb_dzeta+dcc_dbb*d2bb_drsdzeta
2297          d2cc_dzeta2=d2cc_dbb2*dbb_dzeta*dbb_dzeta+dcc_dbb*d2bb_dzeta2
2298 
2299 !        From cc to aa
2300          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
2301          aa=coeff_aa*cc
2302          daa_drs=coeff_aa*dcc_drs
2303          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
2304          d2aa_drs2=coeff_aa*d2cc_drs2
2305          d2aa_drsdzeta=-two*daa_drs*phi_logder+coeff_aa*d2cc_drsdzeta
2306          d2aa_dzeta2=aa*(-two*phi_zeta_inv*phipp_zeta+six*phi_logder*phi_logder)+&
2307 &         coeff_aa*(-four*dcc_dzeta*phi_logder+d2cc_dzeta2)
2308 
2309 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
2310          grrho2=four*grho2_updn(ipts,1)
2311          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
2312 !        Note that tt is (the t variable of PBE divided by phi) squared
2313          tt=half*grrho2*dtt_dg
2314 
2315 !        Get xx from aa and tt
2316          xx=aa*tt
2317          dxx_drs=daa_drs*tt
2318          dxx_dzeta=daa_dzeta*tt
2319          dxx_dtt=aa
2320          d2xx_drs2=d2aa_drs2*tt
2321          d2xx_drsdzeta=d2aa_drsdzeta*tt
2322          d2xx_drsdtt=daa_drs
2323          d2xx_dttdzeta=daa_dzeta
2324          d2xx_dzeta2=d2aa_dzeta2*tt
2325 
2326 !        From xx to pade
2327          pade_den=one/(one+xx*(one+xx))
2328          pade=(one+xx)*pade_den
2329          dpade_dxx=-xx*(two+xx)*pade_den**2
2330          dpade_drs=dpade_dxx*dxx_drs
2331          dpade_dtt=dpade_dxx*dxx_dtt
2332          dpade_dzeta=dpade_dxx*dxx_dzeta
2333          d2pade_dxx2=two*(-one+xx*xx*(three+xx))*pade_den*pade_den*pade_den
2334          d2pade_drs2=d2pade_dxx2*dxx_drs*dxx_drs+dpade_dxx*d2xx_drs2
2335          d2pade_drsdtt=d2pade_dxx2*dxx_drs*dxx_dtt+dpade_dxx*d2xx_drsdtt
2336          d2pade_drsdzeta=d2pade_dxx2*dxx_drs*dxx_dzeta+dpade_dxx*d2xx_drsdzeta
2337          d2pade_dtt2=d2pade_dxx2*dxx_dtt*dxx_dtt
2338          d2pade_dttdzeta=d2pade_dxx2*dxx_dtt*dxx_dzeta+dpade_dxx*d2xx_dttdzeta
2339          d2pade_dzeta2=d2pade_dxx2*dxx_dzeta*dxx_dzeta+dpade_dxx*d2xx_dzeta2
2340 
2341 !        From pade to qq
2342          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
2343          qq=coeff_qq*pade
2344          dqq_drs=coeff_qq*dpade_drs
2345          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
2346          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
2347          d2qq_drs2=coeff_qq*d2pade_drs2
2348          d2qq_drsdtt=phi_zeta_inv*phi_zeta_inv*(dpade_drs+tt*d2pade_drsdtt)
2349          d2qq_drsdzeta=coeff_qq*(d2pade_drsdzeta-two*dpade_drs*phi_logder)
2350          d2qq_dtt2=phi_zeta_inv*phi_zeta_inv*(two*dpade_dtt+tt*d2pade_dtt2)
2351          d2qq_dttdzeta=phi_zeta_inv*phi_zeta_inv*(dpade_dzeta-two*pade*phi_logder)+&
2352 &         coeff_qq*(d2pade_dttdzeta-two*dpade_dtt*phi_logder)
2353          d2qq_dzeta2=coeff_qq*( d2pade_dzeta2-four*dpade_dzeta*phi_logder &
2354 &         +six*pade*phi_logder*phi_logder            &
2355 &         -two*pade*phi_zeta_inv*phipp_zeta)
2356 
2357 !        From qq to rr
2358          arg_rr=one+beta*gamma_inv*qq
2359          div_rr=one/arg_rr
2360          rr=gamma*log(arg_rr)
2361          drr_dqq=beta*div_rr
2362          drr_drs=drr_dqq*dqq_drs
2363          drr_dtt=drr_dqq*dqq_dtt
2364          drr_dzeta=drr_dqq*dqq_dzeta
2365          d2rr_dqq2=-div_rr**2*beta*beta*gamma_inv
2366          d2rr_drs2=d2rr_dqq2*dqq_drs*dqq_drs+drr_dqq*d2qq_drs2
2367          d2rr_drsdtt=d2rr_dqq2*dqq_drs*dqq_dtt+drr_dqq*d2qq_drsdtt
2368          d2rr_drsdzeta=d2rr_dqq2*dqq_drs*dqq_dzeta+drr_dqq*d2qq_drsdzeta
2369          d2rr_dtt2=d2rr_dqq2*dqq_dtt*dqq_dtt+drr_dqq*d2qq_dtt2
2370          d2rr_dttdzeta=d2rr_dqq2*dqq_dtt*dqq_dzeta+drr_dqq*d2qq_dttdzeta
2371          d2rr_dzeta2=d2rr_dqq2*dqq_dzeta*dqq_dzeta+drr_dqq*d2qq_dzeta2
2372 
2373 !        From rr to hh
2374          hh=phi3_zeta*rr
2375          dhh_drs=phi3_zeta*drr_drs
2376          dhh_dtt=phi3_zeta*drr_dtt
2377          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
2378          d2hh_drs2=phi3_zeta*d2rr_drs2
2379          d2hh_drsdtt=phi3_zeta*d2rr_drsdtt
2380          d2hh_drsdzeta=phi3_zeta*(d2rr_drsdzeta+three*drr_drs*phi_logder)
2381          d2hh_dtt2=phi3_zeta*d2rr_dtt2
2382          d2hh_dttdzeta=phi3_zeta*(d2rr_dttdzeta+three*drr_dtt*phi_logder)
2383          d2hh_dzeta2=phi3_zeta*(six*rr*phi_logder*phi_logder+&
2384 &         six*phi_logder*drr_dzeta+d2rr_dzeta2)  &
2385 &         +three*phi_zeta*phi_zeta*rr*phipp_zeta
2386 
2387 !        The GGA correlation energy is added
2388          exci(ipts)=exci(ipts)+hh
2389 
2390 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
2391 
2392 !        From hh to the derivative of the energy wrt the density
2393          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
2394          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
2395 
2396 !        From hh to the derivative of the energy wrt to the gradient of the
2397 !        density, divided by the gradient of the density
2398 !        (The v3.3 definition includes the division by the norm of the gradient)
2399          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
2400 
2401          d2rhohh_drho2=rhotot_inv*&
2402 &         (-two*ninth*rs*dhh_drs +seven*four*ninth*tt*dhh_dtt &
2403 &         +ninth*rs*rs*d2hh_drs2+zeta*zeta*d2hh_dzeta2+(seven*third*tt)**2*d2hh_dtt2 &
2404 &         +two*third*rs*zeta*d2hh_drsdzeta+two*seven*ninth*rs*tt*d2hh_drsdtt &
2405 &         +two*seven*third*tt*zeta*d2hh_dttdzeta)
2406          d2rhohh_drhodg=dtt_dg*(-four*third*dhh_dtt-third*rs*d2hh_drsdtt &
2407 &         -zeta*d2hh_dttdzeta-seven*third*tt*d2hh_dtt2)
2408 
2409 !        Component 12 : first derivative with respect to the gradient
2410 !        of the density, div by the grad of the density
2411          dvxci(ipts,12)=dvxcdgr(ipts,3)
2412 !        Components 9, 10 and 11 : second derivatives with respect to the spin-density
2413 !        Note that there is already a contribution from LSDA
2414          dvxci(ipts,9)=dvxci(ipts,9)+d2rhohh_drho2+rhotot_inv*           &
2415 &         ( d2hh_dzeta2*(one-two*zeta) &
2416 &         -two*third*rs*d2hh_drsdzeta-14.0_dp*third*tt*d2hh_dttdzeta)
2417          dvxci(ipts,10)=dvxci(ipts,10)+d2rhohh_drho2-rhotot_inv*d2hh_dzeta2
2418          dvxci(ipts,11)=dvxci(ipts,11)+d2rhohh_drho2+rhotot_inv*           &
2419 &         ( d2hh_dzeta2*(one+two*zeta) &
2420 &         +two*third*rs*d2hh_drsdzeta+14.0_dp*third*tt*d2hh_dttdzeta)
2421 !        Components 13 and 14 : second derivatives with respect to spin density
2422 !        and gradient, divided by the gradient
2423          dvxci(ipts,13)=d2rhohh_drhodg+dtt_dg*d2hh_dttdzeta
2424          dvxci(ipts,14)=d2rhohh_drhodg-dtt_dg*d2hh_dttdzeta
2425 !        Component 15 : derivative of the (derivative wrt the gradient div by the grad),
2426 !        divided by the grad
2427          dvxci(ipts,15)=rhotot*d2hh_dtt2*dtt_dg*dtt_dg
2428 
2429 
2430 !        End condition of GGA
2431 
2432 
2433 !        Correlation has been added
2434 !        -----------------------------------------------------------------------------
2435 
2436 !        vxci(ipts,2)=vxci(ipts,1)
2437          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
2438 
2439        end do
2440 
2441      else if ((option==6) .or. (option==7)) then
2442        do ipts=1,npts
2443 
2444          rhotot=rhoarr(ipts)
2445          rhotmot=rhom1_3(ipts)
2446          rhotot_inv=rhotmot*rhotmot*rhotmot
2447          rhotmo6=sqrt(rhotmot)
2448          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
2449 !        -----------------------------------------------------------------------
2450 !        First take care of the exchange part of the functional
2451 
2452          exc=zero
2453          dvxcdgr(ipts,3)=zero
2454 !        loop over the spin
2455          ispden=1
2456          rho   =rho_updn(ipts,ispden)
2457          rhomot=rho_updnm1_3(ipts,ispden)
2458          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
2459 !        Perdew-Burke-Ernzerhof GGA, exchange part
2460          rho_inv=rhomot*rhomot*rhomot
2461          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
2462          ss=grho2_updn(ipts,ispden)*coeffss
2463 
2464          if (option==6) then
2465            divss=exp(-mu_divkappa*ss)
2466            dfxdss= mu*divss
2467            d2fxdss2=-mu*mu_divkappa*divss
2468 
2469            fx    = one+kappa*(one-divss)
2470            ex_gga= ex_lsd*fx
2471            dssdn=-eight*third*ss*rho_inv
2472            dfxdn  = dfxdss*dssdn
2473            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
2474 !          The new definition (v3.3) includes the division by the norm of the gradient
2475            dssdg =two*coeffss
2476            dfxdg=dfxdss*dssdg
2477            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
2478            exc=exc+ex_gga*rho
2479 !          This is the Wu and Cohen modification
2480          else
2481            expss=exp(-ss)
2482            p1_wc=b_wc+(mu-b_wc)*(one-ss)*expss+two*c_wc*ss/(one+c_wc*ss*ss)
2483            p2_wc=d_wc*(ss-two)*expss+two*c_wc/(one+c_wc*ss*ss)-&
2484 &           four*c_wc*c_wc*ss*ss/((one+c_wc*ss*ss)*(one+c_wc*ss*ss))
2485            divss=one/(one+(b_wc*ss+d_wc*ss*expss+log(one+c_wc*ss*ss))/kappa)
2486            dfxdss=p1_wc*divss*divss
2487            d2fxdss2=p2_wc*divss*divss-two*divss*divss*divss*p1_wc*p1_wc/kappa
2488 
2489            fx    = one+kappa*(one-divss)
2490            ex_gga= ex_lsd*fx
2491            dssdn=-eight*third*ss*rho_inv
2492            dfxdn  = dfxdss*dssdn
2493            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
2494 !          The new definition (v3.3) includes the division by the norm of the gradient
2495            dssdg =two*coeffss
2496            dfxdg=dfxdss*dssdg
2497            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
2498            exc=exc+ex_gga*rho
2499          end if
2500 
2501 !        Perdew-Burke-Ernzerhof GGA, exchange part
2502 !        Components 3 or 4
2503          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
2504 !        Components 1 or 2
2505          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
2506          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
2507          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
2508 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
2509 !        Components 5 or 6
2510          d2ssdndg=-eight*third*dssdg*rho_inv
2511          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
2512          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
2513 !        Components 7 or 8
2514          d2fxdg2=d2fxdss2*dssdg**2
2515          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
2516 !        For the time being, treat non-spin-polarized like spin-polarized
2517          dvxci(ipts,2)=dvxci(ipts,1)
2518          dvxci(ipts,4)=dvxci(ipts,3)
2519          dvxci(ipts,6)=dvxci(ipts,5)
2520          dvxci(ipts,8)=dvxci(ipts,7)
2521 !        end of loop over the spin
2522 !        If non spin-polarized, treat spin down contribution now, similar to spin up
2523          exc=exc*2
2524          exci(ipts)=exc*rhotot_inv
2525 !        -----------------------------------------------------------------------------
2526 !        Then takes care of the LSD correlation part of the functional
2527 
2528 
2529          rs=rsfac*rhotmot
2530          sqr_rs=sq_rsfac*rhotmo6
2531          rsm1_2=sq_rsfac_inv*rhoto6
2532 
2533 !        Formulas A6-A8 of PW92LSD
2534          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
2535          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
2536          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
2537          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
2538 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
2539          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
2540          ecrs0=ec0_q0*ec0_log
2541          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
2542          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
2543          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
2544 &         -ec0_q0*ec0_q1pp*ec0_den                        &
2545 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
2546 
2547          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
2548          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
2549          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
2550          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
2551          mac_log=-log( mac_q1*mac_q1*mac_den )
2552          macrs=mac_q0*mac_log
2553          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
2554 
2555          ecrs=ecrs0
2556          decrs_drs=decrs0_drs
2557          decrs_dzeta=0.0_dp
2558          d2ecrs_drs2=d2ecrs0_drs2
2559          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
2560          d2ecrs_drsdzeta=zero
2561          zeta=0.0_dp
2562 
2563 
2564 !        Add LSD correlation functional to GGA exchange functional
2565          exci(ipts)=exci(ipts)+ecrs
2566          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
2567 
2568          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
2569 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
2570          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
2571          dvxci(ipts,9)=d2ecrs_drho2
2572          dvxci(ipts,10)=d2ecrs_drho2
2573          dvxci(ipts,11)=d2ecrs_drho2
2574 
2575 !        -----------------------------------------------------------------------------
2576 !        Eventually add the GGA correlation part of the PBE functional
2577 !        Note : the computation of the potential in the spin-unpolarized
2578 !        case could be optimized much further. Other optimizations are left to do.
2579 
2580          phi_zeta=1.0_dp
2581          phip_zeta=0.0_dp
2582          phi_zeta_inv=1.0_dp
2583          phi_logder=0.0_dp
2584          phi3_zeta=1.0_dp
2585          gamphi3inv=gamma_inv
2586          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
2587 
2588 !        From ec to bb
2589          bb=ecrs*gamphi3inv
2590          dbb_drs=decrs_drs*gamphi3inv
2591          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
2592          d2bb_drs2=d2ecrs_drs2*gamphi3inv
2593          d2bb_drsdzeta=gamphi3inv*(d2ecrs_drsdzeta-three*decrs_drs*phi_logder)
2594          d2bb_dzeta2=gamphi3inv*(d2ecrs_dzeta2-six*decrs_dzeta*phi_logder+&
2595 &         12.0_dp*ecrs*phi_logder*phi_logder-three*ecrs*phi_zeta_inv*phipp_zeta)
2596 
2597 !        From bb to cc
2598          exp_pbe=exp(-bb)
2599          cc=one/(exp_pbe-one)
2600          dcc_dbb=cc*cc*exp_pbe
2601          dcc_drs=dcc_dbb*dbb_drs
2602          dcc_dzeta=dcc_dbb*dbb_dzeta
2603          d2cc_dbb2=cc*cc*exp_pbe*(two*cc*exp_pbe-one)
2604          d2cc_drs2=d2cc_dbb2*dbb_drs*dbb_drs+dcc_dbb*d2bb_drs2
2605          d2cc_drsdzeta=d2cc_dbb2*dbb_drs*dbb_dzeta+dcc_dbb*d2bb_drsdzeta
2606          d2cc_dzeta2=d2cc_dbb2*dbb_dzeta*dbb_dzeta+dcc_dbb*d2bb_dzeta2
2607 
2608 !        From cc to aa
2609          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
2610          aa=coeff_aa*cc
2611          daa_drs=coeff_aa*dcc_drs
2612          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
2613          d2aa_drs2=coeff_aa*d2cc_drs2
2614          d2aa_drsdzeta=-two*daa_drs*phi_logder+coeff_aa*d2cc_drsdzeta
2615          d2aa_dzeta2=aa*(-two*phi_zeta_inv*phipp_zeta+six*phi_logder*phi_logder)+&
2616 &         coeff_aa*(-four*dcc_dzeta*phi_logder+d2cc_dzeta2)
2617 
2618 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
2619          grrho2=four*grho2_updn(ipts,1)
2620          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
2621 !        Note that tt is (the t variable of PBE divided by phi) squared
2622          tt=half*grrho2*dtt_dg
2623 
2624 !        Get xx from aa and tt
2625          xx=aa*tt
2626          dxx_drs=daa_drs*tt
2627          dxx_dzeta=daa_dzeta*tt
2628          dxx_dtt=aa
2629          d2xx_drs2=d2aa_drs2*tt
2630          d2xx_drsdzeta=d2aa_drsdzeta*tt
2631          d2xx_drsdtt=daa_drs
2632          d2xx_dttdzeta=daa_dzeta
2633          d2xx_dzeta2=d2aa_dzeta2*tt
2634 
2635 !        From xx to pade
2636          pade_den=one/(one+xx*(one+xx))
2637          pade=(one+xx)*pade_den
2638          dpade_dxx=-xx*(two+xx)*pade_den**2
2639          dpade_drs=dpade_dxx*dxx_drs
2640          dpade_dtt=dpade_dxx*dxx_dtt
2641          dpade_dzeta=dpade_dxx*dxx_dzeta
2642          d2pade_dxx2=two*(-one+xx*xx*(three+xx))*pade_den*pade_den*pade_den
2643          d2pade_drs2=d2pade_dxx2*dxx_drs*dxx_drs+dpade_dxx*d2xx_drs2
2644          d2pade_drsdtt=d2pade_dxx2*dxx_drs*dxx_dtt+dpade_dxx*d2xx_drsdtt
2645          d2pade_drsdzeta=d2pade_dxx2*dxx_drs*dxx_dzeta+dpade_dxx*d2xx_drsdzeta
2646          d2pade_dtt2=d2pade_dxx2*dxx_dtt*dxx_dtt
2647          d2pade_dttdzeta=d2pade_dxx2*dxx_dtt*dxx_dzeta+dpade_dxx*d2xx_dttdzeta
2648          d2pade_dzeta2=d2pade_dxx2*dxx_dzeta*dxx_dzeta+dpade_dxx*d2xx_dzeta2
2649 
2650 !        From pade to qq
2651          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
2652          qq=coeff_qq*pade
2653          dqq_drs=coeff_qq*dpade_drs
2654          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
2655          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
2656          d2qq_drs2=coeff_qq*d2pade_drs2
2657          d2qq_drsdtt=phi_zeta_inv*phi_zeta_inv*(dpade_drs+tt*d2pade_drsdtt)
2658          d2qq_drsdzeta=coeff_qq*(d2pade_drsdzeta-two*dpade_drs*phi_logder)
2659          d2qq_dtt2=phi_zeta_inv*phi_zeta_inv*(two*dpade_dtt+tt*d2pade_dtt2)
2660          d2qq_dttdzeta=phi_zeta_inv*phi_zeta_inv*(dpade_dzeta-two*pade*phi_logder)+&
2661 &         coeff_qq*(d2pade_dttdzeta-two*dpade_dtt*phi_logder)
2662          d2qq_dzeta2=coeff_qq*( d2pade_dzeta2-four*dpade_dzeta*phi_logder &
2663 &         +six*pade*phi_logder*phi_logder            &
2664 &         -two*pade*phi_zeta_inv*phipp_zeta)
2665 
2666 !        From qq to rr
2667          arg_rr=one+beta*gamma_inv*qq
2668          div_rr=one/arg_rr
2669          rr=gamma*log(arg_rr)
2670          drr_dqq=beta*div_rr
2671          drr_drs=drr_dqq*dqq_drs
2672          drr_dtt=drr_dqq*dqq_dtt
2673          drr_dzeta=drr_dqq*dqq_dzeta
2674          d2rr_dqq2=-div_rr**2*beta*beta*gamma_inv
2675          d2rr_drs2=d2rr_dqq2*dqq_drs*dqq_drs+drr_dqq*d2qq_drs2
2676          d2rr_drsdtt=d2rr_dqq2*dqq_drs*dqq_dtt+drr_dqq*d2qq_drsdtt
2677          d2rr_drsdzeta=d2rr_dqq2*dqq_drs*dqq_dzeta+drr_dqq*d2qq_drsdzeta
2678          d2rr_dtt2=d2rr_dqq2*dqq_dtt*dqq_dtt+drr_dqq*d2qq_dtt2
2679          d2rr_dttdzeta=d2rr_dqq2*dqq_dtt*dqq_dzeta+drr_dqq*d2qq_dttdzeta
2680          d2rr_dzeta2=d2rr_dqq2*dqq_dzeta*dqq_dzeta+drr_dqq*d2qq_dzeta2
2681 
2682 !        From rr to hh
2683          hh=phi3_zeta*rr
2684          dhh_drs=phi3_zeta*drr_drs
2685          dhh_dtt=phi3_zeta*drr_dtt
2686          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
2687          d2hh_drs2=phi3_zeta*d2rr_drs2
2688          d2hh_drsdtt=phi3_zeta*d2rr_drsdtt
2689          d2hh_drsdzeta=phi3_zeta*(d2rr_drsdzeta+three*drr_drs*phi_logder)
2690          d2hh_dtt2=phi3_zeta*d2rr_dtt2
2691          d2hh_dttdzeta=phi3_zeta*(d2rr_dttdzeta+three*drr_dtt*phi_logder)
2692          d2hh_dzeta2=phi3_zeta*(six*rr*phi_logder*phi_logder+&
2693 &         six*phi_logder*drr_dzeta+d2rr_dzeta2)  &
2694 &         +three*phi_zeta*phi_zeta*rr*phipp_zeta
2695 
2696 !        The GGA correlation energy is added
2697          exci(ipts)=exci(ipts)+hh
2698 
2699 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
2700 
2701 !        From hh to the derivative of the energy wrt the density
2702          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
2703          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
2704 
2705 !        From hh to the derivative of the energy wrt to the gradient of the
2706 !        density, divided by the gradient of the density
2707 !        (The v3.3 definition includes the division by the norm of the gradient)
2708          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
2709 
2710          d2rhohh_drho2=rhotot_inv*&
2711 &         (-two*ninth*rs*dhh_drs +seven*four*ninth*tt*dhh_dtt &
2712 &         +ninth*rs*rs*d2hh_drs2+zeta*zeta*d2hh_dzeta2+(seven*third*tt)**2*d2hh_dtt2 &
2713 &         +two*third*rs*zeta*d2hh_drsdzeta+two*seven*ninth*rs*tt*d2hh_drsdtt &
2714 &         +two*seven*third*tt*zeta*d2hh_dttdzeta)
2715          d2rhohh_drhodg=dtt_dg*(-four*third*dhh_dtt-third*rs*d2hh_drsdtt &
2716 &         -zeta*d2hh_dttdzeta-seven*third*tt*d2hh_dtt2)
2717 
2718 !        Component 12 : first derivative with respect to the gradient
2719 !        of the density, div by the grad of the density
2720          dvxci(ipts,12)=dvxcdgr(ipts,3)
2721 !        Components 9, 10 and 11 : second derivatives with respect to the spin-density
2722 !        Note that there is already a contribution from LSDA
2723          dvxci(ipts,9)=dvxci(ipts,9)+d2rhohh_drho2+rhotot_inv*           &
2724 &         ( d2hh_dzeta2*(one-two*zeta) &
2725 &         -two*third*rs*d2hh_drsdzeta-14.0_dp*third*tt*d2hh_dttdzeta)
2726          dvxci(ipts,10)=dvxci(ipts,10)+d2rhohh_drho2-rhotot_inv*d2hh_dzeta2
2727          dvxci(ipts,11)=dvxci(ipts,11)+d2rhohh_drho2+rhotot_inv*           &
2728 &         ( d2hh_dzeta2*(one+two*zeta) &
2729 &         +two*third*rs*d2hh_drsdzeta+14.0_dp*third*tt*d2hh_dttdzeta)
2730 !        Components 13 and 14 : second derivatives with respect to spin density
2731 !        and gradient, divided by the gradient
2732          dvxci(ipts,13)=d2rhohh_drhodg+dtt_dg*d2hh_dttdzeta
2733          dvxci(ipts,14)=d2rhohh_drhodg-dtt_dg*d2hh_dttdzeta
2734 !        Component 15 : derivative of the (derivative wrt the gradient div by the grad),
2735 !        divided by the grad
2736          dvxci(ipts,15)=rhotot*d2hh_dtt2*dtt_dg*dtt_dg
2737 
2738 
2739 !        End condition of GGA
2740 
2741 
2742 !        Correlation has been added
2743 !        -----------------------------------------------------------------------------
2744 
2745 !        vxci(ipts,2)=vxci(ipts,1)
2746          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
2747 
2748        end do
2749 
2750      else if (option==-1) then
2751        do ipts=1,npts
2752 
2753          rhotot=rhoarr(ipts)
2754          rhotmot=rhom1_3(ipts)
2755          rhotot_inv=rhotmot*rhotmot*rhotmot
2756          rhotmo6=sqrt(rhotmot)
2757          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
2758 !        -----------------------------------------------------------------------
2759 !        First take care of the exchange part of the functional
2760 
2761          exc=zero
2762 !        loop over the spin
2763          ispden=1
2764          rho   =rho_updn(ipts,ispden)
2765          rhomot=rho_updnm1_3(ipts,ispden)
2766          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
2767 !        Perdew_Wang 91 LSD
2768          vxci(ipts,ispden)=four_thirds*ex_lsd
2769          exc=exc+ex_lsd*rho
2770 
2771 !        Perdew_Wang 91 LSD
2772          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
2773 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
2774          dvxci(ipts,2)=zero
2775 !        If non-spin-polarized, first component of dvxci is second
2776 !        derivative with respect to TOTAL density.
2777          dvxci(ipts,1)=dvxci(ipts,1)*half
2778 !        Compute the second derivative of vx
2779 !        vx^(2) = -2*vx^(1)/(3*rhotot)
2780 !        end of loop over the spin
2781 !        If non spin-polarized, treat spin down contribution now, similar to spin up
2782          exc=exc*2
2783          exci(ipts)=exc*rhotot_inv
2784 
2785 
2786 !        Correlation has been added
2787 !        -----------------------------------------------------------------------------
2788 
2789 
2790        end do
2791 
2792      else if (option==-2) then
2793        do ipts=1,npts
2794 
2795          rhotot=rhoarr(ipts)
2796          rhotmot=rhom1_3(ipts)
2797          rhotot_inv=rhotmot*rhotmot*rhotmot
2798          rhotmo6=sqrt(rhotmot)
2799          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
2800 !        -----------------------------------------------------------------------
2801 !        First take care of the exchange part of the functional
2802 
2803          exc=zero
2804          dvxcdgr(ipts,3)=zero
2805 !        loop over the spin
2806          ispden=1
2807          rho   =rho_updn(ipts,ispden)
2808          rhomot=rho_updnm1_3(ipts,ispden)
2809          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
2810 !        Perdew-Burke-Ernzerhof GGA, exchange part
2811          rho_inv=rhomot*rhomot*rhomot
2812          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
2813          ss=grho2_updn(ipts,ispden)*coeffss
2814          divss=one/(one+mu_divkappa*ss)
2815          dfxdss= mu*divss*divss
2816          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
2817          fx    = one+kappa*(one-divss)
2818          ex_gga= ex_lsd*fx
2819          dssdn=-eight*third*ss*rho_inv
2820          dfxdn  = dfxdss*dssdn
2821          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
2822 !        The new definition (v3.3) includes the division by the norm of the gradient
2823          dssdg =two*coeffss
2824          dfxdg=dfxdss*dssdg
2825          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
2826          exc=exc+ex_gga*rho
2827 
2828 !        Perdew-Burke-Ernzerhof GGA, exchange part
2829 !        Components 3 or 4
2830          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
2831 !        Components 1 or 2
2832          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
2833          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
2834          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
2835 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
2836 !        Components 5 or 6
2837          d2ssdndg=-eight*third*dssdg*rho_inv
2838          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
2839          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
2840 !        Components 7 or 8
2841          d2fxdg2=d2fxdss2*dssdg**2
2842          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
2843 !        For the time being, treat non-spin-polarized like spin-polarized
2844          dvxci(ipts,2)=dvxci(ipts,1)
2845          dvxci(ipts,4)=dvxci(ipts,3)
2846          dvxci(ipts,6)=dvxci(ipts,5)
2847          dvxci(ipts,8)=dvxci(ipts,7)
2848 !        end of loop over the spin
2849 !        If non spin-polarized, treat spin down contribution now, similar to spin up
2850          exc=exc*2
2851          exci(ipts)=exc*rhotot_inv
2852 !        -----------------------------------------------------------------------------
2853 !        Then takes care of the LSD correlation part of the functional
2854 
2855 
2856 !        Correlation has been added
2857 !        -----------------------------------------------------------------------------
2858 
2859 !        vxci(ipts,2)=vxci(ipts,1)
2860          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
2861 
2862        end do
2863 
2864      else if(option==-4) then
2865 
2866 
2867        do ipts=1,npts
2868 
2869          rhotot=rhoarr(ipts)
2870          rhotmot=rhom1_3(ipts)
2871          rhotot_inv=rhotmot*rhotmot*rhotmot
2872          rhotmo6=sqrt(rhotmot)
2873          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
2874 !        -----------------------------------------------------------------------
2875 !        First take care of the exchange part of the functional
2876 
2877          exc=zero
2878          dvxcdgr(ipts,3)=zero
2879 !        loop over the spin
2880          ispden=1
2881          rho   =rho_updn(ipts,ispden)
2882          rhomot=rho_updnm1_3(ipts,ispden)
2883          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
2884 !        VALENTINO R. COOPER C09x GGA, This is an exchange term proposed
2885 !        to use together with vdw-DF (see above).
2886          rho_inv=rhomot*rhomot*rhomot
2887          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
2888 !        the quarter that is lacking is compensated by the grho2_updn in the
2889 !        next line.
2890          ss=grho2_updn(ipts,ispden)*coeffss
2891          alphs2=alpha_c09*ss
2892          alphmu=alpha_c09*mu_c09
2893          dfxdss= mu_c09*exp(-alphs2)*(one-alphs2)+&
2894 &         kappa*alpha_c09*exp(-alphs2/two)/two
2895          d2fxdss2=-alphmu*exp(-alphs2)*(two-alphs2)-&
2896 &         kappa*(alpha_c09**two)*exp(alphs2/two)/four
2897          fx    = one+mu_c09*ss*exp(-alphs2)+kappa*(one-exp(-alphs2/two))
2898          ex_gga= ex_lsd*fx
2899          dssdn=-eight*third*ss*rho_inv
2900          dfxdn  = dfxdss*dssdn
2901          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
2902 !        The new definition (v3.3) includes the division by the norm of the gradient
2903          dssdg =two*coeffss
2904          dfxdg=dfxdss*dssdg
2905          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
2906          exc=exc+ex_gga*rho
2907 
2908 !        Cooper C09x GGA exchange
2909 !        Components 3 or 4
2910          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
2911 !        Components 1 or 2
2912          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
2913          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
2914          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
2915 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
2916 !        Components 5 or 6
2917          d2ssdndg=-eight*third*dssdg*rho_inv
2918          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
2919          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
2920 !        Components 7 or 8
2921          d2fxdg2=d2fxdss2*dssdg**2
2922          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
2923 !        For the time being, treat non-spin-polarized like spin-polarized
2924          dvxci(ipts,2)=dvxci(ipts,1)
2925          dvxci(ipts,4)=dvxci(ipts,3)
2926          dvxci(ipts,6)=dvxci(ipts,5)
2927          dvxci(ipts,8)=dvxci(ipts,7)
2928 
2929 !        end of loop over the spin
2930 !        If non spin-polarized, treat spin down contribution now, similar to spin up
2931          exc=exc*2
2932          exci(ipts)=exc*rhotot_inv
2933 !        -----------------------------------------------------------------------------
2934 !        Then takes care of the LSD correlation part of the functional
2935 
2936 !        Correlation has been added
2937 !        -----------------------------------------------------------------------------
2938 
2939 !        vxci(ipts,2)=vxci(ipts,1)
2940          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
2941 
2942        end do
2943 
2944 
2945      else if (option==1) then
2946        do ipts=1,npts
2947 
2948          rhotot=rhoarr(ipts)
2949          rhotmot=rhom1_3(ipts)
2950          rhotot_inv=rhotmot*rhotmot*rhotmot
2951          rhotmo6=sqrt(rhotmot)
2952          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
2953 !        -----------------------------------------------------------------------
2954 !        First take care of the exchange part of the functional
2955 
2956          exc=zero
2957 !        loop over the spin
2958          ispden=1
2959          rho   =rho_updn(ipts,ispden)
2960          rhomot=rho_updnm1_3(ipts,ispden)
2961          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
2962 !        Perdew_Wang 91 LSD
2963          vxci(ipts,ispden)=four_thirds*ex_lsd
2964          exc=exc+ex_lsd*rho
2965 
2966 
2967 !        Perdew_Wang 91 LSD
2968          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
2969 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
2970          dvxci(ipts,2)=zero
2971 !        If non-spin-polarized, first component of dvxci is second
2972 !        derivative with respect to TOTAL density.
2973          dvxci(ipts,1)=dvxci(ipts,1)*half
2974 !        Compute the second derivative of vx
2975 !        vx^(2) = -2*vx^(1)/(3*rhotot)
2976 !        end of loop over the spin
2977 !        If non spin-polarized, treat spin down contribution now, similar to spin up
2978          exc=exc*2
2979          exci(ipts)=exc*rhotot_inv
2980 !        -----------------------------------------------------------------------------
2981 !        Then takes care of the LSD correlation part of the functional
2982 
2983          rs=rsfac*rhotmot
2984          sqr_rs=sq_rsfac*rhotmo6
2985          rsm1_2=sq_rsfac_inv*rhoto6
2986 
2987 !        Formulas A6-A8 of PW92LSD
2988          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
2989 
2990          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
2991          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
2992          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
2993 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
2994          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
2995          ecrs0=ec0_q0*ec0_log
2996          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
2997          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
2998          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
2999 &         -ec0_q0*ec0_q1pp*ec0_den                        &
3000 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
3001 
3002          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
3003          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
3004          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
3005          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
3006          mac_log=-log( mac_q1*mac_q1*mac_den )
3007          macrs=mac_q0*mac_log
3008          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
3009 
3010          ecrs=ecrs0
3011          decrs_drs=decrs0_drs
3012          decrs_dzeta=0.0_dp
3013          d2ecrs_drs2=d2ecrs0_drs2
3014          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
3015          d2ecrs_drsdzeta=zero
3016          zeta=0.0_dp
3017 
3018 
3019 !        Add LSD correlation functional to GGA exchange functional
3020          exci(ipts)=exci(ipts)+ecrs
3021          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
3022 
3023          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
3024 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
3025          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
3026 
3027          dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drho2
3028          dvxci(ipts,2)=dvxci(ipts,2)+d2ecrs_drho2-d2ecrs_dzeta2*rhotot_inv
3029 
3030 
3031 !        Correlation has been added
3032 !        -----------------------------------------------------------------------------
3033 
3034        end do
3035 
3036      else if (option==3) then
3037        do ipts=1,npts
3038 
3039          rhotot=rhoarr(ipts)
3040          rhotmot=rhom1_3(ipts)
3041          rhotot_inv=rhotmot*rhotmot*rhotmot
3042          rhotmo6=sqrt(rhotmot)
3043          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
3044 !        -----------------------------------------------------------------------
3045 !        First take care of the exchange part of the functional
3046 
3047          exc=zero
3048 !        loop over the spin
3049          ispden=1
3050          rho   =rho_updn(ipts,ispden)
3051          rhomot=rho_updnm1_3(ipts,ispden)
3052          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
3053 !        Perdew_Wang 91 LSD
3054          vxci(ipts,ispden)=four_thirds*ex_lsd
3055          exc=exc+ex_lsd*rho
3056 
3057 !        Perdew_Wang 91 LSD
3058          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
3059 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
3060          dvxci(ipts,2)=zero
3061 !        If non-spin-polarized, first component of dvxci is second
3062 !        derivative with respect to TOTAL density.
3063          dvxci(ipts,1)=dvxci(ipts,1)*half
3064 !        Compute the second derivative of vx
3065 !        vx^(2) = -2*vx^(1)/(3*rhotot)
3066 !        end of loop over the spin
3067 !        If non spin-polarized, treat spin down contribution now, similar to spin up
3068          exc=exc*2
3069          exci(ipts)=exc*rhotot_inv
3070 !        -----------------------------------------------------------------------------
3071 !        Then takes care of the LSD correlation part of the functional
3072 
3073 
3074          rs=rsfac*rhotmot
3075          sqr_rs=sq_rsfac*rhotmo6
3076          rsm1_2=sq_rsfac_inv*rhoto6
3077 
3078 !        Formulas A6-A8 of PW92LSD
3079          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
3080          sqr_sqr_rs=max(1.e-15_dp,sqrt(sqr_rs))
3081          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs/sqr_sqr_rs)
3082          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+3.5_dp*ec0_b4*rs/sqr_sqr_rs)
3083          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
3084 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
3085          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
3086          ecrs0=ec0_q0*ec0_log
3087          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
3088          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
3089          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
3090 &         -ec0_q0*ec0_q1pp*ec0_den                        &
3091 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
3092 
3093          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
3094          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
3095          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
3096          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
3097          mac_log=-log( mac_q1*mac_q1*mac_den )
3098          macrs=mac_q0*mac_log
3099          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
3100 
3101          ecrs=ecrs0
3102          decrs_drs=decrs0_drs
3103          decrs_dzeta=0.0_dp
3104          d2ecrs_drs2=d2ecrs0_drs2
3105          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
3106          d2ecrs_drsdzeta=zero
3107          zeta=0.0_dp
3108 
3109 
3110 !        Add LSD correlation functional to GGA exchange functional
3111          exci(ipts)=exci(ipts)+ecrs
3112          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
3113 
3114          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
3115 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
3116          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
3117          dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drho2
3118          dvxci(ipts,2)=dvxci(ipts,2)+d2ecrs_drho2-d2ecrs_dzeta2*rhotot_inv
3119 
3120 !        Correlation has been added
3121 !        -----------------------------------------------------------------------------
3122 
3123        end do
3124 
3125      end if
3126 
3127 
3128    else if (order**2>1) then
3129 !    separate cases depending to option
3130      if(option==2 .or. option==5) then
3131        do ipts=1,npts
3132 
3133          rhotot=rhoarr(ipts)
3134          rhotmot=rhom1_3(ipts)
3135          rhotot_inv=rhotmot*rhotmot*rhotmot
3136          rhotmo6=sqrt(rhotmot)
3137          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
3138 !        -----------------------------------------------------------------------
3139 !        First take care of the exchange part of the functional
3140 
3141          exc=zero
3142          dvxcdgr(ipts,3)=zero
3143 !        loop over the spin
3144          ispden=1
3145          rho   =rho_updn(ipts,ispden)
3146          rhomot=rho_updnm1_3(ipts,ispden)
3147          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
3148 !        Perdew-Burke-Ernzerhof GGA, exchange part
3149          rho_inv=rhomot*rhomot*rhomot
3150          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
3151          ss=grho2_updn(ipts,ispden)*coeffss
3152 
3153          divss=one/(one+mu_divkappa*ss)
3154          dfxdss= mu*divss*divss
3155          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
3156 
3157          fx    = one+kappa*(one-divss)
3158          ex_gga= ex_lsd*fx
3159          dssdn=-eight*third*ss*rho_inv
3160          dfxdn  = dfxdss*dssdn
3161          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
3162 !        The new definition (v3.3) includes the division by the norm of the gradient
3163          dssdg =two*coeffss
3164          dfxdg=dfxdss*dssdg
3165          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
3166          exc=exc+ex_gga*rho
3167 
3168 !        Perdew-Burke-Ernzerhof GGA, exchange part
3169 !        Components 3 or 4
3170          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
3171 !        Components 1 or 2
3172          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
3173          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
3174          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
3175 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
3176 !        Components 5 or 6
3177          d2ssdndg=-eight*third*dssdg*rho_inv
3178          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
3179          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
3180 !        Components 7 or 8
3181          d2fxdg2=d2fxdss2*dssdg**2
3182          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
3183 !        For the time being, treat non-spin-polarized like spin-polarized
3184          dvxci(ipts,2)=dvxci(ipts,1)
3185          dvxci(ipts,4)=dvxci(ipts,3)
3186          dvxci(ipts,6)=dvxci(ipts,5)
3187          dvxci(ipts,8)=dvxci(ipts,7)
3188 !        end of loop over the spin
3189 !        If non spin-polarized, treat spin down contribution now, similar to spin up
3190          exc=exc*2
3191          exci(ipts)=exc*rhotot_inv
3192 !        -----------------------------------------------------------------------------
3193 !        Then takes care of the LSD correlation part of the functional
3194 
3195          rs=rsfac*rhotmot
3196          sqr_rs=sq_rsfac*rhotmo6
3197          rsm1_2=sq_rsfac_inv*rhoto6
3198 
3199 !        Formulas A6-A8 of PW92LSD
3200          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
3201          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
3202          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
3203          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
3204 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
3205          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
3206          ecrs0=ec0_q0*ec0_log
3207          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
3208          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
3209          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
3210 &         -ec0_q0*ec0_q1pp*ec0_den                        &
3211 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
3212          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
3213          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
3214          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
3215          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
3216          mac_log=-log( mac_q1*mac_q1*mac_den )
3217          macrs=mac_q0*mac_log
3218          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
3219          ecrs=ecrs0
3220          decrs_drs=decrs0_drs
3221          decrs_dzeta=0.0_dp
3222          d2ecrs_drs2=d2ecrs0_drs2
3223          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
3224          d2ecrs_drsdzeta=zero
3225          zeta=0.0_dp
3226 
3227 !        Add LSD correlation functional to GGA exchange functional
3228          exci(ipts)=exci(ipts)+ecrs
3229          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
3230 
3231          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
3232 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
3233          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
3234          dvxci(ipts,9)=d2ecrs_drho2
3235          dvxci(ipts,10)=d2ecrs_drho2
3236          dvxci(ipts,11)=d2ecrs_drho2
3237 
3238 
3239 !        -----------------------------------------------------------------------------
3240 !        Eventually add the GGA correlation part of the PBE functional
3241 !        Note : the computation of the potential in the spin-unpolarized
3242 !        case could be optimized much further. Other optimizations are left to do.
3243 
3244          phi_zeta=1.0_dp
3245          phip_zeta=0.0_dp
3246          phi_zeta_inv=1.0_dp
3247          phi_logder=0.0_dp
3248          phi3_zeta=1.0_dp
3249          gamphi3inv=gamma_inv
3250          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
3251 
3252 !        From ec to bb
3253          bb=ecrs*gamphi3inv
3254          dbb_drs=decrs_drs*gamphi3inv
3255          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
3256          d2bb_drs2=d2ecrs_drs2*gamphi3inv
3257          d2bb_drsdzeta=gamphi3inv*(d2ecrs_drsdzeta-three*decrs_drs*phi_logder)
3258          d2bb_dzeta2=gamphi3inv*(d2ecrs_dzeta2-six*decrs_dzeta*phi_logder+&
3259 &         12.0_dp*ecrs*phi_logder*phi_logder-three*ecrs*phi_zeta_inv*phipp_zeta)
3260 
3261 !        From bb to cc
3262          exp_pbe=exp(-bb)
3263          cc=one/(exp_pbe-one)
3264          dcc_dbb=cc*cc*exp_pbe
3265          dcc_drs=dcc_dbb*dbb_drs
3266          dcc_dzeta=dcc_dbb*dbb_dzeta
3267          d2cc_dbb2=cc*cc*exp_pbe*(two*cc*exp_pbe-one)
3268          d2cc_drs2=d2cc_dbb2*dbb_drs*dbb_drs+dcc_dbb*d2bb_drs2
3269          d2cc_drsdzeta=d2cc_dbb2*dbb_drs*dbb_dzeta+dcc_dbb*d2bb_drsdzeta
3270          d2cc_dzeta2=d2cc_dbb2*dbb_dzeta*dbb_dzeta+dcc_dbb*d2bb_dzeta2
3271 
3272 !        From cc to aa
3273          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
3274          aa=coeff_aa*cc
3275          daa_drs=coeff_aa*dcc_drs
3276          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
3277          d2aa_drs2=coeff_aa*d2cc_drs2
3278          d2aa_drsdzeta=-two*daa_drs*phi_logder+coeff_aa*d2cc_drsdzeta
3279          d2aa_dzeta2=aa*(-two*phi_zeta_inv*phipp_zeta+six*phi_logder*phi_logder)+&
3280 &         coeff_aa*(-four*dcc_dzeta*phi_logder+d2cc_dzeta2)
3281 
3282 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
3283          grrho2=four*grho2_updn(ipts,1)
3284          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
3285 !        Note that tt is (the t variable of PBE divided by phi) squared
3286          tt=half*grrho2*dtt_dg
3287 
3288 !        Get xx from aa and tt
3289          xx=aa*tt
3290          dxx_drs=daa_drs*tt
3291          dxx_dzeta=daa_dzeta*tt
3292          dxx_dtt=aa
3293          d2xx_drs2=d2aa_drs2*tt
3294          d2xx_drsdzeta=d2aa_drsdzeta*tt
3295          d2xx_drsdtt=daa_drs
3296          d2xx_dttdzeta=daa_dzeta
3297          d2xx_dzeta2=d2aa_dzeta2*tt
3298 
3299 !        From xx to pade
3300          pade_den=one/(one+xx*(one+xx))
3301          pade=(one+xx)*pade_den
3302          dpade_dxx=-xx*(two+xx)*pade_den**2
3303          dpade_drs=dpade_dxx*dxx_drs
3304          dpade_dtt=dpade_dxx*dxx_dtt
3305          dpade_dzeta=dpade_dxx*dxx_dzeta
3306          d2pade_dxx2=two*(-one+xx*xx*(three+xx))*pade_den*pade_den*pade_den
3307          d2pade_drs2=d2pade_dxx2*dxx_drs*dxx_drs+dpade_dxx*d2xx_drs2
3308          d2pade_drsdtt=d2pade_dxx2*dxx_drs*dxx_dtt+dpade_dxx*d2xx_drsdtt
3309          d2pade_drsdzeta=d2pade_dxx2*dxx_drs*dxx_dzeta+dpade_dxx*d2xx_drsdzeta
3310          d2pade_dtt2=d2pade_dxx2*dxx_dtt*dxx_dtt
3311          d2pade_dttdzeta=d2pade_dxx2*dxx_dtt*dxx_dzeta+dpade_dxx*d2xx_dttdzeta
3312          d2pade_dzeta2=d2pade_dxx2*dxx_dzeta*dxx_dzeta+dpade_dxx*d2xx_dzeta2
3313 
3314 
3315 !        From pade to qq
3316          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
3317          qq=coeff_qq*pade
3318          dqq_drs=coeff_qq*dpade_drs
3319          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
3320          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
3321          d2qq_drs2=coeff_qq*d2pade_drs2
3322          d2qq_drsdtt=phi_zeta_inv*phi_zeta_inv*(dpade_drs+tt*d2pade_drsdtt)
3323          d2qq_drsdzeta=coeff_qq*(d2pade_drsdzeta-two*dpade_drs*phi_logder)
3324          d2qq_dtt2=phi_zeta_inv*phi_zeta_inv*(two*dpade_dtt+tt*d2pade_dtt2)
3325          d2qq_dttdzeta=phi_zeta_inv*phi_zeta_inv*(dpade_dzeta-two*pade*phi_logder)+&
3326 &         coeff_qq*(d2pade_dttdzeta-two*dpade_dtt*phi_logder)
3327          d2qq_dzeta2=coeff_qq*( d2pade_dzeta2-four*dpade_dzeta*phi_logder &
3328 &         +six*pade*phi_logder*phi_logder            &
3329 &         -two*pade*phi_zeta_inv*phipp_zeta)
3330 
3331 !        From qq to rr
3332          arg_rr=one+beta*gamma_inv*qq
3333          div_rr=one/arg_rr
3334          rr=gamma*log(arg_rr)
3335          drr_dqq=beta*div_rr
3336          drr_drs=drr_dqq*dqq_drs
3337          drr_dtt=drr_dqq*dqq_dtt
3338          drr_dzeta=drr_dqq*dqq_dzeta
3339          d2rr_dqq2=-div_rr**2*beta*beta*gamma_inv
3340          d2rr_drs2=d2rr_dqq2*dqq_drs*dqq_drs+drr_dqq*d2qq_drs2
3341          d2rr_drsdtt=d2rr_dqq2*dqq_drs*dqq_dtt+drr_dqq*d2qq_drsdtt
3342          d2rr_drsdzeta=d2rr_dqq2*dqq_drs*dqq_dzeta+drr_dqq*d2qq_drsdzeta
3343          d2rr_dtt2=d2rr_dqq2*dqq_dtt*dqq_dtt+drr_dqq*d2qq_dtt2
3344          d2rr_dttdzeta=d2rr_dqq2*dqq_dtt*dqq_dzeta+drr_dqq*d2qq_dttdzeta
3345          d2rr_dzeta2=d2rr_dqq2*dqq_dzeta*dqq_dzeta+drr_dqq*d2qq_dzeta2
3346 
3347 !        From rr to hh
3348          hh=phi3_zeta*rr
3349          dhh_drs=phi3_zeta*drr_drs
3350          dhh_dtt=phi3_zeta*drr_dtt
3351          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
3352          d2hh_drs2=phi3_zeta*d2rr_drs2
3353          d2hh_drsdtt=phi3_zeta*d2rr_drsdtt
3354          d2hh_drsdzeta=phi3_zeta*(d2rr_drsdzeta+three*drr_drs*phi_logder)
3355          d2hh_dtt2=phi3_zeta*d2rr_dtt2
3356          d2hh_dttdzeta=phi3_zeta*(d2rr_dttdzeta+three*drr_dtt*phi_logder)
3357          d2hh_dzeta2=phi3_zeta*(six*rr*phi_logder*phi_logder+&
3358 &         six*phi_logder*drr_dzeta+d2rr_dzeta2)  &
3359 &         +three*phi_zeta*phi_zeta*rr*phipp_zeta
3360 
3361 
3362 !        The GGA correlation energy is added
3363          exci(ipts)=exci(ipts)+hh
3364 
3365 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
3366 
3367 
3368 !        From hh to the derivative of the energy wrt the density
3369          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
3370          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
3371 
3372 !        From hh to the derivative of the energy wrt to the gradient of the
3373 !        density, divided by the gradient of the density
3374 !        (The v3.3 definition includes the division by the norm of the gradient)
3375          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
3376 
3377          d2rhohh_drho2=rhotot_inv*&
3378 &         (-two*ninth*rs*dhh_drs +seven*four*ninth*tt*dhh_dtt &
3379 &         +ninth*rs*rs*d2hh_drs2+zeta*zeta*d2hh_dzeta2+(seven*third*tt)**2*d2hh_dtt2 &
3380 &         +two*third*rs*zeta*d2hh_drsdzeta+two*seven*ninth*rs*tt*d2hh_drsdtt &
3381 &         +two*seven*third*tt*zeta*d2hh_dttdzeta)
3382          d2rhohh_drhodg=dtt_dg*(-four*third*dhh_dtt-third*rs*d2hh_drsdtt &
3383 &         -zeta*d2hh_dttdzeta-seven*third*tt*d2hh_dtt2)
3384 
3385 !        Component 12 : first derivative with respect to the gradient
3386 !        of the density, div by the grad of the density
3387          dvxci(ipts,12)=dvxcdgr(ipts,3)
3388 !        Components 9, 10 and 11 : second derivatives with respect to the spin-density
3389 !        Note that there is already a contribution from LSDA
3390          dvxci(ipts,9)=dvxci(ipts,9)+d2rhohh_drho2+rhotot_inv*           &
3391 &         ( d2hh_dzeta2*(one-two*zeta) &
3392 &         -two*third*rs*d2hh_drsdzeta-14.0_dp*third*tt*d2hh_dttdzeta)
3393          dvxci(ipts,10)=dvxci(ipts,10)+d2rhohh_drho2-rhotot_inv*d2hh_dzeta2
3394          dvxci(ipts,11)=dvxci(ipts,11)+d2rhohh_drho2+rhotot_inv*           &
3395 &         ( d2hh_dzeta2*(one+two*zeta) &
3396 &         +two*third*rs*d2hh_drsdzeta+14.0_dp*third*tt*d2hh_dttdzeta)
3397 !        Components 13 and 14 : second derivatives with respect to spin density
3398 !        and gradient, divided by the gradient
3399          dvxci(ipts,13)=d2rhohh_drhodg+dtt_dg*d2hh_dttdzeta
3400          dvxci(ipts,14)=d2rhohh_drhodg-dtt_dg*d2hh_dttdzeta
3401 !        Component 15 : derivative of the (derivative wrt the gradient div by the grad),
3402 !        divided by the grad
3403          dvxci(ipts,15)=rhotot*d2hh_dtt2*dtt_dg*dtt_dg
3404 
3405 
3406 !        End condition of GGA
3407 
3408 !        Correlation has been added
3409 !        -----------------------------------------------------------------------------
3410 
3411 !        vxci(ipts,2)=vxci(ipts,1)
3412          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
3413 
3414        end do
3415 
3416      else if ((option==6) .or. (option==7)) then
3417        do ipts=1,npts
3418 
3419          rhotot=rhoarr(ipts)
3420          rhotmot=rhom1_3(ipts)
3421          rhotot_inv=rhotmot*rhotmot*rhotmot
3422          rhotmo6=sqrt(rhotmot)
3423          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
3424 !        -----------------------------------------------------------------------
3425 !        First take care of the exchange part of the functional
3426 
3427          exc=zero
3428          dvxcdgr(ipts,3)=zero
3429 !        loop over the spin
3430          ispden=1
3431          rho   =rho_updn(ipts,ispden)
3432          rhomot=rho_updnm1_3(ipts,ispden)
3433          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
3434 !        Perdew-Burke-Ernzerhof GGA, exchange part
3435          rho_inv=rhomot*rhomot*rhomot
3436          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
3437          ss=grho2_updn(ipts,ispden)*coeffss
3438 
3439          if (option==6) then
3440            divss=exp(-mu_divkappa*ss)
3441            dfxdss= mu*divss
3442            d2fxdss2=-mu*mu_divkappa*divss
3443 
3444            fx    = one+kappa*(one-divss)
3445            ex_gga= ex_lsd*fx
3446            dssdn=-eight*third*ss*rho_inv
3447            dfxdn  = dfxdss*dssdn
3448            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
3449 !          The new definition (v3.3) includes the division by the norm of the gradient
3450            dssdg =two*coeffss
3451            dfxdg=dfxdss*dssdg
3452            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
3453            exc=exc+ex_gga*rho
3454 !          This is the Wu and Cohen modification
3455          else
3456            expss=exp(-ss)
3457            p1_wc=b_wc+(mu-b_wc)*(one-ss)*expss+two*c_wc*ss/(one+c_wc*ss*ss)
3458            p2_wc=d_wc*(ss-two)*expss+two*c_wc/(one+c_wc*ss*ss)-&
3459 &           four*c_wc*c_wc*ss*ss/((one+c_wc*ss*ss)*(one+c_wc*ss*ss))
3460            divss=one/(one+(b_wc*ss+d_wc*ss*expss+log(one+c_wc*ss*ss))/kappa)
3461            dfxdss=p1_wc*divss*divss
3462            d2fxdss2=p2_wc*divss*divss-two*divss*divss*divss*p1_wc*p1_wc/kappa
3463 
3464            fx    = one+kappa*(one-divss)
3465            ex_gga= ex_lsd*fx
3466            dssdn=-eight*third*ss*rho_inv
3467            dfxdn  = dfxdss*dssdn
3468            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
3469 !          The new definition (v3.3) includes the division by the norm of the gradient
3470            dssdg =two*coeffss
3471            dfxdg=dfxdss*dssdg
3472            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
3473            exc=exc+ex_gga*rho
3474          end if
3475 
3476 !        Perdew-Burke-Ernzerhof GGA, exchange part
3477 !        Components 3 or 4
3478          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
3479 !        Components 1 or 2
3480          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
3481          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
3482          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
3483 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
3484 !        Components 5 or 6
3485          d2ssdndg=-eight*third*dssdg*rho_inv
3486          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
3487          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
3488 !        Components 7 or 8
3489          d2fxdg2=d2fxdss2*dssdg**2
3490          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
3491 !        For the time being, treat non-spin-polarized like spin-polarized
3492          dvxci(ipts,2)=dvxci(ipts,1)
3493          dvxci(ipts,4)=dvxci(ipts,3)
3494          dvxci(ipts,6)=dvxci(ipts,5)
3495          dvxci(ipts,8)=dvxci(ipts,7)
3496 !        end of loop over the spin
3497 !        If non spin-polarized, treat spin down contribution now, similar to spin up
3498          exc=exc*2
3499          exci(ipts)=exc*rhotot_inv
3500 !        -----------------------------------------------------------------------------
3501 !        Then takes care of the LSD correlation part of the functional
3502 
3503 
3504          rs=rsfac*rhotmot
3505          sqr_rs=sq_rsfac*rhotmo6
3506          rsm1_2=sq_rsfac_inv*rhoto6
3507 
3508 !        Formulas A6-A8 of PW92LSD
3509          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
3510          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
3511          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
3512          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
3513 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
3514          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
3515          ecrs0=ec0_q0*ec0_log
3516          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
3517          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
3518          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
3519 &         -ec0_q0*ec0_q1pp*ec0_den                        &
3520 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
3521          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
3522          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
3523          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
3524          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
3525          mac_log=-log( mac_q1*mac_q1*mac_den )
3526          macrs=mac_q0*mac_log
3527          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
3528          ecrs=ecrs0
3529          decrs_drs=decrs0_drs
3530          decrs_dzeta=0.0_dp
3531          d2ecrs_drs2=d2ecrs0_drs2
3532          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
3533          d2ecrs_drsdzeta=zero
3534          zeta=0.0_dp
3535 
3536 !        Add LSD correlation functional to GGA exchange functional
3537          exci(ipts)=exci(ipts)+ecrs
3538          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
3539 
3540          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
3541 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
3542          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
3543          dvxci(ipts,9)=d2ecrs_drho2
3544          dvxci(ipts,10)=d2ecrs_drho2
3545          dvxci(ipts,11)=d2ecrs_drho2
3546 
3547 
3548 !        -----------------------------------------------------------------------------
3549 !        Eventually add the GGA correlation part of the PBE functional
3550 !        Note : the computation of the potential in the spin-unpolarized
3551 !        case could be optimized much further. Other optimizations are left to do.
3552 
3553          phi_zeta=1.0_dp
3554          phip_zeta=0.0_dp
3555          phi_zeta_inv=1.0_dp
3556          phi_logder=0.0_dp
3557          phi3_zeta=1.0_dp
3558          gamphi3inv=gamma_inv
3559          phipp_zeta=-two*ninth*alpha_zeta*alpha_zeta
3560 
3561 !        From ec to bb
3562          bb=ecrs*gamphi3inv
3563          dbb_drs=decrs_drs*gamphi3inv
3564          dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
3565          d2bb_drs2=d2ecrs_drs2*gamphi3inv
3566          d2bb_drsdzeta=gamphi3inv*(d2ecrs_drsdzeta-three*decrs_drs*phi_logder)
3567          d2bb_dzeta2=gamphi3inv*(d2ecrs_dzeta2-six*decrs_dzeta*phi_logder+&
3568 &         12.0_dp*ecrs*phi_logder*phi_logder-three*ecrs*phi_zeta_inv*phipp_zeta)
3569 
3570 !        From bb to cc
3571          exp_pbe=exp(-bb)
3572          cc=one/(exp_pbe-one)
3573          dcc_dbb=cc*cc*exp_pbe
3574          dcc_drs=dcc_dbb*dbb_drs
3575          dcc_dzeta=dcc_dbb*dbb_dzeta
3576          d2cc_dbb2=cc*cc*exp_pbe*(two*cc*exp_pbe-one)
3577          d2cc_drs2=d2cc_dbb2*dbb_drs*dbb_drs+dcc_dbb*d2bb_drs2
3578          d2cc_drsdzeta=d2cc_dbb2*dbb_drs*dbb_dzeta+dcc_dbb*d2bb_drsdzeta
3579          d2cc_dzeta2=d2cc_dbb2*dbb_dzeta*dbb_dzeta+dcc_dbb*d2bb_dzeta2
3580 
3581 !        From cc to aa
3582          coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
3583          aa=coeff_aa*cc
3584          daa_drs=coeff_aa*dcc_drs
3585          daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
3586          d2aa_drs2=coeff_aa*d2cc_drs2
3587          d2aa_drsdzeta=-two*daa_drs*phi_logder+coeff_aa*d2cc_drsdzeta
3588          d2aa_dzeta2=aa*(-two*phi_zeta_inv*phipp_zeta+six*phi_logder*phi_logder)+&
3589 &         coeff_aa*(-four*dcc_dzeta*phi_logder+d2cc_dzeta2)
3590 
3591 !        Introduce tt : do not assume that the spin-dependent gradients are collinear
3592          grrho2=four*grho2_updn(ipts,1)
3593          dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
3594 !        Note that tt is (the t variable of PBE divided by phi) squared
3595          tt=half*grrho2*dtt_dg
3596 
3597 !        Get xx from aa and tt
3598          xx=aa*tt
3599          dxx_drs=daa_drs*tt
3600          dxx_dzeta=daa_dzeta*tt
3601          dxx_dtt=aa
3602          d2xx_drs2=d2aa_drs2*tt
3603          d2xx_drsdzeta=d2aa_drsdzeta*tt
3604          d2xx_drsdtt=daa_drs
3605          d2xx_dttdzeta=daa_dzeta
3606          d2xx_dzeta2=d2aa_dzeta2*tt
3607 
3608 !        From xx to pade
3609          pade_den=one/(one+xx*(one+xx))
3610          pade=(one+xx)*pade_den
3611          dpade_dxx=-xx*(two+xx)*pade_den**2
3612          dpade_drs=dpade_dxx*dxx_drs
3613          dpade_dtt=dpade_dxx*dxx_dtt
3614          dpade_dzeta=dpade_dxx*dxx_dzeta
3615          d2pade_dxx2=two*(-one+xx*xx*(three+xx))*pade_den*pade_den*pade_den
3616          d2pade_drs2=d2pade_dxx2*dxx_drs*dxx_drs+dpade_dxx*d2xx_drs2
3617          d2pade_drsdtt=d2pade_dxx2*dxx_drs*dxx_dtt+dpade_dxx*d2xx_drsdtt
3618          d2pade_drsdzeta=d2pade_dxx2*dxx_drs*dxx_dzeta+dpade_dxx*d2xx_drsdzeta
3619          d2pade_dtt2=d2pade_dxx2*dxx_dtt*dxx_dtt
3620          d2pade_dttdzeta=d2pade_dxx2*dxx_dtt*dxx_dzeta+dpade_dxx*d2xx_dttdzeta
3621          d2pade_dzeta2=d2pade_dxx2*dxx_dzeta*dxx_dzeta+dpade_dxx*d2xx_dzeta2
3622 
3623 
3624 !        From pade to qq
3625          coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
3626          qq=coeff_qq*pade
3627          dqq_drs=coeff_qq*dpade_drs
3628          dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
3629          dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
3630          d2qq_drs2=coeff_qq*d2pade_drs2
3631          d2qq_drsdtt=phi_zeta_inv*phi_zeta_inv*(dpade_drs+tt*d2pade_drsdtt)
3632          d2qq_drsdzeta=coeff_qq*(d2pade_drsdzeta-two*dpade_drs*phi_logder)
3633          d2qq_dtt2=phi_zeta_inv*phi_zeta_inv*(two*dpade_dtt+tt*d2pade_dtt2)
3634          d2qq_dttdzeta=phi_zeta_inv*phi_zeta_inv*(dpade_dzeta-two*pade*phi_logder)+&
3635 &         coeff_qq*(d2pade_dttdzeta-two*dpade_dtt*phi_logder)
3636          d2qq_dzeta2=coeff_qq*( d2pade_dzeta2-four*dpade_dzeta*phi_logder &
3637 &         +six*pade*phi_logder*phi_logder            &
3638 &         -two*pade*phi_zeta_inv*phipp_zeta)
3639 
3640 !        From qq to rr
3641          arg_rr=one+beta*gamma_inv*qq
3642          div_rr=one/arg_rr
3643          rr=gamma*log(arg_rr)
3644          drr_dqq=beta*div_rr
3645          drr_drs=drr_dqq*dqq_drs
3646          drr_dtt=drr_dqq*dqq_dtt
3647          drr_dzeta=drr_dqq*dqq_dzeta
3648          d2rr_dqq2=-div_rr**2*beta*beta*gamma_inv
3649          d2rr_drs2=d2rr_dqq2*dqq_drs*dqq_drs+drr_dqq*d2qq_drs2
3650          d2rr_drsdtt=d2rr_dqq2*dqq_drs*dqq_dtt+drr_dqq*d2qq_drsdtt
3651          d2rr_drsdzeta=d2rr_dqq2*dqq_drs*dqq_dzeta+drr_dqq*d2qq_drsdzeta
3652          d2rr_dtt2=d2rr_dqq2*dqq_dtt*dqq_dtt+drr_dqq*d2qq_dtt2
3653          d2rr_dttdzeta=d2rr_dqq2*dqq_dtt*dqq_dzeta+drr_dqq*d2qq_dttdzeta
3654          d2rr_dzeta2=d2rr_dqq2*dqq_dzeta*dqq_dzeta+drr_dqq*d2qq_dzeta2
3655 
3656 !        From rr to hh
3657          hh=phi3_zeta*rr
3658          dhh_drs=phi3_zeta*drr_drs
3659          dhh_dtt=phi3_zeta*drr_dtt
3660          dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
3661          d2hh_drs2=phi3_zeta*d2rr_drs2
3662          d2hh_drsdtt=phi3_zeta*d2rr_drsdtt
3663          d2hh_drsdzeta=phi3_zeta*(d2rr_drsdzeta+three*drr_drs*phi_logder)
3664          d2hh_dtt2=phi3_zeta*d2rr_dtt2
3665          d2hh_dttdzeta=phi3_zeta*(d2rr_dttdzeta+three*drr_dtt*phi_logder)
3666          d2hh_dzeta2=phi3_zeta*(six*rr*phi_logder*phi_logder+&
3667 &         six*phi_logder*drr_dzeta+d2rr_dzeta2)  &
3668 &         +three*phi_zeta*phi_zeta*rr*phipp_zeta
3669 
3670 
3671 !        The GGA correlation energy is added
3672          exci(ipts)=exci(ipts)+hh
3673 
3674 !        Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
3675 
3676 
3677 !        From hh to the derivative of the energy wrt the density
3678          drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
3679          vxci(ipts,1)=vxci(ipts,1)+drhohh_drho
3680 
3681 !        From hh to the derivative of the energy wrt to the gradient of the
3682 !        density, divided by the gradient of the density
3683 !        (The v3.3 definition includes the division by the norm of the gradient)
3684          dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
3685 
3686          d2rhohh_drho2=rhotot_inv*&
3687 &         (-two*ninth*rs*dhh_drs +seven*four*ninth*tt*dhh_dtt &
3688 &         +ninth*rs*rs*d2hh_drs2+zeta*zeta*d2hh_dzeta2+(seven*third*tt)**2*d2hh_dtt2 &
3689 &         +two*third*rs*zeta*d2hh_drsdzeta+two*seven*ninth*rs*tt*d2hh_drsdtt &
3690 &         +two*seven*third*tt*zeta*d2hh_dttdzeta)
3691          d2rhohh_drhodg=dtt_dg*(-four*third*dhh_dtt-third*rs*d2hh_drsdtt &
3692 &         -zeta*d2hh_dttdzeta-seven*third*tt*d2hh_dtt2)
3693 
3694 !        Component 12 : first derivative with respect to the gradient
3695 !        of the density, div by the grad of the density
3696          dvxci(ipts,12)=dvxcdgr(ipts,3)
3697 !        Components 9, 10 and 11 : second derivatives with respect to the spin-density
3698 !        Note that there is already a contribution from LSDA
3699          dvxci(ipts,9)=dvxci(ipts,9)+d2rhohh_drho2+rhotot_inv*           &
3700 &         ( d2hh_dzeta2*(one-two*zeta) &
3701 &         -two*third*rs*d2hh_drsdzeta-14.0_dp*third*tt*d2hh_dttdzeta)
3702          dvxci(ipts,10)=dvxci(ipts,10)+d2rhohh_drho2-rhotot_inv*d2hh_dzeta2
3703          dvxci(ipts,11)=dvxci(ipts,11)+d2rhohh_drho2+rhotot_inv*           &
3704 &         ( d2hh_dzeta2*(one+two*zeta) &
3705 &         +two*third*rs*d2hh_drsdzeta+14.0_dp*third*tt*d2hh_dttdzeta)
3706 !        Components 13 and 14 : second derivatives with respect to spin density
3707 !        and gradient, divided by the gradient
3708          dvxci(ipts,13)=d2rhohh_drhodg+dtt_dg*d2hh_dttdzeta
3709          dvxci(ipts,14)=d2rhohh_drhodg-dtt_dg*d2hh_dttdzeta
3710 !        Component 15 : derivative of the (derivative wrt the gradient div by the grad),
3711 !        divided by the grad
3712          dvxci(ipts,15)=rhotot*d2hh_dtt2*dtt_dg*dtt_dg
3713 
3714 
3715 !        End condition of GGA
3716 
3717 !        Correlation has been added
3718 !        -----------------------------------------------------------------------------
3719 
3720 !        vxci(ipts,2)=vxci(ipts,1)
3721          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
3722 
3723        end do
3724 
3725      else if (option==-1) then
3726        do ipts=1,npts
3727 
3728          rhotot=rhoarr(ipts)
3729          rhotmot=rhom1_3(ipts)
3730          rhotot_inv=rhotmot*rhotmot*rhotmot
3731          rhotmo6=sqrt(rhotmot)
3732          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
3733 !        -----------------------------------------------------------------------
3734 !        First take care of the exchange part of the functional
3735 
3736          exc=zero
3737 !        loop over the spin
3738          ispden=1
3739          rho   =rho_updn(ipts,ispden)
3740          rhomot=rho_updnm1_3(ipts,ispden)
3741          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
3742 !        Perdew_Wang 91 LSD
3743          vxci(ipts,ispden)=four_thirds*ex_lsd
3744          exc=exc+ex_lsd*rho
3745 
3746 !        Perdew_Wang 91 LSD
3747          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
3748 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
3749          dvxci(ipts,2)=zero
3750 !        If non-spin-polarized, first component of dvxci is second
3751 !        derivative with respect to TOTAL density.
3752          dvxci(ipts,1)=dvxci(ipts,1)*half
3753 !        Compute the second derivative of vx
3754 !        vx^(2) = -2*vx^(1)/(3*rhotot)
3755 !        end of loop over the spin
3756 !        If non spin-polarized, treat spin down contribution now, similar to spin up
3757          exc=exc*2
3758          exci(ipts)=exc*rhotot_inv
3759 !        -----------------------------------------------------------------------------
3760 !        Then takes care of the LSD correlation part of the functional
3761 
3762        end do
3763 
3764      else if (option==-2) then
3765        do ipts=1,npts
3766 
3767          rhotot=rhoarr(ipts)
3768          rhotmot=rhom1_3(ipts)
3769          rhotot_inv=rhotmot*rhotmot*rhotmot
3770          rhotmo6=sqrt(rhotmot)
3771          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
3772 !        -----------------------------------------------------------------------
3773 !        First take care of the exchange part of the functional
3774 
3775          exc=zero
3776          dvxcdgr(ipts,3)=zero
3777 !        loop over the spin
3778          ispden=1
3779          rho   =rho_updn(ipts,ispden)
3780          rhomot=rho_updnm1_3(ipts,ispden)
3781          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
3782 !        Perdew-Burke-Ernzerhof GGA, exchange part
3783          rho_inv=rhomot*rhomot*rhomot
3784          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
3785          ss=grho2_updn(ipts,ispden)*coeffss
3786          divss=one/(one+mu_divkappa*ss)
3787          dfxdss= mu*divss*divss
3788          d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
3789          fx    = one+kappa*(one-divss)
3790          ex_gga= ex_lsd*fx
3791          dssdn=-eight*third*ss*rho_inv
3792          dfxdn  = dfxdss*dssdn
3793          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
3794 !        The new definition (v3.3) includes the division by the norm of the gradient
3795          dssdg =two*coeffss
3796          dfxdg=dfxdss*dssdg
3797          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
3798          exc=exc+ex_gga*rho
3799 
3800 !        Perdew-Burke-Ernzerhof GGA, exchange part
3801 !        Components 3 or 4
3802          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
3803 !        Components 1 or 2
3804          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
3805          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
3806          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
3807 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
3808 !        Components 5 or 6
3809          d2ssdndg=-eight*third*dssdg*rho_inv
3810          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
3811          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
3812 !        Components 7 or 8
3813          d2fxdg2=d2fxdss2*dssdg**2
3814          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
3815 !        For the time being, treat non-spin-polarized like spin-polarized
3816          dvxci(ipts,2)=dvxci(ipts,1)
3817          dvxci(ipts,4)=dvxci(ipts,3)
3818          dvxci(ipts,6)=dvxci(ipts,5)
3819          dvxci(ipts,8)=dvxci(ipts,7)
3820 !        end of loop over the spin
3821 !        If non spin-polarized, treat spin down contribution now, similar to spin up
3822          exc=exc*2
3823          exci(ipts)=exc*rhotot_inv
3824 !        -----------------------------------------------------------------------------
3825 !        Then takes care of the LSD correlation part of the functional
3826 
3827 !        Correlation has been added
3828 !        -----------------------------------------------------------------------------
3829 
3830 !        vxci(ipts,2)=vxci(ipts,1)
3831          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
3832 
3833        end do
3834 
3835 
3836      else if(option==-4) then
3837 
3838 
3839        do ipts=1,npts
3840 
3841          rhotot=rhoarr(ipts)
3842          rhotmot=rhom1_3(ipts)
3843          rhotot_inv=rhotmot*rhotmot*rhotmot
3844          rhotmo6=sqrt(rhotmot)
3845          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
3846 !        -----------------------------------------------------------------------
3847 !        First take care of the exchange part of the functional
3848 
3849          exc=zero
3850          dvxcdgr(ipts,3)=zero
3851 !        loop over the spin
3852          ispden=1
3853          rho   =rho_updn(ipts,ispden)
3854          rhomot=rho_updnm1_3(ipts,ispden)
3855          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
3856 !        VALENTINO R. COOPER C09x GGA, This is an exchange term proposed
3857 !        to use together with vdw-DF (see above).
3858          rho_inv=rhomot*rhomot*rhomot
3859          coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
3860 !        the quarter that is lacking is compensated by the grho2_updn in the
3861 !        next line.
3862          ss=grho2_updn(ipts,ispden)*coeffss
3863          alphs2=alpha_c09*ss
3864          alphmu=alpha_c09*mu_c09
3865          dfxdss= mu_c09*exp(-alphs2)*(one-alphs2)+&
3866 &         kappa*alpha_c09*exp(-alphs2/two)/two
3867          d2fxdss2=-alphmu*exp(-alphs2)*(two-alphs2)-&
3868 &         kappa*(alpha_c09**two)*exp(alphs2/two)/four
3869          fx    = one+mu_c09*ss*exp(-alphs2)+kappa*(one-exp(-alphs2/two))
3870          ex_gga= ex_lsd*fx
3871          dssdn=-eight*third*ss*rho_inv
3872          dfxdn  = dfxdss*dssdn
3873          vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
3874 !        The new definition (v3.3) includes the division by the norm of the gradient
3875          dssdg =two*coeffss
3876          dfxdg=dfxdss*dssdg
3877          dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
3878          exc=exc+ex_gga*rho
3879 
3880 !        Cooper C09x GGA exchange
3881 !        Components 3 or 4
3882          dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
3883 !        Components 1 or 2
3884          d2ssdn2=-11.0_dp*third*dssdn*rho_inv
3885          d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
3886          dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
3887 &         ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
3888 !        Components 5 or 6
3889          d2ssdndg=-eight*third*dssdg*rho_inv
3890          d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
3891          dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
3892 !        Components 7 or 8
3893          d2fxdg2=d2fxdss2*dssdg**2
3894          dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
3895 !        For the time being, treat non-spin-polarized like spin-polarized
3896          dvxci(ipts,2)=dvxci(ipts,1)
3897          dvxci(ipts,4)=dvxci(ipts,3)
3898          dvxci(ipts,6)=dvxci(ipts,5)
3899          dvxci(ipts,8)=dvxci(ipts,7)
3900 
3901 !        end of loop over the spin
3902 !        If non spin-polarized, treat spin down contribution now, similar to spin up
3903          exc=exc*2
3904          exci(ipts)=exc*rhotot_inv
3905 !        -----------------------------------------------------------------------------
3906 !        Then takes care of the LSD correlation part of the functional
3907 
3908 !        Correlation has been added
3909 !        -----------------------------------------------------------------------------
3910 
3911 !        vxci(ipts,2)=vxci(ipts,1)
3912          dvxcdgr(ipts,2)=dvxcdgr(ipts,1)
3913 
3914        end do
3915 
3916 
3917      else if (option==1) then
3918        do ipts=1,npts
3919 
3920          rhotot=rhoarr(ipts)
3921          rhotmot=rhom1_3(ipts)
3922          rhotot_inv=rhotmot*rhotmot*rhotmot
3923          rhotmo6=sqrt(rhotmot)
3924          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
3925 !        -----------------------------------------------------------------------
3926 !        First take care of the exchange part of the functional
3927 
3928          exc=zero
3929 !        loop over the spin
3930          ispden=1
3931          rho   =rho_updn(ipts,ispden)
3932          rhomot=rho_updnm1_3(ipts,ispden)
3933          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
3934 !        Perdew_Wang 91 LSD
3935          vxci(ipts,ispden)=four_thirds*ex_lsd
3936          exc=exc+ex_lsd*rho
3937 
3938 !        Perdew_Wang 91 LSD
3939          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
3940 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
3941          dvxci(ipts,2)=zero
3942 !        If non-spin-polarized, first component of dvxci is second
3943 !        derivative with respect to TOTAL density.
3944          dvxci(ipts,1)=dvxci(ipts,1)*half
3945 !        Compute the second derivative of vx
3946 !        vx^(2) = -2*vx^(1)/(3*rhotot)
3947 !        end of loop over the spin
3948 !        If non spin-polarized, treat spin down contribution now, similar to spin up
3949          exc=exc*2
3950          exci(ipts)=exc*rhotot_inv
3951 !        -----------------------------------------------------------------------------
3952 !        Then takes care of the LSD correlation part of the functional
3953 
3954 
3955          rs=rsfac*rhotmot
3956          sqr_rs=sq_rsfac*rhotmo6
3957          rsm1_2=sq_rsfac_inv*rhoto6
3958 
3959 !        Formulas A6-A8 of PW92LSD
3960          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
3961          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
3962          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
3963          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
3964 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
3965          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
3966          ecrs0=ec0_q0*ec0_log
3967          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
3968          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
3969          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
3970 &         -ec0_q0*ec0_q1pp*ec0_den                        &
3971 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
3972          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
3973          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
3974          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
3975          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
3976          mac_log=-log( mac_q1*mac_q1*mac_den )
3977          macrs=mac_q0*mac_log
3978          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
3979          ecrs=ecrs0
3980          decrs_drs=decrs0_drs
3981          decrs_dzeta=0.0_dp
3982          d2ecrs_drs2=d2ecrs0_drs2
3983          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
3984          d2ecrs_drsdzeta=zero
3985          zeta=0.0_dp
3986 
3987 !        Add LSD correlation functional to GGA exchange functional
3988          exci(ipts)=exci(ipts)+ecrs
3989          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
3990 
3991          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
3992 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
3993          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
3994 
3995          dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drho2
3996 
3997        end do
3998      else if (option==3) then
3999        do ipts=1,npts
4000 
4001          rhotot=rhoarr(ipts)
4002          rhotmot=rhom1_3(ipts)
4003          rhotot_inv=rhotmot*rhotmot*rhotmot
4004          rhotmo6=sqrt(rhotmot)
4005          rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
4006 !        -----------------------------------------------------------------------
4007 !        First take care of the exchange part of the functional
4008 
4009          exc=zero
4010 !        loop over the spin
4011          ispden=1
4012          rho   =rho_updn(ipts,ispden)
4013          rhomot=rho_updnm1_3(ipts,ispden)
4014          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
4015 !        Perdew_Wang 91 LSD
4016          vxci(ipts,ispden)=four_thirds*ex_lsd
4017          exc=exc+ex_lsd*rho
4018 
4019 !        Perdew_Wang 91 LSD
4020          dvxci(ipts,2*ispden-1)=-four_thirds*third*&
4021 &         threefourth_divpi*sixpi2_1_3*rhomot*rhomot
4022          dvxci(ipts,2)=zero
4023 !        If non-spin-polarized, first component of dvxci is second
4024 !        derivative with respect to TOTAL density.
4025          dvxci(ipts,1)=dvxci(ipts,1)*half
4026 !        Compute the second derivative of vx
4027 !        vx^(2) = -2*vx^(1)/(3*rhotot)
4028 !        end of loop over the spin
4029 !        If non spin-polarized, treat spin down contribution now, similar to spin up
4030          exc=exc*2
4031          exci(ipts)=exc*rhotot_inv
4032 !        -----------------------------------------------------------------------------
4033 !        Then takes care of the LSD correlation part of the functional
4034 
4035 
4036          rs=rsfac*rhotmot
4037          sqr_rs=sq_rsfac*rhotmo6
4038          rsm1_2=sq_rsfac_inv*rhoto6
4039 
4040 !        Formulas A6-A8 of PW92LSD
4041          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
4042          sqr_sqr_rs=max(1.e-15_dp,sqrt(sqr_rs))
4043          ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs/sqr_sqr_rs)
4044          ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+3.5_dp*ec0_b4*rs/sqr_sqr_rs)
4045          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
4046 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
4047          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
4048          ecrs0=ec0_q0*ec0_log
4049          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
4050          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
4051          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
4052 &         -ec0_q0*ec0_q1pp*ec0_den                        &
4053 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
4054          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
4055          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
4056          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
4057          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
4058          mac_log=-log( mac_q1*mac_q1*mac_den )
4059          macrs=mac_q0*mac_log
4060          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
4061          ecrs=ecrs0
4062          decrs_drs=decrs0_drs
4063          decrs_dzeta=0.0_dp
4064          d2ecrs_drs2=d2ecrs0_drs2
4065          d2ecrs_dzeta2=alpha_zeta**2*(-macrs)
4066          d2ecrs_drsdzeta=zero
4067          zeta=0.0_dp
4068 
4069 !        Add LSD correlation functional to GGA exchange functional
4070          exci(ipts)=exci(ipts)+ecrs
4071          vxci(ipts,1)=vxci(ipts,1)+ecrs-rs*third*decrs_drs
4072 
4073          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
4074 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
4075          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
4076          dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drho2
4077 
4078        end do
4079      end if
4080 
4081    end if
4082 
4083 
4084 !  fab: here it starts the spin polarized case
4085 
4086  else if(nspden==2) then
4087 
4088 !  we separate different cases depending on order
4089 
4090    if (order**2<=1) then
4091 
4092      do ipts=1,npts
4093 
4094        rhotot=rhoarr(ipts)
4095        rhotmot=rhom1_3(ipts)
4096        rhotot_inv=rhotmot*rhotmot*rhotmot
4097        rhotmo6=sqrt(rhotmot)
4098        rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
4099 !      -----------------------------------------------------------------------
4100 !      First take care of the exchange part of the functional
4101 
4102        exc=zero
4103        if (present(dvxcdgr)) dvxcdgr(ipts,3)=zero
4104        do ispden=1,nspden
4105          rho   =rho_updn(ipts,ispden)
4106          rhomot=rho_updnm1_3(ipts,ispden)
4107          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
4108          if(option==1 .or. option==-1 .or. option==3)then
4109 !          Perdew_Wang 91 LSD
4110            vxci(ipts,ispden)=four_thirds*ex_lsd
4111            if(present(dvxcdgr)) dvxcdgr(ipts,ispden)=0.0_dp
4112            exc=exc+ex_lsd*rho
4113          else
4114            rho_inv=rhomot*rhomot*rhomot
4115            coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
4116            ss=grho2_updn(ipts,ispden)*coeffss
4117            if(option==7) then ! This is WC
4118              expss=exp(-ss)
4119              p1_wc=b_wc+(mu-b_wc)*(one-ss)*expss+two*c_wc*ss/(one+c_wc*ss*ss)
4120              p2_wc=d_wc*(ss-two)*expss+two*c_wc/(one+c_wc*ss*ss)-&
4121 &             four*c_wc*c_wc*ss*ss/((one+c_wc*ss*ss)*(one+c_wc*ss*ss))
4122              divss=one/(one+(b_wc*ss+d_wc*ss*expss+log(one+c_wc*ss*ss))/kappa)
4123              dfxdss=p1_wc*divss*divss
4124              d2fxdss2=p2_wc*divss*divss-two*divss*divss*divss*p1_wc*p1_wc/kappa
4125            else
4126              if(option/=6)then ! This is Perdew-Burke-Ernzerhof GGA, exchange part
4127                divss=one/(one+mu_divkappa*ss)
4128                dfxdss= mu*divss*divss
4129                d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
4130              else  ! This is RPBE modification
4131                divss=exp(-mu_divkappa*ss)
4132                dfxdss= mu*divss
4133                d2fxdss2=-mu*mu_divkappa*divss
4134              end if
4135            end if
4136            fx    = one+kappa*(one-divss)
4137            ex_gga= ex_lsd*fx
4138            dssdn=-eight*third*ss*rho_inv
4139            dfxdn  = dfxdss*dssdn
4140            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
4141 !          The new definition (v3.3) includes the division by the norm of the gradient
4142            dssdg =two*coeffss
4143            dfxdg=dfxdss*dssdg
4144            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
4145            exc=exc+ex_gga*rho
4146          end if
4147 
4148        end do
4149        exci(ipts)=exc*rhotot_inv
4150        if(exexch_==1)cycle
4151 
4152 !      -----------------------------------------------------------------------------
4153 !      Then takes care of the LSD correlation part of the functional
4154 
4155        if(option>0)then
4156 
4157          rs=rsfac*rhotmot
4158          sqr_rs=sq_rsfac*rhotmo6
4159          rsm1_2=sq_rsfac_inv*rhoto6
4160 
4161 !        Formulas A6-A8 of PW92LSD
4162          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
4163          if(option/=3 .and. option/=4)then
4164            ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
4165            ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
4166          else
4167            sqr_sqr_rs=max(1.e-15_dp,sqrt(sqr_rs))
4168            ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs/sqr_sqr_rs)
4169            ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+3.5_dp*ec0_b4*rs/sqr_sqr_rs)
4170          end if
4171          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
4172 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
4173          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
4174          ecrs0=ec0_q0*ec0_log
4175          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
4176 
4177          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
4178          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
4179          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
4180          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
4181          mac_log=-log( mac_q1*mac_q1*mac_den )
4182          macrs=mac_q0*mac_log
4183          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
4184 
4185          zeta=(rho_updn(ipts,1)-rho_updn(ipts,2))*rhotot_inv
4186          ec1_q0=-2.0_dp*ec1_aa*(1.0_dp+ec1_a1*rs)
4187          if(option/=3 .and. option/=4)then
4188            ec1_q1=2.0_dp*ec1_aa*(ec1_b1*sqr_rs+ec1_b2*rs+ec1_b3*rs*sqr_rs+ec1_b4*rs*rs)
4189            ec1_q1p=ec1_aa*(ec1_b1*rsm1_2+2._dp*ec1_b2+3._dp*ec1_b3*sqr_rs+4._dp*ec1_b4*rs)
4190          else
4191            ec1_q1=2.0_dp*ec1_aa*(ec1_b1*sqr_rs+ec1_b2*rs+ec1_b3*rs*sqr_rs+ec1_b4*rs*rs/sqr_sqr_rs)
4192            ec1_q1p=ec1_aa*(ec1_b1*rsm1_2+2._dp*ec1_b2+3._dp*ec1_b3*sqr_rs+3.5_dp*ec1_b4*rs/sqr_sqr_rs)
4193          end if
4194          ec1_den=1.0_dp/(ec1_q1*ec1_q1+ec1_q1)
4195 !        ec1_log=log( 1.0_dp + 1.0_dp / ec1_q1 )
4196          ec1_log=-log( ec1_q1*ec1_q1*ec1_den )
4197          ecrs1=ec1_q0*ec1_log
4198          decrs1_drs= -2.0_dp*ec1_aa*ec1_a1*ec1_log - ec1_q0*ec1_q1p *ec1_den
4199 
4200 !        alpha_zeta is introduced in order to remove singularities for fully
4201 !        polarized systems.
4202          zetp_1_3=(1.0_dp+zeta*alpha_zeta)*zetpm1_3(ipts)**2
4203          zetm_1_3=(1.0_dp-zeta*alpha_zeta)*zetmm1_3(ipts)**2
4204 
4205          f_zeta=( (1.0_dp+zeta*alpha_zeta2)*zetp_1_3 +                      &
4206 &         (1.0_dp-zeta*alpha_zeta2)*zetm_1_3 - 2.0_dp ) * factf_zeta
4207          fp_zeta=( zetp_1_3 - zetm_1_3 ) * factfp_zeta
4208          zeta4=zeta**4
4209 
4210          gcrs=ecrs1-ecrs0+macrs*fsec_inv
4211 !        ecrs=ecrs0+f_zeta*(-macrs*(1.0_dp-zeta4)*fsec_inv+(ecrs1-ecrs0)*zeta4)
4212          ecrs=ecrs0+f_zeta*(zeta4*gcrs-macrs*fsec_inv)
4213 
4214          dgcrs_drs=decrs1_drs-decrs0_drs+dmacrs_drs*fsec_inv
4215 !        decrs_drs=decrs0_drs+f_zeta*&
4216 !        &        (-dmacrs_drs*(1.0_dp-zeta4)*fsec_inv+(decrs1_drs-decrs0_drs)*zeta4)
4217          decrs_drs=decrs0_drs+f_zeta*(zeta4*dgcrs_drs-dmacrs_drs*fsec_inv)
4218          dfzeta4_dzeta=4.0_dp*zeta**3*f_zeta+fp_zeta*zeta4
4219          decrs_dzeta=dfzeta4_dzeta*gcrs-fp_zeta*macrs*fsec_inv
4220 
4221 !        Add LSD correlation functional to GGA exchange functional
4222          exci(ipts)=exci(ipts)+ecrs
4223          vxcadd=ecrs-rs*third*decrs_drs-zeta*decrs_dzeta
4224          vxci(ipts,1)=vxci(ipts,1)+vxcadd+decrs_dzeta
4225          vxci(ipts,2)=vxci(ipts,2)+vxcadd-decrs_dzeta
4226 
4227 !        -----------------------------------------------------------------------------
4228 !        Eventually add the GGA correlation part of the PBE functional
4229 !        Note : the computation of the potential in the spin-unpolarized
4230 !        case could be optimized much further. Other optimizations are left to do.
4231 
4232          if(option==2 .or. option==5 .or. option==6 .or. option==7)then
4233 !          The definition of phi has been slightly changed, because
4234 !          the original PBE one gives divergent behaviour for fully
4235 !          polarized points
4236 !          zetpm1_3=(1.0_dp+zeta*alpha_zeta)**(-third)
4237 !          zetmm1_3=(1.0_dp-zeta*alpha_zeta)**(-third)
4238            phi_zeta=( zetpm1_3(ipts)*(1.0_dp+zeta*alpha_zeta)+ &
4239 &           zetmm1_3(ipts)*(1.0_dp-zeta*alpha_zeta)   )*0.5_dp
4240            phip_zeta=(zetpm1_3(ipts)-zetmm1_3(ipts))*third*alpha_zeta
4241            phi_zeta_inv=1.0_dp/phi_zeta
4242            phi_logder=phip_zeta*phi_zeta_inv
4243            phi3_zeta=phi_zeta*phi_zeta*phi_zeta
4244            gamphi3inv=gamma_inv*phi_zeta_inv*phi_zeta_inv*phi_zeta_inv
4245 
4246 !          From ec to bb
4247            bb=ecrs*gamphi3inv
4248            dbb_drs=decrs_drs*gamphi3inv
4249            dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
4250 !          From bb to cc
4251            exp_pbe=exp(-bb)
4252            cc=one/(exp_pbe-one)
4253            dcc_dbb=cc*cc*exp_pbe
4254            dcc_drs=dcc_dbb*dbb_drs
4255            dcc_dzeta=dcc_dbb*dbb_dzeta
4256 
4257 !          From cc to aa
4258            coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
4259            aa=coeff_aa*cc
4260            daa_drs=coeff_aa*dcc_drs
4261            daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
4262 !          Introduce tt : do not assume that the spin-dependent gradients are collinear
4263            grrho2=grho2_updn(ipts,3)
4264            dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
4265 !          Note that tt is (the t variable of PBE divided by phi) squared
4266            tt=half*grrho2*dtt_dg
4267 
4268 !          Get xx from aa and tt
4269            xx=aa*tt
4270            dxx_drs=daa_drs*tt
4271            dxx_dzeta=daa_dzeta*tt
4272            dxx_dtt=aa
4273 !          From xx to pade
4274            pade_den=one/(one+xx*(one+xx))
4275            pade=(one+xx)*pade_den
4276            dpade_dxx=-xx*(two+xx)*pade_den**2
4277            dpade_drs=dpade_dxx*dxx_drs
4278            dpade_dtt=dpade_dxx*dxx_dtt
4279            dpade_dzeta=dpade_dxx*dxx_dzeta
4280 
4281 !          From pade to qq
4282            coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
4283            qq=coeff_qq*pade
4284            dqq_drs=coeff_qq*dpade_drs
4285            dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
4286            dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
4287 
4288 !          From qq to rr
4289            arg_rr=one+beta*gamma_inv*qq
4290            div_rr=one/arg_rr
4291            rr=gamma*log(arg_rr)
4292            drr_dqq=beta*div_rr
4293            drr_drs=drr_dqq*dqq_drs
4294            drr_dtt=drr_dqq*dqq_dtt
4295            drr_dzeta=drr_dqq*dqq_dzeta
4296 
4297 !          From rr to hh
4298            hh=phi3_zeta*rr
4299            dhh_drs=phi3_zeta*drr_drs
4300            dhh_dtt=phi3_zeta*drr_dtt
4301            dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
4302 
4303 !          The GGA correlation energy is added
4304            exci(ipts)=exci(ipts)+hh
4305 
4306 !          Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
4307 
4308 !          From hh to the derivative of the energy wrt the density
4309            drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
4310            vxci(ipts,1)=vxci(ipts,1)+drhohh_drho+dhh_dzeta
4311            vxci(ipts,2)=vxci(ipts,2)+drhohh_drho-dhh_dzeta
4312 
4313 
4314 !          From hh to the derivative of the energy wrt to the gradient of the
4315 !          density, divided by the gradient of the density
4316 !          (The v3.3 definition includes the division by the norm of the gradient)
4317            dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
4318 
4319 !          End condition of GGA
4320          end if
4321 
4322        else  ! no correlation
4323 
4324 !        End condition of including correlation, and not only exchange
4325        end if
4326 
4327 !      Correlation has been added
4328 !      -----------------------------------------------------------------------------
4329 
4330      end do
4331 
4332 !    fab: the following is the "else" on order !!!
4333 
4334    else
4335 
4336      do ipts=1,npts
4337 
4338        rhotot=rhoarr(ipts)
4339        rhotmot=rhom1_3(ipts)
4340        rhotot_inv=rhotmot*rhotmot*rhotmot
4341        rhotmo6=sqrt(rhotmot)
4342        rhoto6=rhotot*rhotmot*rhotmot*rhotmo6
4343 !      -----------------------------------------------------------------------
4344 !      First take care of the exchange part of the functional
4345 
4346        exc=zero
4347        if (present(dvxcdgr)) dvxcdgr(ipts,3)=zero
4348        do ispden=1,nspden
4349          rho   =rho_updn(ipts,ispden)
4350          rhomot=rho_updnm1_3(ipts,ispden)
4351          ex_lsd= - threefourth_divpi * sixpi2_1_3*rhomot*rhomot*rho
4352          if(option==1 .or. option==-1 .or. option==3)then
4353 !          Perdew_Wang 91 LSD
4354            vxci(ipts,ispden)=four_thirds*ex_lsd
4355            if(present(dvxcdgr)) dvxcdgr(ipts,ispden)=0.0_dp
4356            exc=exc+ex_lsd*rho
4357          else
4358            rho_inv=rhomot*rhomot*rhomot
4359            coeffss=quarter*sixpi2m1_3*sixpi2m1_3*rho_inv*rho_inv*rhomot*rhomot
4360            ss=grho2_updn(ipts,ispden)*coeffss
4361            if(option==7) then ! This is WC
4362              expss=exp(-ss)
4363              p1_wc=b_wc+(mu-b_wc)*(one-ss)*expss+two*c_wc*ss/(one+c_wc*ss*ss)
4364              p2_wc=d_wc*(ss-two)*expss+two*c_wc/(one+c_wc*ss*ss)-&
4365 &             four*c_wc*c_wc*ss*ss/((one+c_wc*ss*ss)*(one+c_wc*ss*ss))
4366              divss=one/(one+(b_wc*ss+d_wc*ss*expss+log(one+c_wc*ss*ss))/kappa)
4367              dfxdss=p1_wc*divss*divss
4368              d2fxdss2=p2_wc*divss*divss-two*divss*divss*divss*p1_wc*p1_wc/kappa
4369            else
4370              if(option/=6)then ! This is Perdew-Burke-Ernzerhof GGA, exchange part
4371                divss=one/(one+mu_divkappa*ss)
4372                dfxdss= mu*divss*divss
4373                d2fxdss2=-mu*two*mu_divkappa*divss*divss*divss
4374              else  ! This is RPBE modification
4375                divss=exp(-mu_divkappa*ss)
4376                dfxdss= mu*divss
4377                d2fxdss2=-mu*mu_divkappa*divss
4378              end if
4379            end if
4380            fx    = one+kappa*(one-divss)
4381            ex_gga= ex_lsd*fx
4382            dssdn=-eight*third*ss*rho_inv
4383            dfxdn  = dfxdss*dssdn
4384            vxci(ipts,ispden)=ex_lsd*(four_thirds*fx+rho*dfxdn)
4385 !          The new definition (v3.3) includes the division by the norm of the gradient
4386            dssdg =two*coeffss
4387            dfxdg=dfxdss*dssdg
4388            dvxcdgr(ipts,ispden)=ex_lsd*rho*dfxdg
4389            exc=exc+ex_gga*rho
4390          end if
4391 
4392          if(option==1 .or. option==-1 .or. option==3)then
4393 
4394 !          Perdew_Wang 91 LSD
4395            dvxci(ipts,2*ispden-1)=-four_thirds*third*&
4396 &           threefourth_divpi*sixpi2_1_3*rhomot*rhomot
4397            dvxci(ipts,2)=zero
4398            if(order==3)then
4399 !            If non-spin-polarized, first component of dvxci is second
4400 !            derivative with respect to TOTAL density.
4401 !            Compute the second derivative of vx
4402 !            vx^(2) = -2*vx^(1)/(3*rhotot)
4403 
4404 !            fab: third order derivatives of the exchange part in the spin polarized case
4405 
4406              d2vxci(ipts,3*ispden-2) = -2._dp*dvxci(ipts,2*ispden-1)*(rhomot*rhomot*rhomot)/3._dp
4407 
4408 !            mixed thir order derivatives of the exchange energy with respect to rho of
4409 !            different spin polarization are zero
4410              d2vxci(ipts,2)=zero
4411              d2vxci(ipts,3)=zero
4412 
4413            end if
4414 
4415          else
4416 !          Perdew-Burke-Ernzerhof GGA, exchange part
4417 !          Components 3 or 4
4418            dvxci(ipts,2+ispden)=dvxcdgr(ipts,ispden)
4419 !          Components 1 or 2
4420            d2ssdn2=-11.0_dp*third*dssdn*rho_inv
4421            d2fxdn2=d2fxdss2*dssdn**2+dfxdss*d2ssdn2
4422            dvxci(ipts,ispden)=third*rho_inv*vxci(ipts,ispden)+&
4423 &           ex_lsd*(seven*third*dfxdn+rho*d2fxdn2)
4424 !          Components 5 or 6
4425            d2ssdndg=-eight*third*dssdg*rho_inv
4426            d2fxdndg=d2fxdss2*dssdn*dssdg+dfxdss*d2ssdndg
4427            dvxci(ipts,4+ispden)=ex_lsd*(four_thirds*dfxdg+rho*d2fxdndg)
4428 !          Components 7 or 8
4429            d2fxdg2=d2fxdss2*dssdg**2
4430            dvxci(ipts,6+ispden)=ex_lsd*rho*d2fxdg2
4431 !          For the time being, treat non-spin-polarized like spin-polarized
4432          end if
4433        end do
4434        exci(ipts)=exc*rhotot_inv
4435 !      -----------------------------------------------------------------------------
4436 !      Then takes care of the LSD correlation part of the functional
4437 
4438        if(option>0)then
4439 
4440          rs=rsfac*rhotmot
4441          sqr_rs=sq_rsfac*rhotmo6
4442          rsm1_2=sq_rsfac_inv*rhoto6
4443 
4444 !        Formulas A6-A8 of PW92LSD
4445          ec0_q0=-2.0_dp*ec0_aa*(1.0_dp+ec0_a1*rs)
4446          if(option/=3 .and. option/=4)then
4447            ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs)
4448            ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+4._dp*ec0_b4*rs)
4449          else
4450            sqr_sqr_rs=max(1.e-15_dp,sqrt(sqr_rs))
4451            ec0_q1=2.0_dp*ec0_aa*(ec0_b1*sqr_rs+ec0_b2*rs+ec0_b3*rs*sqr_rs+ec0_b4*rs*rs/sqr_sqr_rs)
4452            ec0_q1p=ec0_aa*(ec0_b1*rsm1_2+2._dp*ec0_b2+3._dp*ec0_b3*sqr_rs+3.5_dp*ec0_b4*rs/sqr_sqr_rs)
4453          end if
4454          ec0_den=1.0_dp/(ec0_q1*ec0_q1+ec0_q1)
4455 !        ec0_log=log( 1.0_dp + 1.0_dp / ec0_q1 )
4456          ec0_log=-log( ec0_q1*ec0_q1*ec0_den )
4457          ecrs0=ec0_q0*ec0_log
4458          decrs0_drs= -2.0_dp*ec0_aa*ec0_a1*ec0_log - ec0_q0*ec0_q1p *ec0_den
4459          ec0_q1pp=0.5_dp*ec0_aa*(-ec0_b1*rsm1_2**3+3._dp*ec0_b3*rsm1_2+8._dp*ec0_b4)
4460          d2ecrs0_drs2= 4.0_dp*ec0_aa*ec0_a1*ec0_q1p*ec0_den            &
4461 &         -ec0_q0*ec0_q1pp*ec0_den                        &
4462 &         +ec0_q0*ec0_q1p**2*ec0_den**2*(2._dp*ec0_q1+1.0_dp)
4463          if (order==3) then
4464            ec0_q1ppp = 0.75_dp*ec0_aa*(rsm1_2**5)*(ec0_b1-ec0_b3*rs)
4465            ec0_f1 = 1._dp/(ec0_q1*ec0_q1*(1._dp + ec0_q1))
4466            ec0_f2 = 1._dp/(ec0_q1*(1+ec0_q1))
4467            d3ecrs0_drs3 = 6._dp*ec0_q1p*ec0_f1*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + &
4468 &           ec0_q0*ec0_q1pp) - &
4469 &           ec0_f2*(-6._dp*ec0_aa*ec0_a1*ec0_q1pp + ec0_q0*ec0_q1ppp + &
4470 &           ec0_f2*(3._dp*ec0_q1p*(-2._dp*ec0_aa*ec0_a1*ec0_q1p + ec0_q0*ec0_q1pp) + &
4471 &           ec0_f2*2._dp*ec0_q0*(ec0_q1p**3)*(1._dp + 3._dp*ec0_q1*(1._dp + ec0_q1))))
4472          end if
4473 
4474          mac_q0=-2.0_dp*mac_aa*(1.0_dp+mac_a1*rs)
4475          mac_q1=2.0_dp*mac_aa*(mac_b1*sqr_rs+mac_b2*rs+mac_b3*rs*sqr_rs+mac_b4*rs*rs)
4476          mac_q1p=mac_aa*(mac_b1*rsm1_2+2._dp*mac_b2+3._dp*mac_b3*sqr_rs+4._dp*mac_b4*rs)
4477          mac_den=1.0_dp/(mac_q1*mac_q1+mac_q1)
4478          mac_log=-log( mac_q1*mac_q1*mac_den )
4479          macrs=mac_q0*mac_log
4480          dmacrs_drs= -2.0_dp*mac_aa*mac_a1*mac_log - mac_q0*mac_q1p*mac_den
4481 
4482          zeta=(rho_updn(ipts,1)-rho_updn(ipts,2))*rhotot_inv
4483          ec1_q0=-2.0_dp*ec1_aa*(1.0_dp+ec1_a1*rs)
4484          if(option/=3 .and. option/=4)then
4485            ec1_q1=2.0_dp*ec1_aa*(ec1_b1*sqr_rs+ec1_b2*rs+ec1_b3*rs*sqr_rs+ec1_b4*rs*rs)
4486            ec1_q1p=ec1_aa*(ec1_b1*rsm1_2+2._dp*ec1_b2+3._dp*ec1_b3*sqr_rs+4._dp*ec1_b4*rs)
4487          else
4488            ec1_q1=2.0_dp*ec1_aa*(ec1_b1*sqr_rs+ec1_b2*rs+ec1_b3*rs*sqr_rs+ec1_b4*rs*rs/sqr_sqr_rs)
4489            ec1_q1p=ec1_aa*(ec1_b1*rsm1_2+2._dp*ec1_b2+3._dp*ec1_b3*sqr_rs+3.5_dp*ec1_b4*rs/sqr_sqr_rs)
4490          end if
4491          ec1_den=1.0_dp/(ec1_q1*ec1_q1+ec1_q1)
4492 !        ec1_log=log( 1.0_dp + 1.0_dp / ec1_q1 )
4493          ec1_log=-log( ec1_q1*ec1_q1*ec1_den )
4494          ecrs1=ec1_q0*ec1_log
4495          decrs1_drs= -2.0_dp*ec1_aa*ec1_a1*ec1_log - ec1_q0*ec1_q1p *ec1_den
4496 
4497 !        alpha_zeta is introduced in order to remove singularities for fully
4498 !        polarized systems.
4499          zetp_1_3=(1.0_dp+zeta*alpha_zeta)*zetpm1_3(ipts)**2
4500          zetm_1_3=(1.0_dp-zeta*alpha_zeta)*zetmm1_3(ipts)**2
4501 
4502          f_zeta=( (1.0_dp+zeta*alpha_zeta2)*zetp_1_3 +                      &
4503 &         (1.0_dp-zeta*alpha_zeta2)*zetm_1_3 - 2.0_dp ) * factf_zeta
4504          fp_zeta=( zetp_1_3 - zetm_1_3 ) * factfp_zeta
4505          zeta4=zeta**4
4506 
4507          gcrs=ecrs1-ecrs0+macrs*fsec_inv
4508 !        ecrs=ecrs0+f_zeta*(-macrs*(1.0_dp-zeta4)*fsec_inv+(ecrs1-ecrs0)*zeta4)
4509          ecrs=ecrs0+f_zeta*(zeta4*gcrs-macrs*fsec_inv)
4510 
4511          dgcrs_drs=decrs1_drs-decrs0_drs+dmacrs_drs*fsec_inv
4512 !        decrs_drs=decrs0_drs+f_zeta*&
4513 !        &        (-dmacrs_drs*(1.0_dp-zeta4)*fsec_inv+(decrs1_drs-decrs0_drs)*zeta4)
4514          decrs_drs=decrs0_drs+f_zeta*(zeta4*dgcrs_drs-dmacrs_drs*fsec_inv)
4515          dfzeta4_dzeta=4.0_dp*zeta**3*f_zeta+fp_zeta*zeta4
4516          decrs_dzeta=dfzeta4_dzeta*gcrs-fp_zeta*macrs*fsec_inv
4517 
4518          ec1_q1pp=0.5_dp*ec1_aa*(-ec1_b1*rsm1_2**3+3._dp*ec1_b3*rsm1_2+8._dp*ec1_b4)
4519 
4520          d2ecrs1_drs2= 4.0_dp*ec1_aa*ec1_a1*ec1_q1p*ec1_den            &
4521 &         -ec1_q0*ec1_q1pp*ec1_den                        &
4522 &         +ec1_q0*ec1_q1p**2*ec1_den**2*(2._dp*ec1_q1+1.0_dp)
4523 
4524 
4525          mac_q1pp=0.5_dp*mac_aa*(-mac_b1*rsm1_2**3+3._dp*mac_b3*rsm1_2+8._dp*mac_b4)
4526          d2macrs_drs2= 4.0_dp*mac_aa*mac_a1*mac_q1p*mac_den            &
4527 &         -mac_q0*mac_q1pp*mac_den                        &
4528 &         +mac_q0*mac_q1p**2*mac_den**2*(2._dp*mac_q1+1.0_dp)
4529 
4530          d2gcrs_drs2=d2ecrs1_drs2-d2ecrs0_drs2+d2macrs_drs2*fsec_inv
4531          fpp_zeta=(zetpm1_3(ipts)**2+zetmm1_3(ipts)**2) * factfpp_zeta
4532          d2fzeta4_dzeta2=12.0_dp*zeta**2*f_zeta  &
4533 &         + 8.0_dp*zeta**3*fp_zeta &
4534 &         +       zeta4  *fpp_zeta
4535 
4536          d2ecrs_drs2=d2ecrs0_drs2+&
4537 &         f_zeta*(zeta4*d2gcrs_drs2-d2macrs_drs2*fsec_inv)
4538          d2ecrs_drsdzeta=dfzeta4_dzeta*dgcrs_drs-fp_zeta*dmacrs_drs*fsec_inv
4539          d2ecrs_dzeta2=d2fzeta4_dzeta2*gcrs-fpp_zeta*macrs*fsec_inv
4540 
4541 !        End condition of abs(order)>1
4542 !        Add LSD correlation functional to GGA exchange functional
4543          exci(ipts)=exci(ipts)+ecrs
4544          vxcadd=ecrs-rs*third*decrs_drs-zeta*decrs_dzeta
4545 !        decrs_drup=vxcadd+decrs_dzeta
4546 !        decrs_drdn=vxcadd-decrs_dzeta
4547          vxci(ipts,1)=vxci(ipts,1)+vxcadd+decrs_dzeta
4548          vxci(ipts,2)=vxci(ipts,2)+vxcadd-decrs_dzeta
4549 
4550 
4551 
4552          dvcrs_drs=third*(2._dp*decrs_drs-rs*d2ecrs_drs2)
4553 !        And d(vxc)/d(rho)=(-rs/(3*rho))*d(vxc)/d(rs)
4554          d2ecrs_drho2= -rs**4*(four_pi*third)*third*dvcrs_drs
4555          d2ecrs_drup2=d2ecrs_drho2+&
4556 &         two*(-third*rs*d2ecrs_drsdzeta)*(1._dp-zeta)*rhotot_inv+ &
4557 &         d2ecrs_dzeta2*(1._dp-zeta)**2*rhotot_inv
4558          d2ecrs_drdndrup=d2ecrs_drho2+&
4559 &         2.0_dp*(-third*rs*d2ecrs_drsdzeta)*(-zeta)*rhotot_inv+ &
4560 &         d2ecrs_dzeta2*(1._dp-zeta)*(-1._dp-zeta)*rhotot_inv
4561          d2ecrs_drdn2=d2ecrs_drho2+&
4562 &         2.0_dp*(-third*rs*d2ecrs_drsdzeta)*(-1._dp-zeta)*rhotot_inv+ &
4563 &         d2ecrs_dzeta2*(-1._dp-zeta)**2*rhotot_inv
4564 
4565 
4566 
4567 
4568          if (order==3) then
4569 
4570 
4571 !          fab : INGREDIENTS NEEDED
4572 
4573            a1fa=-third*(threefourth_divpi**(third))*((rhotot)**(-4._dp/3._dp))
4574            a2fa=(1._dp-zeta)/rhotot
4575            b2fa=(-2._dp/3._dp)*((threefourth_divpi)**(third))*((7._dp/3._dp)*(-1._dp+zeta)/((rhotot)**(7._dp/3._dp)))
4576            b1fa=a2fa
4577            c2fa=((1._dp-zeta)**2)*(-3._dp*(1._dp/((rhotot)**2)))
4578            c1fa=((1._dp-zeta)*(1._dp-zeta))/rhotot
4579            e2fa=(2._dp/3._dp)*((threefourth_divpi)**(third))*((1._dp-(7._dp/3._dp)*zeta)/((rhotot)**(7._dp/3._dp)))
4580            e1fa=-zeta/rhotot
4581            f2fa=(2._dp*zeta)*(1._dp/((rhotot)**2))-(3._dp*zeta*zeta)*(1._dp/((rhotot)**2))+1._dp/(((rhotot)**2))
4582            f1fa=(zeta*zeta-1._dp)/rhotot
4583            g1fa=a1fa
4584            g2fa=(-1._dp-zeta)/rhotot
4585            h2fa=(2._dp/3._dp)*((threefourth_divpi)**(third))*((-1._dp-(7._dp/3._dp)*zeta)/((rhotot)**(7._dp/3._dp)))
4586            h1fa=e1fa
4587            i2fa=((-2._dp*zeta)-(3*zeta*zeta)+1)/((rhotot)**2)
4588            i1fa=f1fa
4589            m2fa=(-2._dp/3._dp)*((threefourth_divpi)**(third))*(((7._dp/3._dp)*(zeta+1._dp))/((rhotot)**(7._dp/3._dp)))
4590            m1fa=g2fa
4591            n2fa=(-3._dp*(1._dp+zeta)*(1._dp+zeta))/(rhotot*rhotot)
4592            n1fa=((-1._dp-zeta)**2)/rhotot
4593 
4594 
4595 !          TERMS APPEARING IN THE THIRD ORDER DERIVATIVES
4596 !          terms appearing in the third order derivatives of the spin polarized
4597 !          correlation energy with respect to spin densities
4598 
4599 
4600 !          ec1_q0p=-2.0_dp*ec1_aa*ec1_a1
4601 !          ec1_q1ppp=(3._dp/4._dp)*ec1_aa*(ec1_b1*(rsm1_2**5)-ec1_b3*(rsm1_2**3))
4602 !          This must be erroneous ...
4603 !          d3ecrs1_drs3=(ec1_q1pp*(4._dp*ec1_aa*ec1_a1-ec1_q0p)-ec1_q0*ec1_q1ppp)*ec1_den+ &
4604 !          &           ((-ec1_q1p**2)*(4._dp*ec1_aa*ec1_a1)+ec1_q0*ec1_q1pp*ec1_q1+ec1_q0p*(ec1_q1p**2)+ &
4605 !          &           2._dp*ec1_q0*ec1_q1p*ec1_q1pp)*(ec1_den**2)*(2._dp*ec1_q1+1._dp)- &
4606 !          &           (2._dp*ec1_q0*(ec1_q1p**3)*((2._dp*ec1_q1+1._dp)**2))*(ec1_den**3)+  &
4607 !          &           (2._dp*ec1_q0*(ec1_q1p**3))*(ec1_den**2)
4608 
4609            ec1_q1ppp = 0.75_dp*ec1_aa*(rsm1_2**5)*(ec1_b1-ec1_b3*rs)
4610            ec1_f1 = 1._dp/(ec1_q1*ec1_q1*(1._dp + ec1_q1))
4611            ec1_f2 = 1._dp/(ec1_q1*(1+ec1_q1))
4612            d3ecrs1_drs3 = 6._dp*ec1_q1p*ec1_f1*(-2._dp*ec1_aa*ec1_a1*ec1_q1p + &
4613 &           ec1_q0*ec1_q1pp) - &
4614 &           ec1_f2*(-6._dp*ec1_aa*ec1_a1*ec1_q1pp + ec1_q0*ec1_q1ppp + &
4615 &           ec1_f2*(3._dp*ec1_q1p*(-2._dp*ec1_aa*ec1_a1*ec1_q1p + ec1_q0*ec1_q1pp) + &
4616 &           ec1_f2*2._dp*ec1_q0*(ec1_q1p**3)*(1._dp + 3._dp*ec1_q1*(1._dp + ec1_q1))))
4617 
4618 
4619 !          mac_q0p=-2.0_dp*mac_aa*mac_a1
4620 !          mac_q1ppp=(3._dp/4._dp)*mac_aa*(mac_b1*((rsm1_2)**5)-mac_b3*((rsm1_2)**3))
4621 !          This must be erroneous ...
4622 !          d3macrs_drs3=(mac_q1pp*(4._dp*mac_aa*mac_a1-mac_q0p)-mac_q0*mac_q1ppp)*mac_den+ &
4623 !          &           ((-mac_q1p**2)*(4._dp*mac_aa*mac_a1)+mac_q0*mac_q1pp*mac_q1+mac_q0p*(mac_q1p**2)+ &
4624 !          &           2._dp*mac_q0*mac_q1p*mac_q1pp)*(mac_den**2)*(2._dp*mac_q1+1._dp)- &
4625 !          &           (2._dp*mac_q0*(mac_q1p**3)*((2._dp*mac_q1+1._dp)**2))*(mac_den**3)+  &
4626 !          &           (2._dp*mac_q0*(mac_q1p**3))*(mac_den**2)
4627 
4628            mac_q1ppp = 0.75_dp*mac_aa*(rsm1_2**5)*(mac_b1-mac_b3*rs)
4629            mac_f1 = 1._dp/(mac_q1*mac_q1*(1._dp + mac_q1))
4630            mac_f2 = 1._dp/(mac_q1*(1+mac_q1))
4631            d3macrs_drs3 = 6._dp*mac_q1p*mac_f1*(-2._dp*mac_aa*mac_a1*mac_q1p + &
4632 &           mac_q0*mac_q1pp) - &
4633 &           mac_f2*(-6._dp*mac_aa*mac_a1*mac_q1pp + mac_q0*mac_q1ppp + &
4634 &           mac_f2*(3._dp*mac_q1p*(-2._dp*mac_aa*mac_a1*mac_q1p + mac_q0*mac_q1pp) + &
4635 &           mac_f2*2._dp*mac_q0*(mac_q1p**3)*(1._dp + 3._dp*mac_q1*(1._dp + mac_q1))))
4636 
4637            d3gcrs_drs3=d3ecrs1_drs3-d3ecrs0_drs3+d3macrs_drs3*fsec_inv
4638            d3ecrs_drs3=d3ecrs0_drs3+f_zeta*(zeta4*d3gcrs_drs3-d3macrs_drs3*fsec_inv)
4639            factfppp_zeta=-two*third*factfpp_zeta*alpha_zeta2
4640            fppp_zeta=factfppp_zeta*(((zetpm1_3(ipts))**5)-((zetmm1_3(ipts))**5))
4641 
4642            d3ecrs_dzeta3=(24._dp*zeta*f_zeta+36._dp*(zeta**2)*fp_zeta+  &
4643 &           12._dp*(zeta**3)*fpp_zeta+(zeta**4)*fppp_zeta)*gcrs+   &
4644 &           fppp_zeta*(-macrs)*fsec_inv
4645 
4646            d3ecrs_drs2dzeta=dfzeta4_dzeta*(d2gcrs_drs2)+   &
4647 &           fp_zeta*(-d2macrs_drs2)*fsec_inv
4648 
4649            d3ecrs_dzeta2drs=d2fzeta4_dzeta2*dgcrs_drs+  &
4650 &           fpp_zeta*(-dmacrs_drs)*fsec_inv
4651 
4652 
4653 
4654 
4655 !          ***************** all this part has been commented following the suggestion by xavier
4656 
4657 !          The following is the calculations of only a part of the third order derivatives.
4658 !          the term d3ecrs_drho3 is the part which remains in the
4659 !          non spin polarized limit
4660 !          THIS CODING IS CORRECT, but the alternative one is also correct ...
4661 
4662 !          d3ecrs_drho3=(128._dp*pi*pi/243._dp)*(rs**7)*(decrs_drs)- &
4663 !          &           (48._dp*pi*pi/243._dp)*(rs**8)*d2ecrs_drs2- &
4664 !          &           (16._dp*pi*pi/243._dp)*(rs**9)*d3ecrs_drs3
4665 
4666 !          d3ecrs_drhoupdrho2=d3ecrs_drho3+ &
4667 !          &           a2fa*((-8._dp*pi/27._dp)*(rs**4)*d2ecrs_drsdzeta+ &
4668 !          &           (4._dp*pi/27._dp)*(rs**5)*d3ecrs_drs2dzeta)
4669 
4670 
4671 !          d3ecrs_drhodndrho2=d3ecrs_drho3+ &
4672 !          &           g2fa*((-8._dp*pi/27._dp)*(rs**4)*d2ecrs_drsdzeta+ &
4673 !          &           (4._dp*pi/27._dp)*(rs**5)*d3ecrs_drs2dzeta)
4674 
4675 
4676 !          third order derivatives of the exchange-correlation part
4677 
4678 !          d3ecrs_drup3=d3ecrs_drhoupdrho2+b2fa*d2ecrs_drsdzeta-  &
4679 !          &           (2._dp/3._dp)*rs*b1fa*a1fa*d3ecrs_drs2dzeta- &
4680 !          &           (2._dp/3._dp)*rs*b1fa*b1fa*d3ecrs_dzeta2drs+ &
4681 !          &           c2fa*d2ecrs_dzeta2+c1fa*(a1fa*d3ecrs_dzeta2drs+a2fa*d3ecrs_dzeta3)
4682 
4683 
4684 !          d3ecrs_drup2drdn=d3ecrs_drhoupdrho2+e2fa*d2ecrs_drsdzeta-  &
4685 !          &           (2._dp/3._dp)*rs*e1fa*a1fa*d3ecrs_drs2dzeta- &
4686 !          &           (2._dp/3._dp)*rs*e1fa*b1fa*d3ecrs_dzeta2drs+ &
4687 !          &           f2fa*d2ecrs_dzeta2+f1fa*(a1fa*d3ecrs_dzeta2drs+a2fa*d3ecrs_dzeta3)
4688 
4689 !          d3ecrs_drupdrdn2=d3ecrs_drhodndrho2+h2fa*d2ecrs_drsdzeta-  &
4690 !          &           (2._dp/3._dp)*rs*h1fa*g1fa*d3ecrs_drs2dzeta- &
4691 !          &           (2._dp/3._dp)*rs*h1fa*g2fa*d3ecrs_dzeta2drs+ &
4692 !          &           i2fa*d2ecrs_dzeta2+i1fa*(g1fa*d3ecrs_dzeta2drs+g2fa*d3ecrs_dzeta3)
4693 
4694 
4695 !          d3ecrs_drdn3=d3ecrs_drhodndrho2+m2fa*d2ecrs_drsdzeta-  &
4696 !          &           (2._dp/3._dp)*rs*m1fa*g1fa*d3ecrs_drs2dzeta- &
4697 !          &           (2._dp/3._dp)*rs*m1fa*g2fa*d3ecrs_dzeta2drs+ &
4698 !          &           n2fa*d2ecrs_dzeta2+n1fa*(g1fa*d3ecrs_dzeta2drs+g2fa*d3ecrs_dzeta3)
4699 
4700 
4701 !          ********* suggested by xavier  (now corrected, XG100524)
4702 
4703            sp1_up3=three-three*zeta
4704            sp1_up2dn=one-three*zeta
4705            sp1_updn2=-one-three*zeta
4706            sp1_dn3=-three-three*zeta
4707 
4708            sp2_up3=three-six*zeta+three*zeta*zeta
4709            sp2_up2dn=-one-two*zeta+three*zeta*zeta
4710            sp2_updn2=-one+two*zeta+three*zeta*zeta
4711            sp2_dn3=three+six*zeta+three*zeta*zeta
4712 
4713            sp3_up3=(one-zeta)**3
4714            sp3_up2dn=-one+zeta+zeta**2-zeta**3
4715            sp3_updn2=one+zeta-zeta**2-zeta**3
4716            sp3_dn3=-(one+zeta)**3
4717 
4718            d3ecrs_sp0=(eight*rs*decrs_drs-three*rs*rs*d2ecrs_drs2-rs*rs*rs*d3ecrs_drs3)/(rhotot*rhotot*27.0_dp)
4719            d3ecrs_sp1=(four*rs*d2ecrs_drsdzeta+rs*rs*d3ecrs_drs2dzeta)/(rhotot*rhotot*nine)
4720            d3ecrs_sp2=(-three*d2ecrs_dzeta2-rs*d3ecrs_dzeta2drs)/(rhotot*rhotot*three)
4721            d3ecrs_sp3=d3ecrs_dzeta3/(rhotot*rhotot)
4722 
4723            d3ecrs_drup3=d3ecrs_sp0+d3ecrs_sp1*sp1_up3+d3ecrs_sp2*sp2_up3+d3ecrs_sp3*sp3_up3
4724            d3ecrs_drup2drdn=d3ecrs_sp0+d3ecrs_sp1*sp1_up2dn+d3ecrs_sp2*sp2_up2dn+d3ecrs_sp3*sp3_up2dn
4725            d3ecrs_drupdrdn2=d3ecrs_sp0+d3ecrs_sp1*sp1_updn2+d3ecrs_sp2*sp2_updn2+d3ecrs_sp3*sp3_updn2
4726            d3ecrs_drdn3=d3ecrs_sp0+d3ecrs_sp1*sp1_dn3+d3ecrs_sp2*sp2_dn3+d3ecrs_sp3*sp3_dn3
4727 
4728 !          **** end of the section suggested by xavier...
4729 
4730 
4731 !          fab: this is the end of the if over order==3
4732 
4733 
4734          end if
4735 
4736 !        fab: I think the following lines are wrong..indeed we are in the case option>0, so option cannot be -1
4737 !        I comment them and I put the if only for option 1 and 3
4738 
4739 !        if(option==1 .or. option==-1 .or. option==3)then
4740 !        dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drup2
4741 !        dvxci(ipts,2)=dvxci(ipts,2)+d2ecrs_drdndrup
4742 !        dvxci(ipts,3)=dvxci(ipts,3)+d2ecrs_drdn2
4743 
4744 
4745 !        fab: however, here the thing seems a bit strange...option=3 doesn't seem to be completely implemented
4746 !        (the second derivatives of ec1_q and ec0_q are the derived only in correspondance of the first derivative in the case option !=3 and !=4)
4747 !        so..I think that here the case "or option==3 should be cancelled
4748 
4749          if(option==1 .or. option==3)then
4750            dvxci(ipts,1)=dvxci(ipts,1)+d2ecrs_drup2
4751            dvxci(ipts,2)=dvxci(ipts,2)+d2ecrs_drdndrup
4752            dvxci(ipts,3)=dvxci(ipts,3)+d2ecrs_drdn2
4753 
4754            if(order==3) then
4755 !            third order derivatives of the spin polarized exchange+correlation energy
4756              d2vxci(ipts,1)=d2vxci(ipts,1)+d3ecrs_drup3
4757              d2vxci(ipts,2)=d2vxci(ipts,2)+d3ecrs_drup2drdn
4758              d2vxci(ipts,3)=d2vxci(ipts,3)+d3ecrs_drupdrdn2
4759              d2vxci(ipts,4)=d2vxci(ipts,4)+d3ecrs_drdn3
4760 !            DEBUG
4761 !            wecrsz(ipts,1)=ecrs
4762 !            wecrsz(ipts,1)=ecrs*rhotot
4763 !            wecrsz(ipts,2)=rs
4764 !            wecrsz(ipts,2)=rho_updn(ipts,1)
4765 !            wecrsz(ipts,3)=zeta
4766 !            wecrsz(ipts,3)=rho_updn(ipts,2)
4767 !            wecrsz(ipts,5)=ecrs0
4768 !            wecrsz(ipts,6)=gcrs
4769 !            wecrsz(ipts,7)=macrs
4770 !            wecrsz(ipts,8)=ecrs1
4771 !            d1wecrsz(ipts,1)=decrs_drs
4772 !            d1wecrsz(ipts,1)=decrs_drup
4773 !            d1wecrsz(ipts,2)=decrs_dzeta
4774 !            d1wecrsz(ipts,2)=decrs_drdn
4775 !            d1wecrsz(ipts,5)=decrs0_drs
4776 !            d1wecrsz(ipts,6)=dgcrs_drs
4777 !            d1wecrsz(ipts,7)=dmacrs_drs
4778 !            d1wecrsz(ipts,8)=decrs1_drs
4779 !            d2wecrsz(ipts,1)=d2ecrs_drs2
4780 !            d2wecrsz(ipts,1)=d2ecrs_drup2
4781 !            d2wecrsz(ipts,2)=d2ecrs_drsdzeta
4782 !            d2wecrsz(ipts,2)=d2ecrs_drdndrup
4783 !            d2wecrsz(ipts,3)=d2ecrs_dzeta2
4784 !            d2wecrsz(ipts,3)=d2ecrs_drdn2
4785 !            d2wecrsz(ipts,5)=d2ecrs0_drs2
4786 !            d2wecrsz(ipts,6)=d2gcrs_drs2
4787 !            d2wecrsz(ipts,7)=d2macrs_drs2
4788 !            d2wecrsz(ipts,8)=d2ecrs1_drs2
4789 !            d3wecrsz(ipts,1)=d3ecrs_drs3
4790 !            d3wecrsz(ipts,1)=d3ecrs_drup3
4791 !            d3wecrsz(ipts,2)=d3ecrs_drs2dzeta
4792 !            d3wecrsz(ipts,2)=d3ecrs_drup2drdn
4793 !            d3wecrsz(ipts,3)=d3ecrs_dzeta2drs
4794 !            d3wecrsz(ipts,3)=d3ecrs_drupdrdn2
4795 !            d3wecrsz(ipts,4)=d3ecrs_dzeta3
4796 !            d3wecrsz(ipts,4)=d3ecrs_drdn3
4797 !            d3wecrsz(ipts,5)=d3ecrs0_drs3
4798 !            d3wecrsz(ipts,6)=d3gcrs_drs3
4799 !            d3wecrsz(ipts,7)=d3macrs_drs3
4800 !            d3wecrsz(ipts,8)=d3ecrs1_drs3
4801 !            ENDDEBUG
4802            end if
4803 
4804 
4805          else
4806            dvxci(ipts,9)=d2ecrs_drup2
4807            dvxci(ipts,10)=d2ecrs_drdndrup
4808            dvxci(ipts,11)=d2ecrs_drdn2
4809          end if
4810 
4811 
4812 
4813 
4814 !        -----------------------------------------------------------------------------
4815 !        Eventually add the GGA correlation part of the PBE functional
4816 !        Note : the computation of the potential in the spin-unpolarized
4817 !        case could be optimized much further. Other optimizations are left to do.
4818 
4819          if(option==2 .or. option==5 .or. option==6 .or. option==7)then
4820 !          The definition of phi has been slightly changed, because
4821 !          the original PBE one gives divergent behaviour for fully
4822 !          polarized points
4823 !          zetpm1_3=(1.0_dp+zeta*alpha_zeta)**(-third)
4824 !          zetmm1_3=(1.0_dp-zeta*alpha_zeta)**(-third)
4825            phi_zeta=( zetpm1_3(ipts)*(1.0_dp+zeta*alpha_zeta)+ &
4826 &           zetmm1_3(ipts)*(1.0_dp-zeta*alpha_zeta)   )*0.5_dp
4827            phip_zeta=(zetpm1_3(ipts)-zetmm1_3(ipts))*third*alpha_zeta
4828            phi_zeta_inv=1.0_dp/phi_zeta
4829            phi_logder=phip_zeta*phi_zeta_inv
4830            phi3_zeta=phi_zeta*phi_zeta*phi_zeta
4831            gamphi3inv=gamma_inv*phi_zeta_inv*phi_zeta_inv*phi_zeta_inv
4832            phipp_zeta=-alpha_zeta*alpha_zeta*ninth*&
4833 &           (zetpm1_3(ipts)*zetpm1_3(ipts)*zetpm1_3(ipts)*zetpm1_3(ipts) + &
4834 &           zetmm1_3(ipts)*zetmm1_3(ipts)*zetmm1_3(ipts)*zetmm1_3(ipts)  )
4835 
4836 !          From ec to bb
4837            bb=ecrs*gamphi3inv
4838            dbb_drs=decrs_drs*gamphi3inv
4839            dbb_dzeta=gamphi3inv*(decrs_dzeta-three*ecrs*phi_logder)
4840            d2bb_drs2=d2ecrs_drs2*gamphi3inv
4841            d2bb_drsdzeta=gamphi3inv*(d2ecrs_drsdzeta-three*decrs_drs*phi_logder)
4842            d2bb_dzeta2=gamphi3inv*(d2ecrs_dzeta2-six*decrs_dzeta*phi_logder+&
4843 &           12.0_dp*ecrs*phi_logder*phi_logder-three*ecrs*phi_zeta_inv*phipp_zeta)
4844 
4845 !          From bb to cc
4846            exp_pbe=exp(-bb)
4847            cc=one/(exp_pbe-one)
4848            dcc_dbb=cc*cc*exp_pbe
4849            dcc_drs=dcc_dbb*dbb_drs
4850            dcc_dzeta=dcc_dbb*dbb_dzeta
4851            d2cc_dbb2=cc*cc*exp_pbe*(two*cc*exp_pbe-one)
4852            d2cc_drs2=d2cc_dbb2*dbb_drs*dbb_drs+dcc_dbb*d2bb_drs2
4853            d2cc_drsdzeta=d2cc_dbb2*dbb_drs*dbb_dzeta+dcc_dbb*d2bb_drsdzeta
4854            d2cc_dzeta2=d2cc_dbb2*dbb_dzeta*dbb_dzeta+dcc_dbb*d2bb_dzeta2
4855 
4856 !          From cc to aa
4857            coeff_aa=beta*gamma_inv*phi_zeta_inv*phi_zeta_inv
4858            aa=coeff_aa*cc
4859            daa_drs=coeff_aa*dcc_drs
4860            daa_dzeta=-two*aa*phi_logder+coeff_aa*dcc_dzeta
4861            d2aa_drs2=coeff_aa*d2cc_drs2
4862            d2aa_drsdzeta=-two*daa_drs*phi_logder+coeff_aa*d2cc_drsdzeta
4863            d2aa_dzeta2=aa*(-two*phi_zeta_inv*phipp_zeta+six*phi_logder*phi_logder)+&
4864 &           coeff_aa*(-four*dcc_dzeta*phi_logder+d2cc_dzeta2)
4865 
4866 !          Introduce tt : do not assume that the spin-dependent gradients are collinear
4867            grrho2=grho2_updn(ipts,3)
4868            dtt_dg=two*rhotot_inv*rhotot_inv*rhotmot*coeff_tt
4869 !          Note that tt is (the t variable of PBE divided by phi) squared
4870            tt=half*grrho2*dtt_dg
4871 
4872 !          Get xx from aa and tt
4873            xx=aa*tt
4874            dxx_drs=daa_drs*tt
4875            dxx_dzeta=daa_dzeta*tt
4876            dxx_dtt=aa
4877            d2xx_drs2=d2aa_drs2*tt
4878            d2xx_drsdzeta=d2aa_drsdzeta*tt
4879            d2xx_drsdtt=daa_drs
4880            d2xx_dttdzeta=daa_dzeta
4881            d2xx_dzeta2=d2aa_dzeta2*tt
4882 
4883 !          From xx to pade
4884            pade_den=one/(one+xx*(one+xx))
4885            pade=(one+xx)*pade_den
4886            dpade_dxx=-xx*(two+xx)*pade_den**2
4887            dpade_drs=dpade_dxx*dxx_drs
4888            dpade_dtt=dpade_dxx*dxx_dtt
4889            dpade_dzeta=dpade_dxx*dxx_dzeta
4890            d2pade_dxx2=two*(-one+xx*xx*(three+xx))*pade_den*pade_den*pade_den
4891            d2pade_drs2=d2pade_dxx2*dxx_drs*dxx_drs+dpade_dxx*d2xx_drs2
4892            d2pade_drsdtt=d2pade_dxx2*dxx_drs*dxx_dtt+dpade_dxx*d2xx_drsdtt
4893            d2pade_drsdzeta=d2pade_dxx2*dxx_drs*dxx_dzeta+dpade_dxx*d2xx_drsdzeta
4894            d2pade_dtt2=d2pade_dxx2*dxx_dtt*dxx_dtt
4895            d2pade_dttdzeta=d2pade_dxx2*dxx_dtt*dxx_dzeta+dpade_dxx*d2xx_dttdzeta
4896            d2pade_dzeta2=d2pade_dxx2*dxx_dzeta*dxx_dzeta+dpade_dxx*d2xx_dzeta2
4897 
4898 !          From pade to qq
4899            coeff_qq=tt*phi_zeta_inv*phi_zeta_inv
4900            qq=coeff_qq*pade
4901            dqq_drs=coeff_qq*dpade_drs
4902            dqq_dtt=pade*phi_zeta_inv*phi_zeta_inv+coeff_qq*dpade_dtt
4903            dqq_dzeta=coeff_qq*(dpade_dzeta-two*pade*phi_logder)
4904            d2qq_drs2=coeff_qq*d2pade_drs2
4905            d2qq_drsdtt=phi_zeta_inv*phi_zeta_inv*(dpade_drs+tt*d2pade_drsdtt)
4906            d2qq_drsdzeta=coeff_qq*(d2pade_drsdzeta-two*dpade_drs*phi_logder)
4907            d2qq_dtt2=phi_zeta_inv*phi_zeta_inv*(two*dpade_dtt+tt*d2pade_dtt2)
4908            d2qq_dttdzeta=phi_zeta_inv*phi_zeta_inv*(dpade_dzeta-two*pade*phi_logder)+&
4909 &           coeff_qq*(d2pade_dttdzeta-two*dpade_dtt*phi_logder)
4910            d2qq_dzeta2=coeff_qq*( d2pade_dzeta2-four*dpade_dzeta*phi_logder &
4911 &           +six*pade*phi_logder*phi_logder            &
4912 &           -two*pade*phi_zeta_inv*phipp_zeta)
4913 
4914 !          From qq to rr
4915            arg_rr=one+beta*gamma_inv*qq
4916            div_rr=one/arg_rr
4917            rr=gamma*log(arg_rr)
4918            drr_dqq=beta*div_rr
4919            drr_drs=drr_dqq*dqq_drs
4920            drr_dtt=drr_dqq*dqq_dtt
4921            drr_dzeta=drr_dqq*dqq_dzeta
4922            d2rr_dqq2=-div_rr**2*beta*beta*gamma_inv
4923            d2rr_drs2=d2rr_dqq2*dqq_drs*dqq_drs+drr_dqq*d2qq_drs2
4924            d2rr_drsdtt=d2rr_dqq2*dqq_drs*dqq_dtt+drr_dqq*d2qq_drsdtt
4925            d2rr_drsdzeta=d2rr_dqq2*dqq_drs*dqq_dzeta+drr_dqq*d2qq_drsdzeta
4926            d2rr_dtt2=d2rr_dqq2*dqq_dtt*dqq_dtt+drr_dqq*d2qq_dtt2
4927            d2rr_dttdzeta=d2rr_dqq2*dqq_dtt*dqq_dzeta+drr_dqq*d2qq_dttdzeta
4928            d2rr_dzeta2=d2rr_dqq2*dqq_dzeta*dqq_dzeta+drr_dqq*d2qq_dzeta2
4929 
4930 !          From rr to hh
4931            hh=phi3_zeta*rr
4932            dhh_drs=phi3_zeta*drr_drs
4933            dhh_dtt=phi3_zeta*drr_dtt
4934            dhh_dzeta=phi3_zeta*(drr_dzeta+three*rr*phi_logder)
4935            d2hh_drs2=phi3_zeta*d2rr_drs2
4936            d2hh_drsdtt=phi3_zeta*d2rr_drsdtt
4937            d2hh_drsdzeta=phi3_zeta*(d2rr_drsdzeta+three*drr_drs*phi_logder)
4938            d2hh_dtt2=phi3_zeta*d2rr_dtt2
4939            d2hh_dttdzeta=phi3_zeta*(d2rr_dttdzeta+three*drr_dtt*phi_logder)
4940            d2hh_dzeta2=phi3_zeta*(six*rr*phi_logder*phi_logder+&
4941 &           six*phi_logder*drr_dzeta+d2rr_dzeta2)  &
4942 &           +three*phi_zeta*phi_zeta*rr*phipp_zeta
4943 
4944 !          The GGA correlation energy is added
4945            exci(ipts)=exci(ipts)+hh
4946 
4947 !          Change of variables : from (rs,zeta,tt) to (rhoup,rhodn,grrho)
4948 
4949 
4950 !          From hh to the derivative of the energy wrt the density
4951            drhohh_drho=hh-third*rs*dhh_drs-zeta*dhh_dzeta-seven*third*tt*dhh_dtt
4952            vxci(ipts,1)=vxci(ipts,1)+drhohh_drho+dhh_dzeta
4953            vxci(ipts,2)=vxci(ipts,2)+drhohh_drho-dhh_dzeta
4954 
4955 
4956 !          From hh to the derivative of the energy wrt to the gradient of the
4957 !          density, divided by the gradient of the density
4958 !          (The v3.3 definition includes the division by the norm of the gradient)
4959            dvxcdgr(ipts,3)=rhotot*dtt_dg*dhh_dtt
4960 
4961            d2rhohh_drho2=rhotot_inv*&
4962 &           (-two*ninth*rs*dhh_drs +seven*four*ninth*tt*dhh_dtt &
4963 &           +ninth*rs*rs*d2hh_drs2+zeta*zeta*d2hh_dzeta2+(seven*third*tt)**2*d2hh_dtt2 &
4964 &           +two*third*rs*zeta*d2hh_drsdzeta+two*seven*ninth*rs*tt*d2hh_drsdtt &
4965 &           +two*seven*third*tt*zeta*d2hh_dttdzeta)
4966            d2rhohh_drhodg=dtt_dg*(-four*third*dhh_dtt-third*rs*d2hh_drsdtt &
4967 &           -zeta*d2hh_dttdzeta-seven*third*tt*d2hh_dtt2)
4968 
4969 !          Component 12 : first derivative with respect to the gradient
4970 !          of the density, div by the grad of the density
4971            dvxci(ipts,12)=dvxcdgr(ipts,3)
4972 !          Components 9, 10 and 11 : second derivatives with respect to the spin-density
4973 !          Note that there is already a contribution from LSDA
4974            dvxci(ipts,9)=dvxci(ipts,9)+d2rhohh_drho2+rhotot_inv*           &
4975 &           ( d2hh_dzeta2*(one-two*zeta) &
4976 &           -two*third*rs*d2hh_drsdzeta-14.0_dp*third*tt*d2hh_dttdzeta)
4977            dvxci(ipts,10)=dvxci(ipts,10)+d2rhohh_drho2-rhotot_inv*d2hh_dzeta2
4978            dvxci(ipts,11)=dvxci(ipts,11)+d2rhohh_drho2+rhotot_inv*           &
4979 &           ( d2hh_dzeta2*(one+two*zeta) &
4980 &           +two*third*rs*d2hh_drsdzeta+14.0_dp*third*tt*d2hh_dttdzeta)
4981 !          Components 13 and 14 : second derivatives with respect to spin density
4982 !          and gradient, divided by the gradient
4983            dvxci(ipts,13)=d2rhohh_drhodg+dtt_dg*d2hh_dttdzeta
4984            dvxci(ipts,14)=d2rhohh_drhodg-dtt_dg*d2hh_dttdzeta
4985 !          Component 15 : derivative of the (derivative wrt the gradient div by the grad),
4986 !          divided by the grad
4987            dvxci(ipts,15)=rhotot*d2hh_dtt2*dtt_dg*dtt_dg
4988 
4989 !          End condition of GGA
4990          end if
4991 
4992        else  ! no correlation
4993 
4994          if(ndvxci > 8)then
4995 !          Must zero the correlation part of the xc kernel
4996            dvxci(:,9:15)=zero
4997          end if
4998 
4999 !        End condition of including correlation, and not only exchange
5000        end if
5001 
5002 !      Correlation has been added
5003 !      -----------------------------------------------------------------------------
5004 
5005      end do
5006    end if
5007 
5008 !  fab: this should be the else on nspden
5009  else
5010 !  Disallowed value for nspden
5011    write(message, '(a,a,a,i12,a)' )&
5012 &   '  Argument nspden must be 1 or 2; ',ch10,&
5013 &   '  Value provided as argument was ',nspden,'.'
5014    ABI_BUG(message)
5015  end if
5016 
5017 !DEBUG
5018 !Finite-difference debugging, do not take away
5019 !if(debug/=0)then
5020 !do ipts=1,5,5
5021 
5022 !rho=rho_updn(ipts,1)+rho_updn(ipts,2)
5023 !write(std_out,'(a,i5,a,2es16.8)' ) ' Point number',ipts,' with rhoup,rhodn=',rho_updn(ipts,1),rho_updn(ipts,2)
5024 !write(std_out,'(a,i5,a,2es16.8)' ) ' Point number',ipts+1,' with rhoup,rhodn=',rho_updn(ipts+1,1),rho_updn(ipts+1,2)
5025 !write(std_out,'(a,i5,a,2es16.8)' ) ' Point number',ipts+2,' with rhoup,rhodn=',rho_updn(ipts+2,1),rho_updn(ipts+2,2)
5026 !write(std_out,'(a,i5,a,2es16.8)' ) ' Point number',ipts+3,' with rhoup,rhodn=',rho_updn(ipts+3,1),rho_updn(ipts+3,2)
5027 !write(std_out,'(a,i5,a,2es16.8)' ) ' Point number',ipts+4,' with rhoup,rhodn=',rho_updn(ipts+4,1),rho_updn(ipts+4,2)
5028 
5029 !! For rho
5030 !if(debug==1)then
5031 !write(std_out,'(a)' )' Direct values :'
5032 !write(std_out,'(3es16.8)' )exci(ipts)*rho,vxci(ipts,1),vxci(ipts,2)
5033 !else
5034 !!  For grho2
5035 !write(std_out,'(4es16.8)' )exci(ipts)*rho,dvxcdgr(ipts,1),&
5036 !&  dvxcdgr(ipts,2),dvxcdgr(ipts,3)
5037 !end if
5038 
5039 !write(std_out,'(4es16.8)' )dvxci(ipts,1:3)  ! For LDA
5040 !write(std_out,'(a)' )'     3rd-order '
5041 !write(std_out,'(4es16.8)' )d2vxci(ipts,1:4)  ! For LDA
5042 
5043 !write(std_out,'(4es16.8)' )dvxci(ipts,1:4)  ! For exchange
5044 !write(std_out,'(4es16.8)' )dvxci(ipts,5:8)  ! For exchange
5045 !write(std_out,'(4es16.8)' )dvxci(ipts,9:12) ! For correlation
5046 !write(std_out,'(4es16.8)' )dvxci(ipts,13:15) ! For correlation
5047 
5048 !if(debug==1)then
5049 !!  For rho
5050 !write(std_out,'(a)' )' Finite-difference values :'
5051 !write(std_out,'(3es16.8)' )exci(ipts)*rho,&
5052 !&      ( exci(ipts+1)*(rho+delta) - exci(ipts+2)*(rho-delta) )/2._dp/delta,&
5053 !&      ( exci(ipts+3)*(rho+delta) - exci(ipts+4)*(rho-delta) )/2._dp/delta
5054 !write(std_out,'(3es16.8)' )&
5055 !&    ( vxci(ipts+1,1) - vxci(ipts+2,1) )/2._dp/delta,&
5056 !&    ( vxci(ipts+3,1) - vxci(ipts+4,1) )/2._dp/delta,&
5057 !&    ( vxci(ipts+3,2) - vxci(ipts+4,2) )/2._dp/delta
5058 !!This is for order 3
5059 !write(std_out,'(a)' )'     3rd-order by two methods, giving components 1, 2, 3, then on the next line 2, 3, 4 '
5060 !write(std_out,'(3es16.8)' )&
5061 !&    ( dvxci(ipts+1,1) - dvxci(ipts+2,1) )/2._dp/delta,&
5062 !&    ( dvxci(ipts+1,2) - dvxci(ipts+2,2) )/2._dp/delta,&
5063 !&    ( dvxci(ipts+1,3) - dvxci(ipts+2,3) )/2._dp/delta
5064 !write(std_out,'(3es16.8)' )&
5065 !&    ( dvxci(ipts+3,1) - dvxci(ipts+4,1) )/2._dp/delta,&
5066 !&    ( dvxci(ipts+3,2) - dvxci(ipts+4,2) )/2._dp/delta,&
5067 !&    ( dvxci(ipts+3,3) - dvxci(ipts+4,3) )/2._dp/delta
5068 
5069 !write(std_out,*)
5070 !write(std_out,*)' Now for ecrs and derivatives '
5071 !write(std_out,'(a)' )' ecrs, rs, zeta ='
5072 !write(std_out,'(3es16.8)' )wecrsz(ipts,1:3)
5073 !write(std_out,'(3es16.8)' )wecrsz(ipts+1,1:3)
5074 !write(std_out,'(3es16.8)' )wecrsz(ipts+2,1:3)
5075 !write(std_out,'(3es16.8)' )wecrsz(ipts+3,1:3)
5076 !write(std_out,'(3es16.8)' )wecrsz(ipts+4,1:3)
5077 !write(std_out,'(a)' )' ecrs derivatives :'
5078 !write(std_out,'(3es16.8)' )d1wecrsz(ipts,1:2)
5079 !write(std_out,'(3es16.8)' )d2wecrsz(ipts,1:3)
5080 !write(std_out,'(4es16.8)' )d3wecrsz(ipts,1:4)
5081 !write(std_out,'(a)' )' Finite-differences :'
5082 !write(std_out,'(3es16.8)' )&
5083 !&    ( wecrsz(ipts+1,1) - wecrsz(ipts+2,1) )/( wecrsz(ipts+1,2) - wecrsz(ipts+2,2) ),&
5084 !&    ( wecrsz(ipts+3,1) - wecrsz(ipts+4,1) )/( wecrsz(ipts+3,3) - wecrsz(ipts+4,3) )
5085 !write(std_out,'(3es16.8)' )&
5086 !&    ( d1wecrsz(ipts+1,1) - d1wecrsz(ipts+2,1) )/( wecrsz(ipts+1,2) - wecrsz(ipts+2,2) ),&
5087 !&    ( d1wecrsz(ipts+1,2) - d1wecrsz(ipts+2,2) )/( wecrsz(ipts+1,2) - wecrsz(ipts+2,2) )
5088 !write(std_out,'(3es16.8)' )&
5089 !&    ( d1wecrsz(ipts+3,1) - d1wecrsz(ipts+4,1) )/( wecrsz(ipts+3,3) - wecrsz(ipts+4,3) ),&
5090 !&    ( d1wecrsz(ipts+3,2) - d1wecrsz(ipts+4,2) )/( wecrsz(ipts+3,3) - wecrsz(ipts+4,3) )
5091 !write(std_out,'(a)' )' Finite-differences, 3rd order :'
5092 !write(std_out,'(3es16.8)' )&
5093 !&    ( d2wecrsz(ipts+1,1) - d2wecrsz(ipts+2,1) )/( wecrsz(ipts+1,2) - wecrsz(ipts+2,2) ),&
5094 !&    ( d2wecrsz(ipts+1,2) - d2wecrsz(ipts+2,2) )/( wecrsz(ipts+1,2) - wecrsz(ipts+2,2) ),&
5095 !&    ( d2wecrsz(ipts+1,3) - d2wecrsz(ipts+2,3) )/( wecrsz(ipts+1,2) - wecrsz(ipts+2,2) )
5096 !write(std_out,'(3es16.8)' )&
5097 !&    ( d2wecrsz(ipts+3,1) - d2wecrsz(ipts+4,1) )/( wecrsz(ipts+3,3) - wecrsz(ipts+4,3) ),&
5098 !&    ( d2wecrsz(ipts+3,2) - d2wecrsz(ipts+4,2) )/( wecrsz(ipts+3,3) - wecrsz(ipts+4,3) ),&
5099 !&    ( d2wecrsz(ipts+3,3) - d2wecrsz(ipts+4,3) )/( wecrsz(ipts+3,3) - wecrsz(ipts+4,3) )
5100 
5101 
5102 !!This is for GGA
5103 
5104 !write(std_out,'(4es16.8)' )&
5105 !&    ( dvxcdgr(ipts+1,1) - dvxcdgr(ipts+2,1) )/2._dp/delta,&
5106 !&    ( dvxcdgr(ipts+3,2) - dvxcdgr(ipts+4,2) )/2._dp/delta,&
5107 !&    ( dvxcdgr(ipts+1,3) - dvxcdgr(ipts+2,3) )/2._dp/delta,&
5108 !&    ( dvxcdgr(ipts+3,3) - dvxcdgr(ipts+4,3) )/2._dp/delta
5109 !else
5110 !!  For grho2  (should distinguish exchange and correlation ...)
5111 !grr=sqrt(grho2_updn(ipts,1)) ! Analysis of exchange
5112 !grr=sqrt(grho2_updn(ipts,3)) ! Analysis of correlation
5113 !write(std_out,'(3es16.8)' )exci(ipts)*rho,&
5114 !&      ( exci(ipts+1)*rho - exci(ipts+2)*rho )/2._dp/delta/grr,&
5115 !&      ( exci(ipts+3)*rho - exci(ipts+4)*rho )/2._dp/delta/grr
5116 !write(std_out,'(3es16.8)' )&
5117 !&    ( vxci(ipts+1,1) - vxci(ipts+2,1) )/2._dp/delta/grr,&
5118 !&    ( vxci(ipts+3,1) - vxci(ipts+4,1) )/2._dp/delta/grr,&
5119 !&    ( vxci(ipts+3,2) - vxci(ipts+4,2) )/2._dp/delta/grr
5120 !write(std_out,'(4es16.8)' )&
5121 !&    ( dvxcdgr(ipts+1,1) - dvxcdgr(ipts+2,1) )/2._dp/delta/grr,&
5122 !&    ( dvxcdgr(ipts+3,2) - dvxcdgr(ipts+4,2) )/2._dp/delta/grr,&
5123 !&    ( dvxcdgr(ipts+1,3) - dvxcdgr(ipts+2,3) )/2._dp/delta/grr,&
5124 !&    ( dvxcdgr(ipts+3,3) - dvxcdgr(ipts+4,3) )/2._dp/delta/grr
5125 !end if
5126 !end do
5127 !stop
5128 !end if
5129 !ENDDEBUG
5130 
5131  ABI_FREE(rhoarr)
5132  ABI_FREE(rhom1_3)
5133  ABI_FREE(rho_updnm1_3)
5134  ABI_FREE(zetm)
5135  ABI_FREE(zetmm1_3)
5136  ABI_FREE(zetp)
5137  ABI_FREE(zetpm1_3)
5138 
5139 !DEBUG
5140 !deallocate(wecrsz,d1wecrsz,d2wecrsz,d3wecrsz)
5141 !ENDDEBUG
5142 
5143 !DEBUG
5144 !write(std_out,*)' xcpbe : exit'
5145 !write(std_out,*)' nspden=',nspden
5146 !if(order==2)stop
5147 !ENDDEBUG
5148 
5149 end subroutine xcpbe