Scine::Sparrow  5.0.0
Library for fast and agile quantum chemical calculations with semiempirical methods.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
TwoElectronMatrix.h
Go to the documentation of this file.
1 
8 #ifndef SPARROW_TWOELECTRONMATRIX_H
9 #define SPARROW_TWOELECTRONMATRIX_H
10 
12 #include <Utils/Typenames.h>
13 #include <Eigen/Core>
14 
15 namespace Scine {
16 
17 namespace Utils {
18 enum class ElementType : unsigned;
19 class DensityMatrix;
20 class AtomsOrbitalsIndexes;
21 enum class Derivative;
22 } // namespace Utils
23 
24 namespace Sparrow {
25 
26 namespace nddo {
27 class OneCenterIntegralContainer;
28 class OneCenterTwoElectronIntegrals;
29 class TwoCenterIntegralContainer;
30 class ElementParameters;
31 class AtomicParameters;
32 
33 namespace multipole {
34 class Global2c2eMatrix;
35 }
36 
43  public:
44  TwoElectronMatrix(const Utils::ElementTypeCollection& elements, const Utils::DensityMatrix& densityMatrix,
45  const OneCenterIntegralContainer& oneCIntegrals, const TwoCenterIntegralContainer& twoCIntegrals,
46  const ElementParameters& elementPar, const Utils::AtomsOrbitalsIndexes& aoIndexes);
47  void initialize();
48 
49  void calculate(bool spinPolarized);
50  void calculateBlocks();
51  void calculateSameAtomBlock(int startIndex, int nAOs, Utils::ElementType el, Eigen::MatrixXd& G,
52  Eigen::MatrixXd& GAlpha, Eigen::MatrixXd& GBeta);
53  void calculateDifferentAtomsBlock(int startA, int startB, int nAOsA, int nAOsB, const multipole::Global2c2eMatrix& m,
54  Eigen::MatrixXd& G, Eigen::MatrixXd& GAlpha, Eigen::MatrixXd& GBeta);
55  template<Utils::Derivative O>
56  void addDerivatives(Utils::AutomaticDifferentiation::DerivativeContainerType<O>& derivativeContainer) const;
57  const Eigen::MatrixXd& operator()() const {
58  return G_;
59  }
60  const Eigen::MatrixXd& getMatrix() const {
61  return G_;
62  }
63  const Eigen::MatrixXd& getAlpha() const {
64  return GAlpha_;
65  }
66  const Eigen::MatrixXd& getBeta() const {
67  return GBeta_;
68  }
69 
73  const OneCenterIntegralContainer& getOneCenterIntegrals() const;
77  const TwoCenterIntegralContainer& getTwoCenterIntegrals() const;
78 
79  private:
80  template<Utils::Derivative O>
81  void addDerivativesForBlock(Utils::AutomaticDifferentiation::DerivativeContainerType<O>& derivativeContainer, int a, int b,
82  int startA, int startB, int nAOsA, int nAOsB, const multipole::Global2c2eMatrix& m) const;
83 
84  bool spinPolarized_;
85  const Eigen::MatrixXd &P, &PAlpha_, &PBeta_;
86  const OneCenterIntegralContainer& oneCenterIntegrals;
87  const TwoCenterIntegralContainer& twoCenterIntegrals;
88  const ElementParameters& elementParameters;
89  const Utils::AtomsOrbitalsIndexes& aoIndexes_;
90 
91  Eigen::MatrixXd G_, GAlpha_, GBeta_;
92  const Utils::ElementTypeCollection& elementTypes_;
93  int nAOs_;
94  int nAtoms_;
95 };
96 
97 } // namespace nddo
98 
99 } // namespace Sparrow
100 } // namespace Scine
101 #endif // SPARROW_TWOELECTRONMATRIX_H
Definition: Global2c2eMatrix.h:31
Definition: ElementParameters.h:26
Definition: oneCenterIntegralContainer.h:32
Class to generate the two-electron matrix G for semi-empirical methods. This class is parallelized wi...
Definition: TwoElectronMatrix.h:42
ConceptualDftContainer calculate(double energy, const Eigen::VectorXd &atomicCharges, double energyPlus, const Eigen::VectorXd &atomicChargesPlus, double energyMinus, const Eigen::VectorXd &atomicChargesMinus)
This class contains smart pointers to two-center two-electron matrices for different atoms...
Definition: TwoCenterIntegralContainer.h:33