Molassembler  1.0.0
Molecule graph and conformer library
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
NewtonRaphson.h
Go to the documentation of this file.
1 
8 #ifndef INCLUDE_TEMPLE_NEWTON_RAPHSON_H
9 #define INCLUDE_TEMPLE_NEWTON_RAPHSON_H
10 
11 #include <Eigen/Core>
12 #include <Eigen/Eigenvalues>
13 #include <unsupported/Eigen/NumericalDiff>
14 #include "Molassembler/Temple/STL17.h"
15 #include <iostream>
16 
17 namespace Scine {
18 namespace Molassembler {
19 namespace Temple {
20 
30 template<typename FloatType = double>
31 struct NewtonRaphson {
32  using MatrixType = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic>;
33  using VectorType = Eigen::Matrix<FloatType, Eigen::Dynamic, 1>;
34 
38  unsigned iterations;
40  FloatType value;
42  VectorType gradient;
43  };
44 
45  template<
46  typename UpdateFunction,
47  typename Checker
48  > OptimizationReturnType minimize(
49  Eigen::Ref<VectorType> parameters,
50  UpdateFunction&& function,
51  Checker&& check
52  ) {
53  const unsigned P = parameters.size();
54  FloatType value;
55  VectorType gradients = Eigen::VectorXd::Zero(P);
56  MatrixType hessian (P, P);
57  hessian.setZero();
58 
59  VectorType numericalGradient(P);
60  {
61  constexpr FloatType h = 1e-4;
62  FloatType a, b;
63  VectorType diffParameters = parameters;
64  for(unsigned i = 0; i < P; ++i) {
65  diffParameters(i) += h;
66  function(diffParameters, b, gradients, hessian);
67  diffParameters(i) -= 2 * h;
68  function(diffParameters, a, gradients, hessian);
69  diffParameters(i) = parameters(i);
70  numericalGradient(i) = (b - a) / (2 * h);
71  }
72  }
73 
74  function(parameters, value, gradients, hessian);
75 
76  if(numericalGradient.norm() > 1e-4) {
77  if(!numericalGradient.isApprox(gradients, 1e-8)) {
78  std::cout << "MISMATCH: Numerical gradient is " << numericalGradient.transpose() << "\nfunction gradient is : " << gradients.transpose() << "\n";
79  } else {
80  std::cout << "match: norm is " << gradients.norm() << " and diff norm is " << (numericalGradient - gradients).norm() << "\n";
81  }
82  }
83 
84  MatrixType numericalHessian(P, P);
85  {
86  constexpr FloatType h = 1e-4;
87  FloatType a, b, c, d;
88  VectorType diffParameters = parameters;
89  for(unsigned i = 0; i < P; ++i) {
90  // Diagonal second derivative from forward and backward derivative
91  function(diffParameters, b, gradients, hessian);
92  diffParameters(i) += h;
93  function(diffParameters, c, gradients, hessian);
94  diffParameters(i) -= 2 * h;
95  function(diffParameters, a, gradients, hessian);
96  diffParameters(i) = parameters(i);
97  numericalHessian(i, i) = (a - 2 * b + c) / (h * h);
98 
99  // Cross terms from two central derivatives
100  for(unsigned j = i + 1; j < P; ++j) {
101  // (i) + h, (j) + h -> a
102  diffParameters(i) += h;
103  diffParameters(j) += h;
104  function(diffParameters, a, gradients, hessian);
105  // (i) + h, (j) - h - > b
106  diffParameters(j) -= 2 * h;
107  function(diffParameters, b, gradients, hessian);
108  // (i) - h, (j) - h -> d
109  diffParameters(i) -= 2 * h;
110  function(diffParameters, d, gradients, hessian);
111  // (i) - h, (j) + h -> c
112  diffParameters(j) += 2 * h;
113  function(diffParameters, c, gradients, hessian);
114 
115  numericalHessian(i, j) = (a - b - c + d) / (4 * h * h);
116  numericalHessian(j, i) = numericalHessian(i, j);
117  }
118  }
119  }
120 
121  function(parameters, value, gradients, hessian);
122 
123  if(!numericalHessian.isApprox(hessian), 1e-2) {
124  std::cout << "MISMATCH: Hessian is\n" << hessian << "\nNumerical hessian is\n" << numericalHessian << "\n";
125  } else {
126  std::cout << "Hessian matches\n";
127  }
128 
129  unsigned iteration = 0;
130  while(
131  check.shouldContinue(
132  iteration,
133  Stl17::as_const(value),
134  Stl17::as_const(gradients)
135  )
136  ) {
137  // Solve HΔx = g, then apply x_(n+1) = x_n - Δx
138  Eigen::JacobiSVD<MatrixType> svd(
139  hessian,
140  Eigen::ComputeFullV | Eigen::ComputeFullU
141  );
142  svd.setThreshold(svdThreshold);
143  VectorType steps = svd.solve(gradients);
144 
145  // Take a step
146  FloatType rms = std::sqrt(steps.squaredNorm() / P);
147  if (rms > trustRadius) {
148  steps *= trustRadius / rms;
149  }
150  parameters -= steps;
151 
152  function(parameters, value, gradients, hessian);
153  ++iteration;
154  }
155 
156  return {
157  iteration,
158  value,
159  std::move(gradients)
160  };
161  }
162 
164  FloatType svdThreshold = 1.0e-12;
166  FloatType trustRadius = 0.5;
167 };
168 
169 } // namespace Temple
170 } // namespace Molassembler
171 } // namespace Scine
172 
173 #endif
FloatType value
Final function value.
Definition: NewtonRaphson.h:40
unsigned iterations
Number of iterations.
Definition: NewtonRaphson.h:38
Type returned from an optimization.
Definition: NewtonRaphson.h:36
VectorType gradient
Final gradient.
Definition: NewtonRaphson.h:42
FloatType trustRadius
The maximum RMS of a taken step.
Definition: NewtonRaphson.h:166
FloatType svdThreshold
The SVD threshold for the decomposition of the Hessian.
Definition: NewtonRaphson.h:164
A very basic newton-raphson minimizer.
Definition: NewtonRaphson.h:31