TABLE OF CONTENTS


ABINIT/m_xcdata [ Modules ]

[ Top ] [ Modules ]

NAME

  m_xcdata

FUNCTION

  This module provides the definition of
  the xcdata_type used to drive the computation of the XC energy, potential, kernel, etc.

COPYRIGHT

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

NOTES

SOURCE

19 #if defined HAVE_CONFIG_H
20 #include "config.h"
21 #endif
22 
23 #include "abi_common.h"
24 
25 module m_xcdata
26 
27  use defs_basis
28  use m_errors
29  use libxc_functionals
30  use m_dtset, only : dataset_type
31  use m_drivexc, only : size_dvxc
32 
33  implicit none
34 
35  private

m_xcdata/get_auxc_ixc [ Functions ]

[ Top ] [ m_xcdata ] [ Functions ]

NAME

  get_auxc_ixc

FUNCTION

  Returns the ixc of an auxiliary XC functional to be used instead of the input ixc
  For most of the functionals, there is no need of an auxiliary functional, in which case auxc_ixc=0
  For hybrid functionals, on the contrary, some speedup can be achieved by using such an auxiliary functional
  Note that this XC functional intend to replace the whole ixc functional. Generally speakin, it should be
  mistaken for the GGA part of the hybrid functional (that is for exchange only, actually).

  At present, always return ixc=1, but this might change in the future ...

INPUTS

  ixc= index of exchange-correlation functional

OUTPUT

  auxc_ixc= 0 if no need of an auxiliary functional, otherwise, returns the ixc of an auxiliary functional.

SIDE EFFECTS

SOURCE

305 subroutine get_auxc_ixc(auxc_ixc,ixc)
306 
307 !Arguments ------------------------------------
308 !scalars
309  integer, intent(in) :: ixc
310  integer, intent(out) :: auxc_ixc
311 
312 !Local variables-------------------------------
313  integer :: usefock,xclevel
314 !integer :: gga_id(2)
315 
316 ! *************************************************************************
317 
318  auxc_ixc=11
319 !Native hybrid functionals from ABINIT
320  if (ixc==40.or.ixc==41.or.ixc==42) then
321    auxc_ixc = 11
322 !Hybrid functionals from libxc
323  else if (ixc<0) then
324    call get_xclevel(ixc,xclevel,usefock)
325    if(usefock==1)then
326      auxc_ixc=11
327 !    if (libxc_functionals_gga_from_hybrid(hybrid_id=ixc,gga_id=gga_id)) then
328 !      auxc_ixc=-gga_id(1)*1000-gga_id(2)
329 !    endif
330    end if
331  end if
332 
333 end subroutine get_auxc_ixc
334 
335 end module m_xcdata

m_xcdata/get_xclevel [ Functions ]

[ Top ] [ m_xcdata ] [ Functions ]

NAME

  get_xclevel

FUNCTION

  Compute xclevel.

INPUTS

  ixc= index of exchange-correlation functional

OUTPUT

  [usefock = 1 if the XC functional needs the Fock operator]
  xclevel= 0 if no XC functional except possibly Fock; 1 if LDA; 2 if GGA ; 3 for TDDFT kernel tests

SIDE EFFECTS

SOURCE

