21#include <boost/format.hpp>
39 string spintype = options.
get(
"spin").
as<std::string>();
40 if (spintype ==
"all") {
43 }
else if (spintype ==
"triplet") {
46 }
else if (spintype ==
"singlet") {
50 throw std::runtime_error(
52 "Choice % for type not known. Available singlet,triplet,all") %
77 const std::string& name,
78 const Eigen::MatrixXd& mat,
83 for (
Index i = 0; i < mat.rows(); i++) {
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();
89 mat_prop.
setAttribute(
"row_" + std::to_string(i), row_str);
100 Property& coupling_summary = summary.
add(
"coupling",
"");
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());
131 string algorithm =
"j_diag";
133 algorithm =
"j_pert";
137 auto WriteSpinChannel = [&](
const std::string& spin_name,
138 const Eigen::MatrixXd& J_dimer,
139 const Eigen::MatrixXd& S_dimer,
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,
146 Property& spin_summary = bsecoupling.
add(spin_name,
"");
150 for (
Index stateA = 0; stateA <
levA_; ++stateA) {
152 for (
Index stateB = 0; stateB <
levB_; ++stateB) {
164 Property& mono_prop = spin_summary.
add(
"monomer_energies",
"");
166 for (
Index i = 0; i < monA_energies.size(); i++) {
168 QMState st(spin_type, i,
false);
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());
178 for (
Index i = 0; i < monB_energies.size(); i++) {
180 QMState st(spin_type, i,
false);
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());
193 Property& diag_prop = spin_summary.
add(
"diagnostics",
"");
194 diag_prop.
setAttribute(
"xi", (format(
"%1$1.4f") % diag.xi).str());
196 "pt_rm_discrepancy_eV",
201 diag.downfolding_safe ?
"true" :
"false");
208 Index ct = J_dimer.rows() - bse_exc;
210 Property& tb_prop = spin_summary.
add(
"tb_matrices",
"");
219 J_dimer.topLeftCorner(bse_exc, bse_exc),
222 S_dimer.topLeftCorner(bse_exc, bse_exc));
226 J_dimer.topRightCorner(bse_exc, ct),
229 S_dimer.topRightCorner(bse_exc, ct));
231 J_dimer.bottomRightCorner(ct, ct),
234 S_dimer.bottomRightCorner(ct, ct));
251 const std::vector<Eigen::Vector3d> empty_dipoles;
255 empty_dipoles, empty_dipoles,
261 Index methodindex)
const {
267 Index methodindex)
const {
275 const Eigen::MatrixXd& A_AB,
276 const Eigen::MatrixXd& B_AB)
const {
279 Index bseAB_size = bseAB_vtotal * bseAB_ctotal;
280 Eigen::MatrixXd CTstates = Eigen::MatrixXd::Zero(bseAB_size, noAB + noBA);
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_);
287 const Eigen::MatrixXd A_occ_occ = A_occ.topRows(bseAB_vtotal);
288 const Eigen::MatrixXd B_unocc_unocc = B_unocc.bottomRows(bseAB_ctotal);
292 <<
TimeStamp() <<
" Setting up CT-states" << flush;
294#pragma omp parallel for
295 for (
Index a_occ = 0; a_occ <
occA_; a_occ++) {
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);
305 <<
TimeStamp() <<
" " << noBA <<
" CT states A+B- created" << flush;
307 const Eigen::MatrixXd A_unocc_unocc = A_unocc.bottomRows(bseAB_ctotal);
308 const Eigen::MatrixXd B_occ_occ = B_occ.topRows(bseAB_vtotal);
310#pragma omp parallel for
311 for (
Index b_occ = 0; b_occ <
occB_; b_occ++) {
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);
321 <<
TimeStamp() <<
" " << noBA <<
" CT states A-B+ created" << flush;
326 const Eigen::MatrixXd& BSE_Coeffs,
const Eigen::MatrixXd& X_AB,
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());
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());
349 }
else if (value > 0) {
359 <<
TimeStamp() <<
" Calculating exciton couplings" << flush;
367 if ((basisA == 0) || (basisB == 0)) {
368 throw std::runtime_error(
"Basis set size is not stored in monomers");
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;
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;
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;
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;
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"
416 levA_ = bseA_singlet_exc;
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"
424 levB_ = bseB_singlet_exc;
428 if (
levA_ > bseA_triplet_exc) {
431 <<
" Number of Frenkel states you want is greater than stored for "
432 "molecule A. Setting to max number available"
434 levA_ = bseA_triplet_exc;
436 if (
levB_ > bseB_triplet_exc) {
439 <<
" Number of Frenkel states you want is greater than stored for "
440 "molecule B. Setting to max number available"
442 levB_ = bseB_triplet_exc;
448 <<
" Number of unoccupied orbitals in molecule A for CT creation "
449 "exceeds number of KS-orbitals in BSE"
456 <<
" Number of unoccupied orbitals in molecule B for CT creation "
457 "exceeds number of KS-orbitals in BSE"
464 <<
" Number of occupied orbitals in molecule A for CT creation "
465 "exceeds number of KS-orbitals in BSE"
472 <<
" Number of occupied orbitals in molecule B for CT creation "
473 "exceeds number of KS-orbitals in BSE"
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;
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;
498 Eigen::MatrixXd MOsA =
500 Eigen::MatrixXd MOsB =
502 Eigen::MatrixXd MOsAB =
506 <<
TimeStamp() <<
" Calculating overlap matrix for basisset: "
511 <<
TimeStamp() <<
" Projecting monomers onto dimer orbitals" << flush;
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) {
519 <<
"Projection of orbitals of monomer A on dimer is insufficient,mag="
520 << mag_A.minCoeff() << flush;
522 Eigen::VectorXd mag_B = B_AB.colwise().squaredNorm();
523 if (mag_B.any() < 0.95) {
526 <<
"Projection of orbitals of monomer B on dimer is insufficient,mag="
527 << mag_B.minCoeff() << flush;
539 Eigen::MatrixXd Hqp = qpcoeff *
559 <<
TimeStamp() <<
" Evaluating singlets" << flush;
561 <<
TimeStamp() <<
" Setup Hamiltonian" << flush;
573 Eigen::MatrixXd FE_AB = Eigen::MatrixXd::Zero(bseAB_size,
levA_ +
levB_);
574 const Eigen::MatrixXd bseA =
577 bseA, A_AB, bseA_vtotal, bseA_ctotal, bseAB_vtotal, bseAB_ctotal);
578 const Eigen::MatrixXd bseB =
581 bseB, B_AB, bseB_vtotal, bseB_ctotal, bseAB_vtotal, bseAB_ctotal);
583 bseA_vtotal, bseB_vtotal, bseAB_vtotal, bseAB_ctotal, A_AB, B_AB);
588 <<
TimeStamp() <<
" calculated singlet couplings " << flush;
592 <<
TimeStamp() <<
" Evaluating triplets" << flush;
602 Eigen::MatrixXd FE_AB = Eigen::MatrixXd::Zero(bseAB_size,
levA_ +
levB_);
603 const Eigen::MatrixXd bseA =
606 bseA, A_AB, bseA_vtotal, bseA_ctotal, bseAB_vtotal, bseAB_ctotal);
607 const Eigen::MatrixXd bseB =
610 bseB, B_AB, bseB_vtotal, bseB_ctotal, bseAB_vtotal, bseAB_ctotal);
612 bseA_vtotal, bseB_vtotal, bseAB_vtotal, bseAB_ctotal, A_AB, B_AB);
617 <<
TimeStamp() <<
" calculated triplet couplings " << flush;
620 <<
TimeStamp() <<
" Done with exciton couplings" << flush;
646 Eigen::MatrixXd& CTStates)
const {
647 Index ct = CTStates.cols();
649 Index bseAB_size = FE_AB.rows();
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;
666 <<
" CT states into projection (no pre-orthogonalization)" << flush;
669 Eigen::MatrixXd projection(bseAB_size, bse_exc + ct);
671 <<
TimeStamp() <<
" merging projections into one vector " << flush;
672 projection.leftCols(bse_exc) = FE_AB;
675 projection.rightCols(ct) = CTStates;
677 CTStates.resize(0, 0);
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;
699 J_dimer_out = projection.transpose() * (
H * projection).eval();
702 <<
TimeStamp() <<
" Setting up overlap matrix size "
703 << projection.cols() << flush;
704 S_dimer_out = projection.transpose() * projection;
706 projection.resize(0, 0);
709 <<
"---------------------------------------" << flush;
715 <<
"---------------------------------------" << flush;
717 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(S_dimer_out);
719 <<
TimeStamp() <<
" Smallest value of dimer overlapmatrix is "
720 << es.eigenvalues()(0) << flush;
722 Eigen::MatrixXd Sm1 = es.operatorInverseSqrt();
723 Eigen::MatrixXd J_ortho = Sm1 * J_dimer_out * Sm1;
726 <<
"---------------------------------------" << flush;
732 <<
"---------------------------------------" << flush;
748 const Eigen::MatrixXd& J_dimer,
const Eigen::MatrixXd& J_pert,
749 const Eigen::MatrixXd& J_diag)
const {
751 Index ct = J_dimer.rows() - bse_exc;
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));
762 diag.
xi = std::max(diag.
xi, coupling / dE);
765 diag.
xi = std::numeric_limits<double>::infinity();
775 double diff = std::abs(J_pert(i, jd) - J_diag(i, jd));
793template <
class BSE_OPERATOR>
795 Eigen::MatrixXd& FE_AB, Eigen::MatrixXd& CTStates,
BSE_OPERATOR H,
796 Eigen::MatrixXd& J_dimer_out, Eigen::MatrixXd& S_dimer_out,
800 Eigen::MatrixXd J_ortho =
803 std::array<Eigen::MatrixXd, 2> J;
805 <<
TimeStamp() <<
" Running Perturbation algorithm" << flush;
808 <<
TimeStamp() <<
" Running Projection algorithm" << flush;
814 <<
TimeStamp() <<
" Diagnostics: xi=" << diag_out.
xi
820 <<
"---------------------------------------" << flush;
826 <<
"---------------------------------------" << flush;
839 const Eigen::MatrixXd& J_dimer)
const {
841 Index ct = J_dimer.rows() - bse_exc;
842 Eigen::MatrixXd J_result = J_dimer;
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();
853 << J_dimer.topLeftCorner(bse_exc, bse_exc) << flush;
857 J_result = transformation.transpose() * J_dimer * transformation;
859 <<
"---------------------------------------" << flush;
863 <<
"---------------------------------------" << flush;
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++) {
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
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
893 J += 0.5 * J_result(k, stateA) * J_result(k, stateBd) *
894 (1.0 / (Ea - Eab) + 1.0 / (Eb - Eab));
896 Jmatrix(stateA, stateBd) = J;
897 Jmatrix(stateBd, stateA) = J;
914 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(J_dimer);
916 <<
"---------------------------------------" << flush;
922 <<
"---------------------------------------" << flush;
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++) {
929 <<
TimeStamp() <<
" Calculating coupling between exciton A"
930 << stateA + 1 <<
" and exciton B" << stateB + 1 << flush;
932 std::array<int, 2> indexes;
933 std::array<int, 2> signs;
936 es.eigenvectors().row(stateA).cwiseAbs().maxCoeff(&indexes[0]);
937 es.eigenvectors().row(stateBd).cwiseAbs().maxCoeff(&indexes[1]);
938 if (indexes[0] == indexes[1]) {
940 Eigen::RowVectorXd stateamplitudes =
941 es.eigenvectors().row(stateBd).cwiseAbs();
942 stateamplitudes[indexes[1]] = 0.0;
943 stateamplitudes.maxCoeff(&indexes[1]);
946 signs[0] =
GetSign(es.eigenvectors()(stateA, indexes[0]));
947 signs[1] =
GetSign(es.eigenvectors()(stateBd, indexes[1]));
950 <<
TimeStamp() <<
" Order is: [Initial state n->nth eigenvalue]"
953 <<
"->" << indexes[0] + 1 <<
" ";
955 <<
"->" << indexes[1] + 1 <<
" " << flush;
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);
968 Tmat.colwise().normalize();
970 if (Tmat.determinant() < 0) {
972 <<
" Reduced state matrix is not in a right handed basis, "
973 "multiplying second eigenvector by -1 "
979 <<
"---------------------------------------" << flush;
984 Eigen::Matrix2d S_small = Tmat * Tmat.transpose();
988 Eigen::SelfAdjointEigenSolver<Eigen::Matrix2d> ss(S_small);
989 Eigen::Matrix2d sm1 = ss.operatorInverseSqrt();
990 Emat = sm1 * Emat * sm1;
993 <<
TimeStamp() <<
" Smallest value of dimer overlapmatrix is "
994 << ss.eigenvalues()(0) << flush;
1004 <<
"---------------------------------------" << flush;
1007 Eigen::Matrix2d J_small = Tmat * Emat * Tmat.transpose();
1011 Jmat(stateA, stateBd) = J_small(0, 1);
1012 Jmat(stateBd, stateA) = J_small(1, 0);
Container to hold Basisfunctions for all atoms.
Index AOBasisSize() const
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
Diagnostics diag_triplet_
bool output_perturbation_
std::string Identify() const
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
Diagnostics diag_singlet_
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
SingletOperator_TDA getSingletOperator_TDA() const
void configure(const options &opt, const Eigen::VectorXd &RPAEnergies, const Eigen::MatrixXd &Hqp_in)
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.
Index getBSEvmax() const
Return the highest valence orbital included in BSE.
const tools::EigenSystem & BSETriplets() const
Return read-only access to triplet BSE eigenpairs.
bool GetFlagUseHqpOffdiag() const
Index getBSEcmin() const
Return the lowest conduction orbital included in BSE.
const tools::EigenSystem & QPdiag() const
Return read-only access to the diagonalized quasiparticle representation.
Index getHomo() const
Return the default HOMO index used by spin-agnostic callers.
const tools::EigenSystem & BSESinglets() const
Return read-only access to singlet BSE eigenpairs.
Index getBSEcmax() const
Return the highest conduction orbital included in BSE.
const AOBasis & getAuxBasis() const
Return the auxiliary AO basis, throwing if it has not been initialized.
Index getBSEvmin() const
Return the lowest valence orbital included in BSE.
Index getRPAmax() const
Return the upper RPA orbital index.
const std::vector< Eigen::Vector3d > & TransitionDipoles() const
Return the stored transition dipole moments.
Index getBasisSetSize() const
Return the number of AO basis functions in the DFT basis.
const Eigen::VectorXd & RPAInputEnergies() const
Return read-only access to the RPA input energies.
Index getGWAmin() const
Return the lower GW orbital index.
const tools::EigenSystem & MOs() const
Return read-only access to alpha/restricted molecular orbitals.
Index getGWAmax() const
Return the upper GW orbital index.
const std::string & getDFTbasisName() const
Return the DFT basis-set name.
Index getRPAmin() const
Return the lower RPA orbital index.
const AOBasis & getDftBasis() const
Return the DFT AO basis, throwing if it has not been initialized.
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
std::string ToString() const
const QMStateType & Type() const
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().
#define XTP_LOG(level, log)
Charge transport classes.
int GetSign(double value)
Provides a means for comparing floating point numbers.
Diagnostic quantities for assessing whether CT downfolding is safe.