Scine::Readuct  6.0.0
This is the SCINE module ReaDuct.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
RecurringProfileCalculator.h
Go to the documentation of this file.
1 
8 #ifndef READUCT_ELEMENTARYSTEPOPTIMIZATION_RECURRINGPROFILECALCULATOR_H
9 #define READUCT_ELEMENTARYSTEPOPTIMIZATION_RECURRINGPROFILECALCULATOR_H
10 
15 #include <vector>
16 
17 // TODO: Improve documentation
18 
19 namespace Scine {
20 namespace Utils {
21 class State;
22 } // namespace Utils
23 namespace Readuct {
24 
25 namespace ElementaryStepOptimization {
26 class ProfileEnergies;
27 
33 class SCINE_DLLEXPORT RecurringProfileCalculator {
34  public:
35  explicit RecurringProfileCalculator(Core::Calculator& calculator, int numberEquidistantPoints);
36 
37  int pointCount() const;
38  const double& deltaU() const;
39  const std::vector<double>& getCoordinates() const;
40  const std::vector<double>& getEnergies() const;
41  const PointSequence& pointSequence() const;
42  const EnergiesAndGradientsAlongSpline& valuesAlongSpline() const;
43  ProfileEnergies getProfileEnergies() const;
44 
45  void calculateEnergies(const Utils::BSplines::BSpline& spline);
46  void calculateEnergies(const Utils::BSplines::BSpline& spline, std::vector<double>& energies);
47 
48  void calculateEnergiesAndGradients(const Utils::BSplines::BSpline& spline);
49  void calculateEnergiesAndGradients(const Utils::BSplines::BSpline& spline, std::vector<double>& energies,
50  std::vector<Utils::GradientCollection>& gradients);
51 
52  void calculateUpToSecondDerivative(const Utils::BSplines::BSpline& spline, std::vector<double>& energies,
53  std::vector<Utils::GradientCollection>& gradients,
54  std::vector<Eigen::MatrixXd>& hessians);
55 
56  void calculateEnergyAndGradients(const Utils::BSplines::BSpline& spline, double u, double& energy,
57  Utils::GradientCollection& gradients);
58 
59  private:
60  void injectDensity(int index);
61  void saveDensity(int index);
62 
63  Core::Calculator& calculator_;
64  BSplineProfileCalculator profileCalculator_;
65 
66  double deltaU_;
67  EnergiesAndGradientsAlongSpline valuesAlongSpline_;
68  std::vector<std::shared_ptr<Core::State>> densities_;
69 };
70 
71 inline int RecurringProfileCalculator::pointCount() const {
72  return pointSequence().count();
73 }
74 
75 inline const PointSequence& RecurringProfileCalculator::pointSequence() const {
76  return valuesAlongSpline_.uValues;
77 }
78 
79 inline const EnergiesAndGradientsAlongSpline& RecurringProfileCalculator::valuesAlongSpline() const {
80  return valuesAlongSpline_;
81 }
82 
83 } // namespace ElementaryStepOptimization
84 
85 } // namespace Readuct
86 } // namespace Scine
87 #endif // ELEMENTARYSTEPOPTIMIZATION_RECURRINGPROFILECALCULATOR_H