librpa.h File Reference

LibRPA: librpa.h File Reference
LibRPA
librpa.h File Reference

C interface of LibRPA. More...

#include <mpi.h>
Include dependency graph for librpa.h:

Go to the source code of this file.

Classes

struct  LibRPAParams
 C struct to handle input parameters for LibRPA calculation. More...
 

Functions

void initialize_librpa_environment (MPI_Comm comm_global_in, int is_fortran_comm, int redirect_stdout, const char *output_filename)
 Initialize the environment of LibRPA calculation. More...
 
void finalize_librpa_environment ()
 Finalize the environment of LibRPA calculation. More...
 
void set_dimension (int nspins, int nkpts, int nstates, int nbasis, int natoms)
 Set dimension parameters of the system. More...
 
void set_wg_ekb_efermi (int nspins, int nkpts, int nstates, double *wg, double *ekb, double efermi)
 Set eigenvalues, occupation number and fermi level. More...
 
void set_ao_basis_wfc (int ispin, int ik, double *wfc_real, double *wfc_imag)
 Set wave function expansion of AO basis for a particular spin channel and k-point. More...
 
void set_latvec_and_G (double lat_mat[9], double G_mat[9])
 Set the real-space and reciprocal-space lattice vectors. More...
 
void set_kgrids_kvec_tot (int nk1, int nk2, int nk3, double *kvecs)
 Set k-mesh grids. More...
 
void set_ibz2bz_index_and_weight (const int nk_irk, const int *ibz2bz_index, const double *wk_irk)
 Set the mapping of irreducible k-points to the whole k-points set and their weights. More...
 
void set_ao_basis_aux (int I, int J, int nbasis_i, int nbasis_j, int naux_mu, int *R, double *Cs_in, int insert_index_only)
 Insert the atom index and set the local RI triple coefficients. More...
 
void set_aux_bare_coulomb_k_atom_pair (int ik, int I, int J, int naux_mu, int naux_nu, double *Vq_real_in, double *Vq_imag_in)
 Set the atom-pair block of bare Coulomb matrix in auxiliary basis. More...
 
void set_aux_cut_coulomb_k_atom_pair (int ik, int I, int J, int naux_mu, int naux_nu, double *Vq_real_in, double *Vq_imag_in)
 Set the atom-pair block of truncated (cut) Coulomb matrix in auxiliary basis. More...
 
void set_aux_bare_coulomb_k_2D_block (int ik, int max_naux, int mu_begin, int mu_end, int nu_begin, int nu_end, double *Vq_real_in, double *Vq_imag_in)
 Set the bare Coulomb matrix in auxiliary basis under BLACS 2D block format. More...
 
void set_aux_cut_coulomb_k_2D_block (int ik, int max_naux, int mu_begin, int mu_end, int nu_begin, int nu_end, double *Vq_real_in, double *Vq_imag_in)
 Set the truncated (cut) Coulomb matrix in auxiliary basis under BLACS 2D block format. More...
 
void set_librpa_params (LibRPAParams *params)
 Set the control parameters of LibRPA. More...
 
void get_default_librpa_params (LibRPAParams *params)
 Obtain the default controlling parameters of LibRPA. More...
 
void get_frequency_grids (int ngrid, double *freqeuncy_grids)
 Obtain the frequency grid points. More...
 
void run_librpa_main ()
 
void get_rpa_correlation_energy (double *rpa_corr, double *rpa_corr_irk_contrib)
 compute and return the RPA correlation energy More...
 
void compute_exx_orbital_energy (int i_state_low, int i_state_high, int n_kpoints_task, const int *i_kpoints_task, double *exx)
 Compute the exact exchange (EXX) energy for states at specified k-points. More...
 

Detailed Description

C interface of LibRPA.

Author
LibRPA developers
Date
2024-06-28

Function Documentation

◆ compute_exx_orbital_energy()

void compute_exx_orbital_energy ( int  i_state_low,
int  i_state_high,
int  n_kpoints_task,
const int *  i_kpoints_task,
double *  exx 
)

Compute the exact exchange (EXX) energy for states at specified k-points.

Parameters
[in]i_state_lowThe lowest index of state (included) to compute
[in]i_state_highThe highest index of state (excluded) to compute
[in]n_kpoints_taskThe number of k-points to return in the called process. When equal to 0, an empty vector will be returned. When less than 0, all k-points will be computed. Otherwise, the states at k-points whose indices are stored in i_kpoints_task will be computed.
[in]i_kpoints_taskThe indices of k-points to compute EXX energy.
[out]exxExchange energy. It should have length of at least n_spins * n_kpoints_task * (i_state_high - i_state_low). When n_kpoints_task < 0, the length should be at least n_spins * n_kpoints * (i_state_high - i_state_low).
Here is the call graph for this function:

◆ finalize_librpa_environment()

void finalize_librpa_environment ( )

Finalize the environment of LibRPA calculation.

This should be the last LibRPA API function to call.

