8 #ifndef INCLUDE_TEMPLE_OPTIMIZATION_TRUST_REGION_NEWTON_H
9 #define INCLUDE_TEMPLE_OPTIMIZATION_TRUST_REGION_NEWTON_H
14 #include "boost/math/tools/roots.hpp"
15 #include <Eigen/Eigenvalues>
20 namespace Molassembler {
24 template<
typename Derived>
25 auto resolve(
const Eigen::MatrixBase<Derived>& matrix) ->
typename Eigen::MatrixBase<Derived>::Scalar {
26 assert(matrix.rows() == 1 && matrix.cols() == 1);
30 template<
typename FloatType>
31 struct TROFloatingPointTolerances {};
34 struct TROFloatingPointTolerances<double> {
35 static constexpr
unsigned rootBitAccuracy = 20;
39 struct TROFloatingPointTolerances<float> {
40 static constexpr
unsigned rootBitAccuracy = 10;
50 template<
typename FloatType =
double>
52 using MatrixType = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic>;
53 using VectorType = Eigen::Matrix<FloatType, Eigen::Dynamic, 1>;
66 typename UpdateFunction,
69 Eigen::Ref<VectorType> parameters,
70 UpdateFunction&&
function,
74 constexpr FloatType trustRadiusLimit = 4.0;
75 static_assert(trustRadiusLimit > 0,
"Trust radius is not positive!");
78 constexpr FloatType agreementContractionUpperBound = 0.25;
79 constexpr FloatType agreementExpansionLowerBound = 0.75;
82 constexpr FloatType agreementEta = 0.125;
84 0.0 <= agreementEta && agreementEta < 0.25,
85 "Eta not within bounds"
90 FloatType trustRadius = 1.0;
91 FloatType modelAgreement = 0;
92 unsigned iterations = 0;
94 step.initialize(parameters,
function, trustRadius);
100 Stl17::as_const(step.values.current),
101 Stl17::as_const(step.gradients.current)
105 modelAgreement = step.modelAgreement();
108 if(modelAgreement < agreementContractionUpperBound) {
110 trustRadius = agreementContractionUpperBound * trustRadius;
112 modelAgreement > agreementExpansionLowerBound
113 && std::fabs(step.direction.norm() - trustRadius) < 1e-5
116 trustRadius = std::min(2 * trustRadius, trustRadiusLimit);
120 if(modelAgreement > agreementEta) {
122 step.accept(
function, trustRadius);
125 step.prepare(
function, trustRadius);
132 std::move(step.gradients.current)
150 VectorType direction;
152 template<
typename UpdateFunction>
154 Eigen::Ref<VectorType> initialParameters,
155 UpdateFunction&&
function,
156 const FloatType trustRadius
158 parameters.current = std::move(initialParameters);
159 const unsigned P = parameters.current.size();
160 parameters.proposed.resize(P);
161 gradients.current.resize(P);
162 gradients.proposed.resize(P);
163 hessians.current.resize(P, P);
164 hessians.proposed.resize(P, P);
166 function(parameters.current, values.current, gradients.current, hessians.current);
167 prepare(
function, trustRadius);
184 const FloatType trustRadius,
188 using Vector2Type = Eigen::Matrix<FloatType, 2, 1>;
189 using Matrix2Type = Eigen::Matrix<FloatType, 2, 2>;
192 const MatrixType Binverse = B.inverse();
194 const Vector2Type gTilde {g.squaredNorm(), g.dot(Binverse * g)};
196 BTilde << g.dot(B * g), gTilde(0),
197 gTilde(0), gTilde(1);
199 const Vector2Type uStar = - BTilde.inverse() * gTilde;
202 BBar << gTilde(0), gTilde(1),
203 gTilde(1), (Binverse * g).squaredNorm();
206 const FloatType J = FloatType {0.5} * uStar.transpose() * BBar * uStar;
207 const FloatType trustRadiusSquare = trustRadius * trustRadius;
209 if(J <= trustRadiusSquare) {
210 return uStar(0) * g + uStar(1) * Binverse * g;
214 boost::uintmax_t iterations = 1000;
215 auto root_result = boost::math::tools::toms748_solve(
216 [&](
const FloatType lambda) -> FloatType {
217 const Vector2Type intermediate = (BTilde + lambda * BBar).inverse() * gTilde;
218 return FloatType {0.5} * intermediate.transpose() * BBar * intermediate - trustRadiusSquare;
220 std::nextafter(FloatType {0.0}, FloatType {1.0}),
222 boost::math::tools::eps_tolerance<FloatType> {
223 Detail::TROFloatingPointTolerances<FloatType>::rootBitAccuracy
228 if(iterations >= 1000) {
229 throw std::logic_error(
"Could not find lambda to satisfy equation");
232 const FloatType lambda = (root_result.first + root_result.second) / 2;
235 const Vector2Type uMin = -(BTilde + lambda * BBar).inverse() * gTilde;
236 return uMin(0) * g + uMin(1) * Binverse * g;
239 VectorType determineDirection(
const FloatType trustRadius)
const {
250 VectorType eigenvalues = hessians.current.template selfadjointView<Eigen::Lower>().eigenvalues();
255 if((eigenvalues.array() < FloatType {0.0}).any()) {
257 const FloatType mostNegativeEigenvalue = eigenvalues(0);
262 const FloatType alpha = FloatType {-1.5} * mostNegativeEigenvalue;
263 const unsigned P = parameters.current.size();
264 const Eigen::MatrixXd modifiedHessian = hessians.current + alpha * Eigen::MatrixXd::Identity(P, P);
279 const FloatType gradientNorm = gradients.current.norm();
280 const FloatType
tau = [&]() -> FloatType {
281 const FloatType cauchyIntermediate = gradients.current.transpose() * hessians.current * gradients.current;
282 if(cauchyIntermediate <= 0) {
287 std::pow(gradientNorm, 3) / (trustRadius * cauchyIntermediate),
292 return - (tau * trustRadius / gradientNorm) * gradients.current;
295 template<
typename UpdateFunction>
296 void prepare(UpdateFunction&&
function,
const FloatType trustRadius) {
297 direction = determineDirection(trustRadius);
298 parameters.proposed.noalias() = parameters.current + direction;
299 function(parameters.proposed, values.proposed, gradients.proposed, hessians.proposed);
302 FloatType modelAgreement()
const {
303 const FloatType predictedValue = (
305 + Detail::resolve(gradients.current.transpose() * direction)
306 + Detail::resolve(FloatType {0.5} * direction.transpose() * hessians.current * direction)
310 (values.current - values.proposed)
311 / (values.current - predictedValue)
315 template<
typename UpdateFunction>
316 void accept(UpdateFunction&&
function,
const FloatType trustRadius) {
317 parameters.propagate();
319 gradients.propagate();
320 hessians.propagate();
322 prepare(
function, trustRadius);
Sylvester's criterion for positive definiteness and positive semi-definiteness.
VectorType subspaceMinimize(const FloatType trustRadius, const VectorType &g, const MatrixType &B) const
Find the minimizer of the trust region subspace.
Definition: TrustRegion.h:183
unsigned iterations
Number of iterations.
Definition: TrustRegion.h:58
Trust region Newton-Raphson minimizer with 2D subspace minimization.
Definition: TrustRegion.h:51
Type returned from an optimization.
Definition: TrustRegion.h:56
bool positiveDefinite(const Eigen::MatrixBase< Derived > &matrix)
Determine whether a matrix is positive definite.
Definition: SylvestersCriterion.h:26
FloatType value
Final function value.
Definition: TrustRegion.h:60
Data types common to optimizers.
VectorType gradient
Final gradient.
Definition: TrustRegion.h:62
double tau(const std::vector< double > &angles)
Calculates the tau value for four and five angle symmetries.
Definition: TauCriteria.h:63
Algorithms available in the C++17 STL.
Stores all parameter, value, gradient and hessian data of two points.
Definition: TrustRegion.h:144