8 #ifndef INCLUDE_TEMPLE_NEWTON_RAPHSON_H
9 #define INCLUDE_TEMPLE_NEWTON_RAPHSON_H
12 #include <Eigen/Eigenvalues>
13 #include <unsupported/Eigen/NumericalDiff>
17 namespace Molassembler {
29 template<
typename FloatType =
double>
31 using MatrixType = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic>;
32 using VectorType = Eigen::Matrix<FloatType, Eigen::Dynamic, 1>;
45 typename UpdateFunction,
48 Eigen::Ref<VectorType> parameters,
49 UpdateFunction&&
function,
52 const unsigned P = parameters.size();
54 VectorType gradients = Eigen::VectorXd::Zero(P);
55 MatrixType hessian (P, P);
58 VectorType numericalGradient(P);
60 constexpr FloatType h = 1e-4;
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);
74 function(parameters, value, gradients, hessian);
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";
80 std::cout <<
"match: norm is " << gradients.norm() <<
" and diff norm is " << (numericalGradient - gradients).norm() <<
"\n";
84 MatrixType numericalHessian(P, P);
86 constexpr FloatType h = 1e-4;
91 VectorType diffParameters = parameters;
92 for(
unsigned i = 0; i < P; ++i) {
94 function(diffParameters, b, gradients, hessian);
95 diffParameters(i) += h;
96 function(diffParameters, c, gradients, hessian);
97 diffParameters(i) -= 2 * h;
98 function(diffParameters, a, gradients, hessian);
99 diffParameters(i) = parameters(i);
100 numericalHessian(i, i) = (a - 2 * b + c) / (h * h);
103 for(
unsigned j = i + 1; j < P; ++j) {
105 diffParameters(i) += h;
106 diffParameters(j) += h;
107 function(diffParameters, a, gradients, hessian);
109 diffParameters(j) -= 2 * h;
110 function(diffParameters, b, gradients, hessian);
112 diffParameters(i) -= 2 * h;
113 function(diffParameters, d, gradients, hessian);
115 diffParameters(j) += 2 * h;
116 function(diffParameters, c, gradients, hessian);
118 numericalHessian(i, j) = (a - b - c + d) / (4 * h * h);
119 numericalHessian(j, i) = numericalHessian(i, j);
124 function(parameters, value, gradients, hessian);
126 if(!numericalHessian.isApprox(hessian), 1e-2) {
127 std::cout <<
"MISMATCH: Hessian is\n" << hessian <<
"\nNumerical hessian is\n" << numericalHessian <<
"\n";
129 std::cout <<
"Hessian matches\n";
132 unsigned iteration = 0;
134 check.shouldContinue(
136 std::as_const(value),
137 std::as_const(gradients)
141 Eigen::JacobiSVD<MatrixType> svd(
143 Eigen::ComputeFullV | Eigen::ComputeFullU
146 VectorType steps = svd.solve(gradients);
149 FloatType rms = std::sqrt(steps.squaredNorm() / P);
155 function(parameters, value, gradients, hessian);
FloatType value
Final function value.
Definition: NewtonRaphson.h:39
unsigned iterations
Number of iterations.
Definition: NewtonRaphson.h:37
Type returned from an optimization.
Definition: NewtonRaphson.h:35
VectorType gradient
Final gradient.
Definition: NewtonRaphson.h:41
FloatType trustRadius
The maximum RMS of a taken step.
Definition: NewtonRaphson.h:169
FloatType svdThreshold
The SVD threshold for the decomposition of the Hessian.
Definition: NewtonRaphson.h:167
A very basic newton-raphson minimizer.
Definition: NewtonRaphson.h:30