8 #ifndef INCLUDE_TEMPLE_OPTIMIZATION_COMMON_H
9 #define INCLUDE_TEMPLE_OPTIMIZATION_COMMON_H
14 namespace Molassembler {
18 namespace Optimization {
20 template<
typename FloatType>
23 std::is_floating_point<FloatType>::value,
24 "This struct is not intended for anything else but floats"
32 return std::fabs(proposed - current);
35 void propagate() { current = proposed; }
38 template<
typename EigenType>
43 typename EigenType::Scalar deltaNorm()
const {
44 return (proposed - current).norm();
54 std::swap(current, proposed);
58 template<
typename FloatType,
typename VectorType,
typename UpdateFunction>
60 UpdateFunction
function;
65 const VectorType& parameters,
67 Eigen::Ref<VectorType> gradient
69 function(parameters, value, gradient);
75 template<
typename VectorType,
typename UpdateFunction>
76 auto negateFunction(UpdateFunction&& fn) {
77 using FloatType =
typename VectorType::Scalar;
81 template<
typename Function,
typename FloatType>
82 auto numericalGradient(
84 const Eigen::Matrix<FloatType, Eigen::Dynamic, 1>& parameters
86 using VectorType = Eigen::Matrix<FloatType, Eigen::Dynamic, 1>;
87 const unsigned P = parameters.size();
88 VectorType numericalGradient(P);
89 VectorType diffParameters = parameters;
92 constexpr FloatType h = 1e-4;
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);
103 return numericalGradient;
106 template<
typename Function,
typename FloatType>
107 auto numericalHessian(
109 const Eigen::Matrix<FloatType, Eigen::Dynamic, 1>& parameters
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;
118 constexpr FloatType h = 1e-4;
119 FloatType a, b, c, d;
120 for(
unsigned i = 0; i < P; ++i) {
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);
131 for(
unsigned j = i + 1; j < P; ++j) {
133 diffParameters(i) += h;
134 diffParameters(j) += h;
135 a =
function(diffParameters);
137 diffParameters(j) -= 2 * h;
138 b =
function(diffParameters);
140 diffParameters(i) -= 2 * h;
141 d =
function(diffParameters);
143 diffParameters(j) += 2 * h;
144 c =
function(diffParameters);
146 numericalHessian(i, j) = (a - b - c + d) / (4 * h * h);
147 numericalHessian(j, i) = numericalHessian(i, j);
151 return numericalHessian;
FloatType absDelta() const
Returns signed difference: proposed - current.
Definition: Common.h:31