votca 2024.2-dev
Loading...
Searching...
No Matches
populationanalysis.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// Local VOTCA includes
22
23namespace votca {
24namespace xtp {
25
26template <bool T>
28 const Orbitals& orbitals, const QMState& state) const {
29
30 AOBasis basis = orbitals.getDftBasis();
31 Eigen::MatrixXd dmat = orbitals.DensityMatrixFull(state);
32 AOOverlap overlap;
33 overlap.Fill(basis);
34 Eigen::VectorXd charges = -CalcElecChargeperAtom(dmat, overlap, basis);
35 if (!state.isTransition()) {
36 charges += CalcNucChargeperAtom(orbitals.QMAtoms());
37 }
38
39 StaticSegment seg =
40 StaticSegment(orbitals.QMAtoms().getType(), orbitals.QMAtoms().getId());
41 for (Index i = 0; i < orbitals.QMAtoms().size(); ++i) {
42 seg.push_back(StaticSite(orbitals.QMAtoms()[i], charges(i)));
43 }
44 return seg;
45}
46
47template <bool T>
49 std::vector<QMFragment<BSE_Population> >& frags, const Orbitals& orbitals,
50 QMStateType type) const {
51 if (!type.isExciton()) {
52 throw std::runtime_error(
53 "CalcChargeperFragment: QmStateType must be an exciton");
54 }
55 AOBasis basis = orbitals.getDftBasis();
56 AOOverlap overlap;
57 overlap.Fill(basis);
58 Eigen::VectorXd nuccharges = CalcNucChargeperAtom(orbitals.QMAtoms());
59 Eigen::MatrixXd dmatgs = orbitals.DensityMatrixGroundState();
60 if (orbitals.getCalculationType() != "NoEmbedding") {
61 dmatgs += orbitals.getInactiveDensity();
62 }
63 Eigen::VectorXd gscharges =
64 nuccharges - CalcElecChargeperAtom(dmatgs, overlap, basis);
65 Index numofstates = orbitals.NumberofStates(type);
66 for (auto& frag : frags) {
67 frag.value().Initialize(numofstates);
68 frag.value().Gs = frag.ExtractFromVector(gscharges);
69 }
70 for (Index i_state = 0; i_state < numofstates; i_state++) {
71 QMState state(type, i_state, false);
72 std::array<Eigen::MatrixXd, 2> dmat_ex =
73 orbitals.DensityMatrixExcitedState(state);
74 Eigen::VectorXd atom_h = CalcElecChargeperAtom(dmat_ex[0], overlap, basis);
75 Eigen::VectorXd atom_e = -CalcElecChargeperAtom(dmat_ex[1], overlap, basis);
76 for (auto& frag : frags) {
77 frag.value().E(i_state) = frag.ExtractFromVector(atom_e);
78 frag.value().H(i_state) = frag.ExtractFromVector(atom_h);
79 }
80 }
81}
82
83template <bool T>
85 std::vector<QMFragment<BSE_Population> >& frags, const Orbitals& orbitals,
86 const Eigen::MatrixXd& dmat) const {
87
88 AOBasis basis = orbitals.getDftBasis();
89 AOOverlap overlap;
90 overlap.Fill(basis);
91 // Eigen::VectorXd nuccharges = CalcNucChargeperAtom(orbitals.QMAtoms());
92 // Eigen::MatrixXd dmatgs = orbitals.DensityMatrixGroundState();
93 Eigen::VectorXd charges = CalcElecChargeperAtom(dmat, overlap, basis);
94
95 // Index numofstates = orbitals.NumberofStates(type);
96 for (auto& frag : frags) {
97 frag.value().Initialize(1);
98 frag.value().Gs = frag.ExtractFromVector(charges);
99 }
100 // for (Index i_state = 0; i_state < numofstates; i_state++) {
101 // QMState state(type, i_state, false);
102 // std::array<Eigen::MatrixXd, 2> dmat_ex =
103 // orbitals.DensityMatrixExcitedState(state);
104 // Eigen::VectorXd atom_h = CalcElecChargeperAtom(dmat_ex[0], overlap,
105 // basis); Eigen::VectorXd atom_e = -CalcElecChargeperAtom(dmat_ex[1],
106 // overlap, basis);
107 // for (auto& frag : frags) {
108 // frag.value().E(i_state) = frag.ExtractFromVector(atom_e);
109 // frag.value().H(i_state) = frag.ExtractFromVector(atom_h);
110 // }
111}
112
113template <bool T>
115 const QMMolecule& mol) const {
116 Eigen::VectorXd result = Eigen::VectorXd::Zero(mol.size());
117 for (Index i = 0; i < mol.size(); i++) {
118 result(i) = double(mol[i].getNuccharge());
119 }
120 return result;
121}
122
123template <bool T>
125 const Eigen::MatrixXd& dmat, AOOverlap& overlap,
126 const AOBasis& basis) const {
127 Eigen::MatrixXd prodmat;
128 if (T) {
129 Eigen::MatrixXd Smsqrt = overlap.Sqrt();
130 prodmat = Smsqrt * dmat * Smsqrt;
131 } else {
132 prodmat = dmat * overlap.Matrix();
133 }
134 std::vector<Index> funcperatom = basis.getFuncPerAtom();
135 Index noofatoms = Index(funcperatom.size());
136 Eigen::VectorXd charges = Eigen::VectorXd::Zero(noofatoms);
137 Index start = 0;
138 for (Index i = 0; i < charges.size(); ++i) {
139 Index nofunc = funcperatom[i];
140 charges(i) = prodmat.diagonal().segment(start, nofunc).sum();
141 start += nofunc;
142 }
143 return charges;
144}
145
146template class Populationanalysis<true>;
147template class Populationanalysis<false>;
148
149} // namespace xtp
150} // namespace votca
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
const std::vector< Index > & getFuncPerAtom() const
Definition aobasis.h:72
void Fill(const AOBasis &aobasis) final
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:52
Eigen::MatrixXd Sqrt()
Definition aomatrix.cc:45
const std::string & getType() const
void push_back(const T &atom)
container for molecular orbitals
Definition orbitals.h:46
std::string getCalculationType() const
Definition orbitals.h:159
Eigen::MatrixXd DensityMatrixFull(const QMState &state) const
Definition orbitals.cc:140
const Eigen::MatrixXd & getInactiveDensity() const
Definition orbitals.h:424
Eigen::MatrixXd DensityMatrixGroundState() const
Definition orbitals.cc:183
const QMMolecule & QMAtoms() const
Definition orbitals.h:172
std::array< Eigen::MatrixXd, 2 > DensityMatrixExcitedState(const QMState &state) const
Definition orbitals.cc:302
Index NumberofStates(QMStateType type) const
Definition orbitals.h:135
const AOBasis & getDftBasis() const
Definition orbitals.h:208
void CalcChargeperFragment(std::vector< QMFragment< BSE_Population > > &frags, const Orbitals &orbitals, QMStateType type) const
StaticSegment CalcChargeperAtom(const Orbitals &orbitals, const QMState &state) const
void CalcChargeperFragmentTransition(std::vector< QMFragment< BSE_Population > > &frags, const Orbitals &orbitals, const Eigen::MatrixXd &dmat) const
Eigen::VectorXd CalcElecChargeperAtom(const Eigen::MatrixXd &dmat, AOOverlap &overlap, const AOBasis &basis) const
Eigen::VectorXd CalcNucChargeperAtom(const QMMolecule &mol) const
bool isExciton() const
Definition qmstate.h:71
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
bool isTransition() const
Definition qmstate.h:153
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
ClassicalSegment< StaticSite > StaticSegment
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26