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
OneElectronMatrix.h
Go to the documentation of this file.
1 
8 #ifndef SPARROW_ONEELECTRONMATRIX_H
9 #define SPARROW_ONEELECTRONMATRIX_H
10 
12 #include <Utils/Typenames.h>
13 #include <Eigen/Core>
14 
15 namespace Scine {
16 
17 namespace Utils {
18 class MatrixWithDerivatives;
19 class AtomsOrbitalsIndexes;
20 } // namespace Utils
21 
22 namespace Sparrow {
23 
24 namespace nddo {
25 class AtomicParameters;
26 class ElementParameters;
27 class TwoCenterIntegralContainer;
28 
33  public:
36  const Eigen::MatrixXd& densityMatrix, const TwoCenterIntegralContainer& twoCIntegrals,
37  const ElementParameters& elementPar, const Utils::AtomsOrbitalsIndexes& aoIndexes);
39  void initialize();
45  void calculateSameAtomBlock(int a, int startIndex, int nAOs);
49  void calculateDifferentAtomsBlock(int startRow, int startCol, const AtomicParameters& pA, const AtomicParameters& pB,
52  template<Utils::Derivative O>
53  void addDerivatives(Utils::AutomaticDifferentiation::DerivativeContainerType<O>& derivativeContainer,
54  const Utils::MatrixWithDerivatives& S) const;
56  const Eigen::MatrixXd& operator()() const {
57  return H_;
58  }
59  const Eigen::MatrixXd& getMatrix() const {
60  return H_;
61  }
62 
63  private:
64  template<Utils::Derivative O>
65  void addDerivativesContribution1(Utils::AutomaticDifferentiation::DerivativeContainerType<O>& derivativeContainer,
66  int a, int startIndex, int nAOs) const;
67  template<Utils::Derivative O>
68  void addDerivativesContribution2(Utils::AutomaticDifferentiation::DerivativeContainerType<O>& derivativeContainer,
69  int a, int b, int indexA, int indexB, int nAOsA, int nAOsB,
70  const Utils::MatrixWithDerivatives& S) const;
71 
72  const Eigen::MatrixXd& P;
73  const TwoCenterIntegralContainer& twoCenterIntegrals;
74  const ElementParameters& elementParameters;
75  const Utils::AtomsOrbitalsIndexes& aoIndexes_;
76 
77  int nAOs_ = 0;
78  int nAtoms_ = 0;
79  Eigen::MatrixXd H_;
80  const Utils::ElementTypeCollection& elementTypes_;
81  const Utils::PositionCollection& positions_;
82 };
83 
84 } // namespace nddo
85 
86 } // namespace Sparrow
87 } // namespace Scine
88 #endif // SPARROW_ONEELECTRONMATRIX_H
void calculateDifferentAtomsBlocks(const Utils::MatrixWithDerivatives &S)
Calculates all the blocks on different atom pairs.
Definition: OneElectronMatrix.cpp:111
Definition: ElementParameters.h:26
void initialize()
Initializes one-electron matrix H and the used data structures.
Definition: OneElectronMatrix.cpp:34
void calculateDifferentAtomsBlock(int startRow, int startCol, const AtomicParameters &pA, const AtomicParameters &pB, const Utils::MatrixWithDerivatives &S)
Calculates a block between a specific atom pair.
Definition: OneElectronMatrix.cpp:125
This class generates the one-electron matrix H for semi-empirical methods.
Definition: OneElectronMatrix.h:32
void calculate(const Utils::MatrixWithDerivatives &S)
Fills the one-electron matrix H by reading and calculating the required integrals.
Definition: OneElectronMatrix.cpp:43
const Eigen::MatrixXd & operator()() const
Getter for the one-electron matrix H.
Definition: OneElectronMatrix.h:56
void calculateSameAtomBlock(int a, int startIndex, int nAOs)
Calculates a specific block on an atom.
Definition: OneElectronMatrix.cpp:61
void calculateSameAtomBlocks()
Calculates all the blocks on the same atoms.
Definition: OneElectronMatrix.cpp:52
Definition: AtomicParameters.h:29
This class contains smart pointers to two-center two-electron matrices for different atoms...
Definition: TwoCenterIntegralContainer.h:33
void addDerivatives(Utils::AutomaticDifferentiation::DerivativeContainerType< O > &derivativeContainer, const Utils::MatrixWithDerivatives &S) const
Calculates the derivative contribution up to the order.
OneElectronMatrix(const Utils::ElementTypeCollection &elements, const Utils::PositionCollection &positions, const Eigen::MatrixXd &densityMatrix, const TwoCenterIntegralContainer &twoCIntegrals, const ElementParameters &elementPar, const Utils::AtomsOrbitalsIndexes &aoIndexes)
Constructor.
Definition: OneElectronMatrix.cpp:23