librpa.h Source File
LibRPA
|
librpa.h
Go to the documentation of this file.
136 void set_wg_ekb_efermi(int nspins, int nkpts, int nstates, double* wg, double* ekb, double efermi);
176 void set_ibz2bz_index_and_weight(const int nk_irk, const int* ibz2bz_index, const double* wk_irk);
192 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);
207 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);
222 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);
238 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);
254 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);
316 // void compute_gw_quasiparticle_energy_kgrid(int n_qp_state_low, int n_qp_state_high, const double *vxc, double _Complex *screened_coul);
void set_latvec_and_G(double lat_mat[9], double G_mat[9])
Set the real-space and reciprocal-space lattice vectors.
Definition: librpa.cpp:121
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.
Definition: librpa.cpp:189
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.
Definition: librpa.cpp:106
void finalize_librpa_environment()
Finalize the environment of LibRPA calculation.
Definition: librpa.cpp:64
void set_wg_ekb_efermi(int nspins, int nkpts, int nstates, double *wg, double *ekb, double efermi)
Set eigenvalues, occupation number and fermi level.
Definition: librpa.cpp:85
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.
Definition: librpa.cpp:417
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.
Definition: librpa.cpp:41
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.
Definition: librpa.cpp:329
void set_librpa_params(LibRPAParams *params)
Set the control parameters of LibRPA.
Definition: librpa.cpp:335
void get_rpa_correlation_energy(double *rpa_corr, double *rpa_corr_irk_contrib)
compute and return the RPA correlation energy
Definition: librpa.cpp:404
void set_kgrids_kvec_tot(int nk1, int nk2, int nk3, double *kvecs)
Set k-mesh grids.
Definition: librpa.cpp:163
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.
Definition: librpa.cpp:289
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.
Definition: librpa.cpp:323
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.
Definition: librpa.cpp:283
void get_frequency_grids(int ngrid, double *freqeuncy_grids)
Obtain the frequency grid points.
Definition: librpa.cpp:401
void get_default_librpa_params(LibRPAParams *params)
Obtain the default controlling parameters of LibRPA.
Definition: librpa.cpp:360
void set_dimension(int nspins, int nkpts, int nstates, int nbasis, int natoms)
Set dimension parameters of the system.
Definition: librpa.cpp:70
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.
Definition: librpa.cpp:208
double libri_exx_threshold_D
Threshold of real-space density matrices to compute exact exchange using LibRI.
Definition: librpa.h:79
double libri_exx_threshold_CSM
Threshold of Cauchy-Schwarz filtering to compute exact exchange using LibRI.
Definition: librpa.h:73
char output_file[100]
File path for LibRPA output redirection.
Definition: librpa.h:28
double libri_gw_threshold_W
Threshold of correlation screened Coulomb matrix to compute GW correlation self-energy using LibRI.
Definition: librpa.h:91
double sqrt_coulomb_threshold
Threshold of eigenvalues to perform square root of Coulomb matrices.
Definition: librpa.h:64
double libri_exx_threshold_V
Threshold of real-space Coulomb matrices to compute exact exchange using LibRI.
Definition: librpa.h:82
double cs_threshold
Threshold for real-space LRI triple coefficients.
Definition: librpa.h:58
double libri_gw_threshold_G
Threshold of real-space Green's function to compute GW correlation self-energy using LibRI.
Definition: librpa.h:88
double libri_exx_threshold_C
Threshold of real-space LRI triple coefficients to compute exact exchange using LibRI.
Definition: librpa.h:76
int use_scalapack_ecrpa
Flag to control whether to use ScaLAPACK to compute RPA correlation energy.
Definition: librpa.h:51
char output_dir[100]
Directory path for LibRPA debug output.
Definition: librpa.h:31
double libri_gw_threshold_C
Threshold of real-space LRI triple coefficients to compute GW correlation self-energy using LibRI.
Definition: librpa.h:85
char DFT_software[20]
Choice of hosting DFT software.
Definition: librpa.h:40
double libri_chi0_threshold_C
Threshold of real-space LRI triple coefficients to compute response function using LibRI.
Definition: librpa.h:67
char parallel_routing[20]
Scheme of parallelization.
Definition: librpa.h:34
double libri_chi0_threshold_G
Threshold of real-space Green's function to compute response function using LibRI.
Definition: librpa.h:70
double gf_R_threshold
Threshold for real-space Green's function.
Definition: librpa.h:55
double vq_threshold
Threshold for real-space Coulomb matrices.
Definition: librpa.h:61
Generated by