228 subroutine get_xclevel(ixc,xclevel,usefock)
229 
230 !Arguments ------------------------------------
231 !scalars
232  integer, intent(in) :: ixc
233  integer, intent(out) :: xclevel
234  integer, intent(out), optional :: usefock
235 
236 !Local variables-------------------------------
237  integer :: ii,isiz,jj
238  character(len=500) :: message
239 
240 ! *************************************************************************
241 
242  xclevel=0 ; if(present(usefock)) usefock=0
243  if( ( 1<=ixc .and. ixc<=10).or.(30<=ixc .and. ixc<=39).or.(ixc==50) )xclevel=1 ! LDA
244  if( (11<=ixc .and. ixc<=19).or.(23<=ixc .and. ixc<=29).or. ixc==1402000)xclevel=2 ! GGA
245  if( 20<=ixc .and. ixc<=22 )xclevel=3 ! ixc for TDDFT kernel tests
246  if(present(usefock))then
247    if( ixc>=40 .and. ixc<=42 )usefock=1 ! Hartree-Fock or internal hybrid functionals
248  endif
249  if( ixc>=31 .and. ixc<=35)xclevel=2 ! ixc for internal fake mGGA
250  if( ixc>=41 .and. ixc<=42)xclevel=2 ! ixc for internal hybrids using GGA
251  if (ixc<0) then                     ! libXC: metaGGA and hybrid functionals
252    xclevel=1
253    do isiz=1,2
254 !    ixc has ABINIT sign convention
255 !    ii has Libxc sign convention
256      if (isiz==1) ii=-ixc/1000
257      if (isiz==2) ii=-ixc-ii*1000
258      if (ii<=0) cycle
259      jj=libxc_functionals_family_from_id(ii)
260      if (jj==XC_FAMILY_GGA    .or.jj==XC_FAMILY_MGGA) xclevel=2
261      if (jj==XC_FAMILY_HYB_GGA.or.jj==XC_FAMILY_HYB_MGGA) xclevel=2
262      if (present(usefock)) then
263        if (libxc_functionals_is_hybrid_from_id(ii)) usefock=1
264        if (usefock==1) then
265          if (.not.libxc_functionals_gga_from_hybrid(hybrid_id=ii)) then
266            write(message, '(a,i8,3a,2i8,2a)' )&
267 &           'ixc=',ixc,' (libXC hybrid functional) is presently not allowed.',ch10,&
268 &           'ii,jj=',ii,jj,ch10,&
269 &           'Action: try another hybrid functional.'
270            ABI_ERROR(message)
271          end if
272        end if
273      end if
274    end do
275  end if
276 
277 end subroutine get_xclevel

m_xcdata/xcdata_init [ Functions ]

[ Top ] [ m_xcdata ] [ Functions ]

NAME

  xcdata_init

FUNCTION

  Init the structure. Mostly copy input variables, except compute and usefock and xclevel.

INPUTS

  [dtset = the dataset from which the other input variables are taken, if they are not present]
  [auxc_ixc = possibly the index of the auxiliary xc functional, otherwise 0.]
  [hyb_mixing = parameter for mixing Fock exchange in native PBEx functionals]
  [intxc = 1 if the XC functional has to be interpolated on a more refined mesh than the FFT one]
  [ixc= index of exchange-correlation functional]
  [nelect = Number of electrons in the cell (for Fermi-Amaldi only)]
  [tphysel = Physical temperature (for temperature-dependent functional)]
  [vdw_xc = Choice of van-der-Waals density functional]
  [xc_denpos = density positivity value]
  [xc_taupos = kinetic energy density positivity value (mGGA)]

OUTPUT

  xcdata <type(xcdata_type)>= the data to calculate exchange-correlation are initialized

SIDE EFFECTS

SOURCE

145 subroutine xcdata_init(xcdata,auxc_ixc,dtset,hyb_mixing,intxc,ixc,nelect,nspden,tphysel,&
146 &                      vdw_xc,xc_denpos,xc_taupos)
147 
148 !Arguments ------------------------------------
149 !scalars
150  integer, intent(in),optional :: auxc_ixc,intxc,ixc,nspden,vdw_xc
151  real(dp),intent(in),optional :: hyb_mixing,nelect,tphysel,xc_denpos,xc_taupos
152  type(dataset_type), intent(in),optional :: dtset
153  type(xcdata_type), intent(out) :: xcdata
154 !Local variables-------------------------------
155  integer :: nspden_updn
156  character(len=500) :: msg
157 
158 ! *************************************************************************
159 
160  if(present(dtset))then
161    xcdata%auxc_ixc=dtset%auxc_ixc
162    xcdata%intxc=dtset%intxc
163    xcdata%ixc=dtset%ixc
164    xcdata%nspden=dtset%nspden
165    xcdata%vdw_xc=dtset%vdw_xc
166 
167    xcdata%hyb_mixing=abs(dtset%hyb_mixing) ! Warning : the absolute value is needed, because of the singular way
168                                            ! to define the default values for this input variable.
169    xcdata%nelect=dtset%nelect
170    xcdata%tphysel=merge(dtset%tphysel,dtset%tsmear,dtset%tphysel>tol8.and.dtset%occopt/=3.and.dtset%occopt/=9)
171 
172    xcdata%xc_denpos=dtset%xc_denpos
173    xcdata%xc_taupos=dtset%xc_taupos
174 
175  else
176    if(.not.(present(auxc_ixc).and.present(intxc).and.present(ixc).and.&
177 &           present(vdw_xc).and.present(hyb_mixing).and.&
178 &           present(nelect).and.present(nspden).and.&
179 &           present(tphysel).and.present(xc_denpos).and.present(xc_taupos)))then
180      msg='If dtset is not provided, all the other optional arguments must be provided, which is not the case!'
181      ABI_BUG(msg)
182    endif
183  endif
184 
185  if(present(auxc_ixc))  xcdata%auxc_ixc=auxc_ixc
186  if(present(intxc))     xcdata%intxc=intxc
187  if(present(ixc))       xcdata%ixc=ixc
188  if(present(nspden))    xcdata%nspden=nspden
189  if(present(vdw_xc))    xcdata%vdw_xc=vdw_xc
190 
191  if(present(hyb_mixing))xcdata%hyb_mixing=hyb_mixing
192  if(present(nelect))    xcdata%nelect=nelect
193  if(present(tphysel))   xcdata%tphysel=tphysel
194  if(present(xc_denpos)) xcdata%xc_denpos=xc_denpos
195  if(present(xc_taupos)) xcdata%xc_taupos=xc_taupos
196 
197 !Compute xclevel
198  call get_xclevel(xcdata%ixc,xcdata%xclevel,usefock=xcdata%usefock)
199 
200 !Compute usegradient,uselaplacian,usekden
201  nspden_updn=min(xcdata%nspden,2)
202  call size_dvxc(xcdata%ixc,1,nspden_updn,usegradient=xcdata%usegradient,&
203 &               uselaplacian=xcdata%uselaplacian,usekden=xcdata%usekden)
204 
205 end subroutine xcdata_init

