votca 2026-dev
Loading...
Searching...
No Matches
dftcoupling.h
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#pragma once
21#ifndef VOTCA_XTP_DFTCOUPLING_H
22#define VOTCA_XTP_DFTCOUPLING_H
23
24// Local VOTCA includes
25#include "couplingbase.h"
26
27namespace votca {
28namespace xtp {
29
37
38class DFTcoupling : public CouplingBase {
39 public:
40 std::string Identify() const { return "dftcoupling"; }
41
42 void CalculateCouplings(const Orbitals& orbitalsA, const Orbitals& orbitalsB,
43 const Orbitals& orbitalsAB) override;
44
45 void Initialize(tools::Property&) override;
46
47 void Addoutput(tools::Property& type_summary, const Orbitals& orbitalsA,
48 const Orbitals& orbitalsB) const override;
49
50 private:
51 void WriteToProperty(tools::Property& type_summary, const Orbitals& orbitalsA,
52 const Orbitals& orbitalsB, Index a, Index b) const;
53
54 double getCouplingElement(Index levelA, Index levelB,
55 const Orbitals& orbitalsA,
56 const Orbitals& orbitalsB) const;
57
58 std::pair<Index, Index> DetermineRangeOfStates(const Orbitals& orbital,
59 Index numberofstates) const;
60
70 static void WriteMatrixToProperty(tools::Property& prop,
71 const std::string& name,
72 const Eigen::MatrixXd& mat,
73 double conversion = 1.0);
74
75 // --- existing: Lowdin-orthogonalized effective coupling matrix [Hrt] ---
76 // JAB(i,j) where i in [0, levelsA) indexes fragment A MOs and
77 // j in [levelsA, levelsA+levelsB) indexes fragment B MOs.
78 // Both hole and electron blocks share this matrix, partitioned by range.
79 Eigen::MatrixXd JAB;
80
81 // --- new: raw (pre-Lowdin) TB matrices [Hrt for H, dimensionless for S] ---
82 // Stored separately for holes and electrons since they are computed from
83 // different MO ranges and assembled into the TB Hamiltonian independently.
84 //
85 // Block structure of each matrix (hole or electron):
86 // H = [ H_AA H_AB ] S = [ S_AA S_AB ]
87 // [ H_BA H_BB ] [ S_BA S_BB ]
88 // where A = fragment A MOs, B = fragment B MOs for that carrier type.
89 // H is in Hartree; conversion to eV happens in Addoutput.
90 Eigen::MatrixXd JAB_dimer_hole_; // pre-Lowdin H, hole block [Hrt]
91 Eigen::MatrixXd S_AxB_hole_; // overlap, hole block
92 Eigen::MatrixXd JAB_dimer_elec_; // pre-Lowdin H, electron block [Hrt]
93 Eigen::MatrixXd S_AxB_elec_; // overlap, electron block
94
95 // --- new: monomer MO energies [eV] ---
96 // Isolated monomer quasiparticle energies for the relevant MO ranges.
97 // These provide pairwise-consistent TB site energies: the same monomer
98 // calculation is used regardless of which dimer partner is processed,
99 // unlike the dimer Hamiltonian diagonal which depends on the partner.
100 // Environmental corrections (e.g. from polarisation embedding) should
101 // be added on top of these values externally.
102 //
103 // KS: Kohn-Sham DFT eigenvalues, always available.
104 // QP: perturbative GW quasiparticle energies, present only when a GW
105 // calculation was run on the monomer. Stored as empty vector otherwise.
106 // Indexed identically to the KS vectors — same MO range, same ordering.
107 Eigen::VectorXd moEnergiesA_hole_KS_; // fragment A hole MOs, KS [eV]
108 Eigen::VectorXd moEnergiesB_hole_KS_;
109 Eigen::VectorXd moEnergiesA_elec_KS_; // fragment A electron MOs, KS [eV]
110 Eigen::VectorXd moEnergiesB_elec_KS_;
111
112 Eigen::VectorXd moEnergiesA_hole_QP_; // fragment A hole MOs, QPpert [eV]
113 Eigen::VectorXd moEnergiesB_hole_QP_; // empty if GW not available
114 Eigen::VectorXd moEnergiesA_elec_QP_;
115 Eigen::VectorXd moEnergiesB_elec_QP_;
116
117 // --- new: basis diagnostics ---
118 // Minimum eigenvalue of the overlap matrix S_AxB for each carrier type.
119 // Small values indicate near-linear dependence in the MO projection basis,
120 // which would make the Lowdin orthogonalization ill-conditioned.
123
124 double degeneracy_ = 0.0;
127 // When true, write monomer MO energies, raw H/S TB matrices, and
128 // diagnostics to the XML output. Set false (default) for KMC/rate
129 // workflows that only need scalar effective couplings.
130 bool output_tb_ = false;
131
132 std::pair<Index, Index> Range_orbA;
133 std::pair<Index, Index> Range_orbB;
134};
135
136} // namespace xtp
137} // namespace votca
138
139#endif // VOTCA_XTP_DFTCOUPLING_H
class to manage program options with xml serialization functionality
Definition property.h:55
Base Class to derive DFT and BSE coupling from.
Evaluates electronic coupling elements.
Definition dftcoupling.h:38
void WriteToProperty(tools::Property &type_summary, const Orbitals &orbitalsA, const Orbitals &orbitalsB, Index a, Index b) const
Eigen::MatrixXd JAB_dimer_elec_
Definition dftcoupling.h:92
Eigen::VectorXd moEnergiesA_hole_QP_
Eigen::MatrixXd JAB_dimer_hole_
Definition dftcoupling.h:90
std::pair< Index, Index > Range_orbA
std::pair< Index, Index > DetermineRangeOfStates(const Orbitals &orbital, Index numberofstates) const
void Initialize(tools::Property &) override
Eigen::VectorXd moEnergiesB_elec_KS_
double getCouplingElement(Index levelA, Index levelB, const Orbitals &orbitalsA, const Orbitals &orbitalsB) const
Eigen::MatrixXd S_AxB_hole_
Definition dftcoupling.h:91
Eigen::VectorXd moEnergiesB_elec_QP_
void Addoutput(tools::Property &type_summary, const Orbitals &orbitalsA, const Orbitals &orbitalsB) const override
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::VectorXd moEnergiesB_hole_KS_
Eigen::VectorXd moEnergiesA_elec_KS_
Eigen::VectorXd moEnergiesA_hole_KS_
Eigen::VectorXd moEnergiesA_elec_QP_
std::pair< Index, Index > Range_orbB
void CalculateCouplings(const Orbitals &orbitalsA, const Orbitals &orbitalsB, const Orbitals &orbitalsAB) override
Eigen::VectorXd moEnergiesB_hole_QP_
Eigen::MatrixXd S_AxB_elec_
Definition dftcoupling.h:93
std::string Identify() const
Definition dftcoupling.h:40
Eigen::MatrixXd JAB
Definition dftcoupling.h:79
Container for molecular orbitals and derived one-particle data.
Definition orbitals.h:47
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26