Molassembler  1.0.0
Molecule graph and conformer library
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Common.h
Go to the documentation of this file.
1 
8 #ifndef INCLUDE_TEMPLE_OPTIMIZATION_COMMON_H
9 #define INCLUDE_TEMPLE_OPTIMIZATION_COMMON_H
10 
11 #include <Eigen/Core>
12 
13 namespace Scine {
14 namespace Molassembler {
15 namespace Temple {
16 
18 namespace Optimization {
19 
20 template<typename FloatType>
22  static_assert(
23  std::is_floating_point<FloatType>::value,
24  "This struct is not intended for anything else but floats"
25  );
26 
27  FloatType current;
28  FloatType proposed;
29 
31  FloatType absDelta() const {
32  return std::fabs(proposed - current);
33  }
34 
35  void propagate() { current = proposed; }
36 };
37 
38 template<typename EigenType>
40  EigenType current;
41  EigenType proposed;
42 
43  typename EigenType::Scalar deltaNorm() const {
44  return (proposed - current).norm();
45  }
46 
47  void propagate() {
48  /* Swapping instead of moving proposed into current is three moves instead
49  * of one, but avoids allocation of new memory (after moving from proposed,
50  * we need to resize proposed to have the same size). Assuming moves are
51  * cheap (trade ownership of allocated memory) and allocation could be
52  * costly, this might be cheaper.
53  */
54  std::swap(current, proposed);
55  }
56 };
57 
58 template<typename FloatType, typename VectorType, typename UpdateFunction>
60  UpdateFunction function;
61 
62  NegateFunction(UpdateFunction&& fn) : function(fn) {}
63 
64  void operator() (
65  const VectorType& parameters,
66  FloatType& value,
67  Eigen::Ref<VectorType> gradient
68  ) {
69  function(parameters, value, gradient);
70  value *= -1;
71  gradient *= -1;
72  }
73 };
74 
75 template<typename VectorType, typename UpdateFunction>
76 auto negateFunction(UpdateFunction&& fn) {
77  using FloatType = typename VectorType::Scalar;
78  return NegateFunction<FloatType, VectorType, UpdateFunction>(std::forward<UpdateFunction>(fn));
79 }
80 
81 template<typename Function, typename FloatType>
82 auto numericalGradient(
83  Function&& function,
84  const Eigen::Matrix<FloatType, Eigen::Dynamic, 1>& parameters
85 ) {
86  using VectorType = Eigen::Matrix<FloatType, Eigen::Dynamic, 1>;
87  const unsigned P = parameters.size();
88  VectorType numericalGradient(P);
89  VectorType diffParameters = parameters;
90 
91  // Central differences for all individual parameters
92  constexpr FloatType h = 1e-4;
93  FloatType a, b;
94  for(unsigned i = 0; i < P; ++i) {
95  diffParameters(i) += h;
96  b = function(diffParameters);
97  diffParameters(i) -= 2 * h;
98  a = function(diffParameters);
99  diffParameters(i) = parameters(i);
100  numericalGradient(i) = (b - a) / (2 * h);
101  }
102 
103  return numericalGradient;
104 }
105 
106 template<typename Function, typename FloatType>
107 auto numericalHessian(
108  Function&& function,
109  const Eigen::Matrix<FloatType, Eigen::Dynamic, 1>& parameters
110 ) {
111  using VectorType = Eigen::Matrix<FloatType, Eigen::Dynamic, 1>;
112  using MatrixType = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic>;
113  const unsigned P = parameters.size();
114  MatrixType numericalHessian(P, P);
115  VectorType diffParameters = parameters;
116 
117  // Form the hessian using finite differences
118  constexpr FloatType h = 1e-4;
119  FloatType a, b, c, d;
120  for(unsigned i = 0; i < P; ++i) {
121  // Diagonal second derivative from forward and backward derivative
122  b = function(diffParameters);
123  diffParameters(i) += h;
124  c = function(diffParameters);
125  diffParameters(i) -= 2 * h;
126  a = function(diffParameters);
127  diffParameters(i) = parameters(i);
128  numericalHessian(i, i) = (a - 2 * b + c) / (h * h);
129 
130  // Cross terms from two central derivatives
131  for(unsigned j = i + 1; j < P; ++j) {
132  // (i) + h, (j) + h -> a
133  diffParameters(i) += h;
134  diffParameters(j) += h;
135  a = function(diffParameters);
136  // (i) + h, (j) - h - > b
137  diffParameters(j) -= 2 * h;
138  b = function(diffParameters);
139  // (i) - h, (j) - h -> d
140  diffParameters(i) -= 2 * h;
141  d = function(diffParameters);
142  // (i) - h, (j) + h -> c
143  diffParameters(j) += 2 * h;
144  c = function(diffParameters);
145 
146  numericalHessian(i, j) = (a - b - c + d) / (4 * h * h);
147  numericalHessian(j, i) = numericalHessian(i, j);
148  }
149  }
150 
151  return numericalHessian;
152 }
153 
154 } // namespace Optimization
155 } // namespace Temple
156 } // namespace Molassembler
157 } // namespace Scine
158 
159 #endif
FloatType absDelta() const
Returns signed difference: proposed - current.
Definition: Common.h:31