m_xcdata/xcdata_type [ Types ]

[ Top ] [ m_xcdata ] [ Types ]

NAME

  xcdata_type

FUNCTION

   This object stores the input variables (and derived parameters) needed to compute the exchange-correlation functional,
   not simply to define it.

NOTES

SOURCE

 51  type, public :: xcdata_type
 52 
 53 ! Integer scalars
 54 
 55   integer :: auxc_ixc
 56     ! Choice of auxiliary exchange-correlation functional. See input variable documentation
 57     ! If 0, there is no auxiliary xc functional, the one corresponding to ixc has to be used.
 58 
 59   integer :: intxc
 60     ! 1 if the XC functional has to be interpolated on a more refined mesh than the FFT one.
 61     ! 0 stick to the original FFT mesh
 62 
 63   integer :: ixc
 64     ! Choice of exchange-correlation functional. See input variable documentation
 65 
 66   integer :: nspden
 67     ! Number of spin components of the density
 68 
 69   integer :: usefock
 70     ! 1 if the XC functional includes a (possibly screened) Fock contribution
 71 
 72   integer :: usegradient
 73     ! 1 if the XC functional depends on the density gradient
 74 
 75   integer :: uselaplacian
 76     ! 1 if the XC functional depends on the density laplacian
 77 
 78   integer :: usekden
 79     ! 1 if the XC functional depends on the kinetic energy density
 80 
 81   integer :: vdw_xc
 82     ! Choice of van-der-Waals density functional. See input variable documentation
 83 
 84   integer :: xclevel
 85     ! Determined from ixc
 86     ! 0 if no XC functional
 87     ! 1 if LDA-type XC functional
 88     ! 2 if GGA-type XC functional
 89     ! 3 if for TDDFT kernel
 90 
 91 ! Real scalars
 92 
 93   real(dp) :: hyb_mixing
 94     ! Parameter for mixing Fock exchange in native PBEx functionals
 95 
 96   real(dp) :: nelect
 97     ! Number of electrons in the cell (for Fermi-Amaldi only)
 98 
 99   real(dp) :: tphysel
100     ! Physical temperature (for temperature-dependent functional)
101 
102   real(dp) :: xc_denpos
103     ! density positivity value
104 
105   real(dp) :: xc_taupos
106     ! kinetic energy density positivity value (mGGA)
107 
108  end type xcdata_type
109 
110 !----------------------------------------------------------------------
111 
112  public :: xcdata_init                ! Initialize the object.
113  public :: get_xclevel                ! Get the xclevel from ixc (as well as usefock)
114  public :: get_auxc_ixc               ! Get the auxiliary xc functional (if it exists)
115 
116 contains