Exchange-correlation integrator#
For performing exchange-correlation integration, GauXC provides the GauXC::XCIntegrator interface.
The GauXC::XCIntegratorFactory creates integrators based on user-defined settings encapsulated in GauXC::IntegratorSettingsXC, GauXC::IntegratorSettingsKS, and GauXC::IntegratorSettingsEXC_GRAD.
-
template<typename MatrixType>
class XCIntegrator# High-level interface for exchange-correlation and exact exchange integration.
XCIntegrator provides methods to evaluate XC energies, potentials, gradients, and exact exchange contributions using numerical quadrature on a molecular grid.
- Template Parameters:
MatrixType – Matrix type for density and Fock-like matrices.
Public Types
-
using matrix_type = MatrixType#
Matrix type alias.
-
using value_type = typename matrix_type::value_type#
Scalar type alias.
-
using basisset_type = BasisSet<value_type>#
Basis set type alias.
-
using exc_vxc_type_rks = std::tuple<value_type, matrix_type>#
Return type for RKS exc+vxc: (energy, Vxc).
-
using exc_vxc_type_uks = std::tuple<value_type, matrix_type, matrix_type>#
Return type for UKS exc+vxc: (energy, Vxc_alpha, Vxc_beta).
-
using exc_vxc_type_gks = std::tuple<value_type, matrix_type, matrix_type, matrix_type, matrix_type>#
Return type for GKS exc+vxc: (energy, Vxc_scalar, Vxc_x, Vxc_y, Vxc_z).
-
using exc_grad_type = std::vector<value_type>#
Return type for XC gradient: vector of gradient components.
-
using exx_type = matrix_type#
Return type for exact exchange matrix.
-
using fxc_contraction_type_rks = matrix_type#
Return type for RKS Fxc contraction.
-
using fxc_contraction_type_uks = std::tuple<matrix_type, matrix_type>#
Return type for UKS Fxc contraction: (Fxc_alpha, Fxc_beta).
-
using dd_psi_type = std::vector<value_type>#
Return type for dd-Psi evaluation.
-
using dd_psi_potential_type = matrix_type#
Return type for dd-Psi potential.
Public Functions
-
XCIntegrator() = default#
Default constructor.
-
~XCIntegrator() noexcept#
Destructor.
-
XCIntegrator(std::unique_ptr<pimpl_type> &&pimpl)#
Construct from a pre-built implementation.
- Parameters:
pimpl – Unique pointer to implementation instance.
-
XCIntegrator(const XCIntegrator&) = delete#
Deleted copy constructor.
-
XCIntegrator(XCIntegrator&&) noexcept#
Move constructor.
-
value_type integrate_den(const MatrixType &P)#
Integrate the electron density.
- Parameters:
P – Density matrix.
- Returns:
Integrated number of electrons.
-
value_type eval_exc(const MatrixType &P, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate the XC energy (RKS).
- Parameters:
P – Density matrix.
opts – Integration settings.
- Returns:
XC energy.
-
value_type eval_exc(const MatrixType &Palpha, const MatrixType &Pbeta, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate the XC energy (UKS).
- Parameters:
Palpha – Alpha density matrix.
Pbeta – Beta density matrix.
opts – Integration settings.
- Returns:
XC energy.
-
value_type eval_exc(const MatrixType &Pscalar, const MatrixType &Px, const MatrixType &Py, const MatrixType &Pz, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate the XC energy (GKS).
- Parameters:
Pscalar – Scalar density matrix.
Px – X-component spin density matrix.
Py – Y-component spin density matrix.
Pz – Z-component spin density matrix.
opts – Integration settings.
- Returns:
XC energy.
-
exc_vxc_type_rks eval_exc_vxc(const MatrixType &P, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate XC energy and potential (RKS).
- Parameters:
P – Density matrix.
opts – Integration settings.
- Returns:
Tuple of (energy, Vxc).
-
exc_vxc_type_uks eval_exc_vxc(const MatrixType &Palpha, const MatrixType &Pbeta, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate XC energy and potential (UKS).
- Parameters:
Palpha – Alpha density matrix.
Pbeta – Beta density matrix.
opts – Integration settings.
- Returns:
Tuple of (energy, Vxc_alpha, Vxc_beta).
-
exc_vxc_type_gks eval_exc_vxc(const MatrixType &Pscalar, const MatrixType &Px, const MatrixType &Py, const MatrixType &Pz, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate XC energy and potential (GKS).
- Parameters:
Pscalar – Scalar density matrix.
Px – X-component spin density matrix.
Py – Y-component spin density matrix.
Pz – Z-component spin density matrix.
opts – Integration settings.
- Returns:
Tuple of (energy, Vxc_scalar, Vxc_x, Vxc_y, Vxc_z).
-
exc_grad_type eval_exc_grad(const MatrixType &P, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate XC gradient (RKS).
- Parameters:
P – Density matrix.
opts – Integration settings.
- Returns:
Vector of gradient components (3*natoms).
-
exc_grad_type eval_exc_grad(const MatrixType &Palpha, const MatrixType &Pbeta, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate XC gradient (UKS).
- Parameters:
Palpha – Alpha density matrix.
Pbeta – Beta density matrix.
opts – Integration settings.
- Returns:
Vector of gradient components (3*natoms).
-
exx_type eval_exx(const MatrixType &P, const IntegratorSettingsEXX &opts = IntegratorSettingsEXX{})#
Evaluate exact exchange matrix (sn-LinK).
- Parameters:
P – Density matrix.
opts – EXX integration settings.
- Returns:
Exact exchange matrix K.
-
fxc_contraction_type_rks eval_fxc_contraction(const MatrixType &P, const MatrixType &X, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate Fxc contraction (RKS).
- Parameters:
P – Density matrix.
X – Trial vector for contraction.
opts – Integration settings.
- Returns:
Contracted Fxc matrix.
-
fxc_contraction_type_uks eval_fxc_contraction(const MatrixType &Palpha, const MatrixType &Pbeta, const MatrixType &Xalpha, const MatrixType &Xbeta, const IntegratorSettingsXC &opts = IntegratorSettingsXC{})#
Evaluate Fxc contraction (UKS).
- Parameters:
Palpha – Alpha density matrix.
Pbeta – Beta density matrix.
Xalpha – Alpha trial vector.
Xbeta – Beta trial vector.
opts – Integration settings.
- Returns:
Tuple of (Fxc_alpha, Fxc_beta).
-
dd_psi_type eval_dd_psi(const MatrixType &P, unsigned deriv)#
Evaluate dd-Psi diagnostic.
- Parameters:
P – Density matrix.
deriv – Derivative order.
- Returns:
Vector of dd-Psi values.
-
dd_psi_potential_type eval_dd_psi_potential(const MatrixType &P, unsigned deriv)#
Evaluate dd-Psi potential.
- Parameters:
P – Density matrix.
deriv – Derivative order.
- Returns:
dd-Psi potential matrix.
-
const util::Timer &get_timings() const#
Get the internal timing tracker.
- Returns:
Const reference to timer.
-
const LoadBalancer &load_balancer() const#
Get the associated load balancer (const).
- Returns:
Const reference to LoadBalancer.
-
LoadBalancer &load_balancer()#
Get the associated load balancer (mutable).
- Returns:
Reference to LoadBalancer.
-
template<typename MatrixType>
class XCIntegratorFactory# Factory to generate XCIntegrator Instances.
Public Functions
-
inline XCIntegratorFactory(ExecutionSpace ex, std::string integrator_input_type, std::string integrator_kernel_name, std::string local_work_kernel_name, std::string reduction_kernel_name, LocalWorkSettings settings = LocalWorkSettings())#
Construct an XCIntegratorFactory instance
- Parameters:
ex – [in] Execution space for the XCIntegrator instance
integrator_input_type – [in] Input type for XC integration (e.g. “Replicated”)
integrator_kernel_name – [in] Name of Integraion scaffold kernel to load (e.g. “Reference” or “Default”)
local_work_kerenl_name – [in] Name of LWD to load (e.g. “Reference” or “Default”)
setting – [in] Settings to pass to LWD (not currently used)
Generate XCIntegrator instance
- Parameters:
func – [in] XC functional
lb – [in] Preconstructed Load Balancer instance
-
inline XCIntegratorFactory(ExecutionSpace ex, std::string integrator_input_type, std::string integrator_kernel_name, std::string local_work_kernel_name, std::string reduction_kernel_name, LocalWorkSettings settings = LocalWorkSettings())#
-
struct IntegratorSettingsXC#
Base settings for exchange-correlation (XC) integration.
Subclassed by GauXC::IntegratorSettingsKS
-
struct IntegratorSettingsKS : public GauXC::IntegratorSettingsXC#
Settings for Kohn-Sham DFT integration.
Controls tolerances for generalized Kohn-Sham (GKS) calculations.
Subclassed by GauXC::IntegratorSettingsEXC_GRAD
Public Members
-
double gks_dtol = 1e-12#
Tolerance for GKS density evaluation.
-
double gks_dtol = 1e-12#
-
struct IntegratorSettingsEXC_GRAD : public GauXC::IntegratorSettingsKS#
Settings for XC gradient evaluation.
Controls whether grid weight derivatives are included in the gradient computation or if only the Hellmann-Feynman contribution is used.
Public Members
-
bool include_weight_derivatives = true#
Include grid weight contribution and translational invariance.
-
bool include_weight_derivatives = true#
-
struct IntegratorSettingsEXX#
Base settings for exact exchange (EXX) integration.
Subclassed by GauXC::IntegratorSettingsSNLinK
-
struct IntegratorSettingsSNLinK : public GauXC::IntegratorSettingsEXX#
Settings for seminumerical linear exchange (sn-LinK) integration.
Controls screening and tolerance parameters for sn-LinK calculations.