|
votca 2026-dev
|
Electronic ground-state via Density-Functional Theory. More...
#include <dftengine.h>

Classes | |
| struct | SpinDensity |
Public Member Functions | |
| void | Initialize (tools::Property &options) |
| Read DFT, grid, and SCF settings from the user options tree. | |
| void | setLogger (Logger *pLog) |
| Attach the logger used for SCF progress and diagnostics. | |
| void | setExternalcharges (std::vector< std::unique_ptr< StaticSite > > *externalsites) |
| bool | Evaluate (Orbitals &orb) |
| bool | EvaluateActiveRegion (Orbitals &orb) |
| bool | EvaluateTruncatedActiveRegion (Orbitals &trunc_orb) |
| std::string | getDFTBasisName () const |
| Return the configured AO basis-set name for the DFT calculation. | |
| bool | IsRestrictedOpenShell () const |
| Index | NumberOfRestrictedOccupiedOrbitals () const |
| SpinDensity | BuildSpinDensity (const tools::EigenSystem &MOs) const |
| ConvergenceAcc::options | BuildConvergenceOptions () const |
Private Member Functions | |
| void | Prepare (Orbitals &orb, Index numofelectrons=-1) |
| Vxc_Potential< Vxc_Grid > | SetupVxc (const QMMolecule &mol) |
| Eigen::MatrixXd | OrthogonalizeGuess (const Eigen::MatrixXd &GuessMOs) const |
| Orthonormalize an initial MO guess with respect to the AO overlap matrix. | |
| void | PrintMOs (const Eigen::VectorXd &MOEnergies, Log::Level level) |
| Print a one-spin list of orbital energies and occupations to the logger. | |
| void | PrintMOsUKS (const Eigen::VectorXd &alpha_energies, const Eigen::VectorXd &beta_energies, Log::Level level) const |
| Print separate alpha and beta orbital energies for a UKS calculation. | |
| void | CalcElDipole (const Orbitals &orb) const |
| Evaluate and print the electronic dipole moment from the final density. | |
| std::array< Eigen::MatrixXd, 2 > | CalcERIs_EXX (const Eigen::MatrixXd &MOCoeff, const Eigen::MatrixXd &Dmat, double error) const |
| Eigen::MatrixXd | CalcERIs (const Eigen::MatrixXd &Dmat, double error) const |
| Build the Coulomb matrix contribution from the current density matrix. | |
| void | ConfigOrbfile (Orbitals &orb) |
| Propagate basis-set, XC, and metadata settings into the orbital container. | |
| void | SetupInvariantMatrices () |
| Precompute AO matrices that remain unchanged throughout the SCF procedure. | |
| Eigen::MatrixXd | McWeenyPurification (Eigen::MatrixXd &Dmat_in, AOOverlap &overlap) |
| Mat_p_Energy | SetupH0 (const QMMolecule &mol) const |
| Assemble the one-electron core Hamiltonian for the current molecule. | |
| Mat_p_Energy | IntegrateExternalMultipoles (const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const |
| Mat_p_Energy | IntegrateExternalDensity (const QMMolecule &mol, const Orbitals &extdensity) const |
| Eigen::MatrixXd | IntegrateExternalField (const QMMolecule &mol) const |
| Integrate a homogeneous external electric field into the AO basis. | |
| tools::EigenSystem | IndependentElectronGuess (const Mat_p_Energy &H0) const |
| Generate an initial guess by diagonalizing the core Hamiltonian only. | |
| tools::EigenSystem | ModelPotentialGuess (const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const |
| Eigen::MatrixXd | AtomicGuess (const QMMolecule &mol) const |
| Build an atomic-density based initial guess in the AO basis. | |
| Eigen::VectorXd | BuildEHTOrbitalEnergies (const QMMolecule &mol) const |
| Build orbital energies used in the extended-Hückel starting guess. | |
| Eigen::MatrixXd | BuildEHTHamiltonian (const QMMolecule &mol) const |
| Build the extended-Hückel Hamiltonian for the current molecule. | |
| tools::EigenSystem | ExtendedHuckelGuess (const QMMolecule &mol) const |
| tools::EigenSystem | ExtendedHuckelDFTGuess (const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const |
| Eigen::MatrixXd | RunAtomicDFT_unrestricted (const QMAtom &uniqueAtom) const |
| double | NuclearRepulsion (const QMMolecule &mol) const |
| Compute the classical nucleus-nucleus repulsion energy. | |
| double | ExternalRepulsion (const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const |
| Eigen::MatrixXd | SphericalAverageShells (const Eigen::MatrixXd &dmat, const AOBasis &dftbasis) const |
| void | TruncateBasis (Orbitals &orb, std::vector< Index > &activeatoms, Mat_p_Energy &H0, Eigen::MatrixXd InitialActiveDensityMatrix, Eigen::MatrixXd v_embedding, Eigen::MatrixXd InitialInactiveMOs) |
| void | TruncMOsFullBasis (Orbitals &orb, std::vector< Index > activeatoms, std::vector< Index > numfuncpatom) |
| Eigen::MatrixXd | InsertZeroCols (Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerocols) |
| Eigen::MatrixXd | InsertZeroRows (Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerorows) |
| Insert zero rows into an MO coefficient matrix at the requested position. | |
| bool | EvaluateClosedShell (Orbitals &orb, const Mat_p_Energy &H0, const Vxc_Potential< Vxc_Grid > &vxcpotential) |
| Run the restricted closed-shell SCF loop and store the converged result. | |
| bool | EvaluateUKS (Orbitals &orb, const Mat_p_Energy &H0, const Vxc_Potential< Vxc_Grid > &vxcpotential) |
Friends | |
| class | DFTEngineTestAccess |
Electronic ground-state via Density-Functional Theory.
This class assembles the one- and two-electron matrix contributions needed for self-consistent Kohn-Sham calculations in a Gaussian AO basis. The SCF machinery supports restricted closed-shell calculations and unrestricted spin-polarized calculations, while reusing the same integral and numerical XC infrastructure whenever possible.
Definition at line 54 of file dftengine.h.
|
private |
Build an atomic-density based initial guess in the AO basis.
Definition at line 1077 of file dftengine.cc.
| ConvergenceAcc::options votca::xtp::DFTEngine::BuildConvergenceOptions | ( | ) | const |
Assemble SCF acceleration settings consistent with the current spin treatment and occupation model.
|
private |
Build the extended-Hückel Hamiltonian for the current molecule.
Definition at line 1512 of file dftengine.cc.
|
private |
Build orbital energies used in the extended-Hückel starting guess.
Definition at line 1484 of file dftengine.cc.
| SpinDensity votca::xtp::DFTEngine::BuildSpinDensity | ( | const tools::EigenSystem & | MOs | ) | const |
Construct alpha and beta density matrices from a shared MO coefficient matrix and the current occupation metadata.
|
private |
Evaluate and print the electronic dipole moment from the final density.
Definition at line 228 of file dftengine.cc.
|
private |
Build the Coulomb matrix contribution from the current density matrix.
Definition at line 262 of file dftengine.cc.
|
private |
Build Coulomb and exact-exchange matrix contributions from the current MO coefficients and density.
Definition at line 245 of file dftengine.cc.
|
private |
Propagate basis-set, XC, and metadata settings into the orbital container.
Definition at line 1118 of file dftengine.cc.
| bool votca::xtp::DFTEngine::Evaluate | ( | Orbitals & | orb | ) |
Run a full ground-state DFT calculation and store the results in the orbital container.
Definition at line 303 of file dftengine.cc.
| bool votca::xtp::DFTEngine::EvaluateActiveRegion | ( | Orbitals & | orb | ) |
Run an embedded active-region DFT calculation for the supplied orbital container.
Definition at line 63 of file embeddingengine.cc.
|
private |
Run the restricted closed-shell SCF loop and store the converged result.
Definition at line 321 of file dftengine.cc.
| bool votca::xtp::DFTEngine::EvaluateTruncatedActiveRegion | ( | Orbitals & | trunc_orb | ) |
Run the truncated active-region workflow used for reduced embedded calculations.
Definition at line 353 of file embeddingengine.cc.
|
private |
Run the unrestricted Kohn-Sham SCF loop and store alpha and beta orbitals separately.
Definition at line 516 of file dftengine.cc.
|
private |
Generate an extended-Hückel based guess refined with the one-electron DFT Hamiltonian.
Definition at line 1545 of file dftengine.cc.
|
private |
Generate an initial guess by diagonalizing the extended-Hückel Hamiltonian.
Definition at line 1532 of file dftengine.cc.
|
private |
Compute the classical interaction energy between nuclei and external multipoles.
Definition at line 1378 of file dftengine.cc.
|
inline |
Return the configured AO basis-set name for the DFT calculation.
Definition at line 81 of file dftengine.h.
|
private |
Generate an initial guess by diagonalizing the core Hamiltonian only.
Definition at line 271 of file dftengine.cc.
| void votca::xtp::DFTEngine::Initialize | ( | tools::Property & | options | ) |
Read DFT, grid, and SCF settings from the user options tree.
Self-consistent Kohn-Sham implementation.
The SCF cycle solves F C = S C eps in a Gaussian AO basis. In the restricted branch a single density matrix is iterated, whereas in the UKS branch separate alpha and beta densities are propagated while sharing the same one-electron Hamiltonian, Coulomb term, and AO overlap matrix.
Relative to the earlier restricted implementation, the UKS extension keeps the spin channels separate only where the equations require it: exchange, spin-resolved XC potentials, occupations, and convergence acceleration.
Definition at line 60 of file dftengine.cc.
|
private |
Insert zero columns into an MO coefficient matrix at the requested position.
|
private |
Insert zero rows into an MO coefficient matrix at the requested position.
Definition at line 687 of file embeddingengine.cc.
|
private |
Integrate an external electron density represented by another orbital container.
Definition at line 1434 of file dftengine.cc.
|
private |
Integrate a homogeneous external electric field into the AO basis.
Definition at line 1404 of file dftengine.cc.
|
private |
Integrate the electrostatic potential generated by external multipoles into the AO basis.
Definition at line 1417 of file dftengine.cc.
|
inline |
Report whether the current electron counts define a spin-polarized reference.
Definition at line 98 of file dftengine.h.
|
private |
Apply McWeeny purification to improve the idempotency of a density-matrix guess.
Definition at line 649 of file embeddingengine.cc.
|
private |
Generate an initial guess from a model potential including numerical XC contributions.
Definition at line 282 of file dftengine.cc.
|
private |
Compute the classical nucleus-nucleus repulsion energy.
Definition at line 1331 of file dftengine.cc.
|
inline |
Return the number of spatial orbitals occupied by at least one electron in a restricted open-shell reference.
Definition at line 104 of file dftengine.h.
|
private |
Orthonormalize an initial MO guess with respect to the AO overlap matrix.
Definition at line 1472 of file dftengine.cc.
Initialize basis sets, integral engines, and electron counts before entering the SCF loop.
Definition at line 1189 of file dftengine.cc.
|
private |
Print a one-spin list of orbital energies and occupations to the logger.
Definition at line 148 of file dftengine.cc.
|
private |
Print separate alpha and beta orbital energies for a UKS calculation.
Definition at line 168 of file dftengine.cc.
|
private |
Run an unrestricted atomic reference calculation used in open-shell atomic guesses.
Definition at line 893 of file dftengine.cc.
|
inline |
Provide external static sites whose electrostatic potential enters the Hamiltonian.
Definition at line 64 of file dftengine.h.
|
inline |
Attach the logger used for SCF progress and diagnostics.
Definition at line 60 of file dftengine.h.
|
private |
Assemble the one-electron core Hamiltonian for the current molecule.
Definition at line 745 of file dftengine.cc.
|
private |
Precompute AO matrices that remain unchanged throughout the SCF procedure.
Definition at line 855 of file dftengine.cc.
|
private |
Build the numerical exchange-correlation integration object for the current molecule.
Definition at line 1310 of file dftengine.cc.
|
private |
Average density-matrix elements over functions belonging to the same atomic shell.
Definition at line 1347 of file dftengine.cc.
|
private |
Project the full-system Hamiltonian and densities onto the selected truncated active basis.
Definition at line 496 of file embeddingengine.cc.
|
private |
Expand truncated active-region orbitals back into the full AO basis with zero padding.
Definition at line 670 of file embeddingengine.cc.
|
friend |
Definition at line 117 of file dftengine.h.
|
private |
Definition at line 302 of file dftengine.h.
|
private |
Definition at line 289 of file dftengine.h.
|
private |
Definition at line 298 of file dftengine.h.
|
private |
Definition at line 290 of file dftengine.h.
|
private |
Definition at line 283 of file dftengine.h.
|
private |
Definition at line 246 of file dftengine.h.
|
private |
Definition at line 242 of file dftengine.h.
|
private |
Definition at line 311 of file dftengine.h.
|
private |
Definition at line 266 of file dftengine.h.
|
private |
Definition at line 264 of file dftengine.h.
|
private |
Definition at line 257 of file dftengine.h.
|
private |
Definition at line 245 of file dftengine.h.
|
private |
Definition at line 243 of file dftengine.h.
|
private |
Definition at line 300 of file dftengine.h.
|
private |
Definition at line 247 of file dftengine.h.
|
private |
Definition at line 244 of file dftengine.h.
|
private |
Definition at line 268 of file dftengine.h.
|
private |
Definition at line 271 of file dftengine.h.
|
private |
Definition at line 286 of file dftengine.h.
|
private |
Definition at line 249 of file dftengine.h.
|
private |
Definition at line 254 of file dftengine.h.
|
private |
Definition at line 280 of file dftengine.h.
|
private |
Definition at line 294 of file dftengine.h.
|
private |
Definition at line 259 of file dftengine.h.
|
private |
Definition at line 295 of file dftengine.h.
|
private |
Definition at line 277 of file dftengine.h.
|
private |
Definition at line 287 of file dftengine.h.
|
private |
Definition at line 291 of file dftengine.h.
|
private |
Definition at line 263 of file dftengine.h.
|
private |
Definition at line 306 of file dftengine.h.
|
private |
Definition at line 307 of file dftengine.h.
|
private |
Definition at line 308 of file dftengine.h.
|
private |
Definition at line 309 of file dftengine.h.
|
private |
Definition at line 303 of file dftengine.h.
|
private |
Definition at line 262 of file dftengine.h.
|
private |
Definition at line 279 of file dftengine.h.
|
private |
Definition at line 239 of file dftengine.h.
|
private |
Definition at line 274 of file dftengine.h.
|
private |
Definition at line 251 of file dftengine.h.
|
private |
Definition at line 310 of file dftengine.h.
|
private |
Definition at line 281 of file dftengine.h.
|
private |
Definition at line 299 of file dftengine.h.
|
private |
Definition at line 297 of file dftengine.h.
|
private |
Definition at line 301 of file dftengine.h.
|
private |
Definition at line 296 of file dftengine.h.
|
private |
Definition at line 275 of file dftengine.h.