librpa.h File Reference
LibRPA
|
C interface of LibRPA. More...
#include <mpi.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.
- 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_low The lowest index of state (included) to compute [in] i_state_high The highest index of state (excluded) to compute [in] n_kpoints_task The 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_task The indices of k-points to compute EXX energy. [out] exx Exchange energy. It should have length of at least n_spins
*n_kpoints_task
* (i_state_high
-i_state_low
). Whenn_kpoints_task
< 0, the length should be at leastn_spins
*n_kpoints
* (i_state_high
-i_state_low
).
◆ finalize_librpa_environment()
void finalize_librpa_environment | ( | ) |
Finalize the environment of LibRPA calculation.
This should be the last LibRPA API function to call.
◆ get_default_librpa_params()
void get_default_librpa_params | ( | LibRPAParams * | params | ) |
Obtain the default controlling parameters of LibRPA.
- Parameters
-
[in] params Struct containing the control parameters
◆ get_frequency_grids()
void get_frequency_grids | ( | int | ngrid, |
double * | freqeuncy_grids | ||
) |
Obtain the frequency grid points.
- Parameters
-
[in] params Struct 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_corr RPA correlation energy, in Hartree unit. Complex number represented as double array [2]. [out] rpa_corr_irk_contrib Weighted 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_in Global MPI communicator [in] is_fortran_comm Flag to identify whether the input communicator is Fortran. [in] redirect_stdout Flag to control whether output LibRPA will be printed to a file [in] output_filename Name of file for redirected output. Only used when redirect_stdout is true
◆ 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] I Index of atom, where one basis and the auxiliary basis reside. [in] J Index of atom, where the other basis resides. [in] nbasis_i Number of basis functions centered at atom I [in] nbasis_j Number of basis functions centered at atom J [in] naux_mu Number of auxliary basis functions centered at atom I [in] R Real-space lattice vector where atom J resides, [3] [in] Cs_in Local RI triple coefficients, dimension: [nbasis_i][nbasis_j][naux_mu] [in] insert_index_only Flag 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] ispin index of spin channel [in] ik index of k-point [in] wfc_real real part of wave function, dimension: [nstates][nbasis] [in] wfc_imag imaginary 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] ik Index of k-point of parsed Coulomb matrix [in] max_naux Total number of auxiliary basis functions [in] mu_begin Starting row index [in] mu_end End row index [in] nu_begin Starting column index [in] nu_end End column index [in] Vq_real_in Real part of the Coulomb matrix block, dimension: [mu_begin-mu_end+1][nu_begin-nu_end+1] [in] Vq_imag_in Imaginary 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] ik Index of k-point of parsed Coulomb matrix [in] I Index of atom for basis functions of the row indices [in] J Index of atom for basis functions of the column indices [in] naux_mu Number of auxliary basis functions centered at atom I [in] naux_nu Number of auxliary basis functions centered at atom J [in] Vq_real_in Real part of the Coulomb matrix block, dimension: [naux_mu][naux_nu] [in] Vq_imag_in Imaginary 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] ik Index of k-point of parsed Coulomb matrix [in] max_naux Total number of auxiliary basis functions [in] mu_begin Starting row index [in] mu_end End row index [in] nu_begin Starting column index [in] nu_end End column index [in] Vq_real_in Real part of the Coulomb matrix block, dimension: [mu_begin-mu_end+1][nu_begin-nu_end+1] [in] Vq_imag_in Imaginary 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] ik Index of k-point of parsed Coulomb matrix [in] I Index of atom for basis functions of the row indices [in] J Index of atom for basis functions of the column indices [in] naux_mu Number of auxliary basis functions centered at atom I [in] naux_nu Number of auxliary basis functions centered at atom J [in] Vq_real_in Real part of the truncated Coulomb matrix block, dimension: [naux_mu][naux_nu] [in] Vq_imag_in Imaginary 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] nspins Number of spin channels [in] nkpts Number of k-points, on which the Kohn-Sham eigenvectors are computed [in] nstates Number of states [in] nbasis Total number of AO basis [in] natoms Total number of atoms
◆ 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_irk Number of irreducible k-points [in] ibz2bz_index Mapping of irreducible k-points to the full set, [nk_irk] [in] wk_irk Weights 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] nk1 Number of k-grids along the 1st reciprocal lattice vector [in] nk2 Number of k-grids along the 2nd reciprocal lattice vector [in] nk3 Number of k-grids along the 3rd reciprocal lattice vector [in] kvecs Coordinates 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_mat pointer to array of real-space lattice vectors, in Bohr unit [in] G_mat pointer 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] params Struct 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] nspins Number of spin channels [in] nkpts Number of k-points, on which the Kohn-Sham eigenvectors are computed [in] nstates Number of states [in] wg Unnormalized occupation number, dimension: [nspins][nkpts][nstates] [in] ekb Eigenvalues in Hartree unit, dimension: [nspins][nkpts][nstates] [in] efermi Fermi level in Hartree unit
The array wg
and ekb
must be stored in C-style row-major order.
Generated by