votca 2024.2-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
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
47 protected:
49 Eigen::MatrixXd inv_sqrt_;
50};
51
52class TCMatrix_dft final : public TCMatrix {
53 public:
54 void Fill(const AOBasis& auxbasis, const AOBasis& dftbasis);
55
56 Index size() const { return Index(matrix_.size()); }
57
59
60 const Symmetric_Matrix& operator[](Index i) const { return matrix_[i]; }
61
62 private:
63 std::vector<Symmetric_Matrix> matrix_;
64
65 void FillBlock(std::vector<Eigen::MatrixXd>& block, Index shellindex,
66 const AOBasis& dftbasis, const AOBasis& auxbasis);
67};
68
69class TCMatrix_gwbse final : public TCMatrix {
70 public:
71 // returns one level as a constant reference
72 const Eigen::MatrixXd& operator[](Index i) const { return matrix_[i]; }
73
74 // returns one level as a reference
75 Eigen::MatrixXd& operator[](Index i) { return matrix_[i]; }
76 // returns auxbasissize
77 Index auxsize() const { return auxbasissize_; }
78
79 Index get_mmin() const { return mmin_; }
80
81 Index get_mmax() const { return mmax_; }
82
83 Index get_nmin() const { return nmin_; }
84
85 Index get_nmax() const { return nmax_; }
86
87 Index msize() const { return mtotal_; }
88
89 Index nsize() const { return ntotal_; }
90
91 void Initialize(Index basissize, Index mmin, Index mmax, Index nmin,
92 Index nmax);
93
94 void Fill(const AOBasis& auxbasis, const AOBasis& dftbasis,
95 const Eigen::MatrixXd& dft_orbitals);
96 // Rebuilds ThreeCenterIntegrals, only works if the original basisobjects
97 // still exist
99
100 void MultiplyRightWithAuxMatrix(const Eigen::MatrixXd& matrix);
101
102 private:
103 // store vector of matrices
104 std::vector<Eigen::MatrixXd> matrix_;
105
106 // band summation indices
114
115 const AOBasis* auxbasis_ = nullptr;
116 const AOBasis* dftbasis_ = nullptr;
117 const Eigen::MatrixXd* dft_orbitals_ = nullptr;
118
119 void Fill3cMO(const AOBasis& auxbasis, const AOBasis& dftbasis,
120 const Eigen::MatrixXd& dft_orbitals);
121};
122
123} // namespace xtp
124} // namespace votca
125
126#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:60
Symmetric_Matrix & operator[](Index i)
Definition threecenter.h:58
void Fill(const AOBasis &auxbasis, const AOBasis &dftbasis)
std::vector< Symmetric_Matrix > matrix_
Definition threecenter.h:63
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:72
Eigen::MatrixXd & operator[](Index i)
Definition threecenter.h:75
Eigen::MatrixXd inv_sqrt_
Definition threecenter.h:49
virtual ~TCMatrix()=default
Index Removedfunctions() const
Definition threecenter.h:45
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26