TABLE OF CONTENTS


ABINIT/hybrid_corr [ Functions ]

[ Top ] [ Functions ]

NAME

 hybrid_corr

FUNCTION

 Compute the correction to the XC potential and energy due to the hybridation

COPYRIGHT

  Copyright (C) 2015-2017 ABINIT group (FA)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

  dtset <type(dataset_type)>=all input variables in this dataset
  ixc = choice of exchange-correlation functional.
  mpi_enreg=information about MPI parallelization
  nfft = number of fft grid points.
  ngfft(1:3) = integer fft box dimensions, see getng for ngfft(4:8).
  nspden = number of spin-density components.
  rhor(nfft,nspden) = electron density in real space in electrons/bohr**3
   (total in first half and spin-up in second half if nspden = 2).
  rprimd(3,3) = dimensional primitive translations for real space in Bohr.

OUTPUT

SIDE EFFECTS

  enxc = exchange correlation energy corrected for hybrid functionals
  vxc = exchange correlation potential corrected for hybrid functionals

WARNINGS

 Current restrictions are:
  a - Spin-polarized case not tested.

PARENTS

CHILDREN

      dtset_copy,dtset_free,libxc_functionals_end,libxc_functionals_init
      rhohxc

SOURCE

 44 #if defined HAVE_CONFIG_H
 45 #include "config.h"
 46 #endif
 47 
 48 #include "abi_common.h"
 49 
 50 subroutine hybrid_corr(dtset,ixc,nkxc,mpi_enreg,nfft,ngfft,nspden,rhor,rprimd,hybrid_mixing,vxc,enxc)
 51 
 52  use defs_basis
 53  use m_profiling_abi
 54  use m_errors
 55  use defs_abitypes, only : MPI_type, dataset_type
 56  use m_dtset,       only : dtset_copy, dtset_free
 57 
 58 !This section has been created automatically by the script Abilint (TD).
 59 !Do not modify the following lines by hand.
 60 #undef ABI_FUNC
 61 #define ABI_FUNC 'hybrid_corr'
 62 !End of the abilint section
 63 
 64  implicit none
 65 
 66 !Arguments -------------------------------------------------------------
 67 !scalars
 68  integer,intent(in) :: ixc,nfft,nspden,nkxc
 69  real(dp),intent(in) :: hybrid_mixing
 70  real(dp),intent(inout) :: enxc
 71  type(MPI_type),intent(in) :: mpi_enreg
 72  type(dataset_type),intent(in) :: dtset
 73 !arrays
 74  integer,intent(in) :: ngfft(18)
 75  real(dp),intent(in) :: rhor(nfft,nspden),rprimd(3,3)
 76  real(dp),intent(inout) :: vxc(nfft,nspden)
 77 
 78 !Local variables -------------------------------------------------------
 79 !No improved xc quadrature.
 80 !No core correction.
 81 !Dummy here.
 82 !For debugging purposes (see tests below):
 83 !integer :: i1,i2,i3,k1,n1,n2,n3
 84 !real(dp) :: kx,rho,rhomax,ftest
 85 !scalars
 86  character(len=500) :: message
 87  type(dataset_type) :: dtLocal
 88 !arrays
 89  real(dp) :: dum(0)
 90  real(dp) :: enxc_corr
 91  real(dp),allocatable :: kxcr(:,:)
 92  real(dp),allocatable :: xccc3d(:),vxc_corr(:,:)
 93 
 94 !***********************************************************************
 95 !For debugging purposes (see tests below):
 96 
 97 !ftest(i1,n1,k1) = 0._dp+1._dp*cos(k1*two_pi*float(i1)/float(n1))
 98 
 99 !***********************************************************************
100 
101 !Check input parameters.
102 
103  if (nspden > 2) then
104    message = ' kxc_alda does not work yet for nspden > 2.'
105    MSG_ERROR(message)
106  end if
107 
108 !Copy the input variables from the current dataset to a temporary one
109 !to tune some parameters
110  call dtset_copy(dtLocal, dtset)
111  dtLocal%intxc = 0
112  dtLocal%ixc   = -101
113 
114  ! Reinitialize the libxc module with the overriden values
115  if (dtset%ixc<0) then
116    call libxc_functionals_end()
117  end if
118  if (dtLocal%ixc<0) then
119    call libxc_functionals_init(dtLocal%ixc,dtLocal%nspden)
120  end if
121 
122  call rhohxc(dtLocal,enxc_corr,dum,0,kxcr,mpi_enreg,nfft,dum,dum,0,dum,0,nkxc,0,&
123 & dtset%usepawu==4.or.dtset%usepawu==14,nspden,n3xccc,&
124 & 0,dum,rhor,rprimd,strsxc,1,dum,vxc_corr,dum,xccc3d)
125 
126  vxc(:,:) = vxc(:,:) + hybrid_mixing*vxc_corr(:,:)
127  enxc = enxc + hybrid_mixing*enxc_corr
128 
129  call rhohxc(dtLocal,enxc_corr,dum,0,kxcr,mpi_enreg,dum,dum,dum,0,dum,0,nkxc,0,&
130 & dtset%usepawu==4.or.dtset%usepawu==14,nspden,0,&
131 & 0,dum,rhor,rprimd,strsxc,1,dum,vxc_corr,vxcavg,0)
132 
133  vxc(:,:) = vxc(:,:) - hybrid_mixing*vxc_corr(:,:)
134  enxc = enxc - hybrid_mixing*enxc_corr
135 
136 
137 ! Revert libxc module to the original settings
138  if (dtLocal%ixc<0) then
139    call libxc_functionals_end()
140  end if
141  if (dtset%ixc<0) then
142    call libxc_functionals_init(dtset%ixc,dtset%nspden)
143  end if
144 
145 !Free memory.
146  call dtset_free(dtLocal)
147  ABI_FREE(kxcr)
148  ABI_FREE(xccc3d)
149 
150 end subroutine hybrid_corr