Scine::Readuct  6.0.0
This is the SCINE module ReaDuct.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
CostCombiner.h
Go to the documentation of this file.
1 
8 #ifndef READUCT_ELEMENTARYSTEPOPTIMIZATION_COSTCALCULATORS_COSTCOMBINER_H
9 #define READUCT_ELEMENTARYSTEPOPTIMIZATION_COSTCALCULATORS_COSTCOMBINER_H
10 
12 #include <iostream>
13 
14 // TODO: Improve documentation
15 
16 namespace Scine {
17 namespace Readuct {
18 
19 namespace ElementaryStepOptimization {
20 
21 namespace CostBasedOptimization {
22 
28 template<typename C1, typename C2>
29 class SCINE_DLLEXPORT CostCombiner : public ReactionPathCostCalculator {
30  static_assert(std::is_base_of<ReactionPathCostCalculator, C1>::value,
31  "C1 must be a descendant of ReactionPathCostCalculator");
32  static_assert(std::is_base_of<ReactionPathCostCalculator, C2>::value,
33  "C2 must be a descendant of ReactionPathCostCalculator");
34 
35  public:
38  assert(0.0 <= f && f <= 1.0);
39  contributionFromFirstCalculator_ = f;
40  }
41 
42  double getFirstCalculatorContribution() const {
43  return contributionFromFirstCalculator_;
44  }
45 
46  const ReactionPathCostCalculator& firstCalculator() const {
47  return calculator1_;
48  }
49  ReactionPathCostCalculator& firstCalculator() {
50  return calculator1_;
51  }
52 
53  const ReactionPathCostCalculator& secondCalculator() const {
54  return calculator2_;
55  }
56  ReactionPathCostCalculator& secondCalculator() {
57  return calculator2_;
58  }
59 
60  static constexpr double defaultContributionFromFirstCalculator = 0.5;
61 
62  private:
63  std::unique_ptr<ReactionPathCostCalculator> cloneImpl() const override {
64  return std::make_unique<CostCombiner<C1, C2>>(*this);
65  }
66 
67  bool energiesRequiredImpl() const override {
68  return calculator1_.energiesRequired() || calculator2_.energiesRequired();
69  }
70 
71  void calculateCostImpl(const Utils::BSplines::BSpline& spline, const EnergiesAndGradientsAlongSpline& energyValues) override {
72  calculator1_.calculateCost(spline, energyValues);
73  calculator2_.calculateCost(spline, energyValues);
74  }
75 
76  double getCostImpl() const override {
77  double f1 = contributionFromFirstCalculator_;
78  double f2 = 1.0 - f1;
79 
80  return f1 * calculator1_.getCost() + f2 * calculator2_.getCost();
81  }
82 
83  Eigen::MatrixXd getCostDerivativesImpl() const override {
84  double f1 = contributionFromFirstCalculator_;
85  double f2 = 1.0 - f1;
86 
87  return f1 * calculator1_.getCostDerivatives() + f2 * calculator2_.getCostDerivatives();
88  }
89 
90  double contributionFromFirstCalculator_ = defaultContributionFromFirstCalculator;
91  C1 calculator1_{};
92  C2 calculator2_{};
93 };
94 
95 } // namespace CostBasedOptimization
96 
97 } // namespace ElementaryStepOptimization
98 
99 } // namespace Readuct
100 } // namespace Scine
101 #endif // READUCT_ELEMENTARYSTEPOPTIMIZATION_COSTCALCULATORS_COSTCOMBINER_H
Interface for the cost calculation of reaction paths.
Definition: ReactionPathCostCalculator.h:39
void setFirstCalculatorContribution(double f)
Definition: CostCombiner.h:37
bool energiesRequired() const
Whether energies (and their gradients) are required in the cost calculation.
Definition: ReactionPathCostCalculator.h:89