Here is the call graph for this function:

◆ get_default_librpa_params()

void get_default_librpa_params ( LibRPAParams params)

Obtain the default controlling parameters of LibRPA.

Parameters
[in]paramsStruct containing the control parameters

◆ get_frequency_grids()

void get_frequency_grids ( int  ngrid,
double *  freqeuncy_grids 
)

Obtain the frequency grid points.

Parameters
[in]paramsStruct containing the control parameters

◆ get_rpa_correlation_energy()

void get_rpa_correlation_energy ( double *  rpa_corr,
double *  rpa_corr_irk_contrib 
)

compute and return the RPA correlation energy

Parameters
[out]rpa_corrRPA correlation energy, in Hartree unit. Complex number represented as double array [2].
[out]rpa_corr_irk_contribWeighted contribution to RPA correlation energy from each irreducible k-point. Complex array represented as double array [2*n_irk_points]
Warning
Current implementation will free the internal copy of LRI coefficients array. So this can only be called once.

◆ initialize_librpa_environment()

void initialize_librpa_environment ( MPI_Comm  comm_global_in,
int  is_fortran_comm,
int  redirect_stdout,
const char *  output_filename 
)

Initialize the environment of LibRPA calculation.

Parameters
[in]comm_global_inGlobal MPI communicator
[in]is_fortran_commFlag to identify whether the input communicator is Fortran.
[in]redirect_stdoutFlag to control whether output LibRPA will be printed to a file
[in]output_filenameName of file for redirected output. Only used when redirect_stdout is true
Here is the call graph for this function:

◆ set_ao_basis_aux()

void set_ao_basis_aux ( int  I,
int  J,
int  nbasis_i,
int  nbasis_j,
int  naux_mu,
int *  R,
double *  Cs_in,
int  insert_index_only 
)

Insert the atom index and set the local RI triple coefficients.

Parameters
[in]IIndex of atom, where one basis and the auxiliary basis reside.
[in]JIndex of atom, where the other basis resides.
[in]nbasis_iNumber of basis functions centered at atom I
[in]nbasis_jNumber of basis functions centered at atom J
[in]naux_muNumber of auxliary basis functions centered at atom I
[in]RReal-space lattice vector where atom J resides, [3]
[in]Cs_inLocal RI triple coefficients, dimension: [nbasis_i][nbasis_j][naux_mu]
[in]insert_index_onlyFlag to insert the indices of atoms but skip setting Cs_in

The array Cs_in must be stored in C-style row-major order.

◆ set_ao_basis_wfc()

void set_ao_basis_wfc ( int  ispin,
int  ik,
double *  wfc_real,
double *  wfc_imag 
)

Set wave function expansion of AO basis for a particular spin channel and k-point.

Parameters
[in]ispinindex of spin channel
[in]ikindex of k-point
[in]wfc_realreal part of wave function, dimension: [nstates][nbasis]
[in]wfc_imagimaginary part of wave function, dimension: [nstates][nbasis]

The array wfc_real and wfc_imag must be stored in C-style row-major order.

◆ set_aux_bare_coulomb_k_2D_block()

void set_aux_bare_coulomb_k_2D_block ( int  ik,
int  max_naux,
int  mu_begin,
int  mu_end,
int  nu_begin,
int  nu_end,
double *  Vq_real_in,
double *  Vq_imag_in 
)

Set the bare Coulomb matrix in auxiliary basis under BLACS 2D block format.

Parameters
[in]ikIndex of k-point of parsed Coulomb matrix
[in]max_nauxTotal number of auxiliary basis functions
[in]mu_beginStarting row index
[in]mu_endEnd row index
[in]nu_beginStarting column index
[in]nu_endEnd column index
[in]Vq_real_inReal part of the Coulomb matrix block, dimension: [mu_begin-mu_end+1][nu_begin-nu_end+1]
[in]Vq_imag_inImaginary part of the Coulomb matrix block, dimension: [mu_begin-mu_end+1][nu_begin-nu_end+1]

The array Vq_real_in and Vq_imag_in must be stored in C-style row-major order.

◆ set_aux_bare_coulomb_k_atom_pair()

void set_aux_bare_coulomb_k_atom_pair ( int  ik,
int  I,
int  J,
int  naux_mu,
int  naux_nu,
double *  Vq_real_in,
double *  Vq_imag_in 
)

Set the atom-pair block of bare Coulomb matrix in auxiliary basis.

Parameters
[in]ikIndex of k-point of parsed Coulomb matrix
[in]IIndex of atom for basis functions of the row indices
[in]JIndex of atom for basis functions of the column indices
[in]naux_muNumber of auxliary basis functions centered at atom I
[in]naux_nuNumber of auxliary basis functions centered at atom J
[in]Vq_real_inReal part of the Coulomb matrix block, dimension: [naux_mu][naux_nu]
[in]Vq_imag_inImaginary part of the Coulomb matrix block, dimension: [naux_mu][naux_nu]

