Scine::Readuct  5.0.0
This is the SCINE module ReaDuct.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
ReactionPathCostCalculator.h
Go to the documentation of this file.
1 
8 #ifndef READUCT_ELEMENTARYSTEPOPTIMIZATION_COSTCALCULATORS_REACTIONPATHCOSTCALCULATOR_H
9 #define READUCT_ELEMENTARYSTEPOPTIMIZATION_COSTCALCULATORS_REACTIONPATHCOSTCALCULATOR_H
10 
11 #include <Eigen/Core>
12 #include <memory>
13 
14 namespace Scine {
15 namespace Utils {
16 namespace BSplines {
17 class BSpline;
18 }
19 } // namespace Utils
20 namespace Readuct {
21 namespace ElementaryStepOptimization {
22 struct EnergiesAndGradientsAlongSpline;
23 
24 namespace CostBasedOptimization {
25 
39  public:
40  std::unique_ptr<ReactionPathCostCalculator> clone() const;
41 
45  virtual ~ReactionPathCostCalculator() = default;
46 
52  bool energiesRequired() const;
53 
60  void calculateCost(const Utils::BSplines::BSpline& spline, const EnergiesAndGradientsAlongSpline& energyValues);
61 
67  double getCost() const;
68 
73  Eigen::MatrixXd getCostDerivatives() const;
74 
75  private:
76  virtual std::unique_ptr<ReactionPathCostCalculator> cloneImpl() const = 0;
77  virtual bool energiesRequiredImpl() const = 0;
78  virtual void calculateCostImpl(const Utils::BSplines::BSpline& spline,
79  const EnergiesAndGradientsAlongSpline& energyValues) = 0;
80  virtual double getCostImpl() const = 0;
81  virtual Eigen::MatrixXd getCostDerivativesImpl() const = 0;
82 };
83 
84 inline std::unique_ptr<ReactionPathCostCalculator> ReactionPathCostCalculator::clone() const {
85  return cloneImpl();
86 }
87 
89  return energiesRequiredImpl();
90 }
91 
93  const EnergiesAndGradientsAlongSpline& energyValues) {
94  calculateCostImpl(spline, energyValues);
95 }
96 
97 inline double ReactionPathCostCalculator::getCost() const {
98  return getCostImpl();
99 }
100 
101 inline Eigen::MatrixXd ReactionPathCostCalculator::getCostDerivatives() const {
102  return getCostDerivativesImpl();
103 }
104 
105 } // namespace CostBasedOptimization
106 
107 } // namespace ElementaryStepOptimization
108 
109 } // namespace Readuct
110 } // namespace Scine
111 #endif // ELEMENTARYSTEPOPTIMIZATION_COSTBASEDOPTIMIZATION_REACTIONPATHCOSTCALCULATOR_H
Interface for the cost calculation of reaction paths.
Definition: ReactionPathCostCalculator.h:38
double getCost() const
Get the cost associated with a given spline (after it has been evaluated).
Definition: ReactionPathCostCalculator.h:97
bool energiesRequired() const
Whether energies (and their gradients) are required in the cost calculation.
Definition: ReactionPathCostCalculator.h:88
Eigen::MatrixXd getCostDerivatives() const
Get the cost derivatives associated with a given spline.
Definition: ReactionPathCostCalculator.h:101
void calculateCost(const Utils::BSplines::BSpline &spline, const EnergiesAndGradientsAlongSpline &energyValues)
Evaluate the cost associated with a given spline.
Definition: ReactionPathCostCalculator.h:92