Skip to content
This repository has been archived by the owner on Sep 28, 2021. It is now read-only.

Commit

Permalink
Merge pull request #648 from votca/RPA_correction_refactor
Browse files Browse the repository at this point in the history
refactoring energy corrections outside qpmin:qpmax for RPA
  • Loading branch information
Bjoern Baumeier authored Jan 19, 2021
2 parents 10542c6 + e4464ae commit e549fe0
Show file tree
Hide file tree
Showing 3 changed files with 42 additions and 19 deletions.
2 changes: 2 additions & 0 deletions CHANGELOG.rst
Original file line number Diff line number Diff line change
Expand Up @@ -5,8 +5,10 @@ Version 2022-dev
================

- Updated coordinate precision (#638)
- Refactored energy corrections in RPA outside QPs (#577)
- Made SetupCptTable static (#650)


Version 2021-rc.2 (released XX.01.21)
=====================================

Expand Down
8 changes: 7 additions & 1 deletion include/votca/xtp/rpa.h
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
/*
* Copyright 2009-2020 The VOTCA Development Team
* Copyright 2009-2021 The VOTCA Development Team
* (http://www.votca.org)
*
* Licensed under the Apache License, Version 2.0 (the "License")
Expand Down Expand Up @@ -92,6 +92,12 @@ class RPA {
Eigen::MatrixXd Calculate_H2p_ApB() const;
Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> Diagonalize_H2p_C(
const Eigen::MatrixXd& C) const;

void ShiftUncorrectedEnergies(const Eigen::VectorXd& dftenergies, Index qpmin,
Index gwsize);

double getMaxCorrection(const Eigen::VectorXd& dftenergies, Index min,
Index max) const;
};

} // namespace xtp
Expand Down
51 changes: 33 additions & 18 deletions src/libxtp/gwbse/rpa.cc
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
/*
* Copyright 2009-2020 The VOTCA Development Team
* Copyright 2009-2021 The VOTCA Development Team
* (http://www.votca.org)
*
* Licensed under the Apache License, Version 2.0 (the "License")
Expand Down Expand Up @@ -33,26 +33,41 @@ void RPA::UpdateRPAInputEnergies(const Eigen::VectorXd& dftenergies,
Index rpatotal = _rpamax - _rpamin + 1;
_energies = dftenergies.segment(_rpamin, rpatotal);
Index gwsize = Index(gwaenergies.size());
Index lumo = _homo + 1;

Index qpmax = qpmin + gwsize - 1;
_energies.segment(qpmin - _rpamin, gwsize) = gwaenergies;

Eigen::VectorXd corrections_occ =
_energies.segment(qpmin - _rpamin, lumo - qpmin) -
dftenergies.segment(qpmin - _rpamin, lumo - qpmin);
Eigen::VectorXd corrections_virt =
_energies.segment(lumo - qpmin, gwsize - (lumo - qpmin)) -
dftenergies.segment(lumo - qpmin, gwsize - (lumo - qpmin));
double max_correction_occ = (corrections_occ.cwiseAbs()).maxCoeff();
double max_correction_virt = (corrections_virt.cwiseAbs()).maxCoeff();

Index levelaboveqpmax = _rpamax - qpmax;
Index levelbelowqpmin = qpmin - _rpamin;

_energies.segment(0, levelbelowqpmin).array() -= max_correction_occ;
_energies.segment(qpmax + 1 - _rpamin, levelaboveqpmax).array() +=
max_correction_virt;
ShiftUncorrectedEnergies(dftenergies, qpmin, gwsize);
}

// Shifts energies of levels that are not QP corrected but
// used in the RPA:
// between rpamin and qpmin: by maximum abs of explicit QP corrections
// from qpmin to HOMO
// between qpmax and rpamax: by maximum abs of explicit QP corrections
// from LUMO to qpmax
void RPA::ShiftUncorrectedEnergies(const Eigen::VectorXd& dftenergies,
Index qpmin, Index gwsize) {

Index lumo = _homo + 1;
Index qpmax = qpmin + gwsize - 1;

// get max abs QP corrections for occupied/virtual levels
double max_correction_occ = getMaxCorrection(dftenergies, qpmin, _homo);
double max_correction_virt = getMaxCorrection(dftenergies, lumo, qpmax);

// shift energies
_energies.head(qpmin).array() -= max_correction_occ;
_energies.tail(_rpamax - qpmax).array() += max_correction_virt;
}

double RPA::getMaxCorrection(const Eigen::VectorXd& dftenergies, Index min,
Index max) const {

Index range = max - min + 1;
Eigen::VectorXd corrections =
_energies.segment(min, range) - dftenergies.segment(min - _rpamin, range);

return (corrections.cwiseAbs()).maxCoeff();
}

template <bool imag>
Expand Down

0 comments on commit e549fe0

Please sign in to comment.