votca 2024-dev
Loading...
Searching...
No Matches
populationanalysis.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_POPULATIONANALYSIS_H
22#define VOTCA_XTP_POPULATIONANALYSIS_H
23
24// Local VOTCA includes
25#include "aomatrix.h"
26#include "orbitals.h"
27#include "qmfragment.h"
28
37namespace votca {
38namespace xtp {
39
40template <bool T>
42 public:
44 const QMState& state) const;
45
46 void CalcChargeperFragment(std::vector<QMFragment<BSE_Population> >& frags,
47 const Orbitals& orbitals, QMStateType type) const;
48
50 std::vector<QMFragment<BSE_Population> >& frags, const Orbitals& orbitals,
51 const Eigen::MatrixXd& dmat) const;
52
53 private:
54 Eigen::VectorXd CalcNucChargeperAtom(const QMMolecule& mol) const;
55
56 Eigen::VectorXd CalcElecChargeperAtom(const Eigen::MatrixXd& dmat,
57 AOOverlap& overlap,
58 const AOBasis& basis) const;
59};
60
63
64} // namespace xtp
65} // namespace votca
66
67#endif // VOTCA_XTP_POPULATIONANALYSIS_H
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
container for molecular orbitals
Definition orbitals.h:46
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
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
base class for all analysis tools
Definition basebead.h:33