Belle II Software  release-08-01-10
TrackSetEvaluatorHopfieldNNDEVModule.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/modules/trackSetEvaluatorVXD/TrackSetEvaluatorHopfieldNNDEVModule.h"
10 
11 #include "tracking/trackFindingVXD/trackSetEvaluator/HopfieldNetwork.h"
12 
13 using namespace Belle2;
14 
15 
16 REG_MODULE(TrackSetEvaluatorHopfieldNNDEV);
17 
19 {
20  setDescription("This module expects a container of SpacePointTrackCandidates and an OverlapNetwork\
21  and thenselects a subset of non-overlapping TCs determined using a neural network of Hopfield type.");
22 
24 
25  addParam("NameSpacePointTrackCands", m_PARAMtcArrayName, "Sets the name of expected StoreArray with SpacePointTrackCand in it.",
26  std::string(""));
27  addParam("NameOverlapNetworks", m_PARAMtcNetworkName, "Sets the name of expected StoreArray<OverlapNetwork>.",
28  std::string(""));
29  addParam("minActivityState", m_minActivityState, "Sets the minimal value of activity (Neuron Value) for acceptance.",
30  float(0.7));
31 }
32 
33 
35 {
37  //If no SpacePointTrackCands are available, later algorithms may crash.
39  return;
40  }
42 
43  //Prepare the information for the actual HNN.
44  std::vector<OverlapResolverNodeInfo> overlapResolverNodeInfos;
45  overlapResolverNodeInfos.reserve(m_spacePointTrackCands.getEntries());
46 
47  for (const SpacePointTrackCand& sPTC : m_spacePointTrackCands) {
48  overlapResolverNodeInfos.emplace_back(sPTC.getQualityIndicator(), sPTC.getArrayIndex(),
49  m_overlapNetworks[0]->getOverlapForTrackIndex(sPTC.getArrayIndex()), 1.0);
50  }
51 
52  //Performs the actual HNN.
53  //As the parameter is taken as reference, the values are changed and can be reused below.
54  HopfieldNetwork hopfieldNetwork;
55  unsigned maxIterations = 20;
56  if (hopfieldNetwork.doHopfield(overlapResolverNodeInfos, maxIterations) == maxIterations) {
57  B2INFO("Hopfield Network failed converge.");
59  }
60 
61  //Update tcs and kill those which were rejected by the Hopfield algorithm
62  unsigned int nSurvivors = 0;
63  for (const auto& overlapResolverNodeInfo : overlapResolverNodeInfos) {
64  if (overlapResolverNodeInfo.activityState > m_minActivityState) {
65  nSurvivors++;
66  continue;
67  }
68  m_spacePointTrackCands[overlapResolverNodeInfo.trackIndex]->removeRefereeStatus(SpacePointTrackCand::c_isActive);
69  }
70 
71  //Reporting stuff.
72  if (nSurvivors == 0) {
73  B2WARNING("Hopfield network - had no survivors!");
74  }
75 
76  m_nFinalTCs += nSurvivors;
77 }
78 
79 
81 {
82  if (m_eventCounter == 0) { m_eventCounter++; } // prevents division by zero
83  double invEvents = 1. / m_eventCounter;
84 
85  B2INFO("TrackSetEvaluatorHopfieldNNDEVModule-endRun: " <<
86  " nTCs per event: " << float(m_nTCsTotal)*invEvents <<
87  " nFinalTCs per event: " << float(m_nFinalTCs)*invEvents <<
88  "\n nTCs total: " << m_nTCsTotal <<
89  ", nFinalTCs total: " << m_nFinalTCs <<
90  ", number of times Hopfield did not succeed: " << m_nHopfieldFails);
91 
92  //After having evaluated the counters, we reset them to zero.
93  m_eventCounter = 0;
94  m_nTCsTotal = 0;
95  m_nFinalTCs = 0;
96  m_nHopfieldFails = 0;
97 }
Hopfield Algorithm with number based inputs.
unsigned short doHopfield(std::vector< OverlapResolverNodeInfo > &overlapResolverNodeInfos, unsigned short nIterations=20)
Performance of the actual algorithm.
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Storage for (VXD) SpacePoint-based track candidates.
@ c_isActive
bit 11: SPTC is active (i.e.
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
float m_minActivityState
************************************** Member variables *********************************************...
void event() override final
Applies the Hopfield neural network algorithm at given sets of SpacePointTrackCandidates.
std::string m_PARAMtcNetworkName
sets the name of the StoreObjPtr used for storing a TC network.
std::string m_PARAMtcArrayName
************************************** Module Parameters ********************************************...
unsigned int m_nFinalTCs
number of TCs found for final sets of non-overlapping TCs so far.
unsigned int m_nTCsTotal
total number of TCs evaluated so far.
StoreArray< OverlapNetwork > m_overlapNetworks
access to tcNetwork, which will be produced by this module.
StoreArray< SpacePointTrackCand > m_spacePointTrackCands
the storeArray for SpacePointTrackCands as member, is faster than recreating link for each event.
unsigned int m_nHopfieldFails
counts number of times when Hopfield was not able to clean overlaps.
REG_MODULE(arichBtest)
Register the Module.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
Abstract base class for different kinds of events.