Libxc:manual

From OctopusWiki

Jump to: navigation, search

Contents

An example

Probably the best to explain the usage of libxc is through an example. This snippet is taken from the file testsuite/xc-get_data.c

switch(xc_family_from_id(xc.functional))
  {
  case XC_FAMILY_LDA:
    if(xc.functional == XC_LDA_X)
      xc_lda_x_init(&lda_func, xc.nspin, 3, 0);
    else
      xc_lda_init(&lda_func, xc.functional, xc.nspin);
    xc_lda_vxc(&lda_func, xc.rho, &xc.zk, xc.vrho);
    break;
  case XC_FAMILY_GGA:
    xc_gga_init(&gga_func, xc.functional, xc.nspin);
    xc_gga(&gga_func, xc.rho, xc.sigma, &xc.zk, xc.vrho, xc.vsigma);
    xc_gga_end(&gga_func);
    break;
  default:
    fprintf(stderr, "Functional '%d' not found\n", xc.functional);
    exit(1);
  }

The functionals are divided in families (LDA, GGA, etc.). Given a functional identifier, xc.functional returns the corresponding family. Then the functional is initialized by one of the xc_XXX_init, and evaluated by xc_XXX, which returns the energy per unit particle (zk and the derivatives of the energy per unit volume. Finally the functions xc_XXX_end clean up.

This is an example in Fortran:

program lxctest
   use libxc

   implicit none

   real(8) :: rho, e_c, v_c

   TYPE(xc_func) :: xc_c_func
   TYPE(xc_info) :: xc_c_info

   CALL xc_f90_lda_init(xc_c_func, xc_c_info, XC_LDA_C_VWN, XC_UNPOLARIZED)
   CALL xc_f90_lda_vxc(xc_c_func, rho, e_c, v_c)
   CALL xc_f90_lda_end(xc_c_func)

end program lxctest

The info structure

For each functional there is a structure that holds several pieces of information. The relevant part of this structure for the end user is defined as

typedef struct{
  int   number;   /* indentifier number */
  int   kind;     /* XC_EXCHANGE or XC_CORRELATION */

  char *name;     /* name of the functional, e.g. "PBE" */
  int   family;   /* type of the functional, e.g. XC_FAMILY_GGA */
  char *refs;     /* references                       */

  int   provides; /* what the functional provides, e.g. XC_PROVIDES_EXC | XC_PROVIDES_VXC */

  ...
} xc_func_info_type;

For example, for the Slater exchange functional, this structrue is defined as

const xc_func_info_type func_info_lda_x = {
  XC_LDA_X,
  XC_EXCHANGE,
  "Slater exchange",
  XC_FAMILY_LDA,
  "P.A.M. Dirac, Proceedings of the Cambridge Philosophical Society 26, 376 (1930)\n"
  "F. Bloch, Zeitschrift fuer Physik 57, 545 (1929)",
  XC_PROVIDES_EXC | XC_PROVIDES_VXC | XC_PROVIDES_FXC,
  ...
};

Note that the references are separated by a newline. The user of the library can access the information in the following way:

  xc_gga_type b88;

  xc_gga_init(&b88, XC_GGA_X_B88, XC_UNPOLARIZED);
  printf("The functional '%s' is defined in the reference(s):\n%s",
    b88.info->name, b88.info->refs);
  xc_gga_end(&b88);

Varia

xc_family_from_id

If you have the identifier of a functional (that was read, for example from an input file), and want to find out to which family it belongs to, you can use the routine

int xc_family_from_id(int functional);
input:
  functional: identifier of the functional
returns:
  XC_FAMILY_UNKNOWN: could not find the family
  XC_FAMILY_LDA, XC_FAMILY_GGA, etc.

Fortran bindings

The Fortran 90 bindings can be found in the file libxc.F90 (see the download section). In general, calling libxc from Fortran is as simple as from C. There are a few exceptions though. The most notable is how to access the information about the functional, due to the way that Fortran handles strings. You can take a look at the following example (taken from octopus/src/xc_functl.F90) for help

i = xc_f90_info_kind(functl%info)
select case(i)
case(XC_EXCHANGE)
  write(message(1), '(2x,a)') 'Exchange'
case(XC_CORRELATION)
  write(message(1), '(2x,a)') 'Correlation'
case(XC_EXCHANGE_CORRELATION)
  write(message(1), '(2x,a)') 'Exchange-correlation'
end select

call xc_f90_info_name(functl%info, s1)
select case(functl%family)
  case (XC_FAMILY_LDA);  write(s2,'(a)') "LDA"
  case (XC_FAMILY_GGA);  write(s2,'(a)') "GGA"
  case (XC_FAMILY_MGGA); write(s2,'(a)') "MGGA"
  case (XC_FAMILY_LCA);  write(s2,'(a)') "LCA"
end select
write(message(2), '(4x,4a)') trim(s1), ' (', trim(s2), ')'
call write_info(2, iunit)
      
i = 1; str = 0
call xc_f90_info_ref(functl%info, str, s1)
do while(str >= 0)
  write(message(1), '(4x,a,i1,2a)') '[', i, '] ', trim(s1)
  call write_info(1, iunit)
  call xc_f90_info_ref(functl%info, str, s1)
  i = i + 1
end do

The LDAs

Functionals

These are the LDA functionals the library currently knows about. Exc, Vxc and Fxc implementation audit status.

Exchange
Correlation
Exchange-Correlation

Initialization

A generic LDA functional is initialized by the routine

int xc_lda_init(xc_lda_type *p, int functional, int nspin);
input:
  functional: which functional do we want?
  nspin: either XC_UNPOLARIZED or XC_POLARIZED
output:
  p: structure that holds our LDA functional
returns: 0 (OK) or -1 (ERROR)

However, there are a couple of functionals that require some extra information in the initialization. These are

xc_lda_x_init(xc_lda_type *p, int nspin, int dim, int irel);
input:
  dim: either 2 or 3. The dimensionality of the space
  irel: either XC_NON_RELATIVISTIC or XC_RELATIVISTIC

The exchange functional can be evaluated for systems with either 2 or 3 dimensions. It is also possible to apply a relativistic correction.

void xc_lda_c_xalpha_init(xc_lda_type *p, int nspin, int dim, double alpha);
input:
  alpha: Slater's parameter that multiplies the exchange

This correlation functional, added to the exchange functional, produces a total exchange and correlation functional, Exc, equal to 3/2 * alpha * Ex Setting alpha equal to one gives the *usual* Slater Xalpha functional, whereas alpha equal to 2/3 just leaves the exchange functional unchanged.

Evaluation

This is done by the routines

void xc_lda_exc(xc_lda_type *p, double *rho, double *exc);
void xc_lda_vxc(xc_lda_type *p, double *rho, double *exc, double *vxc);
void xc_lda_fxc(xc_lda_type *p, double *rho, double *fxc);
void xc_lda_kxc(xc_lda_type *p, double *rho, double *kxc);
void xc_lda(xc_lda_type *p, double *rho, double *exc, double *vxc, double *fxc, double *kxc);
input:
  p: structure obtained from calling xc_lda_init
  rho[]: the density
output:
  exc: the energy per unit particle
  vxc[]: first derivative of the energy per unit volume
  fxc[]: second derivative of the energy per unit volume
  kxc[]: third derivative of the energy per unit volume

The derivatives are defined as

{\rm vxc}_\alpha = \frac{d \epsilon}{d n_\alpha}

{\rm fxc}_{\alpha\beta} = \frac{d^2 \epsilon}{d n_\alpha d n_\beta}

where ε is the energy per unit volume.

If the functional was initialized with nspin=XC_UNPOLARIZED, the spin indices should be dropped from the previous expressions. Otherwise, the parameters have dimensions rho[2], vxc[2], and fxc[3]. The components of fxc are

f_{\rm xc}[1] = f_{{\rm xc} \uparrow\uparrow}   \qquad   f_{\rm xc}[2] = f_{{\rm xc} \uparrow\downarrow}   \qquad   f_{\rm xc}[3] = f_{{\rm xc} \downarrow\downarrow}   \qquad

Destruction

Due to its simplicity, the LDA functionals do not require, at the moment, an explicit destruction.

The GGAs

Functionals

These are the GGA functionals we know about. Exc, Vxc and Fxc implementation audit status.

Exchange
Correlation
Exchange-Correlation

Initialization

A generic GGA functional is initialized by the routine

int  xc_gga_init(xc_gga_type *p, int functional, int nspin);
input:
  functional: which functional do we want?
  nspin: either XC_UNPOLARIZED or XC_POLARIZED
output:
  p: structure that holds our GGA functional
returns: 0 (OK) or -1 (ERROR)

Evaluation

This is done by the routines

void xc_gga(xc_gga_type *p, double *rho, double *sigma,
            double *exc, double *vrho, double *vsigma)
input:
  p: structure obtained from calling xc_lda_init
  rho[]: the density
  sigma[]: contracted gradients of the density
output:
  exc: the energy per unit particle
  vrho[]: first partial derivative of the energy per unit volume in terms of the density
  vsigma[]: first partial derivative of the energy per unit volume in terms of sigma

The derivatives are defined as

{\rm vrho}_{\alpha} = \frac{\partial \epsilon}{\partial n_\alpha}   \qquad   {\rm vsigma}_{\alpha} = \frac{\partial \epsilon}{\partial \sigma_\alpha}

where ε is the energy per unit volume.

If the functional was initialized with nspin=XC_UNPOLARIZED, the spin indices should be dropped from the previous expressions. Otherwise, the parameters have dimensions rho[2], vxc[2], sigma[3], and vsigma[3]. The components of sigma are

\sigma[1] = \nabla n_\uparrow \cdot \nabla n_\uparrow   \qquad   \sigma[2] = \nabla n_\uparrow \cdot \nabla n_\downarrow   \qquad   \sigma[3] = \nabla n_\downarrow \cdot \nabla n_\downarrow   \qquad

Destruction

If you no longer need your xc functional, you should destruct it with the routine

void xc_gga_end (xc_gga_type *p);
input:
  p: structure holding the functional to destroy

The Hybrid GGAs

Functionals

These are the hybrid GGA functionals we know about. Exc, Vxc and Fxc implementation audit status.

TODO

  • XC_HYB_GGA_XC_mPW91PW91
  • XC_HYB_GGA_XC_PBE0
  • XC_HYB_GGA_B97_1
  • XC_HYB_GGA_B97_1_REV
  • XC_HYB_GGA_B97_2
  • XC_HYB_GGA_SB98H
  • XC_HYB_GGA_X3LYP

Initialization

A generic hybrid GGA functional is initialized by the routine

int  xc_hyb_gga_init(xc_hyb_gga_type *p, int functional, int nspin);
input:
  functional: which functional do we want?
  nspin: either XC_UNPOLARIZED or XC_POLARIZED
output:
  p: structure that holds our hybrid GGA functional
returns: 0 (OK) or -1 (ERROR)

Evaluation

This is done by the routines (cf. the GGAs)

void xc_hyb_gga(xc_hyb_gga_type *p, double *rho, double *sigma,
               double *exc, double *vrho, double *vsigma)
input:
  p: structure obtained from calling xc_hyb_gga_init
  rho[]: the density
  sigma[]: contracted gradients of the density
output:
  exc: the energy per unit particle
  vrho[]: first partial derivative of the energy per unit volume in terms of the density
  vsigma[]: first partial derivative of the energy per unit volume in terms of sigma


Destruction

If you no longer need your xc functional, you should destruct it with the routine

void xc_hyb_gga_end (xc_hyb_gga_type *p);
input:
  p: structure holding the functional to destroy

The MetaGGAs (experimental)

Functionals

These are the MetaGGA functionals we know about. Exc, Vxc and Fxc implementation audit status.

Exchange
Correlation

The LCAs (not working)

Functionals

These are the LCA functionals we know about. Exc, Vxc and Fxc implementation audit status.

Personal tools