|
| | RPA_UKS (Logger &log, const TCMatrix_gwbse_spin &Mmn) |
| void | configure (Index homo_alpha, Index homo_beta, Index rpamin, Index rpamax) |
| | Configure orbital window and spin-resolved HOMO indices.
|
| double | getEta () const |
| | Small positive broadening used in the real-frequency response.
|
| Eigen::MatrixXd | calculate_epsilon_i (double frequency) const |
| | Dielectric matrix on the imaginary frequency axis.
|
| Eigen::MatrixXd | calculate_epsilon_r (double frequency) const |
| | Dielectric matrix on the real frequency axis.
|
| Eigen::MatrixXd | calculate_epsilon_r (std::complex< double > frequency) const |
| | Real part of the dielectric matrix at complex frequency.
|
| const Eigen::VectorXd & | getRPAInputEnergiesAlpha () const |
| | Access the current alpha RPA input energies used in the response.
|
| const Eigen::VectorXd & | getRPAInputEnergiesBeta () const |
| | Access the current beta RPA input energies used in the response.
|
| void | setRPAInputEnergies (const Eigen::VectorXd &rpaenergies_alpha, const Eigen::VectorXd &rpaenergies_beta) |
| | Set alpha and beta RPA input energies directly.
|
| void | UpdateRPAInputEnergies (const Eigen::VectorXd &dftenergies_alpha, const Eigen::VectorXd &dftenergies_beta, const Eigen::VectorXd &gwaenergies_alpha, const Eigen::VectorXd &gwaenergies_beta, Index qpmin) |
| | Update alpha and beta RPA input energies from DFT and GW energies.
|
| const rpa_eigensolution & | Diagonalize_H2p () const |
| | Diagonalize the unrestricted two-particle Hamiltonian.
|
| Index | getHomoAlpha () const |
| Index | getHomoBeta () const |
| void | GetCachedScreeningModes (const Eigen::VectorXd *&omegas, const std::vector< Eigen::VectorXd > *&modes) const |
|
| template<bool imag> |
| Eigen::MatrixXd | calculate_epsilon (double frequency) const |
| | Internal implementation of epsilon(w).
|
| void | InvalidateH2pCache () |
| Eigen::VectorXd | Calculate_H2p_AmB () const |
| | Construct the diagonal AmB block of the two-particle Hamiltonian.
|
| Eigen::MatrixXd | Calculate_H2p_ApB () const |
| | Construct the ApB block of the unrestricted two-particle Hamiltonian.
|
| Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > | Diagonalize_H2p_C (const Eigen::MatrixXd &C) const |
| | Diagonalize the symmetrized two-particle matrix C.
|
| void | ShiftUncorrectedEnergies (Eigen::VectorXd &energies, const Eigen::VectorXd &dftenergies, Index homo, Index qpmin, Index gwsize) |
| | Shift uncorrected states outside the GW window.
|
| double | getMaxCorrection (const Eigen::VectorXd &energies, const Eigen::VectorXd &dftenergies, Index min, Index max) const |
| | Return the largest absolute GW correction in a given index range.
|
| void | BuildCachedScreeningModes () const |
Unrestricted RPA helper for spin-resolved GW screening.
This class is the unrestricted/open-shell analogue of the existing restricted RPA class. It constructs the independent-particle polarizability and related dielectric quantities from explicit alpha and beta particle-hole channels.
Conceptually, for a collinear unrestricted reference,
chi0(w) = chi0_alpha(w) + chi0_beta(w)
where each spin contribution contains only same-spin particle-hole excitations:
i_alpha -> a_alpha i_beta -> a_beta
The Coulomb interaction is spin-independent, so the resulting dielectric matrix is built from the sum of both spin channels. This is exactly what is needed later for unrestricted GW, where Sigma^alpha and Sigma^beta share the same screened interaction W but differ in exchange and in the external spin-resolved states on which Sigma acts.
Compared to the restricted closed-shell implementation:
- we do not assume one doubly occupied spatial-orbital ladder,
- we do not insert a blanket closed-shell spin-degeneracy factor,
- instead, alpha and beta channels are treated explicitly and summed.
The class also provides a two-particle Hamiltonian construction (Diagonalize_H2p), mainly mirroring the restricted implementation. For the initial unrestricted GW work, the dielectric functions are the more critical part; the H2p diagonalization is included for completeness and future use.
Definition at line 68 of file rpa_uks.h.
| Eigen::MatrixXd votca::xtp::RPA_UKS::Calculate_H2p_ApB |
( |
| ) |
const |
|
private |
Construct the ApB block of the unrestricted two-particle Hamiltonian.
This contains Coulomb coupling between particle-hole excitations. In the unrestricted case, all four spin blocks are present:
alpha-alpha, alpha-beta, beta-alpha, beta-beta
because the Coulomb interaction couples total density fluctuations and is spin independent.
Definition at line 464 of file rpa_uks.cc.
| void votca::xtp::RPA_UKS::configure |
( |
Index | homo_alpha, |
|
|
Index | homo_beta, |
|
|
Index | rpamin, |
|
|
Index | rpamax ) |
|
inline |
Configure orbital window and spin-resolved HOMO indices.
The same global orbital window [rpamin_, rpamax_] is used for both spin channels, but the occupied/virtual split is spin dependent:
alpha: occupied <= homo_alpha_, virtual >= homo_alpha_ + 1 beta : occupied <= homo_beta_, virtual >= homo_beta_ + 1
This keeps the interface close to the restricted implementation while still allowing different alpha and beta occupations.
Definition at line 91 of file rpa_uks.h.