votca 2024.2-dev
Loading...
Searching...
No Matches
embeddingengine.cc
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 *References- (1) A Simple, Exact Density-Functional-Theory
21 *Embedding Scheme. Frederick R. Manby, Martina Stella, Jason D. Goodpaster, and
22 *Thomas F. Miller Journal of Chemical Theory and Computation 2012 8 (8),
23 *2564-2568 DOI: 10.1021/ct300544e
24 *(2) Taylor A. Barnes, Jason D. Goodpaster, Frederick R. Manby, and Thomas F.
25 *Miller III , "Accurate basis set truncation for wavefunction embedding", J.
26 *Chem. Phys. 139, 024103 (2013) https://doi.org/10.1063/1.4811112
27 *(3) Simon J. Bennie, Martina Stella, Thomas F. Miller III, and Frederick R.
28 *Manby , "Accelerating wavefunction in density-functional-theory embedding by
29 *truncating the active basis set", J. Chem. Phys. 143, 024105 (2015)
30 *https://doi.org/10.1063/1.4923367
31 *(3) Projection-Based Wavefunction-in-DFT Embedding
32 *Sebastian J. R. Lee, Matthew Welborn, Frederick R. Manby, and Thomas F. Miller
33 *III Accounts of Chemical Research 2019 52 (5), 1359-1368
34 *DOI: 10.1021/acs.accounts.8b00672
35 */
36
37// Third party includes
38#include <algorithm>
39#include <boost/filesystem.hpp>
40#include <boost/format.hpp>
41
42// VOTCA includes
45
46// Local VOTCA includes
50#include "votca/xtp/aomatrix.h"
53#include "votca/xtp/dftengine.h"
55#include "votca/xtp/logger.h"
56#include "votca/xtp/mmregion.h"
57#include "votca/xtp/orbitals.h"
59
60namespace votca {
61namespace xtp {
62
64
65 // reading in the orbitals of the full DFT calculation
66 tools::EigenSystem embeddingMOs = orb.MOs();
67
68 // constructing the full electron density matrix
69 const Eigen::MatrixXd FullDensityMatrix = orb.DensityMatrixGroundState();
70 XTP_LOG(Log::error, *pLog_) << std::flush;
72 << TimeStamp() << " Starting Projection based dft-in-dft embedding"
73 << std::flush;
74
75 // reading the localized orbitals and update the occupied orbitals in MOs
77 << TimeStamp() << " Passing localized orbitals as the initial guess"
78 << std::flush;
79 Eigen::MatrixXd LMOs = orb.getLMOs();
80
81 embeddingMOs.eigenvectors().leftCols(orb.getNumberOfAlphaElectrons()) = LMOs;
82
83 // determine the active and inactive electron densities
84 std::vector<Index> activeatoms =
86
88 << "Indices of active atoms selected are: " << active_atoms_as_string_
89 << std::flush;
90 ActiveDensityMatrix DMAT_A(orb, activeatoms, active_threshold_);
91 const Eigen::MatrixXd InitialActiveDensityMatrix = DMAT_A.compute_Dmat_A()[0];
92 Eigen::MatrixXd InitialInactiveMOs = DMAT_A.compute_Dmat_A()[2];
93
95 << TimeStamp() << " Active density formation done" << std::flush;
96
97 Eigen::MatrixXd InactiveDensityMatrix =
98 FullDensityMatrix - InitialActiveDensityMatrix;
99
100 orb.setInactiveDensity(InactiveDensityMatrix);
101
102 // AOoverlap to calculate number of active/inactive electrons
103 AOBasis aobasis = orb.getDftBasis();
104 AOOverlap overlap;
105 overlap.Fill(aobasis);
106
107 Index all_electrons = static_cast<Index>(
108 std::round(FullDensityMatrix.cwiseProduct(overlap.Matrix()).sum()));
109 active_electrons_ = static_cast<Index>(std::round(
110 InitialActiveDensityMatrix.cwiseProduct(overlap.Matrix()).sum()));
111 Index inactive_electrons = static_cast<Index>(
112 std::round(InactiveDensityMatrix.cwiseProduct(overlap.Matrix()).sum()));
113
115 << TimeStamp() << " Total electrons: " << all_electrons << std::flush;
117 << TimeStamp() << " Active electrons: " << active_electrons_
118 << std::flush;
120 << TimeStamp() << " Inactive electrons: " << inactive_electrons
121 << std::flush;
122
123 // check for consistency
124 if ((active_electrons_ + inactive_electrons) != all_electrons) {
125 throw std::runtime_error(
126 " Sum of active and inactive electrons does "
127 "not match full number of electrons!");
128 return false;
129 }
130
131 // setup the DFT engine with the active electrons only
133
134 // setup one-electron part of the Hamiltonian
135 Mat_p_Energy H0 = SetupH0(orb.QMAtoms());
136
137 // energy of the one-electron part of the Hamiltonian
138 const double E0_full = FullDensityMatrix.cwiseProduct(H0.matrix()).sum();
139 const double E0_initial_active =
140 InitialActiveDensityMatrix.cwiseProduct(H0.matrix()).sum();
141 E_nuc_ = H0.energy();
142
143 // setup the Vxc integrator
144 Vxc_Potential<Vxc_Grid> vxcpotential = SetupVxc(orb.QMAtoms());
145 // get the constant XC contributions for the embedding
146 const Mat_p_Energy xc_full = vxcpotential.IntegrateVXC(FullDensityMatrix);
147 const Mat_p_Energy xc_initial_active =
148 vxcpotential.IntegrateVXC(InitialActiveDensityMatrix);
149
150 // get the constant Coulomb matrices for the embedding
151 Eigen::MatrixXd J_full =
152 Eigen::MatrixXd::Zero(FullDensityMatrix.rows(), FullDensityMatrix.cols());
153 Eigen::MatrixXd J_initial_active =
154 Eigen::MatrixXd::Zero(FullDensityMatrix.rows(), FullDensityMatrix.cols());
155 Eigen::MatrixXd K_full =
156 Eigen::MatrixXd::Zero(FullDensityMatrix.rows(), FullDensityMatrix.cols());
157 Eigen::MatrixXd K_initial_active =
158 Eigen::MatrixXd::Zero(FullDensityMatrix.rows(), FullDensityMatrix.cols());
159 if (ScaHFX_ > 0) {
160 std::array<Eigen::MatrixXd, 2> JandK_full =
161 CalcERIs_EXX(embeddingMOs.eigenvectors(), FullDensityMatrix, 1e-12);
162 std::array<Eigen::MatrixXd, 2> JandK_initial_active = CalcERIs_EXX(
163 embeddingMOs.eigenvectors(), InitialActiveDensityMatrix, 1e-12);
164 J_full = JandK_full[0];
165 K_full = 0.5 * ScaHFX_ * JandK_full[1];
166 J_initial_active = JandK_initial_active[0];
167 K_initial_active = 0.5 * ScaHFX_ * JandK_initial_active[1];
168 } else {
169 J_full = CalcERIs(FullDensityMatrix, 1e-12);
170 J_initial_active = CalcERIs(InitialActiveDensityMatrix, 1e-12);
171 }
172 // get the Hartree energies
173 const double E_Hartree_full =
174 0.5 * FullDensityMatrix.cwiseProduct(J_full + K_full).sum();
175 const double E_Hartree_initial_active =
176 0.5 * InitialActiveDensityMatrix
177 .cwiseProduct(J_initial_active + K_initial_active)
178 .sum();
179
180 // Energy of the full reference system
181 Total_E_full_ = E0_full + E_Hartree_full + xc_full.energy() + E_nuc_;
182
184 << TimeStamp()
185 << " Reference total energy full electron system: " << Total_E_full_
186 << " Ha" << std::flush;
187
188 // projection parameter, to be made an option
189 Eigen::MatrixXd ProjectionOperator =
190 overlap.Matrix() * InactiveDensityMatrix * overlap.Matrix();
191
192 // XC and Hartree contribution to the external embedding potential/energy
193 const Eigen::MatrixXd v_embedding = J_full + K_full + xc_full.matrix() -
194 J_initial_active - K_initial_active -
195 xc_initial_active.matrix();
196
197 const double constant_embedding_energy = Total_E_full_ - E0_initial_active -
198 E_Hartree_initial_active -
199 xc_initial_active.energy();
200
202 << TimeStamp() << " Constant energy embedding terms: " << std::flush;
204 << TimeStamp() << " 1) E_DFT[full]: " << Total_E_full_ << std::flush;
206 << TimeStamp() << " 2) E_0[initial active]: " << E0_initial_active
207 << std::flush;
209 << TimeStamp() << " 3) E_H[initial active]: " << E_Hartree_initial_active
210 << std::flush;
211 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " 4) E_XC[initial active]: "
212 << xc_initial_active.energy() << std::flush;
214 << TimeStamp() << " Total 1 -2 -3 -4: " << constant_embedding_energy
215 << std::flush;
216
217 if (truncate_) { // Truncation starts here
218 orb.setCalculationType("Truncated");
219 TruncateBasis(orb, activeatoms, H0, InitialActiveDensityMatrix, v_embedding,
220 InitialInactiveMOs);
221 active_and_border_atoms_ = activeatoms;
222 }
223 // SCF loop if you don't truncate active region
224 else {
225 orb.setCalculationType("Embedded_noTrunc");
226 Eigen::MatrixXd ActiveDensityMatrix = InitialActiveDensityMatrix;
227 for (Index this_iter = 0; this_iter < max_iter_; this_iter++) {
228 XTP_LOG(Log::error, *pLog_) << std::flush;
230 << TimeStamp() << " Iteration " << this_iter + 1 << " of "
231 << max_iter_ << std::flush;
232
233 // get the XC contribution for the updated active density matrix
234 Mat_p_Energy xc_active = vxcpotential.IntegrateVXC(ActiveDensityMatrix);
235 double E_xc_active = xc_active.energy();
236
237 // get the Hartree contribution for the updated active density matrix
238 Eigen::MatrixXd J_active = Eigen::MatrixXd::Zero(
240 Eigen::MatrixXd K_active = Eigen::MatrixXd::Zero(
242 if (ScaHFX_ > 0) {
243 std::array<Eigen::MatrixXd, 2> JandK_active = CalcERIs_EXX(
244 embeddingMOs.eigenvectors(), ActiveDensityMatrix, 1e-12);
245 J_active = JandK_active[0];
246 K_active = 0.5 * ScaHFX_ * JandK_active[1];
247
248 } else {
249 J_active = CalcERIs(ActiveDensityMatrix, 1e-12);
250 }
251 double E_Hartree_active =
252 0.5 * ActiveDensityMatrix.cwiseProduct(J_active + K_active).sum();
253
254 // update the active Hamiltonian
255 const Eigen::MatrixXd Level_Shift_Operator =
256 levelshift_ * ProjectionOperator;
257 Eigen::MatrixXd H_active = H0.matrix() + v_embedding +
258 Level_Shift_Operator + xc_active.matrix() +
259 J_active + K_active;
260
261 const double E0_active =
262 ActiveDensityMatrix.cwiseProduct(H0.matrix()).sum();
263 const double E_embedding =
264 (ActiveDensityMatrix - InitialActiveDensityMatrix)
265 .cwiseProduct(v_embedding)
266 .sum();
267 const double E_levelshift =
268 ActiveDensityMatrix.cwiseProduct(Level_Shift_Operator).sum();
269
270 double TotalEnergy = E0_active + E_Hartree_active + E_xc_active +
271 constant_embedding_energy + E_levelshift +
272 E_embedding;
273
275 << TimeStamp() << " Total Energy after embedding: " << TotalEnergy
276 << std::flush;
278 << TimeStamp() << " I) E_0[active]: " << E0_active << std::flush;
280 << TimeStamp() << " II) E_H[active]: " << E_Hartree_active
281 << std::flush;
283 << TimeStamp() << " III) E_xc[active]: " << E_xc_active << std::flush;
285 << TimeStamp() << " IV) E_ext_emb :" << constant_embedding_energy
286 << std::flush;
288 << TimeStamp() << " V) E_corr_emb :" << E_embedding << std::flush;
290 << TimeStamp() << " VI) E_levelshift :" << E_levelshift
291 << std::flush;
293 << TimeStamp() << " Total I + II + III + IV + V +VI :" << TotalEnergy
294 << std::flush;
295
296 // get a new density matrix in the active region
298 << std::endl
299 << "Active electrons in = "
300 << ActiveDensityMatrix.cwiseProduct(overlap.Matrix()).sum()
301 << std::flush;
303 ActiveDensityMatrix, H_active, embeddingMOs, TotalEnergy);
305 << std::endl
306 << "Active electrons out= "
307 << ActiveDensityMatrix.cwiseProduct(overlap.Matrix()).sum()
308 << std::flush;
309
310 PrintMOs(embeddingMOs.eigenvalues(), Log::info);
311
314 << TimeStamp() << " Total embedding energy has converged to "
315 << std::setprecision(9) << conv_accelerator_.getDeltaE()
316 << "[Ha] after " << this_iter + 1
317 << " iterations. DIIS error is converged up to "
318 << conv_accelerator_.getDIIsError() << std::flush;
320 << TimeStamp() << " Final Single Point Energy of embedding "
321 << std::setprecision(12) << TotalEnergy << " Ha" << std::flush;
323 << "The difference between the reference energy and the embedding "
324 "energy is "
325 << Total_E_full_ - TotalEnergy << " Ha" << std::flush;
326
327 if (abs(Total_E_full_ - TotalEnergy) > 1e-3) {
329 << "Warning!! The difference is greater than 1e-03 Ha"
330 << std::flush;
331 }
332
333 PrintMOs(embeddingMOs.eigenvalues(), Log::error);
334 orb.setEmbeddedMOs(embeddingMOs);
336 CalcElDipole(orb);
337 break;
338 } else if (this_iter == max_iter_ - 1) {
340 << TimeStamp()
341 << " DFT calculation for active region has not converged after "
342 << max_iter_
343 << " iterations. Use more iterations or another convergence "
344 "acceleration scheme."
345 << std::flush;
346 return false;
347 }
348 }
349 }
350 return true;
351}
352
354 if (truncate_) {
356 trunc_orb.QMAtoms() = activemol_;
357 Prepare(trunc_orb, active_electrons_);
358 Vxc_Potential<Vxc_Grid> vxcpotential = SetupVxc(trunc_orb.QMAtoms());
359 ConfigOrbfile(trunc_orb);
360 AOBasis aobasis = trunc_orb.getDftBasis();
361 AOOverlap overlap;
362 overlap.Fill(aobasis);
363
364 const double E0_initial_truncated =
365 InitialActiveDmat_trunc_.cwiseProduct(H0_trunc_).sum();
366 const Mat_p_Energy xc_initial_truncated =
368 Eigen::MatrixXd J_initial_truncated = Eigen::MatrixXd::Zero(
370 Eigen::MatrixXd K_initial_truncated = Eigen::MatrixXd::Zero(
372 tools::EigenSystem MOs_trunc;
373 MOs_trunc.eigenvalues() = Eigen::VectorXd::Zero(H0_trunc_.cols());
374 MOs_trunc.eigenvectors() =
375 Eigen::MatrixXd::Zero(H0_trunc_.rows(), H0_trunc_.cols());
376 if (ScaHFX_ > 0) {
377 std::array<Eigen::MatrixXd, 2> JandK_initial_truncated = CalcERIs_EXX(
378 MOs_trunc.eigenvectors(), InitialActiveDmat_trunc_, 1e-12);
379 J_initial_truncated = JandK_initial_truncated[0];
380 K_initial_truncated = 0.5 * ScaHFX_ * JandK_initial_truncated[1];
381 } else {
382 J_initial_truncated = CalcERIs(InitialActiveDmat_trunc_, 1e-12);
383 }
384 // get the Hartree energies
385 const double E_Hartree_initial_truncated =
387 .cwiseProduct(J_initial_truncated + K_initial_truncated)
388 .sum();
389 const double Initial_truncated_energy = E0_initial_truncated +
390 E_Hartree_initial_truncated +
391 xc_initial_truncated.energy();
393 << "Initial truncated energy = " << Initial_truncated_energy
394 << std::flush;
395
396 Eigen::MatrixXd PurifiedActiveDmat_trunc =
398 Eigen::MatrixXd TruncatedDensityMatrix = PurifiedActiveDmat_trunc;
399
400 for (Index this_iter = 0; this_iter < max_iter_; this_iter++) {
401 XTP_LOG(Log::error, *pLog_) << std::flush;
403 << TimeStamp() << " Iteration " << this_iter + 1 << " of "
404 << max_iter_ << std::flush;
405
406 // get the XC contribution for the updated truncated density matrix
407 Mat_p_Energy xc_truncated =
408 vxcpotential.IntegrateVXC(TruncatedDensityMatrix);
409 double E_xc_truncated = xc_truncated.energy();
410
411 // get the Hartree contribution for the updated truncated density matrix
412 Eigen::MatrixXd J_truncated = Eigen::MatrixXd::Zero(
413 TruncatedDensityMatrix.rows(), TruncatedDensityMatrix.cols());
414 Eigen::MatrixXd K_truncated = Eigen::MatrixXd::Zero(
415 TruncatedDensityMatrix.rows(), TruncatedDensityMatrix.cols());
416 if (ScaHFX_ > 0) {
417 std::array<Eigen::MatrixXd, 2> JandK_truncated = CalcERIs_EXX(
418 MOs_trunc.eigenvectors(), TruncatedDensityMatrix, 1e-12);
419 J_truncated = JandK_truncated[0];
420 K_truncated = 0.5 * ScaHFX_ * JandK_truncated[1];
421
422 } else {
423 J_truncated = CalcERIs(TruncatedDensityMatrix, 1e-12);
424 }
425 double E_Hartree_truncated =
426 0.5 *
427 TruncatedDensityMatrix.cwiseProduct(J_truncated + K_truncated).sum();
428 // update the truncated Hamiltonian
429 Eigen::MatrixXd H_truncated =
430 H0_trunc_ + xc_truncated.matrix() + J_truncated + K_truncated;
431
432 double TruncatedTotalEnergy =
433 TruncatedDensityMatrix.cwiseProduct(H0_trunc_).sum() +
434 E_Hartree_truncated + E_xc_truncated + Total_E_full_ -
435 Initial_truncated_energy +
436 ((TruncatedDensityMatrix - InitialActiveDmat_trunc_)
437 .cwiseProduct(v_embedding_trunc_)
438 .sum());
439 // get the new truncated density matrix
441 << "Electrons in = "
442 << TruncatedDensityMatrix.cwiseProduct(overlap.Matrix()).sum()
443 << std::flush;
444 TruncatedDensityMatrix = conv_accelerator_.Iterate(
445 TruncatedDensityMatrix, H_truncated, MOs_trunc, TruncatedTotalEnergy);
447 << "Electrons out = "
448 << TruncatedDensityMatrix.cwiseProduct(overlap.Matrix()).sum()
449 << std::flush;
450
452 << TimeStamp() << " E_trunc (Ha) = "
453 << TruncatedDensityMatrix.cwiseProduct(H0_trunc_).sum() +
454 E_Hartree_truncated + E_xc_truncated
455 << "\n \t \t \t"
456 << " E_fullDFT (Ha) = " << Total_E_full_ << "\n \t \t \t"
457 << " E_initial_trunc(Ha) = " << Initial_truncated_energy
458 << "\n \t \t \t"
459 << " E_embedding_correction(Ha) = "
460 << ((TruncatedDensityMatrix - InitialActiveDmat_trunc_)
461 .cwiseProduct(v_embedding_trunc_)
462 .sum())
463 << std::flush;
464
466 << " Truncated Energy of the system: E_trunc + E_fullDFT - "
467 << "\n \t \t"
468 << " E_initial_trunc + E_embedding_correction = "
469 << TruncatedTotalEnergy << std::flush;
470
471 PrintMOs(MOs_trunc.eigenvalues(), Log::info);
472
475 << TimeStamp() << " Total embedding energy has converged to "
476 << std::setprecision(9) << conv_accelerator_.getDeltaE()
477 << "[Ha] after " << this_iter + 1
478 << " iterations. DIIS error is converged up to "
479 << conv_accelerator_.getDIIsError() << std::flush;
481 << TimeStamp() << " Final Single Point Energy of embedding "
482 << std::setprecision(12) << TruncatedTotalEnergy << " Ha"
483 << std::flush;
484
485 PrintMOs(MOs_trunc.eigenvalues(), Log::error);
486 trunc_orb.setEmbeddedMOs(MOs_trunc);
488 break;
489 }
490 }
492 }
493 return true;
494}
495
496void DFTEngine::TruncateBasis(Orbitals& orb, std::vector<Index>& activeatoms,
497 Mat_p_Energy& H0,
498 Eigen::MatrixXd InitialActiveDensityMatrix,
499 Eigen::MatrixXd v_embedding,
500 Eigen::MatrixXd InitialInactiveMOs) {
501 AOBasis aobasis;
502 aobasis = orb.getDftBasis();
503 AOOverlap overlap;
504 overlap.Fill(aobasis);
505 // Mulliken population per basis function on every active atom
506 Eigen::VectorXd DiagonalofDmatA = InitialActiveDensityMatrix.diagonal();
507 Eigen::VectorXd DiagonalofOverlap = overlap.Matrix().diagonal();
508 Eigen::VectorXd MnP = DiagonalofDmatA.cwiseProduct(DiagonalofOverlap);
509
510 // Get a vector containing the number of basis functions per atom
511 const std::vector<Index>& numfuncpatom = aobasis.getFuncPerAtom();
512 numfuncpatom_ = numfuncpatom;
513 Index numofactivebasisfunction = 0;
514 // Store start indices. Will be used later
515 std::vector<Index> start_indices, start_indices_activemolecule;
516 Index start_idx = 0, start_idx_activemolecule = 0;
517 std::vector<Index> borderatoms;
518 for (Index atom_num = 0; atom_num < orb.QMAtoms().size(); atom_num++) {
519 start_indices.push_back(start_idx);
520 /* Condition for basis fn on an atom to be counted: either in active
521 region or MnP of any function > threshold (border atoms) */
522 bool partOfActive =
523 (std::find(activeatoms.begin(), activeatoms.end(),
524 orb.QMAtoms()[atom_num].getId()) != activeatoms.end());
525 if (partOfActive == true) {
526 activemol_.push_back(orb.QMAtoms()[atom_num]);
527 /* if active append the atom to the new molecule. Also increment active
528 basis function size and start indices */
529 start_indices_activemolecule.push_back(start_idx_activemolecule);
530 start_idx_activemolecule += numfuncpatom[atom_num];
531 numofactivebasisfunction += numfuncpatom[atom_num];
532 } else {
533 std::vector<const AOShell*> inactiveshell =
534 aobasis.getShellsofAtom(orb.QMAtoms()[atom_num].getId());
535 bool loop_break = false;
536 for (const AOShell* shell : inactiveshell) {
537 for (Index shell_fn_no = shell->getStartIndex();
538 shell_fn_no < shell->getStartIndex() + shell->getNumFunc();
539 shell_fn_no++) {
540 if (MnP[shell_fn_no] > truncation_threshold_) {
541 activeatoms.push_back(atom_num);
542 borderatoms.push_back(atom_num); // push this index to border
543 // atoms
544 activemol_.push_back(orb.QMAtoms()[atom_num]);
545 loop_break = true; // if any function in the whole atom satisfies
546 // we break loop to check next atom
547 start_indices_activemolecule.push_back(start_idx_activemolecule);
548 start_idx_activemolecule += numfuncpatom[atom_num];
549 numofactivebasisfunction += numfuncpatom[atom_num];
550 break;
551 }
552 }
553 if (loop_break) break;
554 }
555 }
556 start_idx += numfuncpatom[atom_num];
557 }
558
559 // Sort atoms before you cut the Hamiltonian
560 sort(activeatoms.begin(), activeatoms.end());
561 sort(borderatoms.begin(), borderatoms.end());
562 Eigen::MatrixXd InactiveDensityMatrix = orb.getInactiveDensity();
563 XTP_LOG(Log::error, *pLog_) << std::flush;
565 << "Active + Border Molecule Size = " << activeatoms.size() << "\n \t \t "
566 << "Border Molecule Size = " << borderatoms.size() << std::flush;
567 Eigen::MatrixXd ProjectionOperator =
568 overlap.Matrix() * InactiveDensityMatrix * overlap.Matrix();
569
570 Eigen::MatrixXd H_embedding =
571 H0.matrix() + v_embedding + levelshift_ * ProjectionOperator;
572
573 if (borderatoms.size() != 0) {
574 Eigen::MatrixXd BorderMOs;
575 for (Index lmo_index = 0; lmo_index < InitialInactiveMOs.cols();
576 lmo_index++) {
577 double mullikenpop_lmo_borderatoms = 0;
578 for (Index borderatom : borderatoms) {
579 Index start = start_indices[borderatom];
580 Index size = numfuncpatom[borderatom];
581 mullikenpop_lmo_borderatoms +=
582 (InitialInactiveMOs.col(lmo_index) *
583 InitialInactiveMOs.col(lmo_index).transpose() * overlap.Matrix())
584 .diagonal()
585 .segment(start, size)
586 .sum();
587 }
588 /*If more than half of a MO contributes on a border atom include that in
589 * the Border MOs list. Not made an user option. If you want to fix
590 * further for density leak you can reduce this number*/
591 if (mullikenpop_lmo_borderatoms > 0.25) {
592 BorderMOs.conservativeResize(InitialInactiveMOs.rows(),
593 BorderMOs.cols() + 1);
594 BorderMOs.col(BorderMOs.cols() - 1) = InitialInactiveMOs.col(lmo_index);
595 }
596 }
597
598 Eigen::MatrixXd BorderDmat = 2 * BorderMOs * BorderMOs.transpose();
599 Eigen::MatrixXd BorderProjectionOperator =
600 overlap.Matrix() * BorderDmat * overlap.Matrix();
601 Eigen::MatrixXd DistantProjectionOperator =
602 ProjectionOperator - BorderProjectionOperator;
603 // 1e+2 is the suitable projection value for border MOs: see citations
604 H_embedding +=
605 1e+2 * BorderProjectionOperator +
606 levelshift_ * (DistantProjectionOperator - ProjectionOperator);
607 }
608
609 // from here it is time to truncate Hamiltonian H0
610 H0_trunc_ =
611 Eigen::MatrixXd::Zero(numofactivebasisfunction, numofactivebasisfunction);
613 Eigen::MatrixXd::Zero(numofactivebasisfunction, numofactivebasisfunction);
615 Eigen::MatrixXd::Zero(numofactivebasisfunction, numofactivebasisfunction);
616 // Defined needed matrices with basisset_size * basisset_size
617
618 for (Index activeatom1_idx = 0; activeatom1_idx < Index(activeatoms.size());
619 activeatom1_idx++) {
620 Index activeatom1 = activeatoms[activeatom1_idx];
621 for (Index activeatom2_idx = 0; activeatom2_idx < Index(activeatoms.size());
622 activeatom2_idx++) {
623 Index activeatom2 = activeatoms[activeatom2_idx];
624 // Fill the H0, Dmat and v_embedding at the right indices
625 H0_trunc_.block(start_indices_activemolecule[activeatom1_idx],
626 start_indices_activemolecule[activeatom2_idx],
627 numfuncpatom[activeatom1], numfuncpatom[activeatom2]) =
628 H_embedding.block(
629 start_indices[activeatom1], start_indices[activeatom2],
630 numfuncpatom[activeatom1], numfuncpatom[activeatom2]);
632 start_indices_activemolecule[activeatom1_idx],
633 start_indices_activemolecule[activeatom2_idx],
634 numfuncpatom[activeatom1], numfuncpatom[activeatom2]) =
635 InitialActiveDensityMatrix.block(
636 start_indices[activeatom1], start_indices[activeatom2],
637 numfuncpatom[activeatom1], numfuncpatom[activeatom2]);
638 v_embedding_trunc_.block(start_indices_activemolecule[activeatom1_idx],
639 start_indices_activemolecule[activeatom2_idx],
640 numfuncpatom[activeatom1],
641 numfuncpatom[activeatom2]) =
642 v_embedding.block(
643 start_indices[activeatom1], start_indices[activeatom2],
644 numfuncpatom[activeatom1], numfuncpatom[activeatom2]);
645 }
646 }
647}
648
649Eigen::MatrixXd DFTEngine::McWeenyPurification(Eigen::MatrixXd& Dmat_in,
650 AOOverlap& overlap) {
651 Eigen::MatrixXd Ssqrt = overlap.Sqrt();
652 Eigen::MatrixXd InvSsqrt = overlap.Pseudo_InvSqrt(1e-8);
653 Eigen::MatrixXd ModifiedDmat = 0.5 * Ssqrt * Dmat_in * Ssqrt;
654 for (Index iter = 0; iter < 100; iter++) {
655 Eigen::MatrixXd Dmat_new = (3 * ModifiedDmat * ModifiedDmat) -
656 (2 * ModifiedDmat * ModifiedDmat * ModifiedDmat);
657 double IdempotencyError =
658 ((Dmat_new * Dmat_new - Dmat_new) * (Dmat_new * Dmat_new - Dmat_new))
659 .trace();
661 << "Idempotency Error: " << IdempotencyError << std::flush;
662
663 ModifiedDmat = Dmat_new;
664 if (IdempotencyError < 1e-20) break;
665 }
666 Eigen::MatrixXd Dmat_out = 2 * InvSsqrt * ModifiedDmat * InvSsqrt;
667 return 0.5 * (Dmat_out + Dmat_out.transpose());
668}
669
670void DFTEngine::TruncMOsFullBasis(Orbitals& orb, std::vector<Index> activeatoms,
671 std::vector<Index> numfuncpatom) {
672 Eigen::MatrixXd expandtruncorb = orb.getEmbeddedMOs().eigenvectors();
673 Index start_index = 0;
674 Index numofactualactoms = numfuncpatom.size();
675 for (Index atomindex = 0; atomindex < numofactualactoms; atomindex++) {
676 bool partOfActive = (std::find(activeatoms.begin(), activeatoms.end(),
677 atomindex) != activeatoms.end());
678 if (partOfActive == false) {
679 expandtruncorb =
680 InsertZeroRows(expandtruncorb, start_index, numfuncpatom[atomindex]);
681 }
682 start_index += numfuncpatom[atomindex];
683 }
684 orb.setTruncMOsFullBasis(expandtruncorb);
685}
686
687Eigen::MatrixXd DFTEngine::InsertZeroRows(Eigen::MatrixXd MOsMatrix,
688 Index startidx, Index numofzerorows) {
689 Eigen::MatrixXd FinalMatrix =
690 Eigen::MatrixXd::Zero(MOsMatrix.rows() + numofzerorows, MOsMatrix.cols());
691 FinalMatrix.topRows(startidx) = MOsMatrix.topRows(startidx);
692 FinalMatrix.bottomRows(MOsMatrix.rows() - startidx) =
693 MOsMatrix.bottomRows(MOsMatrix.rows() - startidx);
694 return FinalMatrix;
695}
696} // namespace xtp
697} // namespace votca
const Eigen::VectorXd & eigenvalues() const
Definition eigensystem.h:30
const Eigen::MatrixXd & eigenvectors() const
Definition eigensystem.h:33
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
const std::vector< const AOShell * > getShellsofAtom(Index AtomId) const
Definition aobasis.cc:63
const std::vector< Index > & getFuncPerAtom() const
Definition aobasis.h:72
void Fill(const AOBasis &aobasis) final
Eigen::MatrixXd Pseudo_InvSqrt(double etol)
Definition aomatrix.cc:28
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:52
Eigen::MatrixXd Sqrt()
Definition aomatrix.cc:45
std::array< Eigen::MatrixXd, 3 > compute_Dmat_A()
void push_back(const T &atom)
Eigen::MatrixXd Iterate(const Eigen::MatrixXd &dmat, Eigen::MatrixXd &H, tools::EigenSystem &MOs, double totE)
Eigen::MatrixXd CalcERIs(const Eigen::MatrixXd &Dmat, double error) const
Definition dftengine.cc:174
bool EvaluateTruncatedActiveRegion(Orbitals &trunc_orb)
Eigen::MatrixXd InitialActiveDmat_trunc_
Definition dftengine.h:179
void TruncMOsFullBasis(Orbitals &orb, std::vector< Index > activeatoms, std::vector< Index > numfuncpatom)
void TruncateBasis(Orbitals &orb, std::vector< Index > &activeatoms, Mat_p_Energy &H0, Eigen::MatrixXd InitialActiveDensityMatrix, Eigen::MatrixXd v_embedding, Eigen::MatrixXd InitialInactiveMOs)
bool EvaluateActiveRegion(Orbitals &orb)
void Prepare(Orbitals &orb, Index numofelectrons=-1)
Definition dftengine.cc:801
Eigen::MatrixXd InsertZeroRows(Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerorows)
std::array< Eigen::MatrixXd, 2 > CalcERIs_EXX(const Eigen::MatrixXd &MOCoeff, const Eigen::MatrixXd &Dmat, double error) const
Definition dftengine.cc:159
Eigen::MatrixXd H0_trunc_
Definition dftengine.h:178
void CalcElDipole(const Orbitals &orb) const
Definition dftengine.cc:150
std::string active_atoms_as_string_
Definition dftengine.h:173
void ConfigOrbfile(Orbitals &orb)
Definition dftengine.cc:738
QMMolecule activemol_
Definition dftengine.h:167
void PrintMOs(const Eigen::VectorXd &MOEnergies, Log::Level level)
Definition dftengine.cc:134
std::vector< Index > numfuncpatom_
Definition dftengine.h:187
std::vector< Index > active_and_border_atoms_
Definition dftengine.h:186
Mat_p_Energy SetupH0(const QMMolecule &mol) const
Definition dftengine.cc:369
Eigen::MatrixXd v_embedding_trunc_
Definition dftengine.h:180
Eigen::MatrixXd McWeenyPurification(Eigen::MatrixXd &Dmat_in, AOOverlap &overlap)
Vxc_Potential< Vxc_Grid > SetupVxc(const QMMolecule &mol)
Definition dftengine.cc:881
ConvergenceAcc conv_accelerator_
Definition dftengine.h:150
std::vector< Index > CreateIndexVector(const std::string &Ids) const
Eigen::MatrixXd & matrix()
Definition eigen.h:80
double & energy()
Definition eigen.h:81
container for molecular orbitals
Definition orbitals.h:46
void setCalculationType(std::string CalcType)
Definition orbitals.h:158
void setInactiveDensity(Eigen::MatrixXd inactivedensity)
Definition orbitals.h:425
void setNumofActiveElectrons(const Index active_electrons)
Definition orbitals.h:420
void setTruncMOsFullBasis(const Eigen::MatrixXd &expandedMOs)
Definition orbitals.h:58
Index getNumberOfAlphaElectrons() const
Definition orbitals.h:93
const Eigen::MatrixXd & getLMOs() const
Definition orbitals.h:411
const Eigen::MatrixXd & getInactiveDensity() const
Definition orbitals.h:424
const tools::EigenSystem & MOs() const
Definition orbitals.h:122
Eigen::MatrixXd DensityMatrixGroundState() const
Definition orbitals.cc:183
const QMMolecule & QMAtoms() const
Definition orbitals.h:172
void setEmbeddedMOs(tools::EigenSystem &system)
Definition orbitals.h:54
const tools::EigenSystem & getEmbeddedMOs() const
Definition orbitals.h:56
const AOBasis & getDftBasis() const
Definition orbitals.h:208
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
Mat_p_Energy IntegrateVXC(const Eigen::MatrixXd &density_matrix) const
#define XTP_LOG(level, log)
Definition logger.h:40
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26