votca 2024.2-dev
Loading...
Searching...
No Matches
dftcoupling.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 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 <boost/format.hpp>
22
23// VOTCA includes
25
26// Local VOTCA includes
27#include "votca/xtp/aomatrix.h"
29
30namespace votca {
31namespace xtp {
32
33using boost::format;
34using std::flush;
35
37
38 degeneracy_ = options.get("degeneracy").as<double>();
40 numberofstatesA_ = options.get("levA").as<Index>();
41 numberofstatesB_ = options.get("levB").as<Index>();
42}
43
45 const Orbitals& orbitalsA,
46 const Orbitals& orbitalsB, Index a,
47 Index b) const {
48 double J = getCouplingElement(a, b, orbitalsA, orbitalsB);
49 tools::Property& coupling = type_summary.add("coupling", "");
50 coupling.setAttribute("levelA", a);
51 coupling.setAttribute("levelB", b);
52 coupling.setAttribute("j", (format("%1$1.6e") % J).str());
53}
54
56 const Orbitals& orbitalsA,
57 const Orbitals& orbitalsB) const {
58 tools::Property& dftcoupling = type_summary.add(Identify(), "");
59 dftcoupling.setAttribute("homoA", orbitalsA.getHomo());
60 dftcoupling.setAttribute("homoB", orbitalsB.getHomo());
61 tools::Property& hole_summary = dftcoupling.add("hole", "");
62 // hole hole
63 for (Index a = Range_orbA.first; a <= orbitalsA.getHomo(); ++a) {
64 for (Index b = Range_orbB.first; b <= orbitalsB.getHomo(); ++b) {
65 WriteToProperty(hole_summary, orbitalsA, orbitalsB, a, b);
66 }
67 }
68 tools::Property& electron_summary = dftcoupling.add("electron", "");
69 // electron-electron
70 for (Index a = orbitalsA.getLumo();
71 a <= Range_orbA.first + Range_orbA.second - 1; ++a) {
72 for (Index b = orbitalsB.getLumo();
73 b <= Range_orbB.first + Range_orbB.second - 1; ++b) {
74 WriteToProperty(electron_summary, orbitalsA, orbitalsB, a, b);
75 }
76 }
77 return;
78}
79
81 const Orbitals& orbital, Index numberofstates) const {
82 const Eigen::VectorXd& MOEnergies = orbital.MOs().eigenvalues();
83 if (std::abs(MOEnergies(orbital.getHomo()) - MOEnergies(orbital.getLumo())) <
85 throw std::runtime_error(
86 "Homo Lumo Gap is smaller than degeneracy. "
87 "Either your degeneracy is too large or your Homo and Lumo are "
88 "degenerate");
89 }
90
91 Index minimal = orbital.getHomo() - numberofstates + 1;
92 Index maximal = orbital.getLumo() + numberofstates - 1;
93
94 std::vector<Index> deg_min = orbital.CheckDegeneracy(minimal, degeneracy_);
95 minimal = *std::min_element(deg_min.begin(), deg_min.end());
96
97 std::vector<Index> deg_max = orbital.CheckDegeneracy(maximal, degeneracy_);
98 maximal = *std::max_element(deg_max.begin(), deg_max.end());
99
100 std::pair<Index, Index> result;
101 result.first = minimal; // start
102 result.second = maximal - minimal + 1; // size
103
104 return result;
105}
106
108 const Orbitals& orbitalsA,
109 const Orbitals& orbitalsB) const {
110
111 Index levelsA = Range_orbA.second;
112 if (degeneracy_ != 0) {
113 std::vector<Index> list_levelsA =
114 orbitalsA.CheckDegeneracy(levelA, degeneracy_);
115 std::vector<Index> list_levelsB =
116 orbitalsB.CheckDegeneracy(levelB, degeneracy_);
117
118 double JAB_sq = 0;
119
120 for (Index iA : list_levelsA) {
121 Index indexA = iA - Range_orbA.first;
122 for (Index iB : list_levelsB) {
123 Index indexB = iB - Range_orbB.first + levelsA;
124 double JAB_one_level = JAB(indexA, indexB);
125 JAB_sq += JAB_one_level * JAB_one_level;
126 }
127 }
128 return std::sqrt(JAB_sq /
129 double(list_levelsA.size() * list_levelsB.size())) *
131 } else {
132 Index indexA = levelA - Range_orbA.first;
133 Index indexB = levelB - Range_orbB.first + levelsA;
134 return JAB(indexA, indexB) * tools::conv::hrt2ev;
135 }
136}
137
145 const Orbitals& orbitalsB,
146 const Orbitals& orbitalsAB) {
147
148 XTP_LOG(Log::error, *pLog_) << "Calculating electronic couplings" << flush;
149
150 CheckAtomCoordinates(orbitalsA, orbitalsB, orbitalsAB);
151
152 // constructing the direct product orbA x orbB
153 Index basisA = orbitalsA.getBasisSetSize();
154 Index basisB = orbitalsB.getBasisSetSize();
155
156 if ((basisA == 0) || (basisB == 0)) {
157 throw std::runtime_error("Basis set size is not stored in monomers");
158 }
159
162
163 Index levelsA = Range_orbA.second;
164 Index levelsB = Range_orbB.second;
165
167 << "Levels:Basis A[" << levelsA << ":" << basisA << "]"
168 << " B[" << levelsB << ":" << basisB << "]" << flush;
169
170 if ((levelsA == 0) || (levelsB == 0)) {
171 throw std::runtime_error(
172 "No information about number of occupied/unoccupied levels is stored");
173 }
174
175 // constructing merged orbitals
176 auto MOsA = orbitalsA.MOs().eigenvectors().middleCols(Range_orbA.first,
177 Range_orbA.second);
178 auto MOsB = orbitalsB.MOs().eigenvectors().middleCols(Range_orbB.first,
179 Range_orbB.second);
180
181 XTP_LOG(Log::info, *pLog_) << "Calculating overlap matrix for basisset: "
182 << orbitalsAB.getDFTbasisName() << flush;
183
184 Eigen::MatrixXd overlap =
185 CalculateOverlapMatrix(orbitalsAB) * orbitalsAB.MOs().eigenvectors();
186
188 << "Projecting monomers onto dimer orbitals" << flush;
189 Eigen::MatrixXd A_AB = MOsA.transpose() * overlap.topRows(basisA);
190 Eigen::MatrixXd B_AB = MOsB.transpose() * overlap.bottomRows(basisB);
191 Eigen::VectorXd mag_A = A_AB.rowwise().squaredNorm();
192 if (mag_A.any() < 0.95) {
194 << "\nWarning: "
195 << "Projection of orbitals of monomer A on dimer is insufficient,mag="
196 << mag_A.minCoeff() << flush;
197 }
198 Eigen::VectorXd mag_B = B_AB.rowwise().squaredNorm();
199 if (mag_B.any() < 0.95) {
201 << "\nWarning: "
202 << "Projection of orbitals of monomer B on dimer is insufficient,mag="
203 << mag_B.minCoeff() << flush;
204 }
205
206 Eigen::MatrixXd psi_AxB_dimer_basis(A_AB.rows() + B_AB.rows(), A_AB.cols());
207 psi_AxB_dimer_basis.topRows(A_AB.rows()) = A_AB;
208 psi_AxB_dimer_basis.bottomRows(B_AB.rows()) = B_AB;
209
211 << "Projecting the Fock matrix onto the dimer basis" << flush;
212 Eigen::MatrixXd JAB_dimer = psi_AxB_dimer_basis *
213 orbitalsAB.MOs().eigenvalues().asDiagonal() *
214 psi_AxB_dimer_basis.transpose();
215 XTP_LOG(Log::info, *pLog_) << "Constructing Overlap matrix" << flush;
216 Eigen::MatrixXd S_AxB = psi_AxB_dimer_basis * psi_AxB_dimer_basis.transpose();
217 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(S_AxB);
218 Eigen::MatrixXd Sm1 = es.operatorInverseSqrt();
219 XTP_LOG(Log::info, *pLog_) << "Smallest eigenvalue of overlap matrix is "
220 << es.eigenvalues()(0) << flush;
221 JAB = Sm1 * JAB_dimer * Sm1;
222 XTP_LOG(Log::error, *pLog_) << "Done with electronic couplings" << flush;
223}
224
225} // namespace xtp
226} // 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 & add(const std::string &key, const std::string &value)
add a new property to structure
Definition property.cc:108
Property & get(const std::string &key)
get existing property
Definition property.cc:79
T as() const
return value as type
Definition property.h:283
void setAttribute(const std::string &attribute, const T &value)
set an attribute
Definition property.h:314
Eigen::MatrixXd CalculateOverlapMatrix(const Orbitals &orbitalsAB) const
void CheckAtomCoordinates(const Orbitals &orbitalsA, const Orbitals &orbitalsB, const Orbitals &orbitalsAB) const
void WriteToProperty(tools::Property &type_summary, const Orbitals &orbitalsA, const Orbitals &orbitalsB, Index a, Index b) const
std::pair< Index, Index > Range_orbA
Definition dftcoupling.h:66
std::pair< Index, Index > DetermineRangeOfStates(const Orbitals &orbital, Index numberofstates) const
void Initialize(tools::Property &) override
double getCouplingElement(Index levelA, Index levelB, const Orbitals &orbitalsA, const Orbitals &orbitalsB) const
void Addoutput(tools::Property &type_summary, const Orbitals &orbitalsA, const Orbitals &orbitalsB) const override
std::pair< Index, Index > Range_orbB
Definition dftcoupling.h:67
void CalculateCouplings(const Orbitals &orbitalsA, const Orbitals &orbitalsB, const Orbitals &orbitalsAB) override
evaluates electronic couplings
std::string Identify() const
Definition dftcoupling.h:40
Eigen::MatrixXd JAB
Definition dftcoupling.h:60
container for molecular orbitals
Definition orbitals.h:46
Index getHomo() const
Definition orbitals.h:68
std::vector< Index > CheckDegeneracy(Index level, double energy_difference) const
Definition orbitals.cc:47
Index getBasisSetSize() const
Definition orbitals.h:64
const tools::EigenSystem & MOs() const
Definition orbitals.h:122
const std::string & getDFTbasisName() const
Definition orbitals.h:202
Index getLumo() const
Definition orbitals.h:66
#define XTP_LOG(level, log)
Definition logger.h:40
const double ev2hrt
Definition constants.h:54
const double hrt2ev
Definition constants.h:53
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26