votca 2026-dev
Loading...
Searching...
No Matches
dftengine.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2026 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// Third party includes
21#include <algorithm>
22#include <boost/filesystem.hpp>
23#include <boost/format.hpp>
24
25// VOTCA includes
28
29// Local VOTCA includes
33#include "votca/xtp/aomatrix.h"
36#include "votca/xtp/dftengine.h"
38#include "votca/xtp/logger.h"
39#include "votca/xtp/mmregion.h"
40#include "votca/xtp/orbitals.h"
43
44namespace votca {
45namespace xtp {
46
47namespace {
48
49void CanonicalizeOrbitalPhases(Eigen::MatrixXd& coeffs) {
50 constexpr double tol = 1e-14;
51
52 for (Index col = 0; col < coeffs.cols(); ++col) {
53 Eigen::Index pivot = 0;
54 const double maxabs = coeffs.col(col).cwiseAbs().maxCoeff(&pivot);
55
56 if (maxabs <= tol) {
57 continue;
58 }
59
60 if (coeffs(pivot, col) < 0.0) {
61 coeffs.col(col) *= -1.0;
62 }
63 }
64}
65
66void CanonicalizeOrbitalPhases(tools::EigenSystem& mos) {
67 CanonicalizeOrbitalPhases(mos.eigenvectors());
68}
69
70} // namespace
71
84
86
87 const std::string key_xtpdft = "xtpdft";
88 dftbasis_name_ = options.get(".basisset").as<std::string>();
89
90 if (options.exists(".auxbasisset")) {
91 auxbasis_name_ = options.get(".auxbasisset").as<std::string>();
92 }
93
94 if (!auxbasis_name_.empty()) {
95 screening_eps_ = options.get(key_xtpdft + ".screening_eps").as<double>();
97 options.get(key_xtpdft + ".fock_matrix_reset").as<Index>();
98 }
99 if (options.exists(".ecp")) {
100 ecp_name_ = options.get(".ecp").as<std::string>();
101 }
102
103 if (options.exists(key_xtpdft + ".force_uks_path")) {
104 force_uks_path_ = options.get(key_xtpdft + ".force_uks_path").as<bool>();
105 }
106
107 initial_guess_ = options.get(".initial_guess").as<std::string>();
108
109 grid_name_ = options.get(key_xtpdft + ".integration_grid").as<std::string>();
110 xc_functional_name_ = options.get(".functional").as<std::string>();
111
112 if (options.exists(key_xtpdft + ".externaldensity")) {
115 options.get(key_xtpdft + ".externaldensity.orbfile").as<std::string>();
116 gridquality_ = options.get(key_xtpdft + ".externaldensity.gridquality")
117 .as<std::string>();
118 state_ =
119 options.get(key_xtpdft + ".externaldensity.state").as<std::string>();
120 }
121
122 if (options.exists(".externalfield")) {
124 extfield_ = options.get(".externalfield").as<Eigen::Vector3d>();
125 }
126
127 conv_opt_.Econverged =
128 options.get(key_xtpdft + ".convergence.energy").as<double>();
129 conv_opt_.error_converged =
130 options.get(key_xtpdft + ".convergence.error").as<double>();
131 max_iter_ =
132 options.get(key_xtpdft + ".convergence.max_iterations").as<Index>();
133
134 std::string method =
135 options.get(key_xtpdft + ".convergence.method").as<std::string>();
136 if (method == "DIIS") {
137 conv_opt_.usediis = true;
138 } else if (method == "mixing") {
139 conv_opt_.usediis = false;
140 }
141 if (!conv_opt_.usediis) {
142 conv_opt_.histlength = 1;
143 conv_opt_.maxout = false;
144 }
145 conv_opt_.mixingparameter =
146 options.get(key_xtpdft + ".convergence.mixing").as<double>();
147 conv_opt_.levelshift =
148 options.get(key_xtpdft + ".convergence.levelshift").as<double>();
149 conv_opt_.levelshiftend =
150 options.get(key_xtpdft + ".convergence.levelshift_end").as<double>();
151 conv_opt_.maxout =
152 options.get(key_xtpdft + ".convergence.DIIS_maxout").as<bool>();
153 conv_opt_.histlength =
154 options.get(key_xtpdft + ".convergence.DIIS_length").as<Index>();
155 conv_opt_.diis_start =
156 options.get(key_xtpdft + ".convergence.DIIS_start").as<double>();
157 conv_opt_.adiis_start =
158 options.get(key_xtpdft + ".convergence.ADIIS_start").as<double>();
159
160 if (options.exists(key_xtpdft + ".dft_in_dft.activeatoms")) {
162 options.get(key_xtpdft + ".dft_in_dft.activeatoms").as<std::string>();
164 options.get(key_xtpdft + ".dft_in_dft.threshold").as<double>();
166 options.get(key_xtpdft + ".dft_in_dft.levelshift").as<double>();
167 truncate_ =
168 options.get(key_xtpdft + ".dft_in_dft.truncate_basis").as<bool>();
169 if (truncate_) {
171 options.get(key_xtpdft + ".dft_in_dft.truncation_threshold")
172 .as<double>();
173 }
174 }
175}
176
177void DFTEngine::PrintMOs(const Eigen::VectorXd& MOEnergies, Log::Level level) {
178 XTP_LOG(level, *pLog_) << " Orbital energies: " << std::flush;
179 XTP_LOG(level, *pLog_) << " index occupation energy(Hartree) " << std::flush;
180
181 for (Index i = 0; i < MOEnergies.size(); ++i) {
182 Index occupancy = 0;
183 if (i < num_docc_) {
184 occupancy = 2;
185 } else if (i < num_docc_ + num_socc_alpha_) {
186 occupancy = 1;
187 }
188
189 XTP_LOG(level, *pLog_) << (boost::format(" %1$5d %2$1d %3$+1.10f") %
190 i % occupancy % MOEnergies(i))
191 .str()
192 << std::flush;
193 }
194 return;
195}
196
197void DFTEngine::PrintMOsUKS(const Eigen::VectorXd& alpha_energies,
198 const Eigen::VectorXd& beta_energies,
199 Log::Level level) const {
200 XTP_LOG(level, *pLog_) << " UKS orbital energies:" << std::flush;
201 XTP_LOG(level, *pLog_) << " index occ eps_a(Ha) eps_b(Ha)"
202 << std::flush;
203
204 const Index nrows =
205 std::max<Index>(alpha_energies.size(), beta_energies.size());
206
207 for (Index i = 0; i < nrows; ++i) {
208 const bool occ_a = (i < num_alpha_electrons_);
209 const bool occ_b = (i < num_beta_electrons_);
210
211 std::string occ = "0";
212 if (occ_a && occ_b) {
213 occ = "2";
214 } else if (occ_a) {
215 occ = "a";
216 } else if (occ_b) {
217 occ = "b";
218 }
219
220 std::string eps_a = " -";
221 std::string eps_b = " -";
222
223 if (i < alpha_energies.size()) {
224 eps_a = (boost::format("%+1.10f") % alpha_energies(i)).str();
225 }
226 if (i < beta_energies.size()) {
227 eps_b = (boost::format("%+1.10f") % beta_energies(i)).str();
228 }
229
230 XTP_LOG(level, *pLog_) << (boost::format(
231 " %1$5d %2$1s %3$15s %4$15s") %
232 i % occ % eps_a % eps_b)
233 .str()
234 << std::flush;
235 }
236
237 if (num_alpha_electrons_ > 0 &&
238 num_alpha_electrons_ < alpha_energies.size()) {
239 XTP_LOG(level, *pLog_) << (boost::format(
240 " alpha HOMO-LUMO gap: %+1.10f Ha") %
241 (alpha_energies(num_alpha_electrons_) -
242 alpha_energies(num_alpha_electrons_ - 1)))
243 .str()
244 << std::flush;
245 }
246
247 if (num_beta_electrons_ > 0 && num_beta_electrons_ < beta_energies.size()) {
248 XTP_LOG(level, *pLog_) << (boost::format(
249 " beta HOMO-LUMO gap: %+1.10f Ha") %
250 (beta_energies(num_beta_electrons_) -
251 beta_energies(num_beta_electrons_ - 1)))
252 .str()
253 << std::flush;
254 }
255}
256
257void DFTEngine::CalcElDipole(const Orbitals& orb) const {
258 QMState state = QMState("n");
259 Eigen::Vector3d result = orb.CalcElDipole(state);
261 << TimeStamp() << " Electric Dipole is[e*bohr]:\n\t\t dx=" << result[0]
262 << "\n\t\t dy=" << result[1] << "\n\t\t dz=" << result[2] << std::flush;
263 return;
264}
265
266// Build the Coulomb and exact-exchange contributions generated by the current
267// density matrix. The returned pair is conventionally interpreted as
268//
269// (J[P], -K[P]),
270//
271// so that the hybrid Fock update becomes F = H0 + J[P] + a_x (-K[P]) + V_xc.
272// For RI/3c builds the occupied MO block is supplied when available to avoid an
273// unnecessary reconstruction of exchange intermediates.
274std::array<Eigen::MatrixXd, 2> DFTEngine::CalcERIs_EXX(
275 const Eigen::MatrixXd& MOCoeff, const Eigen::MatrixXd& Dmat,
276 double error) const {
277 if (!auxbasis_name_.empty()) {
278 if (conv_accelerator_.getUseMixing() || MOCoeff.rows() == 0) {
279 return ERIs_.CalculateERIs_EXX_3c(Eigen::MatrixXd::Zero(0, 0), Dmat);
280 } else {
281 Eigen::MatrixXd occblock = MOCoeff.leftCols(num_docc_ + num_socc_alpha_);
282 return ERIs_.CalculateERIs_EXX_3c(occblock, Dmat);
283 }
284 } else {
285 return ERIs_.CalculateERIs_EXX_4c(Dmat, error);
286 }
287}
288
289// Pure Coulomb contribution J[P] from the current AO density matrix. The code
290// dispatches to either RI/3c or conventional 4-center integral evaluation.
291Eigen::MatrixXd DFTEngine::CalcERIs(const Eigen::MatrixXd& Dmat,
292 double error) const {
293 if (!auxbasis_name_.empty()) {
294 return ERIs_.CalculateERIs_3c(Dmat);
295 } else {
296 return ERIs_.CalculateERIs_4c(Dmat, error);
297 }
298}
299
301 const Mat_p_Energy& H0) const {
302 return conv_accelerator_.SolveFockmatrix(H0.matrix());
303}
304
305// Construct a self-consistent model-potential guess by starting from an
306// atomic density P^(0), evaluating
307//
308// F[P^(0)] = H0 + J[P^(0)] + a_x (-K[P^(0)]) + V_xc[P^(0)],
309//
310// and diagonalizing the resulting Fock matrix once.
312 const Mat_p_Energy& H0, const QMMolecule& mol,
313 const Vxc_Potential<Vxc_Grid>& vxcpotential) const {
314 Eigen::MatrixXd Dmat = AtomicGuess(mol);
315 Mat_p_Energy e_vxc = vxcpotential.IntegrateVXC(Dmat);
317 << TimeStamp() << " Filled DFT Vxc matrix " << std::flush;
318
319 Eigen::MatrixXd H = H0.matrix() + e_vxc.matrix();
320
321 if (ScaHFX_ > 0) {
322 std::array<Eigen::MatrixXd, 2> both =
323 CalcERIs_EXX(Eigen::MatrixXd::Zero(0, 0), Dmat, 1e-12);
324 H += both[0];
325 H += ScaHFX_ * both[1];
326 } else {
327 H += CalcERIs(Dmat, 1e-12);
328 }
329 return conv_accelerator_.SolveFockmatrix(H);
330}
331
333 Prepare(orb);
334 Mat_p_Energy H0 = SetupH0(orb.QMAtoms());
335 Vxc_Potential<Vxc_Grid> vxcpotential = SetupVxc(orb.QMAtoms());
336 ConfigOrbfile(orb);
337
341 << TimeStamp()
342 << " Forcing closed-shell singlet through UKS development path."
343 << std::flush;
344 }
345 return EvaluateUKS(orb, H0, vxcpotential);
346 }
347 return EvaluateClosedShell(orb, H0, vxcpotential);
348}
349
350// Restricted SCF loop. The total energy is assembled as
351//
352// E = Tr[P H0] + E_nuc + E_coul + E_xc + E_exx,
353//
354// with P = 2 C_occ C_occ^T. DIIS or mixing updates the density until both
355// the energy change and the commutator error are converged.
357 Orbitals& orb, const Mat_p_Energy& H0,
358 const Vxc_Potential<Vxc_Grid>& vxcpotential) {
359
361 MOs.eigenvalues() = Eigen::VectorXd::Zero(H0.cols());
362 MOs.eigenvectors() = Eigen::MatrixXd::Zero(H0.rows(), H0.cols());
363
364 if (initial_guess_ == "orbfile") {
366 << TimeStamp() << " Reading guess from orbitals object/file"
367 << std::flush;
368 MOs = orb.MOs();
370 } else {
372 << TimeStamp() << " Setup Initial Guess using: " << initial_guess_
373 << std::flush;
374 if (initial_guess_ == "independent") {
375 MOs = IndependentElectronGuess(H0);
376 } else if (initial_guess_ == "atom") {
377 MOs = ModelPotentialGuess(H0, orb.QMAtoms(), vxcpotential);
378 } else if (initial_guess_ == "huckel") {
379 MOs = ExtendedHuckelGuess(orb.QMAtoms());
380 } else if (initial_guess_ == "huckel_dft") {
381 MOs = ExtendedHuckelDFTGuess(H0, orb.QMAtoms(), vxcpotential);
382 } else {
383 throw std::runtime_error("Initial guess method not known/implemented");
384 }
385 }
386
388 conv_accelerator_.DensityMatrixSpinResolved(MOs);
389 Eigen::MatrixXd Dmat = spin_dmat.total();
390
392 << TimeStamp() << " Guess Matrix gives N=" << std::setprecision(9)
393 << Dmat.cwiseProduct(dftAOoverlap_.Matrix()).sum() << " electrons."
394 << std::flush;
395
397 << TimeStamp() << " STARTING SCF cycle" << std::flush;
399 << " ----------------------------------------------"
400 "----------------------------"
401 << std::flush;
402
403 Eigen::MatrixXd J = Eigen::MatrixXd::Zero(Dmat.rows(), Dmat.cols());
404 Eigen::MatrixXd K;
405 if (ScaHFX_ > 0) {
406 K = Eigen::MatrixXd::Zero(Dmat.rows(), Dmat.cols());
407 }
408
409 double start_incremental_F_threshold = 1e-4;
410 if (!auxbasis_name_.empty()) {
411 start_incremental_F_threshold = 0.0; // Disable if RI is used
412 }
413 IncrementalFockBuilder incremental_fock(*pLog_, start_incremental_F_threshold,
415 incremental_fock.Configure(Dmat);
416
417 for (Index this_iter = 0; this_iter < max_iter_; this_iter++) {
418 XTP_LOG(Log::error, *pLog_) << std::flush;
419 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Iteration " << this_iter + 1
420 << " of " << max_iter_ << std::flush;
421
422 Mat_p_Energy e_vxc = vxcpotential.IntegrateVXC(Dmat);
424 << TimeStamp() << " Filled DFT Vxc matrix " << std::flush;
425
426 Eigen::MatrixXd H = H0.matrix() + e_vxc.matrix();
427 double Eone = Dmat.cwiseProduct(H0.matrix()).sum();
428 double Etwo = e_vxc.energy();
429 double exx = 0.0;
430
431 incremental_fock.Start(this_iter, conv_accelerator_.getDIIsError());
432 incremental_fock.resetMatrices(J, K, Dmat);
433 incremental_fock.UpdateCriteria(conv_accelerator_.getDIIsError(),
434 this_iter);
435
436 double integral_error =
437 std::min(conv_accelerator_.getDIIsError() * 1e-5, 1e-5);
438
439 if (ScaHFX_ > 0) {
440 std::array<Eigen::MatrixXd, 2> both = CalcERIs_EXX(
441 MOs.eigenvectors(), incremental_fock.getDmat_diff(), integral_error);
442 J += both[0];
443 H += J;
444 Etwo += 0.5 * Dmat.cwiseProduct(J).sum();
445 K += both[1];
446 H += 0.5 * ScaHFX_ * K;
447 exx = 0.25 * ScaHFX_ * Dmat.cwiseProduct(K).sum();
449 << TimeStamp() << " Filled F+K matrix " << std::flush;
450 } else {
451 J += CalcERIs(incremental_fock.getDmat_diff(), integral_error);
453 << TimeStamp() << " Filled F matrix " << std::flush;
454 H += J;
455 Etwo += 0.5 * Dmat.cwiseProduct(J).sum();
456 }
457
458 Etwo += exx;
459 double totenergy = Eone + H0.energy() + Etwo;
460
461 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " Single particle energy "
462 << std::setprecision(12) << Eone << std::flush;
463 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " Two particle energy "
464 << std::setprecision(12) << Etwo << std::flush;
466 << TimeStamp() << std::setprecision(12) << " Local Exc contribution "
467 << e_vxc.energy() << std::flush;
468 if (ScaHFX_ > 0) {
470 << TimeStamp() << std::setprecision(12)
471 << " Non local Ex contribution " << exx << std::flush;
472 }
474 << TimeStamp() << " Total Energy " << std::setprecision(12) << totenergy
475 << std::flush;
476
477 Dmat = conv_accelerator_.Iterate(Dmat, H, MOs, totenergy);
478 incremental_fock.UpdateDmats(Dmat, conv_accelerator_.getDIIsError(),
479 this_iter);
480
482
483 if (num_docc_ + num_socc_alpha_ > 0 &&
484 num_docc_ + num_socc_alpha_ < MOs.eigenvalues().size()) {
486 << "\t\tGAP "
489 << std::flush;
490 }
491
492 if (conv_accelerator_.isConverged()) {
494 << TimeStamp() << " Total Energy has converged to "
495 << std::setprecision(9) << conv_accelerator_.getDeltaE()
496 << "[Ha] after " << this_iter + 1
497 << " iterations. DIIS error is converged up to "
498 << conv_accelerator_.getDIIsError() << std::flush;
500 << TimeStamp() << " Final Single Point Energy "
501 << std::setprecision(12) << totenergy << " Ha" << std::flush;
502 XTP_LOG(Log::error, *pLog_) << TimeStamp() << std::setprecision(12)
503 << " Final Local Exc contribution "
504 << e_vxc.energy() << " Ha" << std::flush;
505 if (ScaHFX_ > 0) {
506 XTP_LOG(Log::error, *pLog_) << TimeStamp() << std::setprecision(12)
507 << " Final Non Local Ex contribution "
508 << exx << " Ha" << std::flush;
509 }
510
512
513 Index nuclear_charge = 0;
514 for (const QMAtom& atom : orb.QMAtoms()) {
515 nuclear_charge += atom.getNuccharge();
516 }
517
518 orb.setQMEnergy(totenergy);
519 orb.MOs() = MOs;
524 nuclear_charge - numofelectrons_,
526
527 CalcElDipole(orb);
528 return true;
529 } else if (this_iter == max_iter_ - 1) {
531 << TimeStamp() << " DFT calculation has not converged after "
532 << max_iter_
533 << " iterations. Use more iterations or another convergence "
534 "acceleration scheme."
535 << std::flush;
536 return false;
537 }
538 }
539
540 return true;
541}
542
543// Unrestricted SCF loop. The alpha and beta channels are iterated through
544// separate Fock matrices
545//
546// F^alpha = H0 + J[P^alpha + P^beta] + V_xc^alpha + K^alpha
547// F^beta = H0 + J[P^alpha + P^beta] + V_xc^beta + K^beta,
548//
549// while the total energy uses the spin-summed one-electron and Coulomb terms
550// together with spin-resolved XC and exact-exchange contributions.
552 const Vxc_Potential<Vxc_Grid>& vxcpotential) {
553 tools::EigenSystem MOs_alpha;
554 tools::EigenSystem MOs_beta;
555
556 MOs_alpha.eigenvalues() = Eigen::VectorXd::Zero(H0.cols());
557 MOs_alpha.eigenvectors() = Eigen::MatrixXd::Zero(H0.rows(), H0.cols());
558 MOs_beta.eigenvalues() = Eigen::VectorXd::Zero(H0.cols());
559 MOs_beta.eigenvectors() = Eigen::MatrixXd::Zero(H0.rows(), H0.cols());
560
561 UKSConvergenceAcc conv_uks;
562
566
570
571 conv_uks.Configure(opt_alpha, opt_beta);
572 conv_uks.setLogger(pLog_);
573 conv_uks.setOverlap(dftAOoverlap_, 1e-8);
574
575 if (initial_guess_ == "orbfile") {
577 << TimeStamp() << " Reading UKS guess from orbitals object/file"
578 << std::flush;
579
580 MOs_alpha = orb.MOs();
581 MOs_alpha.eigenvectors() = OrthogonalizeGuess(MOs_alpha.eigenvectors());
582
583 if (orb.hasBetaMOs()) {
584 MOs_beta = orb.MOs_beta();
585 MOs_beta.eigenvectors() = OrthogonalizeGuess(MOs_beta.eigenvectors());
586 } else {
588 << TimeStamp()
589 << " Orbital file has no beta MOs, using alpha guess for beta."
590 << std::flush;
591 MOs_beta = MOs_alpha;
592 }
593 } else {
595 << TimeStamp() << " Setup UKS Initial Guess using: " << initial_guess_
596 << std::flush;
597
598 tools::EigenSystem guess;
599 if (initial_guess_ == "independent") {
600 guess = IndependentElectronGuess(H0);
601 } else if (initial_guess_ == "atom") {
602 guess = ModelPotentialGuess(H0, orb.QMAtoms(), vxcpotential);
603 } else if (initial_guess_ == "huckel") {
604 guess = ExtendedHuckelGuess(orb.QMAtoms());
605 } else if (initial_guess_ == "huckel_dft") {
606 guess = ExtendedHuckelDFTGuess(H0, orb.QMAtoms(), vxcpotential);
607 } else {
608 throw std::runtime_error("Initial guess method not known/implemented");
609 }
610
611 MOs_alpha = guess;
612 MOs_beta = guess;
613 }
614
615 // Build the initial spin densities P^alpha and P^beta from the chosen
616 // starting orbitals before entering the coupled UKS iterations.
618 conv_uks.DensityMatrix(MOs_alpha, MOs_beta);
619
621 << TimeStamp() << " UKS guess gives Nalpha="
622 << Dspin.alpha.cwiseProduct(dftAOoverlap_.Matrix()).sum()
623 << " Nbeta=" << Dspin.beta.cwiseProduct(dftAOoverlap_.Matrix()).sum()
624 << " Ntot=" << Dspin.total().cwiseProduct(dftAOoverlap_.Matrix()).sum()
625 << std::flush;
626
628 << TimeStamp() << " STARTING UKS SCF cycle" << std::flush;
630 << " ------------------------------------------------------------"
631 << std::flush;
632
633 for (Index this_iter = 0; this_iter < max_iter_; ++this_iter) {
634 XTP_LOG(Log::error, *pLog_) << std::flush;
635 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Iteration " << this_iter + 1
636 << " of " << max_iter_ << std::flush;
637
638 Eigen::MatrixXd H_alpha = H0.matrix();
639 Eigen::MatrixXd H_beta = H0.matrix();
640
641 // The Coulomb contribution depends only on the total density
642 // P = P^alpha + P^beta, while exchange and XC remain spin resolved.
643 const Eigen::MatrixXd D_total = Dspin.total();
644
645 double E_one = Dspin.alpha.cwiseProduct(H0.matrix()).sum() +
646 Dspin.beta.cwiseProduct(H0.matrix()).sum();
647
648 double E_coul = 0.0;
649 double E_xc = 0.0;
650 double E_exx = 0.0;
651
652 double integral_error = std::min(conv_uks.getDIIsError() * 1e-5, 1e-5);
653
654 if (ScaHFX_ > 0) {
655 std::array<Eigen::MatrixXd, 2> both_alpha = CalcERIs_EXX(
656 Eigen::MatrixXd::Zero(0, 0), Dspin.alpha, integral_error);
657 std::array<Eigen::MatrixXd, 2> both_beta =
658 CalcERIs_EXX(Eigen::MatrixXd::Zero(0, 0), Dspin.beta, integral_error);
659
660 Eigen::MatrixXd J = both_alpha[0] + both_beta[0];
661 Eigen::MatrixXd K_alpha = both_alpha[1];
662 Eigen::MatrixXd K_beta = both_beta[1];
663
664 H_alpha += J + ScaHFX_ * K_alpha;
665 H_beta += J + ScaHFX_ * K_beta;
666
667 E_coul = 0.5 * D_total.cwiseProduct(J).sum();
668 E_exx = 0.5 * ScaHFX_ *
669 (Dspin.alpha.cwiseProduct(K_alpha).sum() +
670 Dspin.beta.cwiseProduct(K_beta).sum());
671 } else {
672 Eigen::MatrixXd J = CalcERIs(D_total, integral_error);
673 H_alpha += J;
674 H_beta += J;
675 E_coul = 0.5 * D_total.cwiseProduct(J).sum();
676 }
677
678 auto vxc = vxcpotential.IntegrateVXCSpin(Dspin.alpha, Dspin.beta);
679 H_alpha += vxc.vxc_alpha;
680 H_beta += vxc.vxc_beta;
681 E_xc = vxc.energy;
682
683 double totenergy = H0.energy() + E_one + E_coul + E_xc + E_exx;
684
685 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " One particle energy "
686 << std::setprecision(12) << E_one << std::flush;
687 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " Coulomb contribution "
688 << std::setprecision(12) << E_coul << std::flush;
689 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " XC contribution "
690 << std::setprecision(12) << E_xc << std::flush;
691 if (ScaHFX_ > 0) {
693 << TimeStamp() << " EXX contribution " << std::setprecision(12)
694 << E_exx << std::flush;
695 }
697 << TimeStamp() << " Total Energy " << std::setprecision(12) << totenergy
698 << std::flush;
699
700 UKSConvergenceAcc::SpinFock Hspin{H_alpha, H_beta};
701 Dspin = conv_uks.Iterate(Dspin, Hspin, MOs_alpha, MOs_beta, totenergy);
703 MOs_beta = MOs_alpha;
704 Dspin.beta = Dspin.alpha;
705 }
706
708 << TimeStamp()
709 << " Nalpha=" << Dspin.alpha.cwiseProduct(dftAOoverlap_.Matrix()).sum()
710 << " Nbeta=" << Dspin.beta.cwiseProduct(dftAOoverlap_.Matrix()).sum()
711 << std::flush;
712
714 << TimeStamp() << " <Sz> = "
716 << std::flush;
717
718 PrintMOsUKS(MOs_alpha.eigenvalues(), MOs_beta.eigenvalues(), Log::info);
719
720 if (conv_uks.isConverged()) {
721 Index nuclear_charge = 0;
722 for (const QMAtom& atom : orb.QMAtoms()) {
723 nuclear_charge += atom.getNuccharge();
724 }
725
726 CanonicalizeOrbitalPhases(MOs_alpha);
727 CanonicalizeOrbitalPhases(MOs_beta);
728
729 orb.setQMEnergy(totenergy);
730 orb.MOs() = MOs_alpha;
731 orb.MOs_beta() = MOs_beta;
737 nuclear_charge - numofelectrons_,
739
741 << TimeStamp() << " UKS converged after " << this_iter + 1
742 << " iterations. Delta E=" << conv_uks.getDeltaE()
743 << " DIIS error=" << conv_uks.getDIIsError() << std::flush;
744
746 << TimeStamp() << " Final Single Point Energy "
747 << std::setprecision(12) << totenergy << " Ha" << std::flush;
749 << TimeStamp() << std::setprecision(12) << " Final XC contribution "
750 << E_xc << " Ha" << std::flush;
751 if (ScaHFX_ > 0) {
753 << TimeStamp() << std::setprecision(12)
754 << " Final EXX contribution " << E_exx << " Ha" << std::flush;
755 }
756
758 << TimeStamp() << " <Sz> = "
760 << std::flush;
761
762 PrintMOsUKS(MOs_alpha.eigenvalues(), MOs_beta.eigenvalues(), Log::error);
763
764 CalcElDipole(orb);
765 return true;
766 }
767
768 if (this_iter == max_iter_ - 1) {
770 << TimeStamp() << " UKS calculation has not converged after "
771 << max_iter_ << " iterations." << std::flush;
772 return false;
773 }
774 }
775
776 return false;
777}
778
779// One-electron core Hamiltonian and its constant energy offset.
780//
781// The matrix part is
782//
783// H0 = T + V_nuc + V_ECP + V_ext,
784//
785// while the scalar energy collects all nucleus-nucleus and nucleus-external
786// interaction terms that do not depend on the electronic density.
788
789 AOKinetic dftAOkinetic;
790
791 dftAOkinetic.Fill(dftbasis_);
793 << TimeStamp() << " Filled DFT Kinetic energy matrix ." << std::flush;
794
795 AOMultipole dftAOESP;
796 dftAOESP.FillPotential(dftbasis_, mol);
798 << TimeStamp() << " Filled DFT nuclear potential matrix." << std::flush;
799
800 Eigen::MatrixXd H0 = dftAOkinetic.Matrix() + dftAOESP.Matrix();
802 << TimeStamp() << " Constructed independent particle hamiltonian "
803 << std::flush;
804 double E0 = NuclearRepulsion(mol);
805 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Nuclear Repulsion Energy is "
806 << std::setprecision(9) << E0 << std::flush;
807
808 if (!ecp_name_.empty()) {
809 AOECP dftAOECP;
810 dftAOECP.FillPotential(dftbasis_, ecp_);
811 H0 += dftAOECP.Matrix();
813 << TimeStamp() << " Filled DFT ECP matrix" << std::flush;
814 }
815
816 if (externalsites_ != nullptr) {
817 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " " << externalsites_->size()
818 << " External sites" << std::flush;
819 bool has_quadrupoles = std::any_of(
820 externalsites_->begin(), externalsites_->end(),
821 [](const std::unique_ptr<StaticSite>& s) { return s->getRank() == 2; });
822 std::string header =
823 " Name Coordinates[a0] charge[e] dipole[e*a0] ";
824 if (has_quadrupoles) {
825 header += " quadrupole[e*a0^2]";
826 }
827 XTP_LOG(Log::error, *pLog_) << header << std::flush;
828 Index limit = 50;
829 Index counter = 0;
830 for (const std::unique_ptr<StaticSite>& site : *externalsites_) {
831 if (counter == limit) {
832 break;
833 }
834 std::string output =
835 (boost::format(" %1$s"
836 " %2$+1.4f %3$+1.4f %4$+1.4f"
837 " %5$+1.4f") %
838 site->getElement() % site->getPos()[0] % site->getPos()[1] %
839 site->getPos()[2] % site->getCharge())
840 .str();
841 const Eigen::Vector3d& dipole = site->getDipole();
842 output += (boost::format(" %1$+1.4f %2$+1.4f %3$+1.4f") % dipole[0] %
843 dipole[1] % dipole[2])
844 .str();
845 if (site->getRank() > 1) {
846 Eigen::VectorXd quadrupole = site->Q().tail<5>();
847 output +=
848 (boost::format(" %1$+1.4f %2$+1.4f %3$+1.4f %4$+1.4f %5$+1.4f") %
849 quadrupole[0] % quadrupole[1] % quadrupole[2] % quadrupole[3] %
850 quadrupole[4])
851 .str();
852 }
853 XTP_LOG(Log::error, *pLog_) << output << std::flush;
854 counter++;
855 }
856 if (counter == limit) {
858 << " ... (" << externalsites_->size() - limit
859 << " sites not displayed)\n"
860 << std::flush;
861 }
862
863 Mat_p_Energy ext_multipoles =
866 << TimeStamp() << " Nuclei-external site interaction energy "
867 << std::setprecision(9) << ext_multipoles.energy() << std::flush;
868 E0 += ext_multipoles.energy();
869 H0 += ext_multipoles.matrix();
870 }
871
873 Orbitals extdensity;
874 extdensity.ReadFromCpt(orbfilename_);
875 Mat_p_Energy extdensity_result = IntegrateExternalDensity(mol, extdensity);
876 E0 += extdensity_result.energy();
878 << TimeStamp() << " Nuclei-external density interaction energy "
879 << std::setprecision(9) << extdensity_result.energy() << std::flush;
880 H0 += extdensity_result.matrix();
881 }
882
884
886 << TimeStamp() << " Integrating external electric field with F[Hrt]="
887 << extfield_.transpose() << std::flush;
888 H0 += IntegrateExternalField(mol);
889 }
890
891 return Mat_p_Energy(E0, H0);
892}
893
894// Precompute SCF-invariant matrices: overlap for the generalized eigenvalue
895// problem and the RI/4c electron-repulsion backend that later yields J[P] and
896// K[P].
900 << TimeStamp() << " Filled DFT Overlap matrix." << std::flush;
901
902 conv_opt_.numberofelectrons = numofelectrons_;
903 conv_opt_.number_alpha_electrons = num_alpha_electrons_;
904 conv_opt_.number_beta_electrons = num_beta_electrons_;
908 conv_accelerator_.Configure(conv_opt_);
909 conv_accelerator_.setLogger(pLog_);
910 conv_accelerator_.setOverlap(dftAOoverlap_, 1e-8);
911 conv_accelerator_.PrintConfigOptions();
912
913 if (!auxbasis_name_.empty()) {
914 // prepare invariant part of electron repulsion integrals
915 ERIs_.Initialize(dftbasis_, auxbasis_);
917 << TimeStamp() << " Inverted AUX Coulomb matrix, removed "
918 << ERIs_.Removedfunctions() << " functions from aux basis"
919 << std::flush;
921 << TimeStamp()
922 << " Setup invariant parts of Electron Repulsion integrals "
923 << std::flush;
924 } else {
926 << TimeStamp() << " Calculating 4c diagonals. " << std::flush;
927 ERIs_.Initialize_4c(dftbasis_);
929 << TimeStamp() << " Calculated 4c diagonals. " << std::flush;
930 }
931
932 return;
933}
934
936 const QMAtom& uniqueAtom) const {
937 bool with_ecp = !ecp_name_.empty();
938 if (uniqueAtom.getElement() == "H" || uniqueAtom.getElement() == "He") {
939 with_ecp = false;
940 }
941
942 QMMolecule atom = QMMolecule("individual_atom", 0);
943 atom.push_back(uniqueAtom);
944
945 BasisSet basisset;
946 basisset.Load(dftbasis_name_);
947 AOBasis dftbasis;
948 dftbasis.Fill(basisset, atom);
949 Vxc_Grid grid;
950 grid.GridSetup(grid_name_, atom, dftbasis);
951 Vxc_Potential<Vxc_Grid> gridIntegration(grid);
952 gridIntegration.setXCfunctional(xc_functional_name_);
953
954 ECPAOBasis ecp;
955 if (with_ecp) {
956 ECPBasisSet ecps;
957 ecps.Load(ecp_name_);
958 ecp.Fill(ecps, atom);
959 }
960
961 Index numofelectrons = uniqueAtom.getNuccharge();
962 Index alpha_e = 0;
963 Index beta_e = 0;
964
965 if ((numofelectrons % 2) != 0) {
966 alpha_e = numofelectrons / 2 + numofelectrons % 2;
967 beta_e = numofelectrons / 2;
968 } else {
969 alpha_e = numofelectrons / 2;
970 beta_e = alpha_e;
971 }
972
973 AOOverlap dftAOoverlap;
974 AOKinetic dftAOkinetic;
975 AOMultipole dftAOESP;
976 AOECP dftAOECP;
977 ERIs ERIs_atom;
978
979 dftAOoverlap.Fill(dftbasis);
980 dftAOkinetic.Fill(dftbasis);
981
982 dftAOESP.FillPotential(dftbasis, atom);
983 ERIs_atom.Initialize_4c(dftbasis);
984
985 ConvergenceAcc Convergence_alpha;
986 ConvergenceAcc Convergence_beta;
989 opt_alpha.histlength = 20;
990 opt_alpha.levelshift = 0.1;
991 opt_alpha.levelshiftend = 0.0;
992 opt_alpha.usediis = true;
993 opt_alpha.adiis_start = 0.0;
994 opt_alpha.diis_start = 0.0;
995 opt_alpha.numberofelectrons = alpha_e;
996
997 ConvergenceAcc::options opt_beta = opt_alpha;
998 opt_beta.numberofelectrons = beta_e;
999
1000 Logger log;
1001 Convergence_alpha.Configure(opt_alpha);
1002 Convergence_alpha.setLogger(&log);
1003 Convergence_alpha.setOverlap(dftAOoverlap, 1e-8);
1004 Convergence_beta.Configure(opt_beta);
1005 Convergence_beta.setLogger(&log);
1006 Convergence_beta.setOverlap(dftAOoverlap, 1e-8);
1007
1008 Eigen::MatrixXd H0 = dftAOkinetic.Matrix() + dftAOESP.Matrix();
1009 if (with_ecp) {
1010 dftAOECP.FillPotential(dftbasis, ecp);
1011 H0 += dftAOECP.Matrix();
1012 }
1013
1014 tools::EigenSystem MOs_alpha = Convergence_alpha.SolveFockmatrix(H0);
1015 Eigen::MatrixXd dftAOdmat_alpha = Convergence_alpha.DensityMatrix(MOs_alpha);
1016
1017 if (uniqueAtom.getElement() == "H") {
1018 return dftAOdmat_alpha;
1019 }
1020
1021 tools::EigenSystem MOs_beta = Convergence_beta.SolveFockmatrix(H0);
1022 Eigen::MatrixXd dftAOdmat_beta = Convergence_beta.DensityMatrix(MOs_beta);
1023
1024 Index maxiter = 80;
1025 for (Index this_iter = 0; this_iter < maxiter; this_iter++) {
1026 Eigen::MatrixXd H_alpha = H0;
1027 Eigen::MatrixXd H_beta = H0;
1028
1029 double E_coul = 0.0;
1030 double E_exx = 0.0;
1031 double E_xc = 0.0;
1032
1033 double integral_error = std::min(1e-5 * 0.5 *
1034 (Convergence_alpha.getDIIsError() +
1035 Convergence_beta.getDIIsError()),
1036 1e-5);
1037
1038 if (ScaHFX_ > 0) {
1039 std::array<Eigen::MatrixXd, 2> both_alpha =
1040 ERIs_atom.CalculateERIs_EXX_4c(dftAOdmat_alpha, integral_error);
1041 std::array<Eigen::MatrixXd, 2> both_beta =
1042 ERIs_atom.CalculateERIs_EXX_4c(dftAOdmat_beta, integral_error);
1043
1044 Eigen::MatrixXd Hartree = both_alpha[0] + both_beta[0];
1045 H_alpha += Hartree + ScaHFX_ * both_alpha[1];
1046 H_beta += Hartree + ScaHFX_ * both_beta[1];
1047
1048 E_coul =
1049 0.5 * (dftAOdmat_alpha + dftAOdmat_beta).cwiseProduct(Hartree).sum();
1050 E_exx = 0.5 * ScaHFX_ *
1051 (both_alpha[1].cwiseProduct(dftAOdmat_alpha).sum() +
1052 both_beta[1].cwiseProduct(dftAOdmat_beta).sum());
1053 } else {
1054 Eigen::MatrixXd Hartree = ERIs_atom.CalculateERIs_4c(
1055 dftAOdmat_alpha + dftAOdmat_beta, integral_error);
1056 H_alpha += Hartree;
1057 H_beta += Hartree;
1058 E_coul =
1059 0.5 * (dftAOdmat_alpha + dftAOdmat_beta).cwiseProduct(Hartree).sum();
1060 }
1061
1062 auto vxc =
1063 gridIntegration.IntegrateVXCSpin(dftAOdmat_alpha, dftAOdmat_beta);
1064 H_alpha += vxc.vxc_alpha;
1065 H_beta += vxc.vxc_beta;
1066 E_xc = vxc.energy;
1067
1068 double E_one_alpha = dftAOdmat_alpha.cwiseProduct(H0).sum();
1069 double E_one_beta = dftAOdmat_beta.cwiseProduct(H0).sum();
1070 double totenergy = E_one_alpha + E_one_beta + E_coul + E_exx + E_xc;
1071
1072 dftAOdmat_alpha = Convergence_alpha.Iterate(dftAOdmat_alpha, H_alpha,
1073 MOs_alpha, totenergy);
1074 dftAOdmat_beta =
1075 Convergence_beta.Iterate(dftAOdmat_beta, H_beta, MOs_beta, totenergy);
1076
1078 << TimeStamp() << " Iter " << this_iter << " of " << maxiter << " Etot "
1079 << totenergy << " diise_a " << Convergence_alpha.getDIIsError()
1080 << " diise_b " << Convergence_beta.getDIIsError() << "\n\t\t a_gap "
1081 << MOs_alpha.eigenvalues()(alpha_e) -
1082 MOs_alpha.eigenvalues()(alpha_e - 1)
1083 << " b_gap "
1084 << MOs_beta.eigenvalues()(beta_e) - MOs_beta.eigenvalues()(beta_e - 1)
1085 << " Nalpha="
1086 << dftAOoverlap.Matrix().cwiseProduct(dftAOdmat_alpha).sum()
1087 << " Nbeta=" << dftAOoverlap.Matrix().cwiseProduct(dftAOdmat_beta).sum()
1088 << std::flush;
1089
1090 bool converged =
1091 Convergence_alpha.isConverged() && Convergence_beta.isConverged();
1092 if (converged || this_iter == maxiter - 1) {
1093 if (converged) {
1095 << TimeStamp() << " Converged after " << this_iter + 1
1096 << " iterations" << std::flush;
1097 } else {
1099 << TimeStamp() << " Not converged after " << this_iter + 1
1100 << " iterations. Unconverged density.\n\t\t\t"
1101 << " DIIsError_alpha=" << Convergence_alpha.getDIIsError()
1102 << " DIIsError_beta=" << Convergence_beta.getDIIsError()
1103 << std::flush;
1104 }
1105 break;
1106 }
1107 }
1108
1109 Eigen::MatrixXd avgmatrix =
1110 SphericalAverageShells(dftAOdmat_alpha + dftAOdmat_beta, dftbasis);
1112 << TimeStamp() << " Atomic density Matrix for " << uniqueAtom.getElement()
1113 << " gives N=" << std::setprecision(9)
1114 << avgmatrix.cwiseProduct(dftAOoverlap.Matrix()).sum() << " electrons."
1115 << std::flush;
1116 return avgmatrix;
1117}
1118
1119Eigen::MatrixXd DFTEngine::AtomicGuess(const QMMolecule& mol) const {
1120
1121 std::vector<std::string> elements = mol.FindUniqueElements();
1123 << TimeStamp() << " Scanning molecule of size " << mol.size()
1124 << " for unique elements" << std::flush;
1125 QMMolecule uniqueelements = QMMolecule("uniqueelements", 0);
1126 for (auto element : elements) {
1127 uniqueelements.push_back(QMAtom(0, element, Eigen::Vector3d::Zero()));
1128 }
1129
1130 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " " << uniqueelements.size()
1131 << " unique elements found" << std::flush;
1132 std::vector<Eigen::MatrixXd> uniqueatom_guesses;
1133 for (QMAtom& unique_atom : uniqueelements) {
1135 << TimeStamp() << " Calculating atom density for "
1136 << unique_atom.getElement() << std::flush;
1137 Eigen::MatrixXd dmat_unrestricted = RunAtomicDFT_unrestricted(unique_atom);
1138 uniqueatom_guesses.push_back(dmat_unrestricted);
1139 }
1140
1141 Eigen::MatrixXd guess =
1142 Eigen::MatrixXd::Zero(dftbasis_.AOBasisSize(), dftbasis_.AOBasisSize());
1143 Index start = 0;
1144 for (const QMAtom& atom : mol) {
1145 Index index = 0;
1146 for (; index < uniqueelements.size(); index++) {
1147 if (atom.getElement() == uniqueelements[index].getElement()) {
1148 break;
1149 }
1150 }
1151 Eigen::MatrixXd& dmat_unrestricted = uniqueatom_guesses[index];
1152 guess.block(start, start, dmat_unrestricted.rows(),
1153 dmat_unrestricted.cols()) = dmat_unrestricted;
1154 start += dmat_unrestricted.rows();
1155 }
1156
1157 return guess;
1158}
1159
1161 if (initial_guess_ == "orbfile") {
1162
1163 if (orb.hasDFTbasisName()) {
1164 if (orb.getDFTbasisName() != dftbasis_name_) {
1165 throw std::runtime_error(
1166 (boost::format("Basisset Name in guess orb file "
1167 "and in dftengine option file differ %1% vs %2%") %
1169 .str());
1170 }
1171 } else {
1173 << TimeStamp()
1174 << " WARNING: "
1175 "Orbital file has no basisset information,"
1176 "using it as a guess might work or not for calculation with "
1177 << dftbasis_name_ << std::flush;
1178 }
1179 }
1180
1181 const Index target_charge = orb.getCharge();
1182 const Index multiplicity = orb.getSpin();
1183
1184 orb.setChargeAndSpin(target_charge, multiplicity);
1187
1189 orb.setXCGrid(grid_name_);
1190 orb.setScaHFX(ScaHFX_);
1191 if (!ecp_name_.empty()) {
1192 orb.setECPName(ecp_name_);
1193 }
1194 if (!auxbasis_name_.empty()) {
1196 }
1197
1198 if (initial_guess_ == "orbfile") {
1199 if (orb.hasECPName() || !ecp_name_.empty()) {
1200 if (orb.getECPName() != ecp_name_) {
1201 throw std::runtime_error(
1202 (boost::format("ECPs in orb file: %1% and options %2% differ") %
1203 orb.getECPName() % ecp_name_)
1204 .str());
1205 }
1206 }
1209 throw std::runtime_error(
1210 (boost::format("Number of electrons in guess orb file "
1211 "and in dftengine differ: "
1212 "alpha %1% vs %2%, beta %3% vs %4%.") %
1215 .str());
1216 }
1217 if (orb.getBasisSetSize() != dftbasis_.AOBasisSize()) {
1218 throw std::runtime_error(
1219 (boost::format("Number of levels in guess orb file: "
1220 "%1% and in dftengine: %2% differ.") %
1221 orb.getBasisSetSize() % dftbasis_.AOBasisSize())
1222 .str());
1223 }
1224 } else {
1227 }
1228 return;
1229}
1230
1231void DFTEngine::Prepare(Orbitals& orb, Index numofelectrons) {
1232 QMMolecule& mol = orb.QMAtoms();
1233
1235 << TimeStamp() << " Using " << OPENMP::getMaxThreads() << " threads"
1236 << std::flush;
1237
1238 if (XTP_HAS_MKL_OVERLOAD()) {
1240 << TimeStamp() << " Using MKL overload for Eigen " << std::flush;
1241 } else {
1243 << TimeStamp()
1244 << " Using native Eigen implementation, no BLAS overload "
1245 << std::flush;
1246 }
1247
1248 XTP_LOG(Log::error, *pLog_) << " Molecule Coordinates [A] " << std::flush;
1249 for (const QMAtom& atom : mol) {
1250 const Eigen::Vector3d pos = atom.getPos() * tools::conv::bohr2ang;
1251 std::string output = (boost::format(" %1$s"
1252 " %2$+1.4f %3$+1.4f %4$+1.4f") %
1253 atom.getElement() % pos[0] % pos[1] % pos[2])
1254 .str();
1255
1256 XTP_LOG(Log::error, *pLog_) << output << std::flush;
1257 }
1258
1260 dftbasis_ = orb.getDftBasis();
1261
1263 << TimeStamp() << " Loaded DFT Basis Set " << dftbasis_name_ << " with "
1264 << dftbasis_.AOBasisSize() << " functions" << std::flush;
1265
1266 if (!auxbasis_name_.empty()) {
1267 BasisSet auxbasisset;
1268 auxbasisset.Load(auxbasis_name_);
1269 auxbasis_.Fill(auxbasisset, mol);
1271 << TimeStamp() << " Loaded AUX Basis Set " << auxbasis_name_ << " with "
1272 << auxbasis_.AOBasisSize() << " functions" << std::flush;
1273 }
1274 if (!ecp_name_.empty()) {
1275 ECPBasisSet ecpbasisset;
1276 ecpbasisset.Load(ecp_name_);
1278 << TimeStamp() << " Loaded ECP library " << ecp_name_ << std::flush;
1279
1280 std::vector<std::string> results = ecp_.Fill(ecpbasisset, mol);
1282 << TimeStamp() << " Filled ECP Basis" << std::flush;
1283 if (results.size() > 0) {
1284 std::string message = "";
1285 for (const std::string& element : results) {
1286 message += " " + element;
1287 }
1289 << TimeStamp() << " Found no ECPs for elements" << message
1290 << std::flush;
1291 }
1292 }
1293
1294 numofelectrons_ = 0;
1297 num_docc_ = 0;
1298 num_socc_alpha_ = 0;
1299
1300 Index nuclear_charge = 0;
1301 for (const QMAtom& atom : mol) {
1302 nuclear_charge += atom.getNuccharge();
1303 }
1304
1305 Index target_charge = orb.getCharge();
1306 Index multiplicity = orb.getSpin();
1307
1308 if (multiplicity < 1) {
1309 throw std::runtime_error("Spin multiplicity must be >= 1.");
1310 }
1311
1312 if (numofelectrons >= 0) {
1313 numofelectrons_ = numofelectrons;
1314 } else {
1315 numofelectrons_ = nuclear_charge - target_charge;
1316 }
1317
1318 Index spin_excess = multiplicity - 1;
1319
1320 if (numofelectrons_ < 0) {
1321 throw std::runtime_error("Computed a negative number of electrons.");
1322 }
1323
1324 if (spin_excess > numofelectrons_) {
1325 throw std::runtime_error(
1326 "Spin multiplicity incompatible with total number of electrons.");
1327 }
1328
1329 if (((numofelectrons_ + spin_excess) % 2) != 0) {
1330 throw std::runtime_error(
1331 "Charge and spin multiplicity imply non-integer alpha/beta "
1332 "occupations.");
1333 }
1334
1335 num_alpha_electrons_ = (numofelectrons_ + spin_excess) / 2;
1336 num_beta_electrons_ = (numofelectrons_ - spin_excess) / 2;
1337
1340
1342 << TimeStamp() << " Total number of electrons: " << numofelectrons_
1343 << " (charge=" << target_charge << ", multiplicity=" << multiplicity
1344 << ", alpha=" << num_alpha_electrons_ << ", beta=" << num_beta_electrons_
1345 << ", docc=" << num_docc_ << ", socc=" << num_socc_alpha_ << ")"
1346 << std::flush;
1347
1349 return;
1350}
1351
1354 if (ScaHFX_ > 0) {
1356 << TimeStamp() << " Using hybrid functional with alpha=" << ScaHFX_
1357 << std::flush;
1358 }
1359 Vxc_Grid grid;
1360 grid.GridSetup(grid_name_, mol, dftbasis_);
1361 Vxc_Potential<Vxc_Grid> vxc(grid);
1364 << TimeStamp() << " Setup numerical integration grid " << grid_name_
1365 << " for vxc functional " << xc_functional_name_ << std::flush;
1367 << "\t\t "
1368 << " with " << grid.getGridSize() << " points"
1369 << " divided into " << grid.getBoxesSize() << " boxes" << std::flush;
1370 return vxc;
1371}
1372
1373double DFTEngine::NuclearRepulsion(const QMMolecule& mol) const {
1374 double E_nucnuc = 0.0;
1375
1376 for (Index i = 0; i < mol.size(); i++) {
1377 const Eigen::Vector3d& r1 = mol[i].getPos();
1378 double charge1 = double(mol[i].getNuccharge());
1379 for (Index j = 0; j < i; j++) {
1380 const Eigen::Vector3d& r2 = mol[j].getPos();
1381 double charge2 = double(mol[j].getNuccharge());
1382 E_nucnuc += charge1 * charge2 / (r1 - r2).norm();
1383 }
1384 }
1385 return E_nucnuc;
1386}
1387
1388// spherically average the density matrix belonging to two shells
1390 const Eigen::MatrixXd& dmat, const AOBasis& dftbasis) const {
1391 Eigen::MatrixXd avdmat = Eigen::MatrixXd::Zero(dmat.rows(), dmat.cols());
1392 for (const AOShell& shellrow : dftbasis) {
1393 Index size_row = shellrow.getNumFunc();
1394 Index start_row = shellrow.getStartIndex();
1395 for (const AOShell& shellcol : dftbasis) {
1396 Index size_col = shellcol.getNumFunc();
1397 Index start_col = shellcol.getStartIndex();
1398 Eigen::MatrixXd shelldmat =
1399 dmat.block(start_row, start_col, size_row, size_col);
1400 if (shellrow.getL() == shellcol.getL()) {
1401 double diagavg = shelldmat.diagonal().sum() / double(shelldmat.rows());
1402 Index offdiagelements =
1403 shelldmat.rows() * shelldmat.cols() - shelldmat.cols();
1404 double offdiagavg = (shelldmat.sum() - shelldmat.diagonal().sum()) /
1405 double(offdiagelements);
1406 avdmat.block(start_row, start_col, size_row, size_col).array() =
1407 offdiagavg;
1408 avdmat.block(start_row, start_col, size_row, size_col)
1409 .diagonal()
1410 .array() = diagavg;
1411 } else {
1412 double avg = shelldmat.sum() / double(shelldmat.size());
1413 avdmat.block(start_row, start_col, size_row, size_col).array() = avg;
1414 }
1415 }
1416 }
1417 return avdmat;
1418}
1419
1421 const QMMolecule& mol,
1422 const std::vector<std::unique_ptr<StaticSite> >& multipoles) const {
1423
1424 if (multipoles.size() == 0) {
1425 return 0;
1426 }
1427
1428 double E_ext = 0;
1429 eeInteractor interactor;
1430 for (const QMAtom& atom : mol) {
1431 StaticSite nucleus = StaticSite(atom, double(atom.getNuccharge()));
1432 for (const std::unique_ptr<StaticSite>& site : *externalsites_) {
1433 if ((site->getPos() - nucleus.getPos()).norm() < 1e-7) {
1435 << " External site sits on nucleus, "
1436 "interaction between them is ignored."
1437 << std::flush;
1438 continue;
1439 }
1440 E_ext += interactor.CalcStaticEnergy_site(*site, nucleus);
1441 }
1442 }
1443 return E_ext;
1444}
1445
1446Eigen::MatrixXd DFTEngine::IntegrateExternalField(const QMMolecule& mol) const {
1447
1448 AODipole dipole;
1449 dipole.setCenter(mol.getPos());
1450 dipole.Fill(dftbasis_);
1451 Eigen::MatrixXd result =
1452 Eigen::MatrixXd::Zero(dipole.Dimension(), dipole.Dimension());
1453 for (Index i = 0; i < 3; i++) {
1454 result -= dipole.Matrix()[i] * extfield_[i];
1455 }
1456 return result;
1457}
1458
1460 const QMMolecule& mol,
1461 const std::vector<std::unique_ptr<StaticSite> >& multipoles) const {
1462
1463 Mat_p_Energy result(dftbasis_.AOBasisSize(), dftbasis_.AOBasisSize());
1464 AOMultipole dftAOESP;
1465
1466 dftAOESP.FillPotential(dftbasis_, multipoles);
1468 << TimeStamp() << " Filled DFT external multipole potential matrix"
1469 << std::flush;
1470 result.matrix() = dftAOESP.Matrix();
1471 result.energy() = ExternalRepulsion(mol, multipoles);
1472
1473 return result;
1474}
1475
1477 const QMMolecule& mol, const Orbitals& extdensity) const {
1478 BasisSet basis;
1479 basis.Load(extdensity.getDFTbasisName());
1480 AOBasis aobasis;
1481 aobasis.Fill(basis, extdensity.QMAtoms());
1482 Vxc_Grid grid;
1483 grid.GridSetup(gridquality_, extdensity.QMAtoms(), aobasis);
1484 DensityIntegration<Vxc_Grid> numint(grid);
1485 Eigen::MatrixXd dmat = extdensity.DensityMatrixFull(state_);
1486
1487 numint.IntegrateDensity(dmat);
1489 << TimeStamp() << " Calculated external density" << std::flush;
1490 Eigen::MatrixXd e_contrib = numint.IntegratePotential(dftbasis_);
1492 << TimeStamp() << " Calculated potential from electron density"
1493 << std::flush;
1494 AOMultipole esp;
1495 esp.FillPotential(dftbasis_, extdensity.QMAtoms());
1496
1497 double nuc_energy = 0.0;
1498 for (const QMAtom& atom : mol) {
1499 nuc_energy +=
1500 numint.IntegratePotential(atom.getPos()) * double(atom.getNuccharge());
1501 for (const QMAtom& extatom : extdensity.QMAtoms()) {
1502 const double dist = (atom.getPos() - extatom.getPos()).norm();
1503 nuc_energy +=
1504 double(atom.getNuccharge()) * double(extatom.getNuccharge()) / dist;
1505 }
1506 }
1508 << TimeStamp() << " Calculated potential from nuclei" << std::flush;
1510 << TimeStamp() << " Electrostatic: " << nuc_energy << std::flush;
1511 return Mat_p_Energy(nuc_energy, e_contrib + esp.Matrix());
1512}
1513
1515 const Eigen::MatrixXd& GuessMOs) const {
1516 Eigen::MatrixXd nonortho =
1517 GuessMOs.transpose() * dftAOoverlap_.Matrix() * GuessMOs;
1518 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(nonortho);
1519 Eigen::MatrixXd result = GuessMOs * es.operatorInverseSqrt();
1520 return result;
1521}
1522
1523/*************************************************************
1524 * Extended Hueckel Theory
1525 ************************************************************/
1527 const QMMolecule& mol) const {
1528
1530
1531 const Index nao = dftbasis_.AOBasisSize();
1532 Eigen::VectorXd eps = Eigen::VectorXd::Zero(nao);
1533
1534 for (const AOShell& shell : dftbasis_) {
1535
1536 int l = static_cast<int>(shell.getL());
1537 Index start = shell.getStartIndex();
1538 Index nfunc = shell.getNumFunc();
1539
1540 const QMAtom& atom = mol[shell.getAtomIndex()];
1541 const std::string& element = atom.getElement();
1542
1543 int used_l = l;
1544 double e = params.GetWithFallback(element, l, &used_l);
1545
1546 for (Index i = 0; i < nfunc; ++i) {
1547 eps(start + i) = e;
1548 }
1549 }
1550
1551 return eps;
1552}
1553
1554Eigen::MatrixXd DFTEngine::BuildEHTHamiltonian(const QMMolecule& mol) const {
1555
1556 const Eigen::MatrixXd& S = dftAOoverlap_.Matrix();
1557 const Index nao = S.rows();
1558 Eigen::VectorXd eps = BuildEHTOrbitalEnergies(mol);
1559 Eigen::MatrixXd H = Eigen::MatrixXd::Zero(nao, nao);
1560 constexpr double K = 1.75;
1561
1562 for (Index mu = 0; mu < nao; ++mu) {
1563 H(mu, mu) = eps(mu);
1564 for (Index nu = 0; nu < mu; ++nu) {
1565 double hij = K * S(mu, nu) * 0.5 * (eps(mu) + eps(nu));
1566 H(mu, nu) = hij;
1567 H(nu, mu) = hij;
1568 }
1569 }
1570
1571 return H;
1572}
1573
1575
1577 << TimeStamp() << " Building Extended Huckel guess" << std::flush;
1578
1579 Eigen::MatrixXd H = BuildEHTHamiltonian(mol);
1580
1582 << TimeStamp() << " Solving EHT generalized eigenproblem" << std::flush;
1583
1584 return conv_accelerator_.SolveFockmatrix(H);
1585}
1586
1588 const Mat_p_Energy& H0, const QMMolecule& mol,
1589 const Vxc_Potential<Vxc_Grid>& vxcpotential) const {
1590
1592
1593 Eigen::MatrixXd Dmat = conv_accelerator_.DensityMatrix(eht);
1594
1595 Mat_p_Energy e_vxc = vxcpotential.IntegrateVXC(Dmat);
1596
1597 Eigen::MatrixXd H = H0.matrix() + e_vxc.matrix();
1598
1599 if (ScaHFX_ > 0) {
1600 std::array<Eigen::MatrixXd, 2> both =
1601 CalcERIs_EXX(Eigen::MatrixXd::Zero(0, 0), Dmat, 1e-12);
1602 H += both[0];
1603 H += ScaHFX_ * both[1];
1604 } else {
1605 H += CalcERIs(Dmat, 1e-12);
1606 }
1607
1608 return conv_accelerator_.SolveFockmatrix(H);
1609}
1610
1611} // namespace xtp
1612} // namespace votca
const Eigen::VectorXd & eigenvalues() const
Definition eigensystem.h:30
const Eigen::MatrixXd & eigenvectors() const
Definition eigensystem.h:33
class to manage program options with xml serialization functionality
Definition property.h:55
Property & get(const std::string &key)
get existing property
Definition property.cc:79
bool exists(const std::string &key) const
check whether property exists
Definition property.cc:122
T as() const
return value as type
Definition property.h:283
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
void Fill(const BasisSet &bs, const QMMolecule &atoms)
Definition aobasis.cc:85
Index Dimension() final
Definition aomatrix.h:91
void setCenter(const Eigen::Vector3d &r)
Definition aomatrix.h:94
const std::array< Eigen::MatrixXd, 3 > & Matrix() const
Definition aomatrix.h:92
void Fill(const AOBasis &aobasis) final
void FillPotential(const AOBasis &aobasis, const ECPAOBasis &ecp)
Definition aoecp.cc:59
void Fill(const AOBasis &aobasis) final
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:41
void FillPotential(const AOBasis &aobasis, const QMMolecule &atoms)
void Fill(const AOBasis &aobasis) final
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:52
const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > & Matrix() const
Definition aopotential.h:39
void push_back(const T &atom)
std::vector< std::string > FindUniqueElements() const
const Eigen::Vector3d & getPos() const
void Load(const std::string &name)
Definition basisset.cc:149
Eigen::MatrixXd DensityMatrix(const tools::EigenSystem &MOs) const
void Configure(const ConvergenceAcc::options &opt)
void setLogger(Logger *log)
Attach the logger used for convergence diagnostics.
Eigen::MatrixXd Iterate(const Eigen::MatrixXd &dmat, Eigen::MatrixXd &H, tools::EigenSystem &MOs, double totE)
double getDIIsError() const
Return the DIIS commutator norm from the latest iteration.
tools::EigenSystem SolveFockmatrix(const Eigen::MatrixXd &H) const
Solve the generalized eigenvalue problem for the current Fock matrix.
void setOverlap(AOOverlap &S, double etol)
Precompute overlap-dependent quantities used when solving the Fock matrix.
Eigen::MatrixXd CalcERIs(const Eigen::MatrixXd &Dmat, double error) const
Build the Coulomb matrix contribution from the current density matrix.
Definition dftengine.cc:291
std::string auxbasis_name_
Definition dftengine.h:245
std::string gridquality_
Definition dftengine.h:283
tools::EigenSystem ModelPotentialGuess(const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const
Definition dftengine.cc:311
Mat_p_Energy IntegrateExternalMultipoles(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
tools::EigenSystem IndependentElectronGuess(const Mat_p_Energy &H0) const
Generate an initial guess by diagonalizing the core Hamiltonian only.
Definition dftengine.cc:300
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:356
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:197
std::string dftbasis_name_
Definition dftengine.h:246
bool EvaluateUKS(Orbitals &orb, const Mat_p_Energy &H0, const Vxc_Potential< Vxc_Grid > &vxcpotential)
Definition dftengine.cc:551
Mat_p_Energy IntegrateExternalDensity(const QMMolecule &mol, const Orbitals &extdensity) const
std::string ecp_name_
Definition dftengine.h:247
std::string grid_name_
Definition dftengine.h:257
Eigen::MatrixXd RunAtomicDFT_unrestricted(const QMAtom &uniqueAtom) const
Definition dftengine.cc:935
void Prepare(Orbitals &orb, Index numofelectrons=-1)
std::string initial_guess_
Definition dftengine.h:262
std::string orbfilename_
Definition dftengine.h:282
Eigen::MatrixXd AtomicGuess(const QMMolecule &mol) const
Build an atomic-density based initial guess in the AO basis.
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:289
std::array< Eigen::MatrixXd, 2 > CalcERIs_EXX(const Eigen::MatrixXd &MOCoeff, const Eigen::MatrixXd &Dmat, double error) const
Definition dftengine.cc:274
void Initialize(tools::Property &options)
Read DFT, grid, and SCF settings from the user options tree.
Definition dftengine.cc:85
std::string xc_functional_name_
Definition dftengine.h:278
void CalcElDipole(const Orbitals &orb) const
Evaluate and print the electronic dipole moment from the final density.
Definition dftengine.cc:257
double ExternalRepulsion(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
ConvergenceAcc::options conv_opt_
Definition dftengine.h:267
void SetupInvariantMatrices()
Precompute AO matrices that remain unchanged throughout the SCF procedure.
Definition dftengine.cc:897
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:292
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.
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:177
AOOverlap dftAOoverlap_
Definition dftengine.h:260
Mat_p_Energy SetupH0(const QMMolecule &mol) const
Assemble the one-electron core Hamiltonian for the current molecule.
Definition dftengine.cc:787
tools::EigenSystem ExtendedHuckelGuess(const QMMolecule &mol) const
Eigen::MatrixXd OrthogonalizeGuess(const Eigen::MatrixXd &GuessMOs) const
Orthonormalize an initial MO guess with respect to the AO overlap matrix.
Eigen::MatrixXd SphericalAverageShells(const Eigen::MatrixXd &dmat, const AOBasis &dftbasis) const
bool Evaluate(Orbitals &orb)
Definition dftengine.cc:332
std::vector< std::unique_ptr< StaticSite > > * externalsites_
Definition dftengine.h:274
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:269
double IntegratePotential(const Eigen::Vector3d &rvector) const
double IntegrateDensity(const Eigen::MatrixXd &density_matrix)
Container to hold ECPs for all atoms.
Definition ecpaobasis.h:43
std::vector< std::string > Fill(const ECPBasisSet &bs, QMMolecule &atoms)
Definition ecpaobasis.cc:69
void Load(const std::string &name)
Takes a density matrix and and an auxiliary basis set and calculates the electron repulsion integrals...
Definition ERIs.h:35
std::array< Eigen::MatrixXd, 2 > CalculateERIs_EXX_4c(const Eigen::MatrixXd &DMAT, double error) const
Definition ERIs.h:61
void Initialize_4c(const AOBasis &dftbasis)
Definition ERIs.cc:32
Eigen::MatrixXd CalculateERIs_4c(const Eigen::MatrixXd &DMAT, double error) const
Definition ERIs.h:56
double GetWithFallback(const std::string &element, int l, int *used_l=nullptr) const
void UpdateDmats(const Eigen::MatrixXd &dmat, double DiisError, Index Iteration)
void Configure(const Eigen::MatrixXd &dmat)
void resetMatrices(Eigen::MatrixXd &J, Eigen::MatrixXd &K, const Eigen::MatrixXd &dmat)
void Start(Index iteration, double DiisError)
void UpdateCriteria(double DiisError, Index Iteration)
const Eigen::MatrixXd & getDmat_diff() const
Logger is used for thread-safe output of messages.
Definition logger.h:164
Eigen::MatrixXd & matrix()
Definition eigen.h:80
double & energy()
Definition eigen.h:81
Index cols() const
Definition eigen.h:79
Index rows() const
Definition eigen.h:78
Container for molecular orbitals and derived one-particle data.
Definition orbitals.h:47
void setScaHFX(double ScaHFX)
Store the fraction of exact exchange associated with the functional.
Definition orbitals.h:422
Index getCharge() const
Return the stored total charge.
Definition orbitals.h:254
Index getSpin() const
Return the stored spin multiplicity.
Definition orbitals.h:252
const tools::EigenSystem & MOs_beta() const
Return read-only access to beta-spin molecular orbitals.
Definition orbitals.h:202
Index getNumberOfAlphaElectrons() const
Return the stored number of alpha electrons.
Definition orbitals.h:159
void setNumberOfAlphaElectrons(Index electrons)
Store the total number of alpha electrons.
Definition orbitals.h:139
Eigen::MatrixXd DensityMatrixFull(const QMState &state) const
Definition orbitals.cc:150
void SetupAuxBasis(std::string aux_basis_name)
Definition orbitals.cc:108
void setNumberOfBetaElectrons(Index electrons)
Store the total number of beta electrons.
Definition orbitals.h:144
void setECPName(const std::string &ECP)
Store the effective core potential label.
Definition orbitals.h:170
void setXCGrid(std::string grid)
Store the numerical XC grid quality label.
Definition orbitals.h:283
void setNumberOfOccupiedLevels(Index occupied_levels)
Definition orbitals.h:122
Index getBasisSetSize() const
Return the number of AO basis functions in the DFT basis.
Definition orbitals.h:72
void setQMEnergy(double qmenergy)
Store the total DFT energy.
Definition orbitals.h:295
bool hasECPName() const
Report whether an effective core potential label has been stored.
Definition orbitals.h:164
const tools::EigenSystem & MOs() const
Return read-only access to alpha/restricted molecular orbitals.
Definition orbitals.h:192
const QMMolecule & QMAtoms() const
Return read-only access to the molecular geometry.
Definition orbitals.h:262
void ReadFromCpt(const std::string &filename)
Read the orbital container from a checkpoint file on disk.
Definition orbitals.cc:1063
void setNumberOfOccupiedLevelsBeta(Index occupied_levels_beta)
Store the number of occupied beta-spin orbitals.
Definition orbitals.h:134
void setChargeAndSpin(Index charge, Index spin)
Definition orbitals.h:246
bool hasDFTbasisName() const
Report whether a DFT basis-set name has been stored.
Definition orbitals.h:299
const std::string & getECPName() const
Return the effective core potential label.
Definition orbitals.h:167
const std::string & getDFTbasisName() const
Return the DFT basis-set name.
Definition orbitals.h:304
bool hasBetaMOs() const
Report whether beta-spin molecular orbitals are available.
Definition orbitals.h:187
void SetupDftBasis(std::string basis_name)
Build and attach the DFT AO basis from the stored molecular geometry.
Definition orbitals.cc:99
Eigen::Vector3d CalcElDipole(const QMState &state) const
Compute the electronic dipole moment associated with a state density.
Definition orbitals.cc:288
Index getNumberOfBetaElectrons() const
Return the stored number of beta electrons.
Definition orbitals.h:161
const AOBasis & getDftBasis() const
Return the DFT AO basis, throwing if it has not been initialized.
Definition orbitals.h:314
void setXCFunctionalName(std::string functionalname)
Definition orbitals.h:276
container for QM atoms
Definition qmatom.h:37
const std::string & getElement() const
Definition qmatom.h:63
Index getNuccharge() const
Definition qmatom.h:69
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:135
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
const Eigen::Vector3d & getPos() const
Definition staticsite.h:80
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
SpinDensity DensityMatrix(const tools::EigenSystem &MOs_alpha, const tools::EigenSystem &MOs_beta) const
void setOverlap(AOOverlap &S, double etol)
void Configure(const options &opt_alpha, const options &opt_beta)
SpinDensity Iterate(const SpinDensity &dmat, SpinFock &H, tools::EigenSystem &MOs_alpha, tools::EigenSystem &MOs_beta, double totE)
Index getGridSize() const
Definition vxc_grid.h:41
void GridSetup(const std::string &type, const QMMolecule &atoms, const AOBasis &basis)
Definition vxc_grid.cc:229
Index getBoxesSize() const
Definition vxc_grid.h:42
Mat_p_Energy IntegrateVXC(const Eigen::MatrixXd &density_matrix) const
void setXCfunctional(const std::string &functional)
static double getExactExchange(const std::string &functional)
SpinResult IntegrateVXCSpin(const Eigen::MatrixXd &dmat_alpha, const Eigen::MatrixXd &dmat_beta) const
Mediates interaction between polar and static sites.
double CalcStaticEnergy_site(const StaticSite &site1, const StaticSite &site2) const
#define XTP_LOG(level, log)
Definition logger.h:40
const double bohr2ang
Definition constants.h:49
Index getMaxThreads()
Definition eigen.h:128
Charge transport classes.
Definition ERIs.h:28
bool XTP_HAS_MKL_OVERLOAD()
Definition eigen.h:41
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
Spin-resolved density matrices returned for open-shell SCF updates.
Eigen::MatrixXd total() const
Return the total density P = P^alpha + P^beta.