Belle II Software development
HopfieldNetwork Class Reference

Hopfield Algorithm with number based inputs. More...

#include <HopfieldNetwork.h>

Public Member Functions

 HopfieldNetwork (float omega=0.5, float T=3.1, float Tmin=0.1, float cmax=0.01)
 Constructor taking parameters for the algorithm.
 
unsigned short doHopfield (std::vector< OverlapResolverNodeInfo > &overlapResolverNodeInfos, unsigned short nIterations=20)
 Performance of the actual algorithm.
 

Private Attributes

float m_omega
 tuning parameter of the hopfield network
 
float m_T
 start temperature of annealing
 
float m_Tmin
 minimal temperature allowed
 
float m_cmax
 maximal change of weights between iterations
 

Detailed Description

Hopfield Algorithm with number based inputs.

This class encapsulates the actual algorithm with pure number inputs (instead of objects of TrackCandidates or whatever). The reason for this is mainly better testability and potentially as well reusability. The development of this class started as a copy from Jakob's implementation for the VXDTF2.

Relevant resources: [1] R. Frühwirth, "Selection of optimal subsets of tracks with a feedback neural network", Comput. Phys. Commun., vol. 78, pp. 23–28, 1993.

Definition at line 33 of file HopfieldNetwork.h.

Constructor & Destructor Documentation

◆ HopfieldNetwork()

HopfieldNetwork ( float  omega = 0.5,
float  T = 3.1,
float  Tmin = 0.1,
float  cmax = 0.01 
)
inline

Constructor taking parameters for the algorithm.

Parameters
omegaShould be between 0 and 1; small values lead to large number of nodes, large ones to large sums of quality indicators.
TTemperature for annealing.
TminMinimal reached temperature in the annealing scheme.
cmaxMaximum change of weights between iterations, so we accept the network as converged.

Definition at line 42 of file HopfieldNetwork.h.

42 :
43 m_omega(omega), m_T(T), m_Tmin(Tmin), m_cmax(cmax)
44 {}
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
float m_Tmin
minimal temperature allowed

Member Function Documentation

◆ doHopfield()

unsigned short doHopfield ( std::vector< OverlapResolverNodeInfo > &  overlapResolverNodeInfos,
unsigned short  nIterations = 20 
)

Performance of the actual algorithm.

HINT FOR SPEED OPTIMIZATION: A lot of time is spend for checking the Logger, if you don't have LOG_NO_B2DEBUG defined. If you have done that, a lot of time is taken by the tanh function and drawing random numbers.

Currently this algorithm can only be used once for each instance, as the algorithm parameter variables are changed during the performance.

See also
OverlapResolverNodeInfo

Definition at line 20 of file HopfieldNetwork.cc.

22{
23 //Start value for neurons if they are compatible.
24 //Each compatible connection activates a node with this value.
25 //As the sum of all the activations shall be less than one, we divide the
26 //activation by the total number of Nodes.
27 //Incompatible Nodes get later minus one, which counteracts all activation,
28 //if the incompatible Node is active.
29 if (overlapResolverNodeInfos.size() < 2) {
30 B2DEBUG(20, "No reason to doHopfield with less than 2 nodes!");
31 return 0;
32 }
33
34 const float compatibilityValue = (1.0 - m_omega) / static_cast<float>(overlapResolverNodeInfos.size() - 1);
35
36 const size_t overlapSize = overlapResolverNodeInfos.size();
37
38 //Weight matrix; knows compatibility between each possible pair of Nodes
39 Eigen::MatrixXd W(overlapSize, overlapSize);
40 //A): Set all elements to compatible:
41 W.fill(compatibilityValue);
42
43 //B): Inform weight matrix elements of incompatible neurons:
44 for (const auto& aTC : overlapResolverNodeInfos) {
45 for (unsigned int overlapIndex : aTC.overlaps) {
46 W(aTC.trackIndex, overlapIndex) = -1.0;
47 }
48 }
49
50
51 // Neuron values
52 Eigen::VectorXd x(overlapSize);
53 // randomize neuron values for first iteration:
54 for (unsigned int i = 0; i < overlapSize; i++) {
55 x(i) = gRandom->Uniform(1.0); // WARNING: original does Un(0;0.1) not Un(0;1)!
56 }
57
58 //Store for results from last round:
59 Eigen::VectorXd xOld(overlapSize);
60
61 //Order of execution for neuron values:
62 std::vector<unsigned short> sequenceVector(overlapSize);
63 //iota fills the vector with 0, 1, 2, ... , (size-1)
64 std::iota(sequenceVector.begin(), sequenceVector.end(), 0);
65
66 //The following block will be evaluated to empty, if LOG_NO_B2DEBUG is defined:
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) + ", ");
71 }
72
73 //Store all values of c for protocolling:
74 std::vector<float> cValues(nIterations);
75 //Store for maximum change of weights between iterations.
76 float c = 1.0;
77
78 //Iterate until change in weights is small:
79 unsigned iIterations = 0;
80
81 float T = m_T;
82
83
84 while (c > m_cmax) {
85 std::shuffle(sequenceVector.begin(), sequenceVector.end(), TRandomWrapper());
86
87 xOld = x;
88
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));
93 }
94
95 T = 0.5 * (T + m_Tmin);
96
97 //Determine maximum change in weights:
98 c = (x - xOld).cwiseAbs().maxCoeff();
99 B2DEBUG(21, "c value is " << c << " at iteration " << iIterations);
100 cValues.at(iIterations) = c;
101
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));
106 }
107 break;
108 }
109 iIterations++;
110 }
111
112 //Copy Node values into the activity state of the OverlapResolverNodeInfo objects:
113 for (unsigned int i = 0; i < overlapSize; i++) {
114 overlapResolverNodeInfos[i].activityState = x(i);
115 }
116
117 return iIterations;
118}
Wrap TRandom to be useable as a uniform random number generator with STL algorithms like std::shuffle...

Member Data Documentation

◆ m_cmax

float m_cmax
private

maximal change of weights between iterations

Definition at line 63 of file HopfieldNetwork.h.

◆ m_omega

float m_omega
private

tuning parameter of the hopfield network

Definition at line 60 of file HopfieldNetwork.h.

◆ m_T

float m_T
private

start temperature of annealing

Definition at line 61 of file HopfieldNetwork.h.

◆ m_Tmin

float m_Tmin
private

minimal temperature allowed

Definition at line 62 of file HopfieldNetwork.h.


The documentation for this class was generated from the following files: