votca 2024.2-dev
Loading...
Searching...
No Matches
aomatrix.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 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_AOMATRIX_H
22#define VOTCA_XTP_AOMATRIX_H
23
24// Local VOTCA includes
25#include "aobasis.h"
26
27namespace votca {
28namespace xtp {
29
30class AOMatrix {
31 public:
32 virtual void Fill(const AOBasis& aobasis) = 0;
33 virtual Index Dimension() = 0;
34};
35
36// derived class for kinetic energy
37class AOKinetic : public AOMatrix {
38 public:
39 void Fill(const AOBasis& aobasis) final;
40 Index Dimension() final { return aomatrix_.rows(); }
41 const Eigen::MatrixXd& Matrix() const { return aomatrix_; }
42
43 private:
44 Eigen::MatrixXd aomatrix_;
45};
46
47// derived class for atomic orbital overlap
48class AOOverlap : public AOMatrix {
49 public:
50 void Fill(const AOBasis& aobasis) final;
51 Index Dimension() final { return aomatrix_.rows(); }
52 const Eigen::MatrixXd& Matrix() const { return aomatrix_; }
53
54 Eigen::MatrixXd singleShellOverlap(const AOShell& shell) const;
56 double SmallestEigenValue() const { return smallestEigenvalue; }
57
58 Eigen::MatrixXd Pseudo_InvSqrt(double etol);
59 Eigen::MatrixXd Sqrt();
60
61 private:
64 Eigen::MatrixXd aomatrix_;
65};
66
67// derived class for atomic orbital Coulomb interaction
68class AOCoulomb : public AOMatrix {
69 public:
70 void Fill(const AOBasis& aobasis) final;
71 Index Dimension() final { return aomatrix_.rows(); }
72 const Eigen::MatrixXd& Matrix() const { return aomatrix_; }
73
74 Eigen::MatrixXd Pseudo_InvSqrt_GWBSE(const AOOverlap& auxoverlap,
75 double etol);
76 Eigen::MatrixXd Pseudo_InvSqrt(double etol);
78
79 private:
80 void computeCoulombIntegrals(const AOBasis& aobasis);
82 Eigen::MatrixXd aomatrix_;
83};
84
85/* derived class for atomic orbital electrical dipole matrices, required for
86 * electrical transition dipoles
87 */
88class AODipole : public AOMatrix {
89 public:
90 void Fill(const AOBasis& aobasis) final;
91 Index Dimension() final { return aomatrix_[0].rows(); }
92 const std::array<Eigen::MatrixXd, 3>& Matrix() const { return aomatrix_; }
93
94 void setCenter(const Eigen::Vector3d& r) {
95 for (Index i = 0; i < 3; i++) {
96 r_[i] = r[i];
97 }
98 } // definition of a center around which the moment should be calculated
99
100 private:
101 std::array<Eigen::MatrixXd, 3> aomatrix_;
102 std::array<libint2::Shell::real_t, 3> r_ = {0, 0, 0};
103};
104
105} // namespace xtp
106} // namespace votca
107
108#endif // VOTCA_XTP_AOMATRIX_H
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Index Dimension() final
Definition aomatrix.h:71
void Fill(const AOBasis &aobasis) final
Eigen::MatrixXd aomatrix_
Definition aomatrix.h:82
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:72
void computeCoulombIntegrals(const AOBasis &aobasis)
Eigen::MatrixXd Pseudo_InvSqrt(double etol)
Definition aomatrix.cc:88
Index Removedfunctions() const
Definition aomatrix.h:77
Eigen::MatrixXd Pseudo_InvSqrt_GWBSE(const AOOverlap &auxoverlap, double etol)
Definition aomatrix.cc:53
Index Dimension() final
Definition aomatrix.h:91
void setCenter(const Eigen::Vector3d &r)
Definition aomatrix.h:94
std::array< Eigen::MatrixXd, 3 > aomatrix_
Definition aomatrix.h:101
std::array< libint2::Shell::real_t, 3 > r_
Definition aomatrix.h:102
const std::array< Eigen::MatrixXd, 3 > & Matrix() const
Definition aomatrix.h:92
void Fill(const AOBasis &aobasis) final
void Fill(const AOBasis &aobasis) final
Index Dimension() final
Definition aomatrix.h:40
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:41
Eigen::MatrixXd aomatrix_
Definition aomatrix.h:44
virtual void Fill(const AOBasis &aobasis)=0
virtual Index Dimension()=0
void Fill(const AOBasis &aobasis) final
Eigen::MatrixXd Pseudo_InvSqrt(double etol)
Definition aomatrix.cc:28
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:52
double SmallestEigenValue() const
Definition aomatrix.h:56
Index Removedfunctions() const
Definition aomatrix.h:55
Eigen::MatrixXd singleShellOverlap(const AOShell &shell) const
Index Dimension() final
Definition aomatrix.h:51
double smallestEigenvalue
Definition aomatrix.h:63
Eigen::MatrixXd Sqrt()
Definition aomatrix.cc:45
Eigen::MatrixXd aomatrix_
Definition aomatrix.h:64
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26