votca 2026-dev
Loading...
Searching...
No Matches
dftengine.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2023 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20#pragma once
21#ifndef VOTCA_XTP_DFTENGINE_H
22#define VOTCA_XTP_DFTENGINE_H
23
24// VOTCA includes
26
27// Local VOTCA includes
28#include "ERIs.h"
29#include "convergenceacc.h"
30#include "uks_convergenceacc.h"
31
32#include "ecpaobasis.h"
33#include "extended_hueckel.h"
34#include "logger.h"
35#include "qmmolecule.h"
36#include "staticsite.h"
37#include "vxc_grid.h"
38#include "vxc_potential.h"
39namespace votca {
40namespace xtp {
41class Orbitals;
42class DFTEngineTestAccess;
43
53
54class DFTEngine {
55 public:
57 void Initialize(tools::Property& options);
58
60 void setLogger(Logger* pLog) { pLog_ = pLog; }
61
65 std::vector<std::unique_ptr<StaticSite> >* externalsites) {
66 externalsites_ = externalsites;
67 }
68
71 bool Evaluate(Orbitals& orb);
72
79
81 std::string getDFTBasisName() const { return dftbasis_name_; };
82
85 struct SpinDensity {
86 Eigen::MatrixXd alpha;
87 Eigen::MatrixXd beta;
88
90 Eigen::MatrixXd total() const { return alpha + beta; }
91
93 Eigen::MatrixXd spin() const { return alpha - beta; }
94 };
95
101
107
111
115
116 private:
118
121 void Prepare(Orbitals& orb, Index numofelectrons = -1);
122
126
128 Eigen::MatrixXd OrthogonalizeGuess(const Eigen::MatrixXd& GuessMOs) const;
130 void PrintMOs(const Eigen::VectorXd& MOEnergies, Log::Level level);
132 void PrintMOsUKS(const Eigen::VectorXd& alpha_energies,
133 const Eigen::VectorXd& beta_energies,
134 Log::Level level) const;
136 void CalcElDipole(const Orbitals& orb) const;
137
140 std::array<Eigen::MatrixXd, 2> CalcERIs_EXX(const Eigen::MatrixXd& MOCoeff,
141 const Eigen::MatrixXd& Dmat,
142 double error) const;
143
145 Eigen::MatrixXd CalcERIs(const Eigen::MatrixXd& Dmat, double error) const;
146
148 void ConfigOrbfile(Orbitals& orb);
153 Eigen::MatrixXd McWeenyPurification(Eigen::MatrixXd& Dmat_in,
154 AOOverlap& overlap);
155
157 Mat_p_Energy SetupH0(const QMMolecule& mol) const;
161 const QMMolecule& mol,
162 const std::vector<std::unique_ptr<StaticSite> >& multipoles) const;
166 const Orbitals& extdensity) const;
167
169 Eigen::MatrixXd IntegrateExternalField(const QMMolecule& mol) const;
170
176 const Mat_p_Energy& H0, const QMMolecule& mol,
177 const Vxc_Potential<Vxc_Grid>& vxcpotential) const;
178
180 Eigen::MatrixXd AtomicGuess(const QMMolecule& mol) const;
181
183 Eigen::VectorXd BuildEHTOrbitalEnergies(const QMMolecule& mol) const;
185 Eigen::MatrixXd BuildEHTHamiltonian(const QMMolecule& mol) const;
192 const Mat_p_Energy& H0, const QMMolecule& mol,
193 const Vxc_Potential<Vxc_Grid>& vxcpotential) const;
196 Eigen::MatrixXd RunAtomicDFT_unrestricted(const QMAtom& uniqueAtom) const;
197
199 double NuclearRepulsion(const QMMolecule& mol) const;
202 double ExternalRepulsion(
203 const QMMolecule& mol,
204 const std::vector<std::unique_ptr<StaticSite> >& multipoles) const;
207 Eigen::MatrixXd SphericalAverageShells(const Eigen::MatrixXd& dmat,
208 const AOBasis& dftbasis) const;
209
212 void TruncateBasis(Orbitals& orb, std::vector<Index>& activeatoms,
213 Mat_p_Energy& H0,
214 Eigen::MatrixXd InitialActiveDensityMatrix,
215 Eigen::MatrixXd v_embedding,
216 Eigen::MatrixXd InitialInactiveMOs);
217
220 void TruncMOsFullBasis(Orbitals& orb, std::vector<Index> activeatoms,
221 std::vector<Index> numfuncpatom);
224 Eigen::MatrixXd InsertZeroCols(Eigen::MatrixXd MOsMatrix, Index startidx,
225 Index numofzerocols);
227 Eigen::MatrixXd InsertZeroRows(Eigen::MatrixXd MOsMatrix, Index startidx,
228 Index numofzerorows);
229
231 bool EvaluateClosedShell(Orbitals& orb, const Mat_p_Energy& H0,
232 const Vxc_Potential<Vxc_Grid>& vxcpotential);
233
236 bool EvaluateUKS(Orbitals& orb, const Mat_p_Energy& H0,
237 const Vxc_Potential<Vxc_Grid>& vxcpotential);
238
240
241 // basis sets
242 std::string auxbasis_name_;
243 std::string dftbasis_name_;
244 std::string ecp_name_;
248
250 // Pre-screening
252
253 // numerical integration Vxc
254 std::string grid_name_;
255
256 // AO Matrices
258
259 std::string initial_guess_;
260
261 // Convergence
265 // DIIS variables
267 // Electron repulsion integrals
269
270 // external charges
271 std::vector<std::unique_ptr<StaticSite> >* externalsites_ = nullptr;
272
273 // exchange and correlation
274 double ScaHFX_;
276
278 // integrate external density
279 std::string orbfilename_;
280 std::string gridquality_;
281 std::string state_;
282
284 QMMolecule("molecule made of atoms participating in Active region", 1);
285
286 Eigen::Vector3d extfield_ = Eigen::Vector3d::Zero();
288
292
293 // truncation
294 Eigen::MatrixXd H0_trunc_;
296 Eigen::MatrixXd v_embedding_trunc_;
300 double E_nuc_;
302 std::vector<Index> active_and_border_atoms_;
303 std::vector<Index> numfuncpatom_;
304
305 // Spin-DFT Extension
312};
313
314} // namespace xtp
315} // namespace votca
316
317#endif // VOTCA_XTP_DFTENGINE_H
class to manage program options with xml serialization functionality
Definition property.h:55
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Electronic ground-state via Density-Functional Theory.
Definition dftengine.h:54
Eigen::MatrixXd CalcERIs(const Eigen::MatrixXd &Dmat, double error) const
Build the Coulomb matrix contribution from the current density matrix.
Definition dftengine.cc:262
void setExternalcharges(std::vector< std::unique_ptr< StaticSite > > *externalsites)
Definition dftengine.h:64
bool EvaluateTruncatedActiveRegion(Orbitals &trunc_orb)
std::string auxbasis_name_
Definition dftengine.h:242
std::string getDFTBasisName() const
Return the configured AO basis-set name for the DFT calculation.
Definition dftengine.h:81
std::string gridquality_
Definition dftengine.h:280
tools::EigenSystem ModelPotentialGuess(const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const
Definition dftengine.cc:282
Mat_p_Energy IntegrateExternalMultipoles(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
friend class DFTEngineTestAccess
Definition dftengine.h:117
tools::EigenSystem IndependentElectronGuess(const Mat_p_Energy &H0) const
Generate an initial guess by diagonalizing the core Hamiltonian only.
Definition dftengine.cc:271
Eigen::MatrixXd InitialActiveDmat_trunc_
Definition dftengine.h:295
void TruncMOsFullBasis(Orbitals &orb, std::vector< Index > activeatoms, std::vector< Index > numfuncpatom)
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.
Definition dftengine.cc:321
void setLogger(Logger *pLog)
Attach the logger used for SCF progress and diagnostics.
Definition dftengine.h:60
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.
Definition dftengine.cc:168
std::string dftbasis_name_
Definition dftengine.h:243
bool EvaluateUKS(Orbitals &orb, const Mat_p_Energy &H0, const Vxc_Potential< Vxc_Grid > &vxcpotential)
Definition dftengine.cc:516
Mat_p_Energy IntegrateExternalDensity(const QMMolecule &mol, const Orbitals &extdensity) const
std::string ecp_name_
Definition dftengine.h:244
std::string grid_name_
Definition dftengine.h:254
bool EvaluateActiveRegion(Orbitals &orb)
Eigen::MatrixXd RunAtomicDFT_unrestricted(const QMAtom &uniqueAtom) const
Definition dftengine.cc:893
void Prepare(Orbitals &orb, Index numofelectrons=-1)
std::string initial_guess_
Definition dftengine.h:259
std::string orbfilename_
Definition dftengine.h:279
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_
Definition dftengine.h:286
std::array< Eigen::MatrixXd, 2 > CalcERIs_EXX(const Eigen::MatrixXd &MOCoeff, const Eigen::MatrixXd &Dmat, double error) const
Definition dftengine.cc:245
Eigen::MatrixXd H0_trunc_
Definition dftengine.h:294
Index NumberOfRestrictedOccupiedOrbitals() const
Definition dftengine.h:104
void Initialize(tools::Property &options)
Read DFT, grid, and SCF settings from the user options tree.
Definition dftengine.cc:60
std::string xc_functional_name_
Definition dftengine.h:275
void CalcElDipole(const Orbitals &orb) const
Evaluate and print the electronic dipole moment from the final density.
Definition dftengine.cc:228
double ExternalRepulsion(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
ConvergenceAcc::options conv_opt_
Definition dftengine.h:264
void SetupInvariantMatrices()
Precompute AO matrices that remain unchanged throughout the SCF procedure.
Definition dftengine.cc:855
tools::EigenSystem ExtendedHuckelDFTGuess(const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const
std::string active_atoms_as_string_
Definition dftengine.h:289
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.
QMMolecule activemol_
Definition dftengine.h:283
void PrintMOs(const Eigen::VectorXd &MOEnergies, Log::Level level)
Print a one-spin list of orbital energies and occupations to the logger.
Definition dftengine.cc:148
std::vector< Index > numfuncpatom_
Definition dftengine.h:303
AOOverlap dftAOoverlap_
Definition dftengine.h:257
std::vector< Index > active_and_border_atoms_
Definition dftengine.h:302
Mat_p_Energy SetupH0(const QMMolecule &mol) const
Assemble the one-electron core Hamiltonian for the current molecule.
Definition dftengine.cc:745
tools::EigenSystem ExtendedHuckelGuess(const QMMolecule &mol) const
bool IsRestrictedOpenShell() const
Definition dftengine.h:98
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_
Definition dftengine.h:296
Eigen::MatrixXd SphericalAverageShells(const Eigen::MatrixXd &dmat, const AOBasis &dftbasis) const
ConvergenceAcc::options BuildConvergenceOptions() const
bool Evaluate(Orbitals &orb)
Definition dftengine.cc:303
SpinDensity BuildSpinDensity(const tools::EigenSystem &MOs) const
std::vector< std::unique_ptr< StaticSite > > * externalsites_
Definition dftengine.h:271
Eigen::MatrixXd McWeenyPurification(Eigen::MatrixXd &Dmat_in, AOOverlap &overlap)
Vxc_Potential< Vxc_Grid > SetupVxc(const QMMolecule &mol)
double NuclearRepulsion(const QMMolecule &mol) const
Compute the classical nucleus-nucleus repulsion energy.
ConvergenceAcc conv_accelerator_
Definition dftengine.h:266
Eigen::MatrixXd InsertZeroCols(Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerocols)
Container to hold ECPs for all atoms.
Definition ecpaobasis.h:43
Takes a density matrix and and an auxiliary basis set and calculates the electron repulsion integrals...
Definition ERIs.h:35
Logger is used for thread-safe output of messages.
Definition logger.h:164
Container for molecular orbitals and derived one-particle data.
Definition orbitals.h:47
container for QM atoms
Definition qmatom.h:37
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Level
be loud and noisy
Definition globals.h:28
Eigen::MatrixXd spin() const
Return the spin density P^alpha - P^beta.
Definition dftengine.h:93
Eigen::MatrixXd total() const
Return the total density P = P^alpha + P^beta.
Definition dftengine.h:90