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
FockMatrix.h
Go to the documentation of this file.
1 
8 #ifndef SPARROW_NDDO_FOCKMATRIX_H
9 #define SPARROW_NDDO_FOCKMATRIX_H
10 
11 #include "OneElectronMatrix.h"
12 #include "TwoElectronMatrix.h"
16 #include <memory>
17 
18 namespace Scine {
19 
20 namespace Utils {
21 class AtomsOrbitalsIndexes;
22 class DensityMatrix;
23 class OverlapCalculator;
24 class ElectronicEnergyCalculator;
25 class AdditiveElectronicContribution;
26 } // namespace Utils
27 
28 namespace Sparrow {
29 
30 namespace nddo {
31 
33  public:
34  FockMatrix(const Utils::ElementTypeCollection& elements, const Utils::PositionCollection& positions,
35  const Utils::DensityMatrix& densityMatrix, const OneCenterIntegralContainer& oneCIntegrals,
36  const ElementParameters& elementPar, const Utils::AtomsOrbitalsIndexes& aoIndexes,
37  const Utils::OverlapCalculator& overlapCalculator, const bool& unrestrictedCalculationRunning);
38  ~FockMatrix() final;
39 
40  void initialize() override;
41  void calculateDensityIndependentPart(Utils::DerivativeOrder order) override;
42  void calculateDensityDependentPart(Utils::DerivativeOrder order) override;
43  void finalize(Utils::DerivativeOrder order) override;
44  Utils::SpinAdaptedMatrix getMatrix() const override;
45  double calculateElectronicEnergy() const override;
46  void addDerivatives(Utils::AutomaticDifferentiation::DerivativeContainerType<Utils::Derivative::First>& derivatives) const override;
47  void addDerivatives(Utils::AutomaticDifferentiation::DerivativeContainerType<Utils::Derivative::SecondAtomic>& derivatives) const override;
48  void addDerivatives(Utils::AutomaticDifferentiation::DerivativeContainerType<Utils::Derivative::SecondFull>& derivatives) const override;
49 
50  const OneElectronMatrix& getOneElectronMatrix() const;
51  const TwoElectronMatrix& getTwoElectronMatrix() const;
52  const std::vector<std::shared_ptr<Utils::AdditiveElectronicContribution>>& getDensityDependentContributions() const;
53  const std::vector<std::shared_ptr<Utils::AdditiveElectronicContribution>>& getDensityIndependentContributions() const;
54 
59  void addDensityDependentElectronicContribution(std::shared_ptr<Utils::AdditiveElectronicContribution> contribution) final;
64  void addDensityIndependentElectronicContribution(std::shared_ptr<Utils::AdditiveElectronicContribution> contribution) final;
65  void clearElectronicContributions();
66  void eraseElectronicContribution(std::shared_ptr<Utils::AdditiveElectronicContribution> contribution);
67 
68  private:
69  template<Utils::Derivative O>
70  void addDerivativesImpl(Utils::AutomaticDifferentiation::DerivativeContainerType<O>& derivatives) const;
71 
72  TwoCenterIntegralContainer twoCenterIntegrals_;
75  const Utils::OverlapCalculator& overlapCalculator_;
76  const bool& unrestrictedCalculationRunning_;
77  std::unique_ptr<Utils::ElectronicEnergyCalculator> electronicEnergyCalculator_;
78  std::vector<std::shared_ptr<Utils::AdditiveElectronicContribution>> densityDependentContributions_,
79  densityIndependentContributions_;
80 };
81 
82 } // namespace nddo
83 
84 } // namespace Sparrow
85 } // namespace Scine
86 #endif // SPARROW_NDDO_FOCKMATRIX_H
Definition: ElementParameters.h:26
Definition: FockMatrix.h:32
void addDensityDependentElectronicContribution(std::shared_ptr< Utils::AdditiveElectronicContribution > contribution) final
Definition: FockMatrix.cpp:156
This class generates the one-electron matrix H for semi-empirical methods.
Definition: OneElectronMatrix.h:32
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
void addDensityIndependentElectronicContribution(std::shared_ptr< Utils::AdditiveElectronicContribution > contribution) final
Definition: FockMatrix.cpp:152
This class contains smart pointers to two-center two-electron matrices for different atoms...
Definition: TwoCenterIntegralContainer.h:33