The array Vq_real_in and Vq_imag_in must be stored in C-style row-major order.

◆ set_aux_cut_coulomb_k_2D_block()

void set_aux_cut_coulomb_k_2D_block ( int  ik,
int  max_naux,
int  mu_begin,
int  mu_end,
int  nu_begin,
int  nu_end,
double *  Vq_real_in,
double *  Vq_imag_in 
)

Set the truncated (cut) Coulomb matrix in auxiliary basis under BLACS 2D block format.

Parameters
[in]ikIndex of k-point of parsed Coulomb matrix
[in]max_nauxTotal number of auxiliary basis functions
[in]mu_beginStarting row index
[in]mu_endEnd row index
[in]nu_beginStarting column index
[in]nu_endEnd column index
[in]Vq_real_inReal part of the Coulomb matrix block, dimension: [mu_begin-mu_end+1][nu_begin-nu_end+1]
[in]Vq_imag_inImaginary part of the Coulomb matrix block, dimension: [mu_begin-mu_end+1][nu_begin-nu_end+1]

The array Vq_real_in and Vq_imag_in must be stored in C-style row-major order.

◆ set_aux_cut_coulomb_k_atom_pair()

void set_aux_cut_coulomb_k_atom_pair ( int  ik,
int  I,
int  J,
int  naux_mu,
int  naux_nu,
double *  Vq_real_in,
double *  Vq_imag_in 
)

Set the atom-pair block of truncated (cut) Coulomb matrix in auxiliary basis.

Parameters
[in]ikIndex of k-point of parsed Coulomb matrix
[in]IIndex of atom for basis functions of the row indices
[in]JIndex of atom for basis functions of the column indices
[in]naux_muNumber of auxliary basis functions centered at atom I
[in]naux_nuNumber of auxliary basis functions centered at atom J
[in]Vq_real_inReal part of the truncated Coulomb matrix block, dimension: [naux_mu][naux_nu]
[in]Vq_imag_inImaginary part of the truncated Coulomb matrix block, dimension: [naux_mu][naux_nu]

The array Vq_real_in and Vq_imag_in must be stored in C-style row-major order.

◆ set_dimension()

void set_dimension ( int  nspins,
int  nkpts,
int  nstates,
int  nbasis,
int  natoms 
)

Set dimension parameters of the system.

Parameters
[in]nspinsNumber of spin channels
[in]nkptsNumber of k-points, on which the Kohn-Sham eigenvectors are computed
[in]nstatesNumber of states
[in]nbasisTotal number of AO basis
[in]natomsTotal number of atoms
Here is the call graph for this function:

◆ set_ibz2bz_index_and_weight()

void set_ibz2bz_index_and_weight ( const int  nk_irk,
const int *  ibz2bz_index,
const double *  wk_irk 
)

Set the mapping of irreducible k-points to the whole k-points set and their weights.

Parameters
[in]nk_irkNumber of irreducible k-points
[in]ibz2bz_indexMapping of irreducible k-points to the full set, [nk_irk]
[in]wk_irkWeights of irreducible k-points, [nk_irk]

◆ set_kgrids_kvec_tot()

void set_kgrids_kvec_tot ( int  nk1,
int  nk2,
int  nk3,
double *  kvecs 
)

Set k-mesh grids.

Parameters
[in]nk1Number of k-grids along the 1st reciprocal lattice vector
[in]nk2Number of k-grids along the 2nd reciprocal lattice vector
[in]nk3Number of k-grids along the 3rd reciprocal lattice vector
[in]kvecsCoordinates of k-mesh vectors, in inverse Bohr unit, dimension: [nk1*nk2*nk3][3]

The array kvecs must be stored in C-style row-major order.

◆ set_latvec_and_G()

void set_latvec_and_G ( double  lat_mat[9],
double  G_mat[9] 
)

Set the real-space and reciprocal-space lattice vectors.

Parameters
[in]lat_matpointer to array of real-space lattice vectors, in Bohr unit
[in]G_matpointer to array of reciprocal-space lattice vectors, in inverse Bohr unit

◆ set_librpa_params()

void set_librpa_params ( LibRPAParams params)

Set the control parameters of LibRPA.

Parameters
[in]paramsStruct containing the control parameters

◆ set_wg_ekb_efermi()

void set_wg_ekb_efermi ( int  nspins,
int  nkpts,
int  nstates,
double *  wg,
double *  ekb,
double  efermi 
)

Set eigenvalues, occupation number and fermi level.

Parameters
[in]nspinsNumber of spin channels
[in]nkptsNumber of k-points, on which the Kohn-Sham eigenvectors are computed
[in]nstatesNumber of states
[in]wgUnnormalized occupation number, dimension: [nspins][nkpts][nstates]
[in]ekbEigenvalues in Hartree unit, dimension: [nspins][nkpts][nstates]
[in]efermiFermi level in Hartree unit

The array wg and ekb must be stored in C-style row-major order.