votca 2026-dev
Loading...
Searching...
No Matches
threecenter.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_THREECENTER_H
22#define VOTCA_XTP_THREECENTER_H
23
24// Local VOTCA includes
25#include "aobasis.h"
26#include "eigen.h"
27#include "symmetric_matrix.h"
28
35
36namespace votca {
37namespace xtp {
38
39// due to different requirements for the data format for DFT and GW we have two
40// different classes TCMatrix_gwbse and TCMatrix_dft which inherit from TCMatrix
41class TCMatrix {
42
43 public:
44 virtual ~TCMatrix() = default;
46 enum class SpinChannel { Alpha, Beta };
47
48 protected:
50 Eigen::MatrixXd inv_sqrt_;
51};
52
53class TCMatrix_dft final : public TCMatrix {
54 public:
55 void Fill(const AOBasis& auxbasis, const AOBasis& dftbasis);
56
57 Index size() const { return Index(matrix_.size()); }
58
60
61 const Symmetric_Matrix& operator[](Index i) const { return matrix_[i]; }
62
63 private:
64 std::vector<Symmetric_Matrix> matrix_;
65
66 void FillBlock(std::vector<Eigen::MatrixXd>& block, Index shellindex,
67 const AOBasis& dftbasis, const AOBasis& auxbasis);
68};
69
70class TCMatrix_gwbse final : public TCMatrix {
71 public:
72 // returns one level as a constant reference
73 const Eigen::MatrixXd& operator[](Index i) const { return matrix_[i]; }
74
75 // returns one level as a reference
76 Eigen::MatrixXd& operator[](Index i) { return matrix_[i]; }
77 // returns auxbasissize
78 Index auxsize() const { return auxbasissize_; }
79
80 Index get_mmin() const { return mmin_; }
81
82 Index get_mmax() const { return mmax_; }
83
84 Index get_nmin() const { return nmin_; }
85
86 Index get_nmax() const { return nmax_; }
87
88 Index msize() const { return mtotal_; }
89
90 Index nsize() const { return ntotal_; }
91
92 void Initialize(Index basissize, Index mmin, Index mmax, Index nmin,
93 Index nmax);
94
95 void Fill(const AOBasis& auxbasis, const AOBasis& dftbasis,
96 const Eigen::MatrixXd& dft_orbitals);
97 // Rebuilds ThreeCenterIntegrals, only works if the original basisobjects
98 // still exist
100
101 void MultiplyRightWithAuxMatrix(const Eigen::MatrixXd& matrix);
102
103 private:
104 // store vector of matrices
105 std::vector<Eigen::MatrixXd> matrix_;
106
107 // band summation indices
115
116 const AOBasis* auxbasis_ = nullptr;
117 const AOBasis* dftbasis_ = nullptr;
118 const Eigen::MatrixXd* dft_orbitals_ = nullptr;
119
120 void Fill3cMO(const AOBasis& auxbasis, const AOBasis& dftbasis,
121 const Eigen::MatrixXd& dft_orbitals);
122};
123
136
137} // namespace xtp
138} // namespace votca
139
140#endif // VOTCA_XTP_THREECENTER_H
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
void FillBlock(std::vector< Eigen::MatrixXd > &block, Index shellindex, const AOBasis &dftbasis, const AOBasis &auxbasis)
const Symmetric_Matrix & operator[](Index i) const
Definition threecenter.h:61
Symmetric_Matrix & operator[](Index i)
Definition threecenter.h:59
void Fill(const AOBasis &auxbasis, const AOBasis &dftbasis)
std::vector< Symmetric_Matrix > matrix_
Definition threecenter.h:64
void Fill(const AOBasis &auxbasis, const AOBasis &dftbasis, const Eigen::MatrixXd &dft_orbitals)
const AOBasis * auxbasis_
const Eigen::MatrixXd * dft_orbitals_
void Fill3cMO(const AOBasis &auxbasis, const AOBasis &dftbasis, const Eigen::MatrixXd &dft_orbitals)
void Initialize(Index basissize, Index mmin, Index mmax, Index nmin, Index nmax)
void MultiplyRightWithAuxMatrix(const Eigen::MatrixXd &matrix)
std::vector< Eigen::MatrixXd > matrix_
const AOBasis * dftbasis_
const Eigen::MatrixXd & operator[](Index i) const
Definition threecenter.h:73
Eigen::MatrixXd & operator[](Index i)
Definition threecenter.h:76
Eigen::MatrixXd inv_sqrt_
Definition threecenter.h:50
virtual ~TCMatrix()=default
Index Removedfunctions() const
Definition threecenter.h:45
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
TCMatrix_gwbse & operator[](TCMatrix::SpinChannel spin)
const TCMatrix_gwbse & operator[](TCMatrix::SpinChannel spin) const