Belle II Software  release-08-01-10
MCRecoTracksMatcherModule.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 #include <tracking/modules/mcMatcher/MCRecoTracksMatcherModule.h>
9 
10 // hit types
11 #include <pxd/dataobjects/PXDCluster.h>
12 #include <svd/dataobjects/SVDCluster.h>
13 #include <cdc/dataobjects/CDCHit.h>
14 #include <mdst/dataobjects/Track.h>
15 #include <mdst/dataobjects/TrackFitResult.h>
16 
17 #include <map>
18 
19 #include <Eigen/Dense>
20 
21 using namespace Belle2;
22 
23 REG_MODULE(MCRecoTracksMatcher);
24 
25 namespace {
26  // small anonymous helper construct making converting a pair of iterators usable
27  // with range based for
28  template <class Iter>
29  struct iter_pair_range : std::pair<Iter, Iter> {
30  explicit iter_pair_range(std::pair<Iter, Iter> const& x)
31  : std::pair<Iter, Iter>(x)
32  {
33  }
34 
35  Iter begin() const
36  {
37  return this->first;
38  }
39 
40  Iter end() const
41  {
42  return this->second;
43  }
44 
45  bool empty() const
46  {
47  return begin() == end();
48  }
49  };
50 
51  template <class Iter>
52  inline iter_pair_range<Iter> as_range(std::pair<Iter, Iter> const& x)
53  {
54  return iter_pair_range<Iter>(x);
55  }
56 
57  using DetId = Const::EDetector;
58  using HitId = int;
59  using RecoTrackId = int;
60  struct WeightedRecoTrackId {
61  operator int() const { return id; }
62  RecoTrackId id;
63  double weight;
64  };
65 
66  using DetHitIdPair = std::pair<DetId, HitId>;
67 
68  struct CompDetHitIdPair {
69  bool operator()(const std::pair<DetId, HitId>& lhs,
70  const std::pair<std::pair<DetId, HitId>, WeightedRecoTrackId>& rhs)
71  {
72  return lhs < rhs.first;
73  }
74 
75  bool operator()(const std::pair<std::pair<DetId, HitId>, WeightedRecoTrackId>& lhs,
76  const std::pair<DetId, HitId>& rhs)
77  {
78  return lhs.first < rhs;
79  }
80  };
81 
82  // anonymous helper function to fill a set or a map with the hit IDs and det IDs (we need both a set or a map in the following).
83  template <class AMapOrSet>
84  void fillIDsFromStoreArray(AMapOrSet& recoTrackID_by_hitID,
85  const StoreArray<RecoTrack>& storedRecoTracks)
86  {
87  RecoTrackId recoTrackId = -1;
88  for (const RecoTrack& recoTrack : storedRecoTracks) {
89  ++recoTrackId;
90  std::vector<std::pair<DetHitIdPair, WeightedRecoTrackId> > hitIDsInTrack;
91  double totalWeight = 0;
92  using OriginTrackFinder = RecoHitInformation::OriginTrackFinder;
93  const OriginTrackFinder c_MCTrackFinderAuxiliaryHit =
94  OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
95 
96  for (const RecoHitInformation::UsedCDCHit* cdcHit : recoTrack.getCDCHitList()) {
97  OriginTrackFinder originFinder = recoTrack.getFoundByTrackFinder(cdcHit);
98  double weight = originFinder == c_MCTrackFinderAuxiliaryHit ? 0 : 1;
99  totalWeight += weight;
100  hitIDsInTrack.push_back({{Const::CDC, cdcHit->getArrayIndex()}, {recoTrackId, weight}});
101  }
102  for (const RecoHitInformation::UsedSVDHit* svdHit : recoTrack.getSVDHitList()) {
103  OriginTrackFinder originFinder = recoTrack.getFoundByTrackFinder(svdHit);
104  double weight = originFinder == c_MCTrackFinderAuxiliaryHit ? 0 : 1;
105  totalWeight += weight;
106  hitIDsInTrack.push_back({{Const::SVD, svdHit->getArrayIndex()}, {recoTrackId, weight}});
107  }
108  for (const RecoHitInformation::UsedPXDHit* pxdHit : recoTrack.getPXDHitList()) {
109  OriginTrackFinder originFinder = recoTrack.getFoundByTrackFinder(pxdHit);
110  double weight = originFinder == c_MCTrackFinderAuxiliaryHit ? 0 : 1;
111  totalWeight += weight;
112  hitIDsInTrack.push_back({{Const::PXD, pxdHit->getArrayIndex()}, {recoTrackId, weight}});
113  }
114 
115  // In case all hits are auxiliary for a track - reset all weights to 1
116  if (totalWeight == 0) {
117  for (std::pair<DetHitIdPair, WeightedRecoTrackId>& recoTrack_for_hitID : hitIDsInTrack) {
118  recoTrack_for_hitID.second.weight = 1;
119  }
120  }
121 
122  // Commit to output
123  typename AMapOrSet::iterator itInsertHint = recoTrackID_by_hitID.end();
124  for (std::pair<DetHitIdPair, WeightedRecoTrackId>& recoTrack_for_hitID : hitIDsInTrack) {
125  itInsertHint = recoTrackID_by_hitID.insert(itInsertHint, recoTrack_for_hitID);
126  }
127  }
128  }
129 }
130 
132  : Module()
133 {
134  setDescription("This module compares reconstructed tracks generated by some pattern recognition "
135  "algorithm for PXD, SVD and/or CDC to ideal Monte Carlo tracks and performs a "
136  "matching from the former to the underlying MCParticles.");
137 
139 
140  //Parameter definition
141  // Inputs
142  addParam("prRecoTracksStoreArrayName",
144  "Name of the collection containing the tracks as generate a patter recognition algorithm to be evaluated ",
145  std::string(""));
146 
147  addParam("mcRecoTracksStoreArrayName",
149  "Name of the collection containing the reference tracks as generate by a Monte-Carlo-Tracker (e.g. MCTrackFinder)",
150  std::string("MCGFTrackCands"));
151 
152  addParam("TracksStoreArrayName",
154  "Name of the Tracks StoreArray to be used when checking fitted tracks.",
155  std::string(""));
156 
157  // Hit content to be evaluated
158  addParam("UsePXDHits",
159  m_usePXDHits,
160  "Set true if PXDHits or PXDClusters should be used in the matching in case they are present",
161  true);
162 
163  addParam("UseSVDHits",
164  m_useSVDHits,
165  "Set true if SVDHits or SVDClusters should be used in the matching in case they are present",
166  true);
167 
168  addParam("UseCDCHits",
169  m_useCDCHits,
170  "Set true if CDCHits should be used in the matching in case they are present",
171  true);
172 
173  addParam("UseOnlyAxialCDCHits",
175  "Set true if only the axial CDCHits should be used",
176  false);
177 
178  addParam("MinimalPurity",
180  "Minimal purity of a PRTrack to be considered matchable to a MCTrack. "
181  "This number encodes how many correct hits are minimally need to compensate for a false hits. "
182  "The default 2.0 / 3.0 suggests that for each background hit can be compensated by two correct hits.",
183  2.0 / 3.0);
184 
185  addParam("MinimalEfficiency",
187  "Minimal efficiency of a MCTrack to be considered matchable to a PRTrack. "
188  "This number encodes which fraction of the true hits must at least be in the reconstructed track. "
189  "The default 0.05 suggests that at least 5% of the true hits should have been picked up.",
190  0.05);
191 
192  addParam("useFittedTracks",
194  "If true, it uses fitted tracks for matching. Note that the charge of the track can be different from\
195  the seed charge (that is provided by the pattern recognition) since the DAF can flip tracks.",
197 }
198 
200 {
201  if (m_MCParticles.isOptional()) {
202  m_mcParticlesPresent = true;
203 
204  // Require both RecoTrack arrays and the MCParticles to be present in the DataStore
208 
209  // Purity relation - for each PRTrack to store the purest MCTrack
211 
212  // Efficiency relation - for each MCTrack to store the most efficient PRTrack
214 
215  // MCParticle matching relation - purity
217 
218  // MCParticle matching relation - efficiency
220 
221  // Announce optional store arrays to the hits or clusters in case they should be used
222  // We make them optional in case of limited detector setup.
223  // PXD
224  if (m_usePXDHits) {
226  }
227 
228  // SVD
229  if (m_useSVDHits) {
231  }
232 
233  // CDC
234  if (m_useCDCHits) {
235  m_CDCHits.isOptional();
236  }
237  }
238 }
239 
241 {
242  // Skip in the case there are no MC particles present.
243  if (not m_mcParticlesPresent) {
244  B2DEBUG(23, "Skipping MC Track Matcher as there are no MC Particles registered in the DataStore.");
245  return;
246  }
247 
248  B2DEBUG(23, "########## MCRecoTracksMatcherModule ############");
249 
250  int nMCRecoTracks = m_MCRecoTracks.getEntries();
251  int nPRRecoTracks = m_PRRecoTracks.getEntries();
252 
253  B2DEBUG(23, "Number patter recognition tracks is " << nPRRecoTracks);
254  B2DEBUG(23, "Number Monte-Carlo tracks is " << nMCRecoTracks);
255 
256  if (not nMCRecoTracks or not nPRRecoTracks) {
257  // Neither no pattern recognition tracks
258  // or no Monte Carlo tracks are present in this event
259  // Cannot perform matching.
260  return;
261  }
262 
263  // ### Build a detector_id hit_id to reco track map for easier lookup later ###
264  std::multimap<DetHitIdPair, WeightedRecoTrackId > mcId_by_hitId;
265  fillIDsFromStoreArray(mcId_by_hitId, m_MCRecoTracks);
266 
267  // Use set instead of multimap to handel to following situation
268  // * One hit may be assigned to multiple tracks should contribute to the efficiency of both tracks
269  // * One hit assigned twice or more to the same track should not contribute to the purity multiple times
270  // The first part is handled well by the multimap. But to enforce that one hit is also only assigned
271  // once to a track we use a set.
272  std::set<std::pair<DetHitIdPair, WeightedRecoTrackId>> prId_by_hitId;
273  fillIDsFromStoreArray(prId_by_hitId, m_PRRecoTracks);
274 
275  // ### Get the number of relevant hits for each detector ###
276  // Since we are mostly dealing with indices in this module, this is all we need from the StoreArray
277  // Silently skip store arrays that are not present in reduced detector setups.
278  std::map<DetId, int> nHits_by_detId;
279 
280  // PXD
281  if (m_usePXDHits) {
282  nHits_by_detId[Const::PXD] = m_PXDClusters.getEntries();
283  }
284 
285  // SVD
286  if (m_useSVDHits) {
287  nHits_by_detId[Const::SVD] = m_SVDClusters.getEntries();
288  }
289 
290  // CDC
291  if (m_useCDCHits) {
292  nHits_by_detId[Const::CDC] = m_CDCHits.getEntries();
293  }
294 
295  //### Build the confusion matrix ###
296 
297  // Reserve enough space for the confusion matrix
298  // The last row is meant for hits not assigned to a mcRecoTrack (aka background hits)
299  Eigen::MatrixXd confusionMatrix = Eigen::MatrixXd::Zero(nPRRecoTracks, nMCRecoTracks + 1);
300  Eigen::MatrixXd weightedConfusionMatrix = Eigen::MatrixXd::Zero(nPRRecoTracks, nMCRecoTracks + 1);
301 
302  // Accumulated the total number of hits/ndf for each Monte-Carlo track separately to avoid double counting,
303  // in case patter recognition tracks share hits.
304  Eigen::RowVectorXd totalNDF_by_mcId = Eigen::RowVectorXd::Zero(nMCRecoTracks + 1);
305  Eigen::RowVectorXd totalWeight_by_mcId = Eigen::RowVectorXd::Zero(nMCRecoTracks + 1);
306 
307  // Accumulated the total number of hits/ndf for each pattern recognition track separately to avoid double counting,
308  // in case Monte-Carlo tracks share hits.
309  Eigen::VectorXd totalNDF_by_prId = Eigen::VectorXd::Zero(nPRRecoTracks);
310 
311  // Column index for the hits not assigned to any MCRecoTrack
312  const int mcBkgId = nMCRecoTracks;
313 
314  // for each detector examine every hit to which mcRecoTrack and prRecoTrack it belongs
315  // if the hit is not part of any mcRecoTrack put the hit in the background column.
316  for (const std::pair<const DetId, NDF>& detId_nHits_pair : nHits_by_detId) {
317 
318  DetId detId = detId_nHits_pair.first;
319  int nHits = detId_nHits_pair.second;
320  NDF ndfForOneHit = m_ndf_by_detId[detId];
321 
322  for (HitId hitId = 0; hitId < nHits; ++hitId) {
323  DetHitIdPair detId_hitId_pair(detId, hitId);
324 
325  if (m_useOnlyAxialCDCHits and detId == Const::CDC) {
326  StoreArray<CDCHit> cdcHits;
327  const CDCHit* cdcHit = cdcHits[hitId];
328  if (cdcHit->getISuperLayer() % 2 != 0) {
329  // Skip stereo hits
330  continue;
331  }
332  }
333 
334  // Seek all Monte Carlo tracks with the given hit
335  const auto mcIds_for_detId_hitId_pair =
336  as_range(mcId_by_hitId.equal_range(detId_hitId_pair));
337 
338  // Seek all pattern recognition tracks with the given hit
339  const auto prIds_for_detId_hitId_pair =
340  as_range(std::equal_range(prId_by_hitId.begin(),
341  prId_by_hitId.end(),
342  detId_hitId_pair,
343  CompDetHitIdPair()));
344 
345  // Assign the hits/ndf to the total ndf vector separately to avoid double counting,
346  // if patter recognition track share hits.
347  if (mcIds_for_detId_hitId_pair.empty()) {
348  // If the hit is not assigned to any mcRecoTrack
349  // The hit is assigned to the background column
350  RecoTrackId mcId = mcBkgId;
351  double mcWeight = 1;
352  totalNDF_by_mcId(mcId) += ndfForOneHit;
353  totalWeight_by_mcId(mcId) += ndfForOneHit * mcWeight;
354  } else {
355  for (const auto& detId_hitId_pair_and_mcId : mcIds_for_detId_hitId_pair) {
356  int mcId = detId_hitId_pair_and_mcId.second;
357  double mcWeight = detId_hitId_pair_and_mcId.second.weight;
358  totalNDF_by_mcId(mcId) += ndfForOneHit;
359  totalWeight_by_mcId(mcId) += ndfForOneHit * mcWeight;
360  }
361  }
362 
363  // Assign the hits/ndf to the total ndf vector separately here to avoid double counting,
364  // if Monte-Carlo track share hits.
365  for (const auto& detId_hitId_pair_and_prId : prIds_for_detId_hitId_pair) {
366  RecoTrackId prId = detId_hitId_pair_and_prId.second;
367  totalNDF_by_prId(prId) += ndfForOneHit;
368  }
369 
370  // Fill the confusion matrix
371  for (const auto& detId_hitId_pair_and_prId : prIds_for_detId_hitId_pair) {
372  int prId = detId_hitId_pair_and_prId.second;
373  if (mcIds_for_detId_hitId_pair.empty()) {
374  int mcId = mcBkgId;
375  double mcWeight = 1;
376  confusionMatrix(prId, mcId) += ndfForOneHit;
377  weightedConfusionMatrix(prId, mcId) += ndfForOneHit * mcWeight;
378  } else {
379  for (const auto& detId_hitId_pair_and_mcId : mcIds_for_detId_hitId_pair) {
380  int mcId = detId_hitId_pair_and_mcId.second;
381  double mcWeight = detId_hitId_pair_and_mcId.second.weight;
382  confusionMatrix(prId, mcId) += ndfForOneHit;
383  weightedConfusionMatrix(prId, mcId) += ndfForOneHit * mcWeight;
384  }
385  }
386  } // end for
387  } // end for hitId
388  } // end for detId
389 
390  B2DEBUG(24, "Confusion matrix of the event : " << std::endl << confusionMatrix);
391  B2DEBUG(24, "Weighted confusion matrix of the event : " << std::endl << weightedConfusionMatrix);
392 
393  B2DEBUG(24, "totalNDF_by_mcId : " << std::endl << totalNDF_by_mcId);
394  B2DEBUG(24, "totalWeight_by_mcId : " << std::endl << totalWeight_by_mcId);
395 
396  B2DEBUG(24, "totalNDF_by_prId : " << std::endl << totalNDF_by_prId);
397 
398  Eigen::MatrixXd purityMatrix = confusionMatrix.array().colwise() / totalNDF_by_prId.array();
399  Eigen::MatrixXd efficiencyMatrix = confusionMatrix.array().rowwise() / totalNDF_by_mcId.array();
400  Eigen::MatrixXd weightedEfficiencyMatrix = weightedConfusionMatrix.array().rowwise() / totalWeight_by_mcId.array();
401 
402  B2DEBUG(23, "Purities");
403  B2DEBUG(23, purityMatrix);
404 
405  B2DEBUG(23, "Efficiencies");
406  B2DEBUG(23, efficiencyMatrix);
407 
408  B2DEBUG(23, "Weighted efficiencies");
409  B2DEBUG(23, weightedEfficiencyMatrix);
410 
411  // ### Building the Monte-Carlo track to highest efficiency patter recognition track relation ###
412  // Weighted efficiency
413  using Efficiency = float;
414  using Purity = float;
415 
416  struct MostWeightEfficientPRId {
417  RecoTrackId id;
418  Efficiency weightedEfficiency;
419  // cppcheck-suppress unusedStructMember
420  Efficiency efficiency;
421  };
422  std::vector<MostWeightEfficientPRId> mostWeightEfficientPRId_by_mcId(nMCRecoTracks);
423  for (RecoTrackId mcId = 0; mcId < nMCRecoTracks; ++mcId) {
424  Eigen::VectorXd efficiencyCol = efficiencyMatrix.col(mcId);
425  Eigen::VectorXd weightedEfficiencyCol = weightedEfficiencyMatrix.col(mcId);
426 
427  RecoTrackId bestPrId = 0;
428  Efficiency bestWeightedEfficiency = weightedEfficiencyCol(0);
429  Efficiency bestEfficiency = efficiencyCol(0);
430  Purity bestPurity = purityMatrix.row(0)(mcId);
431 
432  // Reject efficiency smaller than the minimal one
433  if (bestWeightedEfficiency < m_minimalEfficiency) {
434  bestWeightedEfficiency = 0;
435  }
436 
437  // In case of a tie in the weighted efficiency we use the regular efficiency to break it.
438  for (RecoTrackId prId = 1; prId < nPRRecoTracks; ++prId) {
439  Eigen::RowVectorXd purityRow = purityMatrix.row(prId);
440 
441  Efficiency currentWeightedEfficiency = weightedEfficiencyCol(prId);
442  Efficiency currentEfficiency = efficiencyCol(prId);
443  Purity currentPurity = purityRow(mcId);
444 
445  // Reject efficiency smaller than the minimal one
446  if (currentWeightedEfficiency < m_minimalEfficiency) {
447  currentWeightedEfficiency = 0;
448  }
449 
450  if (std::tie(currentWeightedEfficiency, currentEfficiency, currentPurity) >
451  std::tie(bestWeightedEfficiency, bestEfficiency, bestPurity)) {
452  bestPrId = prId;
453  bestEfficiency = currentEfficiency;
454  bestWeightedEfficiency = currentWeightedEfficiency;
455  bestPurity = currentPurity;
456  }
457  }
458 
459  bestWeightedEfficiency = weightedEfficiencyCol(bestPrId);
460  bestEfficiency = efficiencyCol(bestPrId);
461  mostWeightEfficientPRId_by_mcId[mcId] = {bestPrId, bestWeightedEfficiency, bestEfficiency};
462  }
463 
464  // ### Building the patter recognition track to highest purity Monte-Carlo track relation ###
465  // Unweighted purity
466  struct MostPureMCId {
467  RecoTrackId id;
468  Purity purity;
469  };
470 
471  std::vector<MostPureMCId> mostPureMCId_by_prId(nPRRecoTracks);
472  for (int prId = 0; prId < nPRRecoTracks; ++prId) {
473  Eigen::RowVectorXd purityRow = purityMatrix.row(prId);
474 
475  int mcId;
476  Purity highestPurity = purityRow.maxCoeff(&mcId);
477 
478  mostPureMCId_by_prId[prId] = {mcId, highestPurity};
479  }
480 
481  // Log the Monte-Carlo track to highest efficiency patter recognition track relation
482  // Weighted efficiency
483  {
484  RecoTrackId mcId = -1;
485  B2DEBUG(24, "MCTrack to highest weighted efficiency PRTrack relation");
486  for (const auto& mostWeightEfficientPRId_for_mcId : mostWeightEfficientPRId_by_mcId) {
487  ++mcId;
488  const Efficiency& weightedEfficiency = mostWeightEfficientPRId_for_mcId.weightedEfficiency;
489  const RecoTrackId& prId = mostWeightEfficientPRId_for_mcId.id;
490  B2DEBUG(24,
491  "mcId : " << mcId << " -> prId : " << prId << " with weighted efficiency "
492  << weightedEfficiency);
493  }
494  }
495 
496  // Log the patter recognition track to highest purity Monte-Carlo track relation
497  // Unweighted purity
498  {
499  RecoTrackId prId = -1;
500  B2DEBUG(24, "PRTrack to highest purity MCTrack relation");
501  for (const auto& mostPureMCId_for_prId : mostPureMCId_by_prId) {
502  ++prId;
503  const RecoTrackId& mcId = mostPureMCId_for_prId.id;
504  const Purity& purity = mostPureMCId_for_prId.purity;
505  B2DEBUG(24, "prId : " << prId << " -> mcId : " << mcId << " with purity " << purity);
506  }
507  }
508 
509  // Count the categories
510  int nMatched{}, nWrongCharge{}, nBackground{}, nClones{}, nClonesWrongCharge{}, nGhost{};
511 
512  // ### Classify the patter recognition tracks ###
513  // Means saving the highest purity relation to the data store
514  // + setup the PRTrack to MCParticle relation
515  // + save the set the MatchingStatus
516  for (RecoTrackId prId = 0; prId < nPRRecoTracks; ++prId) {
517  RecoTrack* prRecoTrack = m_PRRecoTracks[prId];
518 
519  const MostPureMCId& mostPureMCId = mostPureMCId_by_prId[prId];
520 
521  const RecoTrackId& mcId = mostPureMCId.id;
522  const Purity& purity = mostPureMCId.purity;
523 
524  // GHOST
525  if (not(purity > 0) or not(purity >= m_minimalPurity)) {
526  prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_ghost);
527 
528  B2DEBUG(23, "Stored PRTrack " << prId << " as ghost because of too low purity");
529  ++nGhost;
530  continue;
531  }
532 
533  // BACKGROUND
534  if (mcId == mcBkgId) {
535  prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_background);
536 
537  B2DEBUG(23, "Stored PRTrack " << prId << " as background because of too low purity.");
538  ++nBackground;
539  continue;
540  }
541 
542  // After the classification for bad purity and background we examine,
543  // whether the highest purity Monte-Carlo track has in turn the highest efficiency
544  // patter recognition track equal to this track.
545  // All extra patter recognition tracks are marked as clones.
546 
547  RecoTrack* mcRecoTrack = m_MCRecoTracks[mcId];
548  MCParticle* mcParticle = mcRecoTrack->getRelated<MCParticle>();
549 
550  B2ASSERT("No relation from MCRecoTrack to MCParticle.", mcParticle);
551 
552  const MostWeightEfficientPRId& mostWeightEfficientPRId_for_mcId =
553  mostWeightEfficientPRId_by_mcId[mcId];
554 
555  const RecoTrackId& mostWeightEfficientPRId = mostWeightEfficientPRId_for_mcId.id;
556  const Efficiency& weightedEfficiency = mostWeightEfficientPRId_for_mcId.weightedEfficiency;
557  // const Efficiency& efficiency = mostWeightEfficientPRId_for_mcId.efficiency;
558 
559  // find the true charge and reconstructed charge
560  const short MCParticleTrackCharge = mcParticle->getCharge() > 0 ? 1 : -1;
561  short foundTrackCharge = prRecoTrack->getChargeSeed();
562  if (m_useFittedTracks) {
563  const RelationVector<Track> fittedTracks = prRecoTrack->getRelationsFrom<Track>(m_TracksStoreArrayName);
564  short nPositiveCharges = 0;
565  short nNegativeCharges = 0;
566 
567 
568 
569  if (fittedTracks.size() > 0) {
570 
571  for (const auto& fittedTrack : fittedTracks) {
572 
573  // catch rare case we track long lived non-chargedStable particles (e.g. Sigma)
574  try {
575  const TrackFitResult* trackFitResult = fittedTrack.getTrackFitResultWithClosestMass(Const::ChargedStable(std::abs(
576  mcParticle->getPDG())));
577  trackFitResult->getChargeSign() > 0 ? nPositiveCharges++ : nNegativeCharges++;
578  } catch (...) {
579  const TrackFitResult* trackFitResult = fittedTrack.getTrackFitResultWithClosestMass(Const::pion);
580  trackFitResult->getChargeSign() > 0 ? nPositiveCharges++ : nNegativeCharges++;
581  }
582  }
583  }
584  if (nPositiveCharges > 0 and nNegativeCharges > 0) {
585  B2DEBUG(23,
586  "There are different charges attributed to the same track, this shouldn't happen. Continue with the majority of positive or negative charges");
587  }
588  // Only use nPositiveCharges and nNegativeCharges to assign a new value to foundTrackCharge if at least one fitted Tracks exists
589  // and at least one of the two values is > 0
590  if (fittedTracks.size() > 0 and (nPositiveCharges > 0 or nNegativeCharges > 0)) {
591  foundTrackCharge = nPositiveCharges > nNegativeCharges ? 1 : -1;
592  }
593  }
594 
595  // Note : The matched category may also contain higher order clones recognisable by their low
596  // absolute efficiency
597 
598  // MATCHED
599  if (prId == mostWeightEfficientPRId) {
600  if (foundTrackCharge != MCParticleTrackCharge) {
601  prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_matchedWrongCharge);
602  ++nWrongCharge;
603  } else {
604  prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_matched);
605  ++nMatched;
606  }
607  // Setup the relation purity relation
608  // regardless of the charge matching
609  prRecoTrack->addRelationTo(mcRecoTrack, purity);
610 
611  // Add the mc matching relation
612  prRecoTrack->addRelationTo(mcParticle, purity);
613 
614  B2DEBUG(23, "Stored PRTrack " << prId << " as matched.");
615  B2DEBUG(23, "MC Match prId " << prId << " to mcPartId " << mcParticle->getArrayIndex());
616  B2DEBUG(23, "Purity rel: prId " << prId << " -> mcId " << mcId << " : " << purity);
617  continue;
618  }
619 
620  // GHOST
621  // Pattern recognition track fails the minimal efficiency requirement to be matched.
622  // We might want to introduce a different classification here, if we see problems
623  // with too many ghosts and want to investigate the specific source of the mismatch.
624  if (not(weightedEfficiency >= m_minimalEfficiency)) {
625  prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_ghost);
626  B2DEBUG(23, "Stored PRTrack " << prId << " as ghost because of too low efficiency.");
627  ++nGhost;
628  continue;
629  }
630 
631  // Final category
632  // CLONE
633  if (foundTrackCharge != MCParticleTrackCharge) {
634  prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_cloneWrongCharge);
635  ++nClonesWrongCharge;
636  } else {
637  prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_clone);
638  ++nClones;
639  }
640  // Setup the relation purity relation
641  // regardless whether the charge is correctly reconstructed or not
642  prRecoTrack->addRelationTo(mcRecoTrack, -purity);
643  prRecoTrack->addRelationTo(mcParticle, -purity);
644  ++nClones;
645 
646  // Add the mc matching relation to the mc particle
647  B2DEBUG(23, "Stored PRTrack " << prId << " as clone.");
648  B2DEBUG(23, "MC Match prId " << prId << " to mcPartId " << mcParticle->getArrayIndex());
649  B2DEBUG(23, "Purity rel: prId " << prId << " -> mcId " << mcId << " : " << -purity);
650  } // end for prId
651 
652 
653 
654  B2DEBUG(23, "Number of matches " << nMatched);
655  B2DEBUG(23, "Number of wrongCharge " << nWrongCharge);
656  B2DEBUG(23, "Number of clones " << nClones);
657  B2DEBUG(23, "Number of clones wrongCharge" << nClonesWrongCharge);
658  B2DEBUG(23, "Number of bkg " << nBackground);
659  B2DEBUG(23, "Number of ghost " << nGhost);
660 
661  // ### Classify the Monte-Carlo tracks ###
662  // Meaning save the highest weighted efficiency relation to the data store.
663  for (RecoTrackId mcId = 0; mcId < nMCRecoTracks; ++mcId) {
664  RecoTrack* mcRecoTrack = m_MCRecoTracks[mcId];
665  MCParticle* mcParticle = mcRecoTrack->getRelated<MCParticle>();
666 
667  const MostWeightEfficientPRId& mostWeightEfficiencyPRId = mostWeightEfficientPRId_by_mcId[mcId];
668 
669  const RecoTrackId& prId = mostWeightEfficiencyPRId.id;
670  const Efficiency& weightedEfficiency = mostWeightEfficiencyPRId.weightedEfficiency;
671  // const Efficiency& efficiency = mostWeightEfficiencyPRId.efficiency;
672 
673  B2ASSERT("Index of pattern recognition tracks out of range.", prId < nPRRecoTracks and prId >= 0);
674 
675  RecoTrack* prRecoTrack = m_PRRecoTracks[prId];
676 
677  const MostPureMCId& mostPureMCId_for_prId = mostPureMCId_by_prId[prId];
678  const RecoTrackId& mostPureMCId = mostPureMCId_for_prId.id;
679 
680  // MATCHED
681  if (mcId == mostPureMCId and
682  (prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_matched or
683  prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_matchedWrongCharge or
684  prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_clone or
685  prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_cloneWrongCharge)) {
686  // Setup the relation with positive weighted efficiency for this case.
687  mcRecoTrack->addRelationTo(prRecoTrack, weightedEfficiency);
688  mcParticle->addRelationTo(prRecoTrack, weightedEfficiency);
689  B2DEBUG(23, "Efficiency rel: mcId " << mcId << " -> prId " << prId << " : " << weightedEfficiency);
690  continue;
691  }
692 
693  // MERGED
694  // This MCTrack has a significant portion of hits in a PRTrack
695  // which in turn better describes a MCTrack different form this.
696  // Setup the relation with negative weighted efficiency for this case.
697  bool isMergedMCRecoTrack =
698  (weightedEfficiency >= m_minimalEfficiency) and
699  (prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_matched or
700  prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_matchedWrongCharge or
701  prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_clone or
702  prRecoTrack->getMatchingStatus() == RecoTrack::MatchingStatus::c_cloneWrongCharge);
703 
704  if (isMergedMCRecoTrack) {
705  mcRecoTrack->addRelationTo(prRecoTrack, -weightedEfficiency);
706  mcParticle->addRelationTo(prRecoTrack, -weightedEfficiency);
707  B2DEBUG(23, "Efficiency rel: mcId " << mcId << " -> prId " << prId << " : " << -weightedEfficiency);
708  continue;
709  }
710 
711  // MISSING
712  // No related pattern recognition track
713  // Do not create a relation.
714  B2DEBUG(23, "mcId " << mcId << " is missing. No relation created.");
715  B2DEBUG(23, "is Primary" << m_MCRecoTracks[mcId]->getRelatedTo<MCParticle>()->isPrimaryParticle());
716  B2DEBUG(23, "best prId " << prId << " with purity " << mostPureMCId_for_prId.purity << " -> " << mostPureMCId);
717  B2DEBUG(23, "MC Total ndf " << totalNDF_by_mcId[mcId]);
718  B2DEBUG(23, "MC Total weight" << totalWeight_by_mcId[mcId]);
719  B2DEBUG(23, "MC Overlap ndf\n " << confusionMatrix.col(mcId).transpose());
720  B2DEBUG(23, "MC Overlap weight\n " << weightedConfusionMatrix.col(mcId).transpose());
721  B2DEBUG(23, "MC Efficiencies for the track\n" << efficiencyMatrix.col(mcId).transpose());
722  B2DEBUG(23, "MC Weighted efficiencies for the track\n" << weightedEfficiencyMatrix.col(mcId).transpose());
723  } // end for mcId
724 
725  B2DEBUG(23, "########## End MCRecoTracksMatcherModule ############");
726 
727 } //end event()
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
unsigned short getISuperLayer() const
Getter for iSuperLayer.
Definition: CDCHit.h:184
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:580
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition: Const.h:42
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
int getArrayIndex() const
Get 0-based index of the particle in the corresponding MCParticle list.
Definition: MCParticle.h:244
float getCharge() const
Return the particle charge defined in TDatabasePDG.
Definition: MCParticle.cc:36
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
bool m_useCDCHits
Parameter : Switch whether CDCHits should be used in the matching.
StoreArray< SVDCluster > m_SVDClusters
StoreArray containing SVDClusters.
void initialize() final
Signal the beginning of the event processing.
std::string m_TracksStoreArrayName
Parameter : Name of the Tracks StoreArray.
int NDF
Descriptive type defintion for a number of degrees of freedom.
bool m_mcParticlesPresent
Flag to indicated whether the Monte Carlo track are on the DataStore.
MCRecoTracksMatcherModule()
Constructor setting up the parameter.
StoreArray< RecoTrack > m_PRRecoTracks
StoreArray containing PR RecoTracks.
double m_minimalPurity
Parameter : Minimal purity of a PRTrack to be considered matchable to a MCTrack.
double m_minimalEfficiency
Parameter : Minimal efficiency for a MCTrack to be considered matchable to a PRTrack.
void event() final
Process the event.
bool m_useFittedTracks
Use fitted tracks for matching.
StoreArray< PXDCluster > m_PXDClusters
StoreArray containing PXDClusters.
bool m_usePXDHits
Parameter : Switch whether PXDHits should be used in the matching.
StoreArray< RecoTrack > m_MCRecoTracks
StoreArray containing MC RecoTracks.
std::string m_prRecoTracksStoreArrayName
Parameter : Name of the RecoTracks StoreArray from pattern recognition.
std::map< int, NDF > m_ndf_by_detId
Map storing the standard number degrees of freedom for a single hit by detector *‍/.
bool m_useOnlyAxialCDCHits
Parameter : Switch whether only axial CDCHits should be used.
std::string m_mcRecoTracksStoreArrayName
Parameter : Name of the RecoTracks StoreArray from MC track finding.
StoreArray< MCParticle > m_MCParticles
StoreArray containing MCParticles.
StoreArray< CDCHit > m_CDCHits
StoreArray containing CDCHits.
bool m_useSVDHits
Parameter : Switch whether SVDHits should be used in the matching.
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
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
OriginTrackFinder
The TrackFinder which has added the hit to the track.
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
void setMatchingStatus(MatchingStatus matchingStatus)
Set the matching status (used by the TrackMatcher module)
Definition: RecoTrack.h:835
MatchingStatus getMatchingStatus() const
Return the matching status set by the TrackMatcher module.
Definition: RecoTrack.h:829
short int getChargeSeed() const
Return the charge seed stored in the reco track. ATTENTION: This is not the fitted charge.
Definition: RecoTrack.h:508
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
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).
RelationVector< FROM > getRelationsFrom(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from another store array to this object.
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
Class that bundles various TrackFitResults.
Definition: Track.h:25
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.