8#ifndef VOTCA_XTP_GW_UKS_H
9#define VOTCA_XTP_GW_UKS_H
12#include <unordered_set>
82 const Eigen::MatrixXd& vxc_alpha,
const Eigen::MatrixXd& vxc_beta,
83 const Eigen::VectorXd& dft_energies_alpha,
84 const Eigen::VectorXd& dft_energies_beta);
109 std::pair<double, double>
operator()(
double frequency)
const {
110 std::pair<double, double> result;
111 result.first =
value(frequency, EvalStage::Other);
112 result.second =
deriv(frequency);
120 if (!insert_result.second) {
121 ++
stats_.sigma_repeat_calls;
123 ++
stats_.sigma_unique_frequencies;
134 double deriv(
double frequency)
const {
145 std::uint64_t key = 0;
146 static_assert(
sizeof(double) ==
sizeof(std::uint64_t),
147 "Unexpected double size");
148 std::memcpy(&key, &x,
sizeof(
double));
154 case EvalStage::Scan:
155 ++
stats_.sigma_scan_calls;
157 case EvalStage::Refine:
158 ++
stats_.sigma_refine_calls;
160 case EvalStage::Derivative:
161 ++
stats_.sigma_derivative_calls;
163 case EvalStage::Other:
165 ++
stats_.sigma_other_calls;
180 const Eigen::MatrixXd&
Vxc(
Spin spin)
const;
183 const Eigen::MatrixXd&
SigmaX(
Spin spin)
const;
184 const Eigen::MatrixXd&
SigmaC(
Spin spin)
const;
196 const Eigen::VectorXd& qp_diag_energies)
const;
197 Eigen::VectorXd
SolveQP(
Spin spin,
const Eigen::VectorXd& frequencies)
const;
200 double frequency0,
Index gw_level,
201 QPStats* stats =
nullptr)
const;
204 Spin spin,
double intercept0,
double frequency0,
Index gw_level,
205 double left_limit,
double right_limit,
bool allow_rejected_return =
true,
206 QPStats* stats =
nullptr)
const;
209 Spin spin,
double intercept0,
double frequency0,
Index gw_level,
210 double left_limit,
double right_limit,
bool allow_rejected_return =
true,
211 QPStats* stats =
nullptr)
const;
214 Spin spin,
double intercept0,
double frequency0,
Index gw_level,
215 double left_limit,
double right_limit,
bool allow_rejected_return =
true,
216 QPStats* stats =
nullptr)
const;
219 double frequency0,
Index gw_level,
220 QPStats* stats =
nullptr)
const;
225 QPStats* stats =
nullptr)
const;
227 bool Converged(
const Eigen::VectorXd& e1,
const Eigen::VectorXd& e2,
228 double epsilon)
const;
231 double lowerbound,
double f_lowerbound,
double upperbound,
232 double f_upperbound,
const QPFunc& f,
double reference)
const;
void CountSigmaStage(EvalStage stage) const
std::pair< double, double > operator()(double frequency) const
std::unordered_set< std::uint64_t > seen_frequencies_
double sigma(double frequency, EvalStage stage=EvalStage::Other) const
const QPStats & GetStats() const
double deriv(double frequency) const
double value(double frequency, EvalStage stage=EvalStage::Other) const
static std::uint64_t FrequencyKey(double x)
QPFunc(Index gw_level, const Sigma_base_UKS &sigma, double offset)
const Sigma_base_UKS & sigma_c_func_
boost::optional< double > SolveQP_Grid_Windowed(Spin spin, double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
TCMatrix_gwbse_spin & Mmn_
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianAlpha() const
bool Converged(const Eigen::VectorXd &e1, const Eigen::VectorXd &e2, double epsilon) const
void PrintQP_Energies(Spin spin, const Eigen::VectorXd &qp_diag_energies) const
Eigen::MatrixXd & SigmaC(Spin spin)
const Eigen::VectorXd & DftEnergies(Spin spin) const
std::unique_ptr< Sigma_base_UKS > sigma_alpha_
boost::optional< double > SolveQP_Linearisation(Spin spin, double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
std::string LevelLabel(Spin spin, Index level) const
Eigen::MatrixXd Sigma_x_alpha_
std::unique_ptr< Sigma_base_UKS > sigma_beta_
boost::optional< double > SolveQP_Grid_Windowed_Adaptive(Spin spin, 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(Spin spin, double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Eigen::MatrixXd Sigma_c_beta_
Eigen::MatrixXd & SigmaX(Spin spin)
Eigen::VectorXd getGWAResultsAlpha() const
boost::optional< QPRootCandidate > RefineQPInterval(double lowerbound, double f_lowerbound, double upperbound, double f_upperbound, const QPFunc &f, double reference) const
const Eigen::MatrixXd & vxc_alpha_
Eigen::VectorXd ScissorShift_DFTlevel(const Eigen::VectorXd &dft_energies, Index homo) const
const Eigen::VectorXd & RPAInputEnergiesAlpha() const
qp_solver::WindowDiagnostics QPWindowDiagnostics
const char * OccupationTag(Spin spin, Index level) const
Sigma_base_UKS & SigmaEvaluator(Spin spin)
Eigen::VectorXd SolveQP(Spin spin, const Eigen::VectorXd &frequencies) const
const Eigen::VectorXd & RPAInputEnergiesBeta() const
Eigen::MatrixXd getHQPBeta() const
void PrintGWA_Energies(Spin spin) const
Eigen::VectorXd getGWAResultsBeta() const
Eigen::MatrixXd getHQPAlpha() const
Eigen::MatrixXd Sigma_x_beta_
const char * SpinName(Spin spin) const
GW_UKS(Logger &log, TCMatrix_gwbse_spin &Mmn, const Eigen::MatrixXd &vxc_alpha, const Eigen::MatrixXd &vxc_beta, const Eigen::VectorXd &dft_energies_alpha, const Eigen::VectorXd &dft_energies_beta)
const Eigen::VectorXd & dft_energies_alpha_
double CalcSpinHomoLumoShift(const Eigen::VectorXd &frequencies, Spin spin) const
boost::optional< double > SolveQP_Grid_Windowed_Dense(Spin spin, double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
const Eigen::MatrixXd & vxc_beta_
const Eigen::VectorXd & dft_energies_beta_
qp_solver::EvalStage EvalStage
qp_solver::RootCandidate QPRootCandidate
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianBeta() const
Eigen::MatrixXd Sigma_c_alpha_
const Eigen::MatrixXd & Vxc(Spin spin) const
void configure(const options &opt)
void CalculateGWPerturbation()
boost::optional< double > SolveQP_FixedPoint(Spin spin, double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Index Homo(Spin spin) const
Logger is used for thread-safe output of messages.
Unrestricted RPA helper for spin-resolved GW screening.
Provides a means for comparing floating point numbers.
Index qp_adaptive_shell_count
double qp_virtual_min_energy
std::string quadrature_scheme
double qp_full_window_half_width
std::string qp_grid_search_mode
std::string sigma_integration
Index g_sc_max_iterations
std::string qp_root_finder
Index gw_sc_max_iterations
double qp_adaptive_shell_width