votca 2026-dev
Loading...
Searching...
No Matches
bsecoupling.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"
28#include "votca/xtp/bse.h"
31
32namespace votca {
33namespace xtp {
34using namespace std;
35using boost::format;
36using namespace tools;
37
39 string spintype = options.get("spin").as<std::string>();
40 if (spintype == "all") {
41 doSinglets_ = true;
42 doTriplets_ = true;
43 } else if (spintype == "triplet") {
44 doTriplets_ = true;
45 doSinglets_ = false;
46 } else if (spintype == "singlet") {
47 doSinglets_ = true;
48 doTriplets_ = false;
49 } else {
50 throw std::runtime_error(
51 (boost::format(
52 "Choice % for type not known. Available singlet,triplet,all") %
53 spintype)
54 .str());
55 }
56 output_perturbation_ = options.get("use_perturbation").as<bool>();
57 // When true, write monomer energies, transition dipoles, diagnostics,
58 // and raw H/S TB matrices to the XML output. Default false for compact
59 // output in KMC/rate workflows that only need scalar couplings.
61 options.ifExistsReturnElseReturnDefault<bool>("output_tb", false);
62
63 levA_ = options.get("moleculeA.states").as<Index>();
64 levB_ = options.get("moleculeB.states").as<Index>();
65 occA_ = options.get("moleculeA.occLevels").as<Index>();
66 occB_ = options.get("moleculeB.occLevels").as<Index>();
67 unoccA_ = options.get("moleculeA.unoccLevels").as<Index>();
68 unoccB_ = options.get("moleculeB.unoccLevels").as<Index>();
69}
70
71// =============================================================================
72// WriteMatrixToProperty
73// Write an Eigen matrix as space-separated rows into an XML property node.
74// Each row is stored as attribute "row_N". Unit conversion is applied if given.
75// =============================================================================
77 const std::string& name,
78 const Eigen::MatrixXd& mat,
79 double conversion) {
80 tools::Property& mat_prop = prop.add(name, "");
81 mat_prop.setAttribute("rows", mat.rows());
82 mat_prop.setAttribute("cols", mat.cols());
83 for (Index i = 0; i < mat.rows(); i++) {
84 std::string row_str;
85 for (Index j = 0; j < mat.cols(); j++) {
86 if (j > 0) row_str += " ";
87 row_str += (format("%1$1.6e") % (mat(i, j) * conversion)).str();
88 }
89 mat_prop.setAttribute("row_" + std::to_string(i), row_str);
90 }
91}
92
93// =============================================================================
94// WriteToProperty
95// Write the pairwise effective couplings for one (stateA, stateB) pair.
96// Backward compatible: attributes unchanged.
97// =============================================================================
98void BSECoupling::WriteToProperty(Property& summary, const QMState& stateA,
99 const QMState& stateB) const {
100 Property& coupling_summary = summary.add("coupling", "");
101 double JAB_pert = 0;
102 double JAB_diag = 0;
103 if (stateA.Type() == QMStateType::Singlet) {
104 JAB_pert =
105 getSingletCouplingElement(stateA.StateIdx(), stateB.StateIdx(), 0);
106 JAB_diag =
107 getSingletCouplingElement(stateA.StateIdx(), stateB.StateIdx(), 1);
108 } else if (stateA.Type() == QMStateType::Triplet) {
109 JAB_pert =
110 getTripletCouplingElement(stateA.StateIdx(), stateB.StateIdx(), 0);
111 JAB_diag =
112 getTripletCouplingElement(stateA.StateIdx(), stateB.StateIdx(), 1);
113 }
114 coupling_summary.setAttribute("stateA", stateA.ToString());
115 coupling_summary.setAttribute("stateB", stateB.ToString());
116 coupling_summary.setAttribute("j_pert", (format("%1$1.6e") % JAB_pert).str());
117 coupling_summary.setAttribute("j_diag", (format("%1$1.6e") % JAB_diag).str());
118}
119
120// =============================================================================
121// Addoutput
122// Write all results to the XML property tree.
123// Extends the existing pairwise coupling output with:
124// - diagnostics node (xi, PT/RM discrepancy, downfolding_safe flag)
125// - tb_matrices node (raw H and S blocks in FE+CT basis, in eV)
126// The existing <coupling> elements are unchanged for backward compatibility.
127// =============================================================================
128void BSECoupling::Addoutput(Property& type_summary, const Orbitals& orbitalsA,
129 const Orbitals& orbitalsB) const {
130 tools::Property& bsecoupling = type_summary.add(Identify(), "");
131 string algorithm = "j_diag";
133 algorithm = "j_pert";
134 }
135
136 // Lambda to write one spin channel — avoids duplicating singlet/triplet logic
137 auto WriteSpinChannel = [&](const std::string& spin_name,
138 const Eigen::MatrixXd& J_dimer,
139 const Eigen::MatrixXd& S_dimer,
140 const Diagnostics& diag,
141 const Eigen::VectorXd& monA_energies,
142 const Eigen::VectorXd& monB_energies,
143 const std::vector<Eigen::Vector3d>& tdipsA,
144 const std::vector<Eigen::Vector3d>& tdipsB,
145 QMStateType spin_type) {
146 Property& spin_summary = bsecoupling.add(spin_name, "");
147 spin_summary.setAttribute("algorithm", algorithm);
148
149 // --- existing pairwise effective couplings (backward compatible) ---
150 for (Index stateA = 0; stateA < levA_; ++stateA) {
151 QMState qmstateA = QMState(spin_type, stateA, false);
152 for (Index stateB = 0; stateB < levB_; ++stateB) {
153 QMState qmstateB = QMState(spin_type, stateB, false);
154 WriteToProperty(spin_summary, qmstateA, qmstateB);
155 }
156 }
157
158 // --- TB output: monomer energies, diagnostics, raw matrices ---
159 // Only written when output_tb_ = true. For KMC/rate workflows that only
160 // need scalar effective couplings, leave output_tb_ = false.
161 if (output_tb_) {
162 // Monomer excitation energies and transition dipoles
163 {
164 Property& mono_prop = spin_summary.add("monomer_energies", "");
165 Property& propA = mono_prop.add("fragmentA", "");
166 for (Index i = 0; i < monA_energies.size(); i++) {
167 Property& e = propA.add("energy", "");
168 QMState st(spin_type, i, false);
169 e.setAttribute("state", st.ToString());
170 e.setAttribute("eV", (format("%1$1.6e") % monA_energies(i)).str());
171 if (i < static_cast<Index>(tdipsA.size())) {
172 e.setAttribute("tdx", (format("%1$1.6e") % tdipsA[i].x()).str());
173 e.setAttribute("tdy", (format("%1$1.6e") % tdipsA[i].y()).str());
174 e.setAttribute("tdz", (format("%1$1.6e") % tdipsA[i].z()).str());
175 }
176 }
177 Property& propB = mono_prop.add("fragmentB", "");
178 for (Index i = 0; i < monB_energies.size(); i++) {
179 Property& e = propB.add("energy", "");
180 QMState st(spin_type, i, false);
181 e.setAttribute("state", st.ToString());
182 e.setAttribute("eV", (format("%1$1.6e") % monB_energies(i)).str());
183 if (i < static_cast<Index>(tdipsB.size())) {
184 e.setAttribute("tdx", (format("%1$1.6e") % tdipsB[i].x()).str());
185 e.setAttribute("tdy", (format("%1$1.6e") % tdipsB[i].y()).str());
186 e.setAttribute("tdz", (format("%1$1.6e") % tdipsB[i].z()).str());
187 }
188 }
189 }
190
191 // Diagnostics
192 {
193 Property& diag_prop = spin_summary.add("diagnostics", "");
194 diag_prop.setAttribute("xi", (format("%1$1.4f") % diag.xi).str());
195 diag_prop.setAttribute(
196 "pt_rm_discrepancy_eV",
197 (format("%1$1.6e") %
198 (diag.pt_rm_discrepancy * votca::tools::conv::hrt2ev))
199 .str());
200 diag_prop.setAttribute("downfolding_safe",
201 diag.downfolding_safe ? "true" : "false");
202 }
203
204 // Raw TB matrices (H in eV, S dimensionless)
205 // Block layout: [ H_FE_FE H_FE_CT ] / [ H_CT_FE H_CT_CT ]
206 {
207 Index bse_exc = levA_ + levB_;
208 Index ct = J_dimer.rows() - bse_exc;
209
210 Property& tb_prop = spin_summary.add("tb_matrices", "");
211 tb_prop.setAttribute("n_FE", bse_exc);
212 tb_prop.setAttribute("n_CT", ct);
213 tb_prop.setAttribute("n_occA", occA_);
214 tb_prop.setAttribute("n_unoccA", unoccA_);
215 tb_prop.setAttribute("n_occB", occB_);
216 tb_prop.setAttribute("n_unoccB", unoccB_);
217
218 WriteMatrixToProperty(tb_prop, "H_FE_FE",
219 J_dimer.topLeftCorner(bse_exc, bse_exc),
221 WriteMatrixToProperty(tb_prop, "S_FE_FE",
222 S_dimer.topLeftCorner(bse_exc, bse_exc));
223
224 if (ct > 0) {
225 WriteMatrixToProperty(tb_prop, "H_FE_CT",
226 J_dimer.topRightCorner(bse_exc, ct),
228 WriteMatrixToProperty(tb_prop, "S_FE_CT",
229 S_dimer.topRightCorner(bse_exc, ct));
230 WriteMatrixToProperty(tb_prop, "H_CT_CT",
231 J_dimer.bottomRightCorner(ct, ct),
233 WriteMatrixToProperty(tb_prop, "S_CT_CT",
234 S_dimer.bottomRightCorner(ct, ct));
235 }
236 }
237 } // output_tb_
238 };
239
240 if (doSinglets_) {
241 WriteSpinChannel(
246 }
247 if (doTriplets_) {
248 // Triplet states have zero electric transition dipole from the singlet
249 // ground state — pass empty vectors so no tdx/tdy/tdz attributes are
250 // written, making clear that oscillator strengths are not meaningful.
251 const std::vector<Eigen::Vector3d> empty_dipoles;
252 WriteSpinChannel(QMStateType(QMStateType::Triplet).ToLongString(),
255 empty_dipoles, empty_dipoles,
257 }
258}
259
261 Index methodindex) const {
262 return JAB_singlet[methodindex](levelA, levelB + levA_) *
264}
265
267 Index methodindex) const {
268 return JAB_triplet[methodindex](levelA, levelB + levA_) *
270}
271
272Eigen::MatrixXd BSECoupling::SetupCTStates(Index bseA_vtotal, Index bseB_vtotal,
273 Index bseAB_vtotal,
274 Index bseAB_ctotal,
275 const Eigen::MatrixXd& A_AB,
276 const Eigen::MatrixXd& B_AB) const {
277 Index noAB = occA_ * unoccB_;
278 Index noBA = unoccA_ * occB_;
279 Index bseAB_size = bseAB_vtotal * bseAB_ctotal;
280 Eigen::MatrixXd CTstates = Eigen::MatrixXd::Zero(bseAB_size, noAB + noBA);
281
282 auto A_occ = A_AB.middleCols(bseA_vtotal - occA_, occA_);
283 auto A_unocc = A_AB.middleCols(bseA_vtotal, unoccA_);
284 auto B_occ = B_AB.middleCols(bseB_vtotal - occB_, occB_);
285 auto B_unocc = B_AB.middleCols(bseB_vtotal, unoccB_);
286
287 const Eigen::MatrixXd A_occ_occ = A_occ.topRows(bseAB_vtotal);
288 const Eigen::MatrixXd B_unocc_unocc = B_unocc.bottomRows(bseAB_ctotal);
289
290 // notation AB is CT states with A+B-, BA is the counterpart
292 << TimeStamp() << " Setting up CT-states" << flush;
293
294#pragma omp parallel for
295 for (Index a_occ = 0; a_occ < occA_; a_occ++) {
296 for (Index b_unocc = 0; b_unocc < unoccB_; b_unocc++) {
297 Index index = a_occ * unoccB_ + b_unocc;
298 Eigen::MatrixXd Coeff =
299 B_unocc_unocc.col(b_unocc) * A_occ_occ.col(a_occ).transpose();
300 CTstates.col(index) =
301 Eigen::Map<Eigen::VectorXd>(Coeff.data(), bseAB_size);
302 }
303 }
305 << TimeStamp() << " " << noBA << " CT states A+B- created" << flush;
306
307 const Eigen::MatrixXd A_unocc_unocc = A_unocc.bottomRows(bseAB_ctotal);
308 const Eigen::MatrixXd B_occ_occ = B_occ.topRows(bseAB_vtotal);
309
310#pragma omp parallel for
311 for (Index b_occ = 0; b_occ < occB_; b_occ++) {
312 for (Index a_unocc = 0; a_unocc < unoccA_; a_unocc++) {
313 Index index = b_occ * unoccA_ + a_unocc + noAB;
314 Eigen::MatrixXd Coeff =
315 A_unocc_unocc.col(a_unocc) * B_occ_occ.col(b_occ).transpose();
316 CTstates.col(index) =
317 Eigen::Map<Eigen::VectorXd>(Coeff.data(), bseAB_size);
318 }
319 }
321 << TimeStamp() << " " << noBA << " CT states A-B+ created" << flush;
322 return CTstates;
323}
324
326 const Eigen::MatrixXd& BSE_Coeffs, const Eigen::MatrixXd& X_AB,
327 Index bseX_vtotal, Index bseX_ctotal, Index bseAB_vtotal,
328 Index bseAB_ctotal) const {
329 Index bseAB_size = bseAB_vtotal * bseAB_ctotal;
330 auto X_occ = X_AB.leftCols(bseX_vtotal);
331 auto X_unocc = X_AB.rightCols(bseX_ctotal);
332 const Eigen::MatrixXd X_occ_occ = X_occ.topRows(bseAB_vtotal);
333 const Eigen::MatrixXd X_unocc_unocc = X_unocc.bottomRows(bseAB_ctotal);
334 Eigen::MatrixXd result = Eigen::MatrixXd::Zero(bseAB_size, BSE_Coeffs.cols());
335 // no pragma here because often we will only have one Coeff
336 for (Index i = 0; i < BSE_Coeffs.cols(); i++) {
337 Eigen::VectorXd coeff = BSE_Coeffs.col(i);
338 Eigen::Map<Eigen::MatrixXd> coeffmatrix =
339 Eigen::Map<Eigen::MatrixXd>(coeff.data(), bseX_ctotal, bseX_vtotal);
340 Eigen::MatrixXd proj = X_unocc_unocc * coeffmatrix * X_occ_occ.transpose();
341 result.col(i) = Eigen::Map<Eigen::VectorXd>(proj.data(), proj.size());
342 }
343 return result;
344}
345
346int GetSign(double value) {
347 if (value < 0) {
348 return -1;
349 } else if (value > 0) {
350 return 1;
351 }
352 return 0;
353}
354
356 const Orbitals& orbitalsB,
357 const Orbitals& orbitalsAB) {
359 << TimeStamp() << " Calculating exciton couplings" << flush;
360 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Using "
361 << OPENMP::getMaxThreads() << " threads" << flush;
362
363 CheckAtomCoordinates(orbitalsA, orbitalsB, orbitalsAB);
364
365 Index basisB = orbitalsB.getBasisSetSize();
366 Index basisA = orbitalsA.getBasisSetSize();
367 if ((basisA == 0) || (basisB == 0)) {
368 throw std::runtime_error("Basis set size is not stored in monomers");
369 }
370
371 // get exciton information of molecule A
372 Index bseA_cmax = orbitalsA.getBSEcmax();
373 Index bseA_cmin = orbitalsA.getBSEcmin();
374 Index bseA_vmax = orbitalsA.getBSEvmax();
375 Index bseA_vmin = orbitalsA.getBSEvmin();
376 Index bseA_vtotal = bseA_vmax - bseA_vmin + 1;
377 Index bseA_ctotal = bseA_cmax - bseA_cmin + 1;
378 Index bseA_total = bseA_vtotal + bseA_ctotal;
379 Index bseA_size = bseA_vtotal * bseA_ctotal;
380 Index bseA_singlet_exc = orbitalsA.BSESinglets().eigenvectors().cols();
381 Index bseA_triplet_exc = orbitalsA.BSETriplets().eigenvectors().cols();
382
384 << TimeStamp() << " molecule A has " << bseA_singlet_exc
385 << " singlet excitons with dimension " << bseA_size << flush;
387 << TimeStamp() << " molecule A has " << bseA_triplet_exc
388 << " triplet excitons with dimension " << bseA_size << flush;
389
390 // get exciton information of molecule B
391 Index bseB_cmax = orbitalsB.getBSEcmax();
392 Index bseB_cmin = orbitalsB.getBSEcmin();
393 Index bseB_vmax = orbitalsB.getBSEvmax();
394 Index bseB_vmin = orbitalsB.getBSEvmin();
395 Index bseB_vtotal = bseB_vmax - bseB_vmin + 1;
396 Index bseB_ctotal = bseB_cmax - bseB_cmin + 1;
397 Index bseB_total = bseB_vtotal + bseB_ctotal;
398 Index bseB_size = bseB_vtotal * bseB_ctotal;
399 Index bseB_singlet_exc = orbitalsB.BSESinglets().eigenvectors().cols();
400 Index bseB_triplet_exc = orbitalsB.BSETriplets().eigenvectors().cols();
401
403 << TimeStamp() << " molecule B has " << bseB_singlet_exc
404 << " singlet excitons with dimension " << bseB_size << flush;
406 << TimeStamp() << " molecule B has " << bseB_triplet_exc
407 << " triplet excitons with dimension " << bseB_size << flush;
408
409 if (doSinglets_) {
410 if (levA_ > bseA_singlet_exc) {
412 << " Number of excitons you want is greater "
413 "than stored for molecule "
414 "A. Setting to max number available"
415 << flush;
416 levA_ = bseA_singlet_exc;
417 }
418 if (levB_ > bseB_singlet_exc) {
420 << " Number of excitons you want is greater "
421 "than stored for molecule "
422 "B. Setting to max number available"
423 << flush;
424 levB_ = bseB_singlet_exc;
425 }
426 }
427 if (doTriplets_) {
428 if (levA_ > bseA_triplet_exc) {
430 << TimeStamp()
431 << " Number of Frenkel states you want is greater than stored for "
432 "molecule A. Setting to max number available"
433 << flush;
434 levA_ = bseA_triplet_exc;
435 }
436 if (levB_ > bseB_triplet_exc) {
438 << TimeStamp()
439 << " Number of Frenkel states you want is greater than stored for "
440 "molecule B. Setting to max number available"
441 << flush;
442 levB_ = bseB_triplet_exc;
443 }
444 }
445 if (unoccA_ > bseA_ctotal || unoccA_ < 0) {
447 << TimeStamp()
448 << " Number of unoccupied orbitals in molecule A for CT creation "
449 "exceeds number of KS-orbitals in BSE"
450 << flush;
451 unoccA_ = bseA_ctotal;
452 }
453 if (unoccB_ > bseB_ctotal || unoccB_ < 0) {
455 << TimeStamp()
456 << " Number of unoccupied orbitals in molecule B for CT creation "
457 "exceeds number of KS-orbitals in BSE"
458 << flush;
459 unoccB_ = bseB_ctotal;
460 }
461 if (occA_ > bseA_vtotal || occA_ < 0) {
463 << TimeStamp()
464 << " Number of occupied orbitals in molecule A for CT creation "
465 "exceeds number of KS-orbitals in BSE"
466 << flush;
467 occA_ = bseA_vtotal;
468 }
469 if (occB_ > bseB_vtotal || occB_ < 0) {
471 << TimeStamp()
472 << " Number of occupied orbitals in molecule B for CT creation "
473 "exceeds number of KS-orbitals in BSE"
474 << flush;
475 occB_ = bseB_vtotal;
476 }
477
478 // get exciton information of pair AB
479 Index bseAB_cmax = orbitalsAB.getBSEcmax();
480 Index bseAB_cmin = orbitalsAB.getBSEcmin();
481 Index bseAB_vmax = orbitalsAB.getBSEvmax();
482 Index bseAB_vmin = orbitalsAB.getBSEvmin();
483 Index bseAB_vtotal = bseAB_vmax - bseAB_vmin + 1;
484 Index bseAB_ctotal = bseAB_cmax - bseAB_cmin + 1;
485 Index bseAB_total = bseAB_vtotal + bseAB_ctotal;
486 Index bseAB_size = bseAB_vtotal * bseAB_ctotal;
487
489 << TimeStamp() << " levels used for BSE of molA: " << bseA_vmin
490 << " to " << bseA_cmax << " total: " << bseA_total << flush;
492 << TimeStamp() << " levels used for BSE of molB: " << bseB_vmin
493 << " to " << bseB_cmax << " total: " << bseB_total << flush;
495 << TimeStamp() << " levels used for BSE of dimer AB: " << bseAB_vmin
496 << " to " << bseAB_cmax << " total: " << bseAB_total << flush;
497
498 Eigen::MatrixXd MOsA =
499 orbitalsA.MOs().eigenvectors().middleCols(bseA_vmin, bseA_total);
500 Eigen::MatrixXd MOsB =
501 orbitalsB.MOs().eigenvectors().middleCols(bseB_vmin, bseB_total);
502 Eigen::MatrixXd MOsAB =
503 orbitalsAB.MOs().eigenvectors().middleCols(bseAB_vmin, bseAB_total);
504
506 << TimeStamp() << " Calculating overlap matrix for basisset: "
507 << orbitalsAB.getDFTbasisName() << flush;
508
509 Eigen::MatrixXd overlap = CalculateOverlapMatrix(orbitalsAB) * MOsAB;
511 << TimeStamp() << " Projecting monomers onto dimer orbitals" << flush;
512
513 Eigen::MatrixXd A_AB = overlap.topRows(basisA).transpose() * MOsA;
514 Eigen::MatrixXd B_AB = overlap.bottomRows(basisB).transpose() * MOsB;
515 Eigen::VectorXd mag_A = A_AB.colwise().squaredNorm();
516 if (mag_A.any() < 0.95) {
518 << "\nWarning: "
519 << "Projection of orbitals of monomer A on dimer is insufficient,mag="
520 << mag_A.minCoeff() << flush;
521 }
522 Eigen::VectorXd mag_B = B_AB.colwise().squaredNorm();
523 if (mag_B.any() < 0.95) {
525 << "\nWarning: "
526 << "Projection of orbitals of monomer B on dimer is insufficient,mag="
527 << mag_B.minCoeff() << flush;
528 }
529
530 AOBasis dftbasis = orbitalsAB.getDftBasis();
531 AOBasis auxbasis = orbitalsAB.getAuxBasis();
532 TCMatrix_gwbse Mmn;
533 Mmn.Initialize(auxbasis.AOBasisSize(), orbitalsAB.getRPAmin(),
534 orbitalsAB.getGWAmax(), orbitalsAB.getRPAmin(),
535 orbitalsAB.getRPAmax());
536 Mmn.Fill(auxbasis, dftbasis, orbitalsAB.MOs().eigenvectors());
537
538 const Eigen::MatrixXd& qpcoeff = orbitalsAB.QPdiag().eigenvectors();
539 Eigen::MatrixXd Hqp = qpcoeff *
540 orbitalsAB.QPdiag().eigenvalues().asDiagonal() *
541 qpcoeff.transpose();
542
543 BSE::options opt;
544 opt.cmax = orbitalsAB.getBSEcmax();
545 opt.homo = orbitalsAB.getHomo();
546 opt.qpmin = orbitalsAB.getGWAmin();
547 opt.qpmax = orbitalsAB.getGWAmax();
548 opt.rpamax = orbitalsAB.getRPAmax();
549 opt.rpamin = orbitalsAB.getRPAmin();
550 opt.useTDA = true;
551 opt.vmin = orbitalsAB.getBSEvmin();
552 opt.use_Hqp_offdiag = orbitalsAB.GetFlagUseHqpOffdiag();
553 BSE bse(*pLog_, Mmn);
554 bse.configure(opt, orbitalsAB.RPAInputEnergies(), Hqp);
555 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Setup BSE operator" << flush;
556
557 if (doSinglets_) {
559 << TimeStamp() << " Evaluating singlets" << flush;
561 << TimeStamp() << " Setup Hamiltonian" << flush;
562
563 // Store isolated monomer singlet excitation energies [eV].
564 // These are pairwise-consistent site energies for TB assembly —
565 // independent of which dimer partner is used.
567 orbitalsA.BSESinglets().eigenvalues().head(levA_) *
570 orbitalsB.BSESinglets().eigenvalues().head(levB_) *
572
573 Eigen::MatrixXd FE_AB = Eigen::MatrixXd::Zero(bseAB_size, levA_ + levB_);
574 const Eigen::MatrixXd bseA =
575 orbitalsA.BSESinglets().eigenvectors().leftCols(levA_);
576 FE_AB.leftCols(levA_) = ProjectFrenkelExcitons(
577 bseA, A_AB, bseA_vtotal, bseA_ctotal, bseAB_vtotal, bseAB_ctotal);
578 const Eigen::MatrixXd bseB =
579 orbitalsB.BSESinglets().eigenvectors().leftCols(levB_);
580 FE_AB.rightCols(levB_) = ProjectFrenkelExcitons(
581 bseB, B_AB, bseB_vtotal, bseB_ctotal, bseAB_vtotal, bseAB_ctotal);
582 Eigen::MatrixXd CTStates = SetupCTStates(
583 bseA_vtotal, bseB_vtotal, bseAB_vtotal, bseAB_ctotal, A_AB, B_AB);
585 ProjectExcitons(FE_AB, CTStates, bse.getSingletOperator_TDA(),
588 << TimeStamp() << " calculated singlet couplings " << flush;
589 }
590 if (doTriplets_) {
592 << TimeStamp() << " Evaluating triplets" << flush;
593
594 // Store isolated monomer triplet excitation energies [eV].
596 orbitalsA.BSETriplets().eigenvalues().head(levA_) *
599 orbitalsB.BSETriplets().eigenvalues().head(levB_) *
601
602 Eigen::MatrixXd FE_AB = Eigen::MatrixXd::Zero(bseAB_size, levA_ + levB_);
603 const Eigen::MatrixXd bseA =
604 orbitalsA.BSETriplets().eigenvectors().leftCols(levA_);
605 FE_AB.leftCols(levA_) = ProjectFrenkelExcitons(
606 bseA, A_AB, bseA_vtotal, bseA_ctotal, bseAB_vtotal, bseAB_ctotal);
607 const Eigen::MatrixXd bseB =
608 orbitalsB.BSETriplets().eigenvectors().leftCols(levB_);
609 FE_AB.rightCols(levB_) = ProjectFrenkelExcitons(
610 bseB, B_AB, bseB_vtotal, bseB_ctotal, bseAB_vtotal, bseAB_ctotal);
611 Eigen::MatrixXd CTStates = SetupCTStates(
612 bseA_vtotal, bseB_vtotal, bseAB_vtotal, bseAB_ctotal, A_AB, B_AB);
614 ProjectExcitons(FE_AB, CTStates, bse.getTripletOperator_TDA(),
617 << TimeStamp() << " calculated triplet couplings " << flush;
618 }
620 << TimeStamp() << " Done with exciton couplings" << flush;
621}
622
623// =============================================================================
624// OrthogonalizeCTs
625// Merge FE and CT projection vectors into a single matrix without any
626// pre-orthogonalization. The non-orthogonality between FE and CT states,
627// and within each block, is handled correctly downstream:
628// - For the reduction method: by the joint Lowdin in CalcJ_dimer.
629// - For TB use of J_dimer/S_dimer: by the generalized eigenvalue problem
630// or SRG.
631//
632// Pre-orthogonalizing CT states against FE states (as done previously) used
633// the incorrect projector P = F*F^T, which is only exact when FE_AB columns
634// are orthonormal. Since projected monomer FE states are generally not
635// orthonormal, this introduced a systematic error — particularly for larger
636// CT manifolds (many occ/unocc pairs) where higher orbital pairs have
637// non-negligible overlap with the FE space. Removing the pre-orthogonalization
638// preserves the physical FE and CT state characters and lets the subsequent
639// linear algebra handle the non-orthogonality correctly and consistently.
640//
641// A linear dependence check is performed on the CT states: if any CT state
642// has norm below threshold after projection onto the dimer basis, it indicates
643// near-linear dependence with other basis states and is flagged.
644// =============================================================================
645Eigen::MatrixXd BSECoupling::OrthogonalizeCTs(Eigen::MatrixXd& FE_AB,
646 Eigen::MatrixXd& CTStates) const {
647 Index ct = CTStates.cols();
648 Index bse_exc = levA_ + levB_;
649 Index bseAB_size = FE_AB.rows();
650
651 // Check for near-linear dependence in CT states via their norms.
652 // A very small norm indicates a CT state that is nearly identical to
653 // another basis state and would cause numerical problems in the Lowdin step.
654 if (ct > 0) {
655 Eigen::VectorXd norm = CTStates.colwise().norm();
656 Index minstateindex = 0;
657 double minnorm = norm.minCoeff(&minstateindex);
658 if (minnorm < 0.95) {
660 << TimeStamp() << " WARNING: CT-state " << minstateindex
661 << " norm is only " << minnorm
662 << " -- near-linear dependence in projection basis" << flush;
663 }
665 << TimeStamp() << " Merging " << ct
666 << " CT states into projection (no pre-orthogonalization)" << flush;
667 }
668
669 Eigen::MatrixXd projection(bseAB_size, bse_exc + ct);
671 << TimeStamp() << " merging projections into one vector " << flush;
672 projection.leftCols(bse_exc) = FE_AB;
673 FE_AB.resize(0, 0);
674 if (ct > 0) {
675 projection.rightCols(ct) = CTStates;
676 }
677 CTStates.resize(0, 0);
678 return projection;
679}
680
681// =============================================================================
682// CalcJ_dimer
683// Form J_dimer and S_dimer from the projection, expose them via output
684// references for TB use, then Lowdin-orthogonalize to produce J_ortho.
685// The block structure of J_dimer and S_dimer is:
686// [ H_FE_FE H_FE_CT ] [ S_FE_FE S_FE_CT ]
687// [ H_CT_FE H_CT_CT ] [ S_CT_FE S_CT_CT ]
688// =============================================================================
689template <class BSE_OPERATOR>
691 Eigen::MatrixXd& projection,
692 Eigen::MatrixXd& J_dimer_out,
693 Eigen::MatrixXd& S_dimer_out) const {
695 << TimeStamp() << " Setting up coupling matrix size "
696 << projection.cols() << flush;
697
698 // this only works for hermitian/symmetric H so only in TDA
699 J_dimer_out = projection.transpose() * (H * projection).eval();
700
702 << TimeStamp() << " Setting up overlap matrix size "
703 << projection.cols() << flush;
704 S_dimer_out = projection.transpose() * projection;
705
706 projection.resize(0, 0);
707
709 << "---------------------------------------" << flush;
710 XTP_LOG(Log::debug, *pLog_) << " J_dimer_[Hrt]" << flush;
711 XTP_LOG(Log::debug, *pLog_) << J_dimer_out << flush;
712 XTP_LOG(Log::debug, *pLog_) << " S_dimer_" << flush;
713 XTP_LOG(Log::debug, *pLog_) << S_dimer_out << flush;
715 << "---------------------------------------" << flush;
716
717 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(S_dimer_out);
719 << TimeStamp() << " Smallest value of dimer overlapmatrix is "
720 << es.eigenvalues()(0) << flush;
721
722 Eigen::MatrixXd Sm1 = es.operatorInverseSqrt();
723 Eigen::MatrixXd J_ortho = Sm1 * J_dimer_out * Sm1;
724
726 << "---------------------------------------" << flush;
727 XTP_LOG(Log::debug, *pLog_) << " J_ortho_[Hrt]" << flush;
728 XTP_LOG(Log::debug, *pLog_) << J_ortho << flush;
729 XTP_LOG(Log::debug, *pLog_) << " S_-1/2" << flush;
730 XTP_LOG(Log::debug, *pLog_) << Sm1 << flush;
732 << "---------------------------------------" << flush;
733
734 return J_ortho;
735}
736
737// =============================================================================
738// ComputeDiagnostics
739// Assess whether CT downfolding is reliable for this pair.
740// xi: max |H_FE_CT| / |E_FE - E_CT| over all FE-CT pairs.
741// Small xi => perturbation theory is safe.
742// pt_rm_discrepancy: max |J_pert - J_diag| over all FE pairs [Hrt].
743// Large discrepancy flags near-resonant CT states.
744// downfolding_safe: true if xi < 0.3 and discrepancy < 1e-4 Hrt.
745// Note: J_dimer here is the raw pre-Lowdin matrix (in Hrt), not J_ortho.
746// =============================================================================
748 const Eigen::MatrixXd& J_dimer, const Eigen::MatrixXd& J_pert,
749 const Eigen::MatrixXd& J_diag) const {
750 Index bse_exc = levA_ + levB_;
751 Index ct = J_dimer.rows() - bse_exc;
752 Diagnostics diag;
753
754 if (ct > 0) {
755 for (Index i = 0; i < bse_exc; i++) {
756 double Ei = J_dimer(i, i);
757 for (Index k = bse_exc; k < bse_exc + ct; k++) {
758 double Ek = J_dimer(k, k);
759 double dE = std::abs(Ei - Ek);
760 double coupling = std::abs(J_dimer(i, k));
761 if (dE > 1e-10) {
762 diag.xi = std::max(diag.xi, coupling / dE);
763 } else {
764 // exact degeneracy — downfolding is meaningless
765 diag.xi = std::numeric_limits<double>::infinity();
766 }
767 }
768 }
769 }
770
771 // PT vs RM discrepancy over all FE pairs
772 for (Index i = 0; i < levA_; i++) {
773 for (Index j = 0; j < levB_; j++) {
774 Index jd = j + levA_;
775 double diff = std::abs(J_pert(i, jd) - J_diag(i, jd));
776 diag.pt_rm_discrepancy = std::max(diag.pt_rm_discrepancy, diff);
777 }
778 }
779
780 diag.downfolding_safe = std::isfinite(diag.xi) && (diag.xi < 0.3) &&
781 (diag.pt_rm_discrepancy < 1e-4);
782
783 return diag;
784}
785
786// =============================================================================
787// ProjectExcitons
788// Top-level driver: orthogonalize CTs, form J_dimer/S_dimer, Lowdin-
789// orthogonalize, run perturbation and reduction methods, compute diagnostics.
790// Raw J_dimer and S_dimer (pre-Lowdin, in Hrt) and diagnostics are returned
791// via output references for storage as member variables.
792// =============================================================================
793template <class BSE_OPERATOR>
794std::array<Eigen::MatrixXd, 2> BSECoupling::ProjectExcitons(
795 Eigen::MatrixXd& FE_AB, Eigen::MatrixXd& CTStates, BSE_OPERATOR H,
796 Eigen::MatrixXd& J_dimer_out, Eigen::MatrixXd& S_dimer_out,
797 Diagnostics& diag_out) const {
798
799 Eigen::MatrixXd projection = OrthogonalizeCTs(FE_AB, CTStates);
800 Eigen::MatrixXd J_ortho =
801 CalcJ_dimer(H, projection, J_dimer_out, S_dimer_out);
802
803 std::array<Eigen::MatrixXd, 2> J;
805 << TimeStamp() << " Running Perturbation algorithm" << flush;
806 J[0] = Perturbation(J_ortho);
808 << TimeStamp() << " Running Projection algorithm" << flush;
809 J[1] = Fulldiag(J_ortho);
810
811 diag_out = ComputeDiagnostics(J_dimer_out, J[0], J[1]);
812
814 << TimeStamp() << " Diagnostics: xi=" << diag_out.xi
815 << " PT/RM discrepancy[Hrt]=" << diag_out.pt_rm_discrepancy
816 << " downfolding_safe=" << (diag_out.downfolding_safe ? "true" : "false")
817 << flush;
818
820 << "---------------------------------------" << flush;
821 XTP_LOG(Log::debug, *pLog_) << "Jeff_pert[Hrt]" << flush;
822 XTP_LOG(Log::debug, *pLog_) << J[0] << flush;
823 XTP_LOG(Log::debug, *pLog_) << "Jeff_diag[Hrt]" << flush;
824 XTP_LOG(Log::debug, *pLog_) << J[1] << flush;
826 << "---------------------------------------" << flush;
827
828 return J;
829}
830
831// =============================================================================
832// Perturbation
833// First-order perturbation theory correction for CT-mediated coupling.
834// The factor of 0.5 in the correction term comes from the symmetric treatment
835// of corrections |delta Phi_A> and |delta Phi_B> in eq 28 of Wehner/Baumeier
836// JCTC 2017 — each contributes half the total CT-mediated coupling.
837// =============================================================================
839 const Eigen::MatrixXd& J_dimer) const {
840 Index bse_exc = levA_ + levB_;
841 Index ct = J_dimer.rows() - bse_exc;
842 Eigen::MatrixXd J_result = J_dimer;
843 if (ct > 0) {
844 Eigen::MatrixXd transformation =
845 Eigen::MatrixXd::Identity(J_dimer.rows(), J_dimer.cols());
846 Eigen::MatrixXd Ct = J_dimer.bottomRightCorner(ct, ct);
847 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(Ct);
848 transformation.bottomRightCorner(ct, ct) = es.eigenvectors();
849 Ct.resize(0, 0);
850
851 XTP_LOG(Log::debug, *pLog_) << "FE state hamiltonian" << flush;
853 << J_dimer.topLeftCorner(bse_exc, bse_exc) << flush;
854 XTP_LOG(Log::debug, *pLog_) << "eigenvalues of CT states" << flush;
855 XTP_LOG(Log::debug, *pLog_) << es.eigenvalues().transpose() << flush;
856
857 J_result = transformation.transpose() * J_dimer * transformation;
859 << "---------------------------------------" << flush;
860 XTP_LOG(Log::debug, *pLog_) << " J_ortho_[Hrt] CT-state diag" << flush;
861 XTP_LOG(Log::debug, *pLog_) << J_result << flush;
863 << "---------------------------------------" << flush;
864 }
865
866 Eigen::MatrixXd Jmatrix = Eigen::MatrixXd::Zero(bse_exc, bse_exc);
867 for (Index stateA = 0; stateA < levA_; stateA++) {
868 double Ea = J_result(stateA, stateA);
869 for (Index stateB = 0; stateB < levB_; stateB++) {
870 Index stateBd = stateB + levA_;
872 << TimeStamp() << " Calculating coupling between exciton A"
873 << stateA + 1 << " and exciton B" << stateB + 1 << flush;
874 double J = J_result(stateA, stateBd);
875 double Eb = J_result(stateBd, stateBd);
876 for (Index k = bse_exc; k < (bse_exc + ct); k++) {
877 double Eab = J_result(k, k);
878 if (std::abs(Eab - Ea) < 0.001) {
880 << TimeStamp() << "Energydifference between state A "
881 << stateA + 1 << "and CT state " << k + 1 << " is " << Eab - Ea
882 << "[Hrt]" << flush;
883 }
884 if (std::abs(Eab - Eb) < 0.001) {
886 << TimeStamp() << "Energydifference between state B "
887 << stateB + 1 << "and CT state " << k + 1 << " is " << Eab - Eb
888 << "[Hrt]" << flush;
889 }
890 // Factor 0.5 from symmetric first-order correction: see eq 28,
891 // Wehner/Baumeier JCTC 2017. Each of |delta Phi_A> and |delta Phi_B>
892 // contributes equally, giving the 1/2 * (1/Delta_A + 1/Delta_B) form.
893 J += 0.5 * J_result(k, stateA) * J_result(k, stateBd) *
894 (1.0 / (Ea - Eab) + 1.0 / (Eb - Eab));
895 }
896 Jmatrix(stateA, stateBd) = J;
897 Jmatrix(stateBd, stateA) = J;
898 }
899 }
900 return Jmatrix;
901}
902
903// =============================================================================
904// Fulldiag (reduction method)
905// Diagonalizes J_ortho (full FE+CT system), then for each (stateA, stateB)
906// pair independently extracts the two eigenvectors most resembling those
907// states, projects to the 2x2 FE subspace, Lowdin-orthogonalizes, and
908// rotates back to recover the effective coupling. This per-pair 2x2 approach
909// is consistent with the reduction method of Wehner/Baumeier JCTC 2017 and
910// is robust at CT/FE near-resonance where perturbation theory diverges.
911// =============================================================================
912Eigen::MatrixXd BSECoupling::Fulldiag(const Eigen::MatrixXd& J_dimer) const {
913 Index bse_exc = levA_ + levB_;
914 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(J_dimer);
916 << "---------------------------------------" << flush;
917 XTP_LOG(Log::debug, *pLog_) << "Eigenvectors of J" << flush;
918 XTP_LOG(Log::debug, *pLog_) << es.eigenvectors() << flush;
919 XTP_LOG(Log::debug, *pLog_) << "J_eigenvalues[Hrt]" << flush;
920 XTP_LOG(Log::debug, *pLog_) << es.eigenvalues() << flush;
922 << "---------------------------------------" << flush;
923
924 Eigen::MatrixXd Jmat = Eigen::MatrixXd::Zero(bse_exc, bse_exc);
925 for (Index stateA = 0; stateA < levA_; stateA++) {
926 for (Index stateB = 0; stateB < levB_; stateB++) {
927 Index stateBd = stateB + levA_;
929 << TimeStamp() << " Calculating coupling between exciton A"
930 << stateA + 1 << " and exciton B" << stateB + 1 << flush;
931
932 std::array<int, 2> indexes;
933 std::array<int, 2> signs;
934
935 // Find the eigenstates closest to stateA and stateBd by maximum overlap
936 es.eigenvectors().row(stateA).cwiseAbs().maxCoeff(&indexes[0]);
937 es.eigenvectors().row(stateBd).cwiseAbs().maxCoeff(&indexes[1]);
938 if (indexes[0] == indexes[1]) {
939 // same eigenvector selected for both: pick next best for stateBd
940 Eigen::RowVectorXd stateamplitudes =
941 es.eigenvectors().row(stateBd).cwiseAbs();
942 stateamplitudes[indexes[1]] = 0.0;
943 stateamplitudes.maxCoeff(&indexes[1]);
944 }
945
946 signs[0] = GetSign(es.eigenvectors()(stateA, indexes[0]));
947 signs[1] = GetSign(es.eigenvectors()(stateBd, indexes[1]));
948
950 << TimeStamp() << " Order is: [Initial state n->nth eigenvalue]"
951 << flush;
952 XTP_LOG(Log::info, *pLog_) << " A" << stateA + 1 << ":" << stateA + 1
953 << "->" << indexes[0] + 1 << " ";
954 XTP_LOG(Log::info, *pLog_) << " B" << stateB + 1 << ":" << stateBd + 1
955 << "->" << indexes[1] + 1 << " " << flush;
956
957 // Build 2x2 transformation matrix from the two selected eigenvectors
958 // (rows stateA and stateBd only) and diagonal eigenvalue matrix
959 Eigen::Matrix2d Emat = Eigen::Matrix2d::Zero();
960 Eigen::Matrix2d Tmat = Eigen::Matrix2d::Zero();
961 for (Index i = 0; i < 2; i++) {
962 Index k = indexes[i];
963 double sign = signs[i];
964 Tmat(0, i) = sign * es.eigenvectors()(stateA, k);
965 Tmat(1, i) = sign * es.eigenvectors()(stateBd, k);
966 Emat(i, i) = es.eigenvalues()(k);
967 }
968 Tmat.colwise().normalize();
969
970 if (Tmat.determinant() < 0) {
972 << " Reduced state matrix is not in a right handed basis, "
973 "multiplying second eigenvector by -1 "
974 << flush;
975 Tmat.col(1) *= -1;
976 }
977
979 << "---------------------------------------" << flush;
980 XTP_LOG(Log::debug, *pLog_) << " T_" << flush;
981 XTP_LOG(Log::debug, *pLog_) << Tmat << flush;
982
983 // Lowdin orthogonalize the 2x2 projected subspace
984 Eigen::Matrix2d S_small = Tmat * Tmat.transpose();
985 XTP_LOG(Log::debug, *pLog_) << "S_small" << flush;
986 XTP_LOG(Log::debug, *pLog_) << S_small << flush;
987
988 Eigen::SelfAdjointEigenSolver<Eigen::Matrix2d> ss(S_small);
989 Eigen::Matrix2d sm1 = ss.operatorInverseSqrt();
990 Emat = sm1 * Emat * sm1;
991
993 << TimeStamp() << " Smallest value of dimer overlapmatrix is "
994 << ss.eigenvalues()(0) << flush;
995 XTP_LOG(Log::debug, *pLog_) << "S-1/2" << flush;
996 XTP_LOG(Log::debug, *pLog_) << sm1 << flush;
997 XTP_LOG(Log::debug, *pLog_) << "E_ortho" << flush;
998 XTP_LOG(Log::debug, *pLog_) << Emat << flush;
999
1000 Tmat = Tmat * sm1;
1001 XTP_LOG(Log::debug, *pLog_) << "T_ortho" << flush;
1002 XTP_LOG(Log::debug, *pLog_) << Tmat << flush;
1004 << "---------------------------------------" << flush;
1005
1006 // H_eff = T * E * T^T; off-diagonal is the effective coupling
1007 Eigen::Matrix2d J_small = Tmat * Emat * Tmat.transpose();
1008 XTP_LOG(Log::debug, *pLog_) << "T_ortho*E_ortho*T_ortho^T" << flush;
1009 XTP_LOG(Log::debug, *pLog_) << J_small << flush;
1010
1011 Jmat(stateA, stateBd) = J_small(0, 1);
1012 Jmat(stateBd, stateA) = J_small(1, 0);
1013 }
1014 }
1015 return Jmat;
1016}
1017
1018} // namespace xtp
1019} // 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
T ifExistsReturnElseReturnDefault(const std::string &key, T defaultvalue) const
Definition property.h:321
void setAttribute(const std::string &attribute, const T &value)
set an attribute
Definition property.h:314
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Index AOBasisSize() const
Definition aobasis.h:46
Eigen::VectorXd monomerA_energies_triplet_
Eigen::MatrixXd J_dimer_triplet_
double getSingletCouplingElement(Index levelA, Index levelB, Index methodindex) const
Eigen::VectorXd monomerB_energies_singlet_
double getTripletCouplingElement(Index levelA, Index levelB, Index methodindex) const
std::string Identify() const
Definition bsecoupling.h:42
void CalculateCouplings(const Orbitals &orbitalsA, const Orbitals &orbitalsB, const Orbitals &orbitalsAB) override
evaluates electronic couplings
void Addoutput(tools::Property &type_summary, const Orbitals &orbitalsA, const Orbitals &orbitalsB) const override
void WriteToProperty(tools::Property &summary, const QMState &stateA, const QMState &stateB) const
Diagnostics ComputeDiagnostics(const Eigen::MatrixXd &J_dimer, const Eigen::MatrixXd &J_pert, const Eigen::MatrixXd &J_diag) const
Compute diagnostics from raw J_dimer and the two effective coupling matrices.
std::array< Eigen::MatrixXd, 2 > JAB_triplet
Eigen::MatrixXd Perturbation(const Eigen::MatrixXd &J_dimer) const
Eigen::VectorXd monomerB_energies_triplet_
std::array< Eigen::MatrixXd, 2 > ProjectExcitons(Eigen::MatrixXd &FE_AB, Eigen::MatrixXd &CTStates, BSE_OPERATOR H, Eigen::MatrixXd &J_dimer_out, Eigen::MatrixXd &S_dimer_out, Diagnostics &diag_out) const
Compute J_pert and J_diag, and store raw J_dimer/S_dimer and diagnostics for later output.
void Initialize(tools::Property &options) override
Eigen::MatrixXd S_dimer_singlet_
Eigen::VectorXd monomerA_energies_singlet_
Eigen::MatrixXd S_dimer_triplet_
std::array< Eigen::MatrixXd, 2 > JAB_singlet
Eigen::MatrixXd ProjectFrenkelExcitons(const Eigen::MatrixXd &BSE_Coeffs, const Eigen::MatrixXd &X_AB, Index bseX_vtotal, Index bseX_ctotal, Index bseAB_vtotal, Index bseAB_ctotal) const
static void WriteMatrixToProperty(tools::Property &prop, const std::string &name, const Eigen::MatrixXd &mat, double conversion=1.0)
Write an Eigen matrix as rows of space-separated values into an XML property node....
Eigen::MatrixXd SetupCTStates(Index bseA_vtotal, Index bseB_vtotal, Index bseAB_vtotal, Index bseAB_ctotal, const Eigen::MatrixXd &A_AB, const Eigen::MatrixXd &B_AB) const
Eigen::MatrixXd Fulldiag(const Eigen::MatrixXd &J_dimer) const
Eigen::MatrixXd CalcJ_dimer(BSE_OPERATOR &H, Eigen::MatrixXd &projection, Eigen::MatrixXd &J_dimer_out, Eigen::MatrixXd &S_dimer_out) const
Form J_dimer and S_dimer from the projection, then Lowdin orthogonalize to produce J_ortho.
Eigen::MatrixXd OrthogonalizeCTs(Eigen::MatrixXd &FE_AB, Eigen::MatrixXd &CTStates) const
Merge FE and CT projection vectors into a single projection matrix without pre-orthogonalization.
Eigen::MatrixXd J_dimer_singlet_
TripletOperator_TDA getTripletOperator_TDA() const
Definition bse.cc:259
SingletOperator_TDA getSingletOperator_TDA() const
Definition bse.cc:252
void configure(const options &opt, const Eigen::VectorXd &RPAEnergies, const Eigen::MatrixXd &Hqp_in)
Definition bse.cc:43
Eigen::MatrixXd CalculateOverlapMatrix(const Orbitals &orbitalsAB) const
void CheckAtomCoordinates(const Orbitals &orbitalsA, const Orbitals &orbitalsB, const Orbitals &orbitalsAB) const
Container for molecular orbitals and derived one-particle data.
Definition orbitals.h:47
Index getBSEvmax() const
Return the highest valence orbital included in BSE.
Definition orbitals.h:410
const tools::EigenSystem & BSETriplets() const
Return read-only access to triplet BSE eigenpairs.
Definition orbitals.h:464
bool GetFlagUseHqpOffdiag() const
Definition orbitals.h:595
Index getBSEcmin() const
Return the lowest conduction orbital included in BSE.
Definition orbitals.h:413
const tools::EigenSystem & QPdiag() const
Return read-only access to the diagonalized quasiparticle representation.
Definition orbitals.h:454
Index getHomo() const
Return the default HOMO index used by spin-agnostic callers.
Definition orbitals.h:107
const tools::EigenSystem & BSESinglets() const
Return read-only access to singlet BSE eigenpairs.
Definition orbitals.h:477
Index getBSEcmax() const
Return the highest conduction orbital included in BSE.
Definition orbitals.h:416
const AOBasis & getAuxBasis() const
Return the auxiliary AO basis, throwing if it has not been initialized.
Definition orbitals.h:325
Index getBSEvmin() const
Return the lowest valence orbital included in BSE.
Definition orbitals.h:407
Index getRPAmax() const
Return the upper RPA orbital index.
Definition orbitals.h:381
const std::vector< Eigen::Vector3d > & TransitionDipoles() const
Return the stored transition dipole moments.
Definition orbitals.h:521
Index getBasisSetSize() const
Return the number of AO basis functions in the DFT basis.
Definition orbitals.h:72
const Eigen::VectorXd & RPAInputEnergies() const
Return read-only access to the RPA input energies.
Definition orbitals.h:429
Index getGWAmin() const
Return the lower GW orbital index.
Definition orbitals.h:361
const tools::EigenSystem & MOs() const
Return read-only access to alpha/restricted molecular orbitals.
Definition orbitals.h:192
Index getGWAmax() const
Return the upper GW orbital index.
Definition orbitals.h:364
const std::string & getDFTbasisName() const
Return the DFT basis-set name.
Definition orbitals.h:304
Index getRPAmin() const
Return the lower RPA orbital index.
Definition orbitals.h:378
const AOBasis & getDftBasis() const
Return the DFT AO basis, throwing if it has not been initialized.
Definition orbitals.h:314
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:135
std::string ToString() const
Definition qmstate.cc:155
Index StateIdx() const
Definition qmstate.h:157
const QMStateType & Type() const
Definition qmstate.h:154
void Fill(const AOBasis &auxbasis, const AOBasis &dftbasis, const Eigen::MatrixXd &dft_orbitals)
void Initialize(Index basissize, Index mmin, Index mmax, Index nmin, Index nmax)
Timestamp returns the current time as a string Example: cout << TimeStamp().
Definition logger.h:224
#define XTP_LOG(level, log)
Definition logger.h:40
const double hrt2ev
Definition constants.h:53
Index getMaxThreads()
Definition eigen.h:128
Charge transport classes.
Definition ERIs.h:28
int GetSign(double value)
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Diagnostic quantities for assessing whether CT downfolding is safe.
Definition bsecoupling.h:69