Belle II Software  release-05-02-19
CDCTriggerMCMatcherModule.cc
1 #include "trg/cdc/modules/mcmatcher/CDCTriggerMCMatcherModule.h"
2 
3 #include <framework/gearbox/Const.h>
4 
5 #include <Eigen/Dense>
6 
7 namespace {
8  //small anonymous helper construct making converting a pair of iterators usable
9  //with range based for
10  template<class Iter>
11  struct iter_pair_range : std::pair<Iter, Iter> {
12  explicit iter_pair_range(std::pair<Iter, Iter> const& x) : std::pair<Iter, Iter>(x) {}
13  Iter begin() const {return this->first;}
14  Iter end() const {return this->second;}
15  };
16 
17  template<class Iter>
18  inline iter_pair_range<Iter> as_range(std::pair<Iter, Iter> const& x)
19  {
20  return iter_pair_range<Iter>(x);
21  }
22 
23  typedef int HitId;
24  typedef int TrackId;
25  typedef float Purity;
26  typedef float Efficiency;
27 }
28 
29 
30 
31 using namespace Belle2;
32 using namespace std;
33 
34 //this line registers the module with the framework and actually makes it available
35 //in steering files or the the module list (basf2 -m).
36 REG_MODULE(CDCTriggerMCMatcher)
37 
39 {
40  setDescription("A module to match CDCTriggerTracks to MCParticles.\n"
41  "Creates an array MCTracks of trackable MCParticles "
42  "and makes relations from MCTracks to CDCTriggerTracks "
43  "and vice-versa.");
44 
45  addParam("MCParticleCollectionName", m_MCParticleCollectionName,
46  "Name of the MCParticle StoreArray to be matched.",
47  string("MCParticles"));
48  addParam("TrgTrackCollectionName", m_TrgTrackCollectionName,
49  "Name of the CDCTriggerTrack StoreArray to be matched.",
50  string("TRGCDC2DFinderTracks"));
51  addParam("MCTrackableCollectionName", m_MCTrackableCollectionName,
52  "Name of a new StoreArray holding MCParticles considered as trackable.",
53  string("MCTracks"));
54  addParam("hitCollectionName", m_hitCollectionName,
55  "Name of the StoreArray of CDCTriggerSegmentHits used for the matching.",
56  string(""));
57  addParam("minAxial", m_minAxial,
58  "Minimum number of axial hits in different super layers.",
59  0);
60  addParam("minStereo", m_minStereo,
61  "Minimum number of stereo hits in different super layers.",
62  0);
63  addParam("axialOnly", m_axialOnly,
64  "Switch to ignore stereo hits (= 2D matching).",
65  false);
66  addParam("minPurity", m_minPurity,
67  "Minimum purity for reconstructed tracks.",
68  0.1);
69  addParam("minEfficiency", m_minEfficiency,
70  "Minimum efficiency for MC tracks.",
71  0.1);
72  addParam("onlyPrimaries", m_onlyPrimaries,
73  "If true, MCTracks are only made from primary particles.",
74  false);
75  addParam("relateClonesAndMerged", m_relateClonesAndMerged,
76  "If true, create relations for clones and merged tracks "
77  "(will get negative weight).",
78  true);
79 }
80 
81 
82 void
84 {
85  m_segmentHits.isRequired(m_hitCollectionName);
86  m_prTracks.isRequired(m_TrgTrackCollectionName);
87  m_mcParticles.isRequired(m_MCParticleCollectionName);
88  m_mcTracks.registerInDataStore(m_MCTrackableCollectionName);
89 
90  m_mcParticles.requireRelationTo(m_segmentHits);
91  m_prTracks.requireRelationTo(m_segmentHits);
92  m_mcParticles.registerRelationTo(m_mcTracks);
93  m_mcParticles.registerRelationTo(m_prTracks);
94  m_mcTracks.registerRelationTo(m_segmentHits);
95  m_mcTracks.registerRelationTo(m_prTracks);
96  m_prTracks.registerRelationTo(m_mcParticles);
97  m_prTracks.registerRelationTo(m_mcTracks);
98 }
99 
100 
101 void
103 {
104  // get all trackable particles
105  for (int imc = 0; imc < m_mcParticles.getEntries(); ++imc) {
106  MCParticle* mcParticle = m_mcParticles[imc];
107 
108  // minimum requirement: seen in CDC, unless monopoles
109  if (!mcParticle->hasSeenInDetector(Const::CDC) || (mcParticle->getCharge() == 0 && abs(mcParticle->getPDG()) != 99666))
110  continue;
111 
112  // reject secondaries
113  if (m_onlyPrimaries && !mcParticle->hasStatus(MCParticle::c_PrimaryParticle))
114  continue;
115 
116  // count super layers with related hits
117  int nAxial = 0;
118  int nStereo = 0;
119  vector<bool> SLcounted;
120  SLcounted.assign(9, false);
121 
123  mcParticle->getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
124  for (unsigned ihit = 0; ihit < relHits.size(); ++ihit) {
125  unsigned iSL = relHits[ihit]->getISuperLayer();
126  if (SLcounted[iSL]) continue;
127  if (iSL % 2) ++nStereo;
128  else ++nAxial;
129  SLcounted[iSL] = true;
130  }
131 
132  if (nAxial < m_minAxial || nStereo < m_minStereo) continue;
133 
134  // copy particle to list of trackable particles
135  MCParticle* mcTrack = m_mcTracks.appendNew(m_mcTracks.getPtr(), *mcParticle);
136  mcParticle->addRelationTo(mcTrack);
137  for (unsigned ihit = 0; ihit < relHits.size(); ++ihit) {
138  mcTrack->addRelationTo(relHits[ihit]);
139  }
140  }
141 
142  // derive relations between trigger tracks and MC tracks,
143  // following the basic ideas of the MCMatcherTracksModule (tracking)
144 
145  B2DEBUG(100, "########## start MC matching ############");
146 
147  int nMCTracks = m_mcTracks.getEntries();
148  int nPRTracks = m_prTracks.getEntries();
149 
150  B2DEBUG(100, "Number of pattern recognition tracks is " << nPRTracks);
151  B2DEBUG(100, "Number of Monte-Carlo tracks is " << nMCTracks);
152 
153  if (not nMCTracks or not nPRTracks) {
154  // Either no pattern recognition tracks
155  // or no Monte Carlo tracks are present in this event
156  // Cannot perform matching.
157  return;
158  }
159 
160  //### Build a hit_id to track map for easier lookup later ###
161  std::map<HitId, TrackId> mcTrackId_by_hitId;
162  {
163  std::map<HitId, TrackId>::iterator itMCInsertHit = mcTrackId_by_hitId.end();
164  TrackId mcTrackId = -1;
165  for (const MCParticle& mcTrack : m_mcTracks) {
166  ++mcTrackId;
168  mcTrack.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
169  for (unsigned iHit = 0; iHit < relHits.size(); ++iHit) {
170  const HitId hitId = relHits[iHit]->getArrayIndex();
171  itMCInsertHit = mcTrackId_by_hitId.insert(itMCInsertHit,
172  make_pair(hitId, mcTrackId));
173  B2DEBUG(250, "hitId " << hitId << " in SL " << relHits[iHit]->getISuperLayer()
174  << ", mcTrackId " << mcTrackId);
175  }
176  }
177  }
178 
179  //### Build a hit_id to track map for easier lookup later ###
180  std::multimap<HitId, TrackId> prTrackId_by_hitId;
181  {
182  std::multimap<HitId, TrackId>::iterator itPRInsertHit = prTrackId_by_hitId.end();
183  TrackId prTrackId = -1;
184  for (const CDCTriggerTrack& prTrack : m_prTracks) {
185  ++prTrackId;
187  prTrack.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
188  for (unsigned int iHit = 0; iHit < relHits.size(); ++iHit) {
189  const HitId hitId = relHits[iHit]->getArrayIndex();
190  itPRInsertHit = prTrackId_by_hitId.insert(itPRInsertHit,
191  make_pair(hitId, prTrackId));
192  B2DEBUG(250, "hitId " << hitId << " in SL " << relHits[iHit]->getISuperLayer()
193  << ", prTrackId " << prTrackId);
194  }
195  }
196  }
197 
198  //### Build the confusion matrix ###
199 
200  // Reserve enough space for the confusion matrix
201  // The last column is meant for hits not assigned to a mcTrack (aka background hits)
202  Eigen::MatrixXi confusionMatrix = Eigen::MatrixXi::Zero(nPRTracks, nMCTracks + 1);
203 
204  // Count total number of hits for each Monte-Carlo track separately
205  // to avoid double counting (in case pattern recognition tracks share hits)
206  // and to count hits not associated to any pattern recognition track.
207  Eigen::RowVectorXi totalHits_by_mcTrackId = Eigen::RowVectorXi::Zero(nMCTracks + 1);
208 
209  // Column index for the hits not assigned to any MCTrackCand
210  const int mcBkgId = nMCTracks;
211 
212  // examine every hit to which mcTrack and prTrack it belongs.
213  // if the hit is not part of any mcTrack put the hit in the background column.
214  for (HitId hitId = 0; hitId < m_segmentHits.getEntries(); ++hitId) {
215  // skip stereo hits
216  if (m_axialOnly && m_segmentHits[hitId]->getISuperLayer() % 2) continue;
217 
218  // First search the unique mcTrackId for the hit.
219  // If the hit is not assigned to any mcTrack the Id is set to the background column.
220  auto it_mcTrackId = mcTrackId_by_hitId.find(hitId);
221  TrackId mcTrackId =
222  (it_mcTrackId == mcTrackId_by_hitId.end()) ? mcBkgId : it_mcTrackId->second;
223 
224  // Assign the hits to the total vector.
225  totalHits_by_mcTrackId(mcTrackId) += 1;
226 
227  // Seek all prTrackCands
228  auto range_prTrackIds = prTrackId_by_hitId.equal_range(hitId);
229 
230  // count for every prTrack that has this hit
231  for (const pair<HitId, TrackId>& hitId_and_prTrackId :
232  as_range(range_prTrackIds)) {
233  TrackId prTrackId = hitId_and_prTrackId.second;
234  B2DEBUG(200, " prTrackId : " << prTrackId << "; mcTrackId : " << mcTrackId);
235  confusionMatrix(prTrackId, mcTrackId) += 1;
236  } //end for prTrackId
237  } //end for hitId
238 
239  Eigen::VectorXi totalHits_by_prTrackId = confusionMatrix.rowwise().sum();
240 
241  B2DEBUG(200, "Confusion matrix of the event : " << endl << confusionMatrix);
242 
243  B2DEBUG(200, "totalHits_by_prTrackId : " << endl << totalHits_by_prTrackId);
244  B2DEBUG(200, "totalHits_by_mcTrackId : " << endl << totalHits_by_mcTrackId);
245 
246  // ### Building the patter recognition track to highest purity Monte-Carlo track relation ###
247  vector<pair<TrackId, Purity>> purestMCTrackId_by_prTrackId(nPRTracks);
248 
249  for (TrackId prTrackId = 0; prTrackId < nPRTracks; ++prTrackId) {
250  Eigen::RowVectorXi prTrackRow = confusionMatrix.row(prTrackId);
251  Eigen::RowVectorXi::Index purestMCTrackId;
252 
253  //Also sets the index of the highest entry in the row vector
254  int highestHits = prTrackRow.maxCoeff(&purestMCTrackId);
255  int totalHits = totalHits_by_prTrackId(prTrackId);
256 
257  Purity highestPurity = Purity(highestHits) / totalHits;
258 
259  purestMCTrackId_by_prTrackId[prTrackId] =
260  pair<TrackId, Purity>(purestMCTrackId, highestPurity);
261  }
262 
263  // Log the pattern recognition track to highest purity Monte-Carlo track
264  // relation to debug output
265  {
266  TrackId prTrackId = -1;
267  B2DEBUG(200, "PRTrack to highest purity MCTrack relation");
268  for (const pair< TrackId, Purity>& purestMCTrackId :
269  purestMCTrackId_by_prTrackId) {
270  ++prTrackId;
271  const Purity& purity = purestMCTrackId.second;
272  const TrackId& mcTrackId = purestMCTrackId.first;
273  B2DEBUG(200, "prTrackId : " << prTrackId << " -> mcTrackId : " << mcTrackId
274  << " with purity " << purity);
275  }
276  }
277 
278  // ### Building the Monte-Carlo track to highest efficiency pattern recognition track relation ###
279  vector<pair<TrackId, Efficiency>> mostEfficientPRTrackId_by_mcTrackId(nMCTracks);
280 
281  for (TrackId mcTrackId = 0; mcTrackId < nMCTracks; ++mcTrackId) {
282  Eigen::VectorXi mcTrackCol = confusionMatrix.col(mcTrackId);
283  Eigen::VectorXi::Index highestHitsPRTrackId;
284 
285  //Also sets the index of the highest entry in the column vector
286  int highestHits = mcTrackCol.maxCoeff(&highestHitsPRTrackId);
287  int totalHits = totalHits_by_mcTrackId(mcTrackId);
288 
289  Efficiency highestEfficiency = Efficiency(highestHits) / totalHits;
290 
291  mostEfficientPRTrackId_by_mcTrackId[mcTrackId] =
292  pair<TrackId, Efficiency> (highestHitsPRTrackId, highestEfficiency);
293  }
294 
295  // Log the Monte-Carlo track to highest efficiency pattern recognition track
296  // relation to debug output
297  {
298  TrackId mcTrackId = -1;
299  B2DEBUG(200, "MCTrack to highest efficiency PRTrack relation");
300  for (const pair<TrackId, Efficiency>& mostEfficientPRTrackId :
301  mostEfficientPRTrackId_by_mcTrackId) {
302  ++mcTrackId;
303  const Efficiency& efficiency = mostEfficientPRTrackId.second;
304  const TrackId& prTrackId = mostEfficientPRTrackId.first;
305  B2DEBUG(200, "mcTrackId : " << mcTrackId << " -> prTrackId : " << prTrackId
306  << " with efficiency " << efficiency);
307  }
308  }
309 
310  // ### Classify the pattern recognition tracks ###
311  for (TrackId prTrackId = 0; prTrackId < nPRTracks; ++prTrackId) {
312  CDCTriggerTrack* prTrack = m_prTracks[prTrackId];
313 
314  const pair<TrackId, Purity>& purestMCTrackId = purestMCTrackId_by_prTrackId[prTrackId];
315  const TrackId& mcTrackId = purestMCTrackId.first;
316  const Purity& purity = purestMCTrackId.second;
317 
318  if (!(purity >= m_minPurity)) {
319  // GHOST
320  B2DEBUG(100, "Classified PRTrack " << prTrackId << " as ghost because of too low purity.");
321  } else if (mcTrackId == mcBkgId) {
322  // BACKGROUND
323  B2DEBUG(100, "Classified PRTrack " << prTrackId << " as background.");
324  } else {
325  // check whether the highest purity Monte-Carlo track has in turn
326  // the highest efficiency pattern recognition track equal to this track.
327  MCParticle* mcTrack = m_mcTracks[mcTrackId];
328 
329  const pair<TrackId, Efficiency>& mostEfficientPRTrackId =
330  mostEfficientPRTrackId_by_mcTrackId[mcTrackId];
331 
332  const TrackId& prTrackIdCompare = mostEfficientPRTrackId.first;
333  const Efficiency& efficiency = mostEfficientPRTrackId.second;
334 
335  if (prTrackId != prTrackIdCompare) {
336  if (efficiency >= m_minEfficiency) {
337  // CLONE
338  if (m_relateClonesAndMerged) {
339  // Add the mc matching relation with negative weight
340  prTrack->addRelationTo(mcTrack, -purity);
341  prTrack->addRelationTo(mcTrack->getRelatedFrom<MCParticle>(m_MCParticleCollectionName), -purity);
342  }
343  B2DEBUG(100, "Classified PRTrack " << prTrackId << " as clone -> mcTrackId "
344  << mcTrackId << " : " << -purity);
345  } else {
346  // GHOST
347  B2DEBUG(100, "Classified PRTrack " << prTrackId
348  << " as ghost because of too low efficiency in purest MCTrack "
349  << "(mcTrackId=" << mcTrackId << ").");
350  }
351  } else {
352  // MATCHED
353  //Add the mc matching relation
354  prTrack->addRelationTo(mcTrack, purity);
355  prTrack->addRelationTo(mcTrack->getRelatedFrom<MCParticle>(m_MCParticleCollectionName), purity);
356  B2DEBUG(100, "Classified PRTrack " << prTrackId << " as match -> mcTrackId "
357  << mcTrackId << " : " << purity);
358  }
359  }
360 
361  }
362 
363  // ### Classify the Monte-Carlo tracks ###
364  for (TrackId mcTrackId = 0; mcTrackId < nMCTracks; ++mcTrackId) {
365  MCParticle* mcTrack = m_mcTracks[mcTrackId];
366 
367  const pair<TrackId, Efficiency>& mostEfficiencyPRTrackId =
368  mostEfficientPRTrackId_by_mcTrackId[mcTrackId];
369  const TrackId& prTrackId = mostEfficiencyPRTrackId.first;
370  const Efficiency& efficiency = mostEfficiencyPRTrackId.second;
371 
372  if (!(efficiency >= m_minEfficiency)) {
373  // MISSING
374  B2DEBUG(100, "Classified MCTrack " << mcTrackId << " as missing.");
375  } else {
376  // check whether the highest efficiency pattern recognition track has in turn
377  // the highest purity Monte-Carlo track equal to this track.
378  CDCTriggerTrack* prTrack = m_prTracks[prTrackId];
379 
380  const pair<TrackId, Purity>& purestMCTrackId = purestMCTrackId_by_prTrackId[prTrackId];
381  const TrackId& mcTrackIdCompare = purestMCTrackId.first;
382 
383  if (mcTrackId != mcTrackIdCompare) {
384  // MERGED
385  if (m_relateClonesAndMerged) {
386  // Add the pr matching relation with negative weight
387  mcTrack->addRelationTo(prTrack, -efficiency);
388  mcTrack->getRelatedFrom<MCParticle>(m_MCParticleCollectionName)->addRelationTo(prTrack, -efficiency);
389  }
390  B2DEBUG(100, "Classifid MCTrack " << mcTrackId << " as merge -> prTrackId "
391  << prTrackId << " : " << -efficiency);
392 
393  } else {
394  // MATCHED
395  // Add the pr matching relation
396  mcTrack->addRelationTo(prTrack, efficiency);
397  mcTrack->getRelatedFrom<MCParticle>(m_MCParticleCollectionName)->addRelationTo(prTrack, efficiency);
398  B2DEBUG(100, "Classified MCTrack " << mcTrackId << " as match -> prTrackId "
399  << prTrackId << " : " << efficiency);
400  }
401  }
402  }
403 
404  B2DEBUG(100, "########## end MC matching ############");
405 }
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
Belle2::MCParticle::getCharge
float getCharge() const
Return the particle charge defined in TDatabasePDG.
Definition: MCParticle.cc:35
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::CDCTriggerTrack
Track created by the CDC trigger.
Definition: CDCTriggerTrack.h:15
Belle2::RelationsInterface::addRelationTo
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
Definition: RelationsObject.h:144
Belle2::MCParticle::hasSeenInDetector
bool hasSeenInDetector(Const::DetectorSet set) const
Return if the seen-in flag for a specific subdetector is set or not.
Definition: MCParticle.h:318
Belle2::CDCTriggerMCMatcherModule
A module to match CDCTriggerTracks to MCParticles.
Definition: CDCTriggerMCMatcherModule.h:22
Belle2::RelationsInterface::getRelationsTo
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
Definition: RelationsObject.h:199
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::MCParticle::hasStatus
bool hasStatus(unsigned short int bitmask) const
Return if specific status bit is set.
Definition: MCParticle.h:140
Belle2::CDCTriggerMCMatcherModule::initialize
virtual void initialize() override
Initialize the module.
Definition: CDCTriggerMCMatcherModule.cc:83
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::MCParticle::getPDG
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:123
Belle2::CDCTriggerMCMatcherModule::event
virtual void event() override
Called once for each event.
Definition: CDCTriggerMCMatcherModule.cc:102
Belle2::MCParticle
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:43
Belle2::CDCTriggerSegmentHit
Combination of several CDCHits to a track segment hit for the trigger.
Definition: CDCTriggerSegmentHit.h:16
Belle2::MCParticle::c_PrimaryParticle
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:58
Belle2::RelationsInterface::getRelatedFrom
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
Definition: RelationsObject.h:265