votca 2026-dev
Loading...
Searching...
No Matches
votca::xtp::GW Class Reference

#include <gw.h>

Collaboration diagram for votca::xtp::GW:

Classes

struct  options
class  QPFunc

Public Member Functions

 GW (Logger &log, TCMatrix_gwbse &Mmn, const Eigen::MatrixXd &vxc, const Eigen::VectorXd &dft_energies)
void configure (const options &opt)
Eigen::VectorXd getGWAResults () const
void CalculateGWPerturbation ()
void CalculateHQP ()
void CalculateQSGW ()
 Run quasiparticle self-consistent GW (QSGW).
const Eigen::MatrixXd & getQSGWRotation () const
const Eigen::VectorXd & getQSGWSeedEnergies () const
 Return the seed (G0W0/evGW) energies that QSGW started from.
void PrintQSGW_Energies (const std::string &seed_label, const Eigen::VectorXd &seed_energies, const Eigen::VectorXd &qsgw_energies) const
 Print a two-column comparison of seed (G0W0 or evGW) vs converged QSGW quasiparticle energies.
void PrintQSGW_Composition (double threshold=0.01) const
 Print the dominant DFT-KS orbital contributions to each converged QSGW quasiparticle state.
Eigen::MatrixXd getHQP () const
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonian () const
void PlotSigma (std::string filename, Index steps, double spacing, std::string states) const
Eigen::VectorXd RPAInputEnergies () const

Private Types

using EvalStage = qp_solver::EvalStage
using QPStats = qp_solver::Stats
using QPRootCandidate = qp_solver::RootCandidate
using QPWindowDiagnostics = qp_solver::WindowDiagnostics

Private Member Functions

double CalcHomoLumoShift (Eigen::VectorXd frequencies) const
Eigen::VectorXd ScissorShift_DFTlevel (const Eigen::VectorXd &dft_energies) const
void PrintQP_Energies (const Eigen::VectorXd &qp_diag_energies) const
void PrintGWA_Energies () const
Eigen::VectorXd SolveQP (const Eigen::VectorXd &frequencies) const
boost::optional< double > SolveQP_Grid (double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
boost::optional< double > SolveQP_Grid_Windowed (double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
boost::optional< double > SolveQP_Grid_Windowed_Adaptive (double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
boost::optional< double > SolveQP_Grid_Windowed_Dense (double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
boost::optional< double > SolveQP_FixedPoint (double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
boost::optional< double > SolveQP_Linearisation (double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
bool Converged (const Eigen::VectorXd &e1, const Eigen::VectorXd &e2, double epsilon) const
boost::optional< QPRootCandidateRefineQPInterval (double lowerbound, double f_lowerbound, double upperbound, double f_upperbound, const QPFunc &f, double reference) const

Private Attributes

Index qptotal_
Eigen::MatrixXd Sigma_x_
Eigen::MatrixXd Sigma_c_
Eigen::MatrixXd qsgw_rotation_
Eigen::VectorXd qsgw_seed_energies_
Eigen::VectorXd qsgw_final_energies_
options opt_
std::unique_ptr< Sigma_basesigma_ = nullptr
Loggerlog_
TCMatrix_gwbseMmn_
const Eigen::MatrixXd & vxc_
const Eigen::VectorXd & dft_energies_
Index gw_sc_iteration_
RPA rpa_

Detailed Description

Definition at line 37 of file gw.h.

Member Typedef Documentation

◆ EvalStage

Definition at line 39 of file gw.h.

◆ QPRootCandidate

Definition at line 41 of file gw.h.

◆ QPStats

Definition at line 40 of file gw.h.

◆ QPWindowDiagnostics

Definition at line 42 of file gw.h.

Constructor & Destructor Documentation

◆ GW()

votca::xtp::GW::GW ( Logger & log,
TCMatrix_gwbse & Mmn,
const Eigen::MatrixXd & vxc,
const Eigen::VectorXd & dft_energies )
inline

Definition at line 45 of file gw.h.

Member Function Documentation

◆ CalcHomoLumoShift()

double votca::xtp::GW::CalcHomoLumoShift ( Eigen::VectorXd frequencies) const
private

Definition at line 60 of file gw.cc.

◆ CalculateGWPerturbation()

void votca::xtp::GW::CalculateGWPerturbation ( )

Definition at line 218 of file gw.cc.

◆ CalculateHQP()

void votca::xtp::GW::CalculateHQP ( )

Definition at line 772 of file gw.cc.

◆ CalculateQSGW()

void votca::xtp::GW::CalculateQSGW ( )

Run quasiparticle self-consistent GW (QSGW).

At each iteration:

  1. Compute the full off-diagonal self-energy matrix tilde_Sigma (symmetrised static approximation) in the current QP basis.
  2. Diagonalise H_QSGW = diag(e_DFT - v_xc) + tilde_Sigma to obtain new QP energies and a rotation matrix U.
  3. Rotate the three-centre integrals Mmn via U so that W and Sigma are computed from the updated QP wavefunctions next iteration. Converges when max|e_QP^{k+1} - e_QP^k| < qsgw_sc_limit.

The accumulated rotation U (DFT MOs -> converged QP wavefunctions) is stored in qsgw_rotation_ and returned via getQSGWRotation(). QPdiag energies and eigenvectors are set to the converged QSGW values.

Definition at line 798 of file gw.cc.

◆ configure()

void votca::xtp::GW::configure ( const options & opt)

Definition at line 35 of file gw.cc.

◆ Converged()

bool votca::xtp::GW::Converged ( const Eigen::VectorXd & e1,
const Eigen::VectorXd & e2,
double epsilon ) const
private

Definition at line 759 of file gw.cc.

◆ DiagonalizeQPHamiltonian()

Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > votca::xtp::GW::DiagonalizeQPHamiltonian ( ) const

Definition at line 73 of file gw.cc.

◆ getGWAResults()

Eigen::VectorXd votca::xtp::GW::getGWAResults ( ) const

Definition at line 312 of file gw.cc.

◆ getHQP()

Eigen::MatrixXd votca::xtp::GW::getHQP ( ) const

Definition at line 67 of file gw.cc.

◆ getQSGWRotation()

const Eigen::MatrixXd & votca::xtp::GW::getQSGWRotation ( ) const
inline

Return the accumulated QSGW rotation matrix U (DFT MOs -> QP wavefunctions)

Definition at line 150 of file gw.h.

◆ getQSGWSeedEnergies()

const Eigen::VectorXd & votca::xtp::GW::getQSGWSeedEnergies ( ) const
inline

Return the seed (G0W0/evGW) energies that QSGW started from.

Definition at line 153 of file gw.h.

◆ PlotSigma()

void votca::xtp::GW::PlotSigma ( std::string filename,
Index steps,
double spacing,
std::string states ) const

Definition at line 1133 of file gw.cc.

◆ PrintGWA_Energies()

void votca::xtp::GW::PrintGWA_Energies ( ) const
private

Definition at line 80 of file gw.cc.

◆ PrintQP_Energies()

void votca::xtp::GW::PrintQP_Energies ( const Eigen::VectorXd & qp_diag_energies) const
private

Definition at line 183 of file gw.cc.

◆ PrintQSGW_Composition()

void votca::xtp::GW::PrintQSGW_Composition ( double threshold = 0.01) const

Print the dominant DFT-KS orbital contributions to each converged QSGW quasiparticle state.

For each QP state n, prints the KS orbitals m with |U_{mn}|^2 > threshold, where U = qsgw_rotation_ (columns = QP states in DFT-KS basis).

Parameters
thresholdMinimum weight to print (default 0.01 = 1%)

Definition at line 113 of file gw.cc.

◆ PrintQSGW_Energies()

void votca::xtp::GW::PrintQSGW_Energies ( const std::string & seed_label,
const Eigen::VectorXd & seed_energies,
const Eigen::VectorXd & qsgw_energies ) const

Print a two-column comparison of seed (G0W0 or evGW) vs converged QSGW quasiparticle energies.

Parameters
seed_labelLabel for the seed column, e.g. "G0W0" or "evGW".
seed_energiesQP energies from the seed calculation (Ha).
qsgw_energiesConverged QSGW eigenvalues (Ha).

Definition at line 156 of file gw.cc.

◆ RefineQPInterval()

boost::optional< QPRootCandidate > votca::xtp::GW::RefineQPInterval ( double lowerbound,
double f_lowerbound,
double upperbound,
double f_upperbound,
const QPFunc & f,
double reference ) const
private

◆ RPAInputEnergies()

Eigen::VectorXd votca::xtp::GW::RPAInputEnergies ( ) const
inline

Definition at line 189 of file gw.h.

◆ ScissorShift_DFTlevel()

Eigen::VectorXd votca::xtp::GW::ScissorShift_DFTlevel ( const Eigen::VectorXd & dft_energies) const
private

Definition at line 210 of file gw.cc.

◆ SolveQP()

Eigen::VectorXd votca::xtp::GW::SolveQP ( const Eigen::VectorXd & frequencies) const
private

Definition at line 323 of file gw.cc.

◆ SolveQP_FixedPoint()

boost::optional< double > votca::xtp::GW::SolveQP_FixedPoint ( double intercept0,
double frequency0,
Index gw_level,
QPStats * stats = nullptr ) const
private

Definition at line 741 of file gw.cc.

◆ SolveQP_Grid()

boost::optional< double > votca::xtp::GW::SolveQP_Grid ( double intercept0,
double frequency0,
Index gw_level,
QPStats * stats = nullptr ) const
private

Definition at line 677 of file gw.cc.

◆ SolveQP_Grid_Windowed()

boost::optional< double > votca::xtp::GW::SolveQP_Grid_Windowed ( double intercept0,
double frequency0,
Index gw_level,
double left_limit,
double right_limit,
bool allow_rejected_return = true,
QPStats * stats = nullptr ) const
private

Definition at line 623 of file gw.cc.

◆ SolveQP_Grid_Windowed_Adaptive()

boost::optional< double > votca::xtp::GW::SolveQP_Grid_Windowed_Adaptive ( double intercept0,
double frequency0,
Index gw_level,
double left_limit,
double right_limit,
bool allow_rejected_return = true,
QPStats * stats = nullptr ) const
private

Definition at line 433 of file gw.cc.

◆ SolveQP_Grid_Windowed_Dense()

boost::optional< double > votca::xtp::GW::SolveQP_Grid_Windowed_Dense ( double intercept0,
double frequency0,
Index gw_level,
double left_limit,
double right_limit,
bool allow_rejected_return = true,
QPStats * stats = nullptr ) const
private

Definition at line 503 of file gw.cc.

◆ SolveQP_Linearisation()

boost::optional< double > votca::xtp::GW::SolveQP_Linearisation ( double intercept0,
double frequency0,
Index gw_level,
QPStats * stats = nullptr ) const
private

Definition at line 412 of file gw.cc.

Member Data Documentation

◆ dft_energies_

const Eigen::VectorXd& votca::xtp::GW::dft_energies_
private

Definition at line 213 of file gw.h.

◆ gw_sc_iteration_

Index votca::xtp::GW::gw_sc_iteration_
private

Definition at line 215 of file gw.h.

◆ log_

Logger& votca::xtp::GW::log_
private

Definition at line 210 of file gw.h.

◆ Mmn_

TCMatrix_gwbse& votca::xtp::GW::Mmn_
private

Definition at line 211 of file gw.h.

◆ opt_

options votca::xtp::GW::opt_
private

Definition at line 207 of file gw.h.

◆ qptotal_

Index votca::xtp::GW::qptotal_
private

Definition at line 194 of file gw.h.

◆ qsgw_final_energies_

Eigen::VectorXd votca::xtp::GW::qsgw_final_energies_
private

Definition at line 205 of file gw.h.

◆ qsgw_rotation_

Eigen::MatrixXd votca::xtp::GW::qsgw_rotation_
private

Definition at line 198 of file gw.h.

◆ qsgw_seed_energies_

Eigen::VectorXd votca::xtp::GW::qsgw_seed_energies_
private

Definition at line 200 of file gw.h.

◆ rpa_

RPA votca::xtp::GW::rpa_
private

Definition at line 217 of file gw.h.

◆ sigma_

std::unique_ptr<Sigma_base> votca::xtp::GW::sigma_ = nullptr
private

Definition at line 209 of file gw.h.

◆ Sigma_c_

Eigen::MatrixXd votca::xtp::GW::Sigma_c_
private

Definition at line 197 of file gw.h.

◆ Sigma_x_

Eigen::MatrixXd votca::xtp::GW::Sigma_x_
private

Definition at line 196 of file gw.h.

◆ vxc_

const Eigen::MatrixXd& votca::xtp::GW::vxc_
private

Definition at line 212 of file gw.h.


The documentation for this class was generated from the following files:
  • xtp/include/votca/xtp/gw.h
  • xtp/src/libxtp/gwbse/gw.cc