votca 2024.2-dev
Loading...
Searching...
No Matches
ERIs.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_ERIS_H
22#define VOTCA_XTP_ERIS_H
23
24// Local VOTCA includes
25#include "threecenter.h"
26
27namespace votca {
28namespace xtp {
29
35class ERIs {
36
37 public:
38 void Initialize(const AOBasis& dftbasis, const AOBasis& auxbasis);
39 void Initialize_4c(const AOBasis& dftbasis);
40
41 Eigen::MatrixXd CalculateERIs_3c(const Eigen::MatrixXd& DMAT) const;
42
43 std::array<Eigen::MatrixXd, 2> CalculateERIs_EXX_3c(
44 const Eigen::MatrixXd& occMos, const Eigen::MatrixXd& DMAT) const {
45 std::array<Eigen::MatrixXd, 2> result;
46 result[0] = CalculateERIs_3c(DMAT);
47 if (occMos.rows() > 0 && occMos.cols() > 0) {
48 assert(occMos.rows() == DMAT.rows() && "occMos.rows()==DMAT.rows()");
49 result[1] = CalculateEXX_mos(occMos);
50 } else {
51 result[1] = CalculateEXX_dmat(DMAT);
52 }
53 return result;
54 }
55
56 Eigen::MatrixXd CalculateERIs_4c(const Eigen::MatrixXd& DMAT,
57 double error) const {
58 return Compute4c<false>(DMAT, error)[0];
59 }
60
61 std::array<Eigen::MatrixXd, 2> CalculateERIs_EXX_4c(
62 const Eigen::MatrixXd& DMAT, double error) const {
63 return Compute4c<true>(DMAT, error);
64 }
65
67
68 static double CalculateEnergy(const Eigen::MatrixXd& DMAT,
69 const Eigen::MatrixXd& matrix_operator) {
70 return matrix_operator.cwiseProduct(DMAT).sum();
71 }
72
73 private:
74 std::vector<libint2::Shell> basis_;
75 std::vector<Index> starts_;
76
77 std::vector<std::vector<Index>> shellpairs_;
78 std::vector<std::vector<libint2::ShellPair>> shellpairdata_;
81
82 Eigen::MatrixXd CalculateEXX_dmat(const Eigen::MatrixXd& DMAT) const;
83 Eigen::MatrixXd CalculateEXX_mos(const Eigen::MatrixXd& occMos) const;
84
85 std::vector<std::vector<libint2::ShellPair>> ComputeShellPairData(
86 const std::vector<libint2::Shell>& basis,
87 const std::vector<std::vector<Index>>& shellpairs) const;
88
89 Eigen::MatrixXd ComputeSchwarzShells(const AOBasis& dftbasis) const;
90 Eigen::MatrixXd ComputeShellBlockNorm(const Eigen::MatrixXd& dmat) const;
91
92 template <bool with_exchange>
93 std::array<Eigen::MatrixXd, 2> Compute4c(const Eigen::MatrixXd& dmat,
94 double error) const;
95
97
98 Eigen::MatrixXd schwarzscreen_; // Square matrix containing <ab|ab> for all
99 // shells
100}; // namespace xtp
101
102} // namespace xtp
103} // namespace votca
104
105#endif // VOTCA_XTP_ERIS_H
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Takes a density matrix and and an auxiliary basis set and calculates the electron repulsion integrals...
Definition ERIs.h:35
Eigen::MatrixXd CalculateERIs_3c(const Eigen::MatrixXd &DMAT) const
Definition ERIs.cc:80
Index Removedfunctions() const
Definition ERIs.h:66
std::vector< std::vector< libint2::ShellPair > > shellpairdata_
Definition ERIs.h:78
std::vector< libint2::Shell > basis_
Definition ERIs.h:74
Eigen::MatrixXd ComputeSchwarzShells(const AOBasis &dftbasis) const
std::vector< std::vector< Index > > shellpairs_
Definition ERIs.h:77
std::array< Eigen::MatrixXd, 2 > CalculateERIs_EXX_4c(const Eigen::MatrixXd &DMAT, double error) const
Definition ERIs.h:61
void Initialize_4c(const AOBasis &dftbasis)
Definition ERIs.cc:32
static double CalculateEnergy(const Eigen::MatrixXd &DMAT, const Eigen::MatrixXd &matrix_operator)
Definition ERIs.h:68
std::array< Eigen::MatrixXd, 2 > Compute4c(const Eigen::MatrixXd &dmat, double error) const
std::array< Eigen::MatrixXd, 2 > CalculateERIs_EXX_3c(const Eigen::MatrixXd &occMos, const Eigen::MatrixXd &DMAT) const
Definition ERIs.h:43
Index maxL_
Definition ERIs.h:80
std::vector< std::vector< libint2::ShellPair > > ComputeShellPairData(const std::vector< libint2::Shell > &basis, const std::vector< std::vector< Index > > &shellpairs) const
Definition ERIs.cc:46
void Initialize(const AOBasis &dftbasis, const AOBasis &auxbasis)
Definition ERIs.cc:27
Eigen::MatrixXd CalculateERIs_4c(const Eigen::MatrixXd &DMAT, double error) const
Definition ERIs.h:56
Eigen::MatrixXd schwarzscreen_
Definition ERIs.h:98
Eigen::MatrixXd ComputeShellBlockNorm(const Eigen::MatrixXd &dmat) const
Definition ERIs.cc:63
std::vector< Index > starts_
Definition ERIs.h:75
Eigen::MatrixXd CalculateEXX_dmat(const Eigen::MatrixXd &DMAT) const
Definition ERIs.cc:98
TCMatrix_dft threecenter_
Definition ERIs.h:96
Index maxnprim_
Definition ERIs.h:79
Eigen::MatrixXd CalculateEXX_mos(const Eigen::MatrixXd &occMos) const
Definition ERIs.cc:112
Index Removedfunctions() const
Definition threecenter.h:45
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26