Belle II Software  release-06-01-15
HopfieldNetwork.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <tracking/trackFindingVXD/trackSetEvaluator/HopfieldNetwork.h>
10 
11 #include <framework/logging/Logger.h>
12 
13 #include <Eigen/Dense>
14 
15 #include <numeric>
16 
17 using namespace Belle2;
18 
20  std::vector<OverlapResolverNodeInfo>& overlapResolverNodeInfos, unsigned short nIterations)
21 {
22  //Start value for neurons if they are compatible.
23  //Each compatible connection activates a node with this value.
24  //As the sum of all the activations shall be less than one, we divide the
25  //activiation by the total number of Nodes.
26  //Incompatible Nodes get later minus one, which counteracts all activation,
27  //if the incompatible Node is active.
28  if (overlapResolverNodeInfos.size() < 2) {
29  B2DEBUG(20, "No reason to doHopfield with less than 2 nodes!");
30  return 0;
31  }
32 
33  const float compatibilityValue = (1.0 - m_omega) / static_cast<float>(overlapResolverNodeInfos.size() - 1);
34 
35  const size_t overlapSize = overlapResolverNodeInfos.size();
36 
37  //Weight matrix; knows compatibility between each possible pair of Nodes
38  Eigen::MatrixXd W(overlapSize, overlapSize);
39  //A): Set all elements to compatible:
40  W.fill(compatibilityValue);
41 
42  //B): Inform weight matrix elements of incompatible neurons:
43  for (const auto& aTC : overlapResolverNodeInfos) {
44  for (unsigned int overlapIndex : aTC.overlaps) {
45  W(aTC.trackIndex, overlapIndex) = -1.0;
46  }
47  }
48 
49 
50  // Neuron values
51  Eigen::VectorXd x(overlapSize);
52  // randomize neuron values for first iteration:
53  for (unsigned int i = 0; i < overlapSize; i++) {
54  x(i) = gRandom->Uniform(1.0); // WARNING: original does Un(0;0.1) not Un(0;1)!
55  }
56 
57  //Store for results from last round:
58  Eigen::VectorXd xOld(overlapSize);
59 
60  //Order of execution for neuron values:
61  std::vector<unsigned short> sequenceVector(overlapSize);
62  //iota fills the vector with 0, 1, 2, ... , (size-1)
63  std::iota(sequenceVector.begin(), sequenceVector.end(), 0);
64 
65  //The following block will be evaluated to empty, if LOG_NO_B2DEBUG is defined:
66  B2DEBUG(100, "sequenceVector with length " << sequenceVector.size());
67  B2DEBUG(100, "Entries are from begin to end:");
68  for (auto && entry : sequenceVector) {
69  B2DEBUG(100, std::to_string(entry) + ", ");
70  }
71 
72  //Store all values of c for protocolling:
73  std::vector<float> cValues(nIterations);
74  //Store for maximum change of weights between iterations.
75  float c = 1.0;
76 
77  //Iterate until change in weights is small:
78  unsigned iIterations = 0;
79 
80  float T = m_T;
81 
82 
83  while (c > m_cmax) {
84  std::shuffle(sequenceVector.begin(), sequenceVector.end(), TRandomWrapper());
85 
86  xOld = x;
87 
88  for (unsigned int i : sequenceVector) {
89  float aTempVal = W.row(i).dot(x);
90  float act = aTempVal + m_omega * overlapResolverNodeInfos[i].qualityIndicator;
91  x(i) = 0.5 * (1. + tanh(act / T));
92  }
93 
94  T = 0.5 * (T + m_Tmin);
95 
96  //Determine maximum change in weights:
97  c = (x - xOld).cwiseAbs().maxCoeff();
98  B2DEBUG(10, "c value is " << c << " at iteration " << iIterations);
99  cValues.at(iIterations) = c;
100 
101  if (iIterations + 1 == nIterations) {
102  B2INFO("Hopfield reached maximum number of iterations without convergence. cValues are:");
103  for (auto && entry : cValues) {
104  B2INFO(std::to_string(entry));
105  }
106  break;
107  }
108  iIterations++;
109  }
110 
111  //Copy Node values into the activity state of the OverlapResolverNodeInfo objects:
112  for (unsigned int i = 0; i < overlapSize; i++) {
113  overlapResolverNodeInfos[i].activityState = x(i);
114  }
115 
116  return iIterations;
117 }
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 std algorithms like std::shuffle...