9#include <tracking/trackFindingVXD/trackSetEvaluator/HopfieldNetwork.h>
11#include <framework/logging/Logger.h>
12#include <framework/utilities/TRandomWrapper.h>
21 std::vector<OverlapResolverNodeInfo>& overlapResolverNodeInfos,
unsigned short nIterations)
29 if (overlapResolverNodeInfos.size() < 2) {
30 B2DEBUG(20,
"No reason to doHopfield with less than 2 nodes!");
34 const float compatibilityValue = (1.0 -
m_omega) /
static_cast<float>(overlapResolverNodeInfos.size() - 1);
36 const size_t overlapSize = overlapResolverNodeInfos.size();
39 Eigen::MatrixXd W(overlapSize, overlapSize);
41 W.fill(compatibilityValue);
44 for (
const auto& aTC : overlapResolverNodeInfos) {
45 for (
unsigned int overlapIndex : aTC.overlaps) {
46 W(aTC.trackIndex, overlapIndex) = -1.0;
52 Eigen::VectorXd x(overlapSize);
54 for (
unsigned int i = 0; i < overlapSize; i++) {
55 x(i) = gRandom->Uniform(1.0);
59 Eigen::VectorXd xOld(overlapSize);
62 std::vector<unsigned short> sequenceVector(overlapSize);
64 std::iota(sequenceVector.begin(), sequenceVector.end(), 0);
67 B2DEBUG(29,
"sequenceVector with length " << sequenceVector.size());
68 B2DEBUG(29,
"Entries are from begin to end:");
69 for (
auto&& entry : sequenceVector) {
70 B2DEBUG(29, std::to_string(entry) +
", ");
74 std::vector<float> cValues(nIterations);
79 unsigned iIterations = 0;
85 std::shuffle(sequenceVector.begin(), sequenceVector.end(),
TRandomWrapper());
89 for (
unsigned int i : sequenceVector) {
90 float aTempVal = W.row(i).dot(x);
91 float act = aTempVal +
m_omega * overlapResolverNodeInfos[i].qualityIndicator;
92 x(i) = 0.5 * (1. + tanh(act / T));
98 c = (x - xOld).cwiseAbs().maxCoeff();
99 B2DEBUG(21,
"c value is " << c <<
" at iteration " << iIterations);
100 cValues.at(iIterations) = c;
102 if (iIterations + 1 == nIterations) {
103 B2INFO(
"Hopfield reached maximum number of iterations without convergence. cValues are:");
104 for (
auto&& entry : cValues) {
105 B2INFO(std::to_string(entry));
113 for (
unsigned int i = 0; i < overlapSize; i++) {
114 overlapResolverNodeInfos[i].activityState = x(i);
float m_omega
tuning parameter of the hopfield network
float m_T
start temperature of annealing
float m_cmax
maximal change of weights between iterations
unsigned short doHopfield(std::vector< OverlapResolverNodeInfo > &overlapResolverNodeInfos, unsigned short nIterations=20)
Performance of the actual algorithm.
float m_Tmin
minimal temperature allowed
Abstract base class for different kinds of events.
Wrap TRandom to be useable as a uniform random number generator with STL algorithms like std::shuffle...