8 #ifndef SPARROW_ONEELECTRONMATRIX_H
9 #define SPARROW_ONEELECTRONMATRIX_H
18 class MatrixWithDerivatives;
19 class AtomsOrbitalsIndexes;
25 class AtomicParameters;
26 class ElementParameters;
27 class TwoCenterIntegralContainer;
52 template<Utils::Derivative O>
53 void addDerivatives(Utils::AutomaticDifferentiation::DerivativeContainerType<O>& derivativeContainer,
59 const Eigen::MatrixXd& getMatrix()
const {
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,
72 const Eigen::MatrixXd& P;
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