21#ifndef VOTCA_XTP_DFTENGINE_H
22#define VOTCA_XTP_DFTENGINE_H
42class DFTEngineTestAccess;
65 std::vector<std::unique_ptr<StaticSite> >* externalsites) {
132 void PrintMOsUKS(
const Eigen::VectorXd& alpha_energies,
133 const Eigen::VectorXd& beta_energies,
140 std::array<Eigen::MatrixXd, 2>
CalcERIs_EXX(
const Eigen::MatrixXd& MOCoeff,
141 const Eigen::MatrixXd& Dmat,
145 Eigen::MatrixXd
CalcERIs(
const Eigen::MatrixXd& Dmat,
double error)
const;
162 const std::vector<std::unique_ptr<StaticSite> >& multipoles)
const;
204 const std::vector<std::unique_ptr<StaticSite> >& multipoles)
const;
208 const AOBasis& dftbasis)
const;
214 Eigen::MatrixXd InitialActiveDensityMatrix,
215 Eigen::MatrixXd v_embedding,
216 Eigen::MatrixXd InitialInactiveMOs);
221 std::vector<Index> numfuncpatom);
225 Index numofzerocols);
228 Index numofzerorows);
284 QMMolecule(
"molecule made of atoms participating in Active region", 1);
Container to hold Basisfunctions for all atoms.
Electronic ground-state via Density-Functional Theory.
Eigen::MatrixXd CalcERIs(const Eigen::MatrixXd &Dmat, double error) const
Build the Coulomb matrix contribution from the current density matrix.
void setExternalcharges(std::vector< std::unique_ptr< StaticSite > > *externalsites)
bool EvaluateTruncatedActiveRegion(Orbitals &trunc_orb)
std::string auxbasis_name_
std::string getDFTBasisName() const
Return the configured AO basis-set name for the DFT calculation.
Index num_alpha_electrons_
tools::EigenSystem ModelPotentialGuess(const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const
Mat_p_Energy IntegrateExternalMultipoles(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
friend class DFTEngineTestAccess
tools::EigenSystem IndependentElectronGuess(const Mat_p_Energy &H0) const
Generate an initial guess by diagonalizing the core Hamiltonian only.
Eigen::MatrixXd InitialActiveDmat_trunc_
void TruncMOsFullBasis(Orbitals &orb, std::vector< Index > activeatoms, std::vector< Index > numfuncpatom)
Index num_beta_electrons_
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.
void setLogger(Logger *pLog)
Attach the logger used for SCF progress and diagnostics.
void TruncateBasis(Orbitals &orb, std::vector< Index > &activeatoms, Mat_p_Energy &H0, Eigen::MatrixXd InitialActiveDensityMatrix, Eigen::MatrixXd v_embedding, Eigen::MatrixXd InitialInactiveMOs)
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.
std::string dftbasis_name_
bool EvaluateUKS(Orbitals &orb, const Mat_p_Energy &H0, const Vxc_Potential< Vxc_Grid > &vxcpotential)
double truncation_threshold_
Mat_p_Energy IntegrateExternalDensity(const QMMolecule &mol, const Orbitals &extdensity) const
bool EvaluateActiveRegion(Orbitals &orb)
Eigen::MatrixXd RunAtomicDFT_unrestricted(const QMAtom &uniqueAtom) const
void Prepare(Orbitals &orb, Index numofelectrons=-1)
std::string initial_guess_
Eigen::MatrixXd AtomicGuess(const QMMolecule &mol) const
Build an atomic-density based initial guess in the AO basis.
Eigen::MatrixXd InsertZeroRows(Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerorows)
Insert zero rows into an MO coefficient matrix at the requested position.
Eigen::MatrixXd BuildEHTHamiltonian(const QMMolecule &mol) const
Build the extended-Hückel Hamiltonian for the current molecule.
Eigen::MatrixXd IntegrateExternalField(const QMMolecule &mol) const
Integrate a homogeneous external electric field into the AO basis.
Eigen::Vector3d extfield_
std::array< Eigen::MatrixXd, 2 > CalcERIs_EXX(const Eigen::MatrixXd &MOCoeff, const Eigen::MatrixXd &Dmat, double error) const
Eigen::MatrixXd H0_trunc_
Index NumberOfRestrictedOccupiedOrbitals() const
void Initialize(tools::Property &options)
Read DFT, grid, and SCF settings from the user options tree.
std::string xc_functional_name_
void CalcElDipole(const Orbitals &orb) const
Evaluate and print the electronic dipole moment from the final density.
double ExternalRepulsion(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
ConvergenceAcc::options conv_opt_
void SetupInvariantMatrices()
Precompute AO matrices that remain unchanged throughout the SCF procedure.
tools::EigenSystem ExtendedHuckelDFTGuess(const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const
std::string active_atoms_as_string_
void ConfigOrbfile(Orbitals &orb)
Propagate basis-set, XC, and metadata settings into the orbital container.
Eigen::VectorXd BuildEHTOrbitalEnergies(const QMMolecule &mol) const
Build orbital energies used in the extended-Hückel starting guess.
bool integrate_ext_density_
void PrintMOs(const Eigen::VectorXd &MOEnergies, Log::Level level)
Print a one-spin list of orbital energies and occupations to the logger.
std::vector< Index > numfuncpatom_
std::vector< Index > active_and_border_atoms_
Mat_p_Energy SetupH0(const QMMolecule &mol) const
Assemble the one-electron core Hamiltonian for the current molecule.
tools::EigenSystem ExtendedHuckelGuess(const QMMolecule &mol) const
bool IsRestrictedOpenShell() const
Eigen::MatrixXd OrthogonalizeGuess(const Eigen::MatrixXd &GuessMOs) const
Orthonormalize an initial MO guess with respect to the AO overlap matrix.
Eigen::MatrixXd v_embedding_trunc_
Eigen::MatrixXd SphericalAverageShells(const Eigen::MatrixXd &dmat, const AOBasis &dftbasis) const
ConvergenceAcc::options BuildConvergenceOptions() const
bool Evaluate(Orbitals &orb)
SpinDensity BuildSpinDensity(const tools::EigenSystem &MOs) const
std::vector< std::unique_ptr< StaticSite > > * externalsites_
Eigen::MatrixXd McWeenyPurification(Eigen::MatrixXd &Dmat_in, AOOverlap &overlap)
Vxc_Potential< Vxc_Grid > SetupVxc(const QMMolecule &mol)
bool integrate_ext_field_
double NuclearRepulsion(const QMMolecule &mol) const
Compute the classical nucleus-nucleus repulsion energy.
ConvergenceAcc conv_accelerator_
Eigen::MatrixXd InsertZeroCols(Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerocols)
Container to hold ECPs for all atoms.
Takes a density matrix and and an auxiliary basis set and calculates the electron repulsion integrals...
Logger is used for thread-safe output of messages.
Container for molecular orbitals and derived one-particle data.
Provides a means for comparing floating point numbers.
Eigen::MatrixXd spin() const
Return the spin density P^alpha - P^beta.
Eigen::MatrixXd total() const
Return the total density P = P^alpha + P^beta.