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