1 #include "trg/cdc/modules/mcmatcher/CDCTriggerMCMatcherModule.h"
3 #include <framework/gearbox/Const.h>
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;}
18 inline iter_pair_range<Iter> as_range(std::pair<Iter, Iter>
const& x)
20 return iter_pair_range<Iter>(x);
26 typedef float Efficiency;
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 "
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.",
54 addParam(
"hitCollectionName", m_hitCollectionName,
55 "Name of the StoreArray of CDCTriggerSegmentHits used for the matching.",
57 addParam(
"minAxial", m_minAxial,
58 "Minimum number of axial hits in different super layers.",
60 addParam(
"minStereo", m_minStereo,
61 "Minimum number of stereo hits in different super layers.",
63 addParam(
"axialOnly", m_axialOnly,
64 "Switch to ignore stereo hits (= 2D matching).",
66 addParam(
"minPurity", m_minPurity,
67 "Minimum purity for reconstructed tracks.",
69 addParam(
"minEfficiency", m_minEfficiency,
70 "Minimum efficiency for MC tracks.",
72 addParam(
"onlyPrimaries", m_onlyPrimaries,
73 "If true, MCTracks are only made from primary particles.",
75 addParam(
"relateClonesAndMerged", m_relateClonesAndMerged,
76 "If true, create relations for clones and merged tracks "
77 "(will get negative weight).",
85 m_segmentHits.isRequired(m_hitCollectionName);
86 m_prTracks.isRequired(m_TrgTrackCollectionName);
87 m_mcParticles.isRequired(m_MCParticleCollectionName);
88 m_mcTracks.registerInDataStore(m_MCTrackableCollectionName);
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);
105 for (
int imc = 0; imc < m_mcParticles.getEntries(); ++imc) {
119 vector<bool> SLcounted;
120 SLcounted.assign(9,
false);
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;
129 SLcounted[iSL] =
true;
132 if (nAxial < m_minAxial || nStereo < m_minStereo)
continue;
135 MCParticle* mcTrack = m_mcTracks.appendNew(m_mcTracks.getPtr(), *mcParticle);
137 for (
unsigned ihit = 0; ihit < relHits.
size(); ++ihit) {
145 B2DEBUG(100,
"########## start MC matching ############");
147 int nMCTracks = m_mcTracks.getEntries();
148 int nPRTracks = m_prTracks.getEntries();
150 B2DEBUG(100,
"Number of pattern recognition tracks is " << nPRTracks);
151 B2DEBUG(100,
"Number of Monte-Carlo tracks is " << nMCTracks);
153 if (not nMCTracks or not nPRTracks) {
161 std::map<HitId, TrackId> mcTrackId_by_hitId;
163 std::map<HitId, TrackId>::iterator itMCInsertHit = mcTrackId_by_hitId.end();
164 TrackId mcTrackId = -1;
165 for (
const MCParticle& mcTrack : m_mcTracks) {
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);
180 std::multimap<HitId, TrackId> prTrackId_by_hitId;
182 std::multimap<HitId, TrackId>::iterator itPRInsertHit = prTrackId_by_hitId.end();
183 TrackId prTrackId = -1;
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);
202 Eigen::MatrixXi confusionMatrix = Eigen::MatrixXi::Zero(nPRTracks, nMCTracks + 1);
207 Eigen::RowVectorXi totalHits_by_mcTrackId = Eigen::RowVectorXi::Zero(nMCTracks + 1);
210 const int mcBkgId = nMCTracks;
214 for (HitId hitId = 0; hitId < m_segmentHits.getEntries(); ++hitId) {
216 if (m_axialOnly && m_segmentHits[hitId]->getISuperLayer() % 2)
continue;
220 auto it_mcTrackId = mcTrackId_by_hitId.find(hitId);
222 (it_mcTrackId == mcTrackId_by_hitId.end()) ? mcBkgId : it_mcTrackId->second;
225 totalHits_by_mcTrackId(mcTrackId) += 1;
228 auto range_prTrackIds = prTrackId_by_hitId.equal_range(hitId);
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;
239 Eigen::VectorXi totalHits_by_prTrackId = confusionMatrix.rowwise().sum();
241 B2DEBUG(200,
"Confusion matrix of the event : " << endl << confusionMatrix);
243 B2DEBUG(200,
"totalHits_by_prTrackId : " << endl << totalHits_by_prTrackId);
244 B2DEBUG(200,
"totalHits_by_mcTrackId : " << endl << totalHits_by_mcTrackId);
247 vector<pair<TrackId, Purity>> purestMCTrackId_by_prTrackId(nPRTracks);
249 for (TrackId prTrackId = 0; prTrackId < nPRTracks; ++prTrackId) {
250 Eigen::RowVectorXi prTrackRow = confusionMatrix.row(prTrackId);
251 Eigen::RowVectorXi::Index purestMCTrackId;
254 int highestHits = prTrackRow.maxCoeff(&purestMCTrackId);
255 int totalHits = totalHits_by_prTrackId(prTrackId);
257 Purity highestPurity = Purity(highestHits) / totalHits;
259 purestMCTrackId_by_prTrackId[prTrackId] =
260 pair<TrackId, Purity>(purestMCTrackId, highestPurity);
266 TrackId prTrackId = -1;
267 B2DEBUG(200,
"PRTrack to highest purity MCTrack relation");
268 for (
const pair< TrackId, Purity>& purestMCTrackId :
269 purestMCTrackId_by_prTrackId) {
271 const Purity& purity = purestMCTrackId.second;
272 const TrackId& mcTrackId = purestMCTrackId.first;
273 B2DEBUG(200,
"prTrackId : " << prTrackId <<
" -> mcTrackId : " << mcTrackId
274 <<
" with purity " << purity);
279 vector<pair<TrackId, Efficiency>> mostEfficientPRTrackId_by_mcTrackId(nMCTracks);
281 for (TrackId mcTrackId = 0; mcTrackId < nMCTracks; ++mcTrackId) {
282 Eigen::VectorXi mcTrackCol = confusionMatrix.col(mcTrackId);
283 Eigen::VectorXi::Index highestHitsPRTrackId;
286 int highestHits = mcTrackCol.maxCoeff(&highestHitsPRTrackId);
287 int totalHits = totalHits_by_mcTrackId(mcTrackId);
289 Efficiency highestEfficiency = Efficiency(highestHits) / totalHits;
291 mostEfficientPRTrackId_by_mcTrackId[mcTrackId] =
292 pair<TrackId, Efficiency> (highestHitsPRTrackId, highestEfficiency);
298 TrackId mcTrackId = -1;
299 B2DEBUG(200,
"MCTrack to highest efficiency PRTrack relation");
300 for (
const pair<TrackId, Efficiency>& mostEfficientPRTrackId :
301 mostEfficientPRTrackId_by_mcTrackId) {
303 const Efficiency& efficiency = mostEfficientPRTrackId.second;
304 const TrackId& prTrackId = mostEfficientPRTrackId.first;
305 B2DEBUG(200,
"mcTrackId : " << mcTrackId <<
" -> prTrackId : " << prTrackId
306 <<
" with efficiency " << efficiency);
311 for (TrackId prTrackId = 0; prTrackId < nPRTracks; ++prTrackId) {
314 const pair<TrackId, Purity>& purestMCTrackId = purestMCTrackId_by_prTrackId[prTrackId];
315 const TrackId& mcTrackId = purestMCTrackId.first;
316 const Purity& purity = purestMCTrackId.second;
318 if (!(purity >= m_minPurity)) {
320 B2DEBUG(100,
"Classified PRTrack " << prTrackId <<
" as ghost because of too low purity.");
321 }
else if (mcTrackId == mcBkgId) {
323 B2DEBUG(100,
"Classified PRTrack " << prTrackId <<
" as background.");
329 const pair<TrackId, Efficiency>& mostEfficientPRTrackId =
330 mostEfficientPRTrackId_by_mcTrackId[mcTrackId];
332 const TrackId& prTrackIdCompare = mostEfficientPRTrackId.first;
333 const Efficiency& efficiency = mostEfficientPRTrackId.second;
335 if (prTrackId != prTrackIdCompare) {
336 if (efficiency >= m_minEfficiency) {
338 if (m_relateClonesAndMerged) {
340 prTrack->addRelationTo(mcTrack, -purity);
343 B2DEBUG(100,
"Classified PRTrack " << prTrackId <<
" as clone -> mcTrackId "
344 << mcTrackId <<
" : " << -purity);
347 B2DEBUG(100,
"Classified PRTrack " << prTrackId
348 <<
" as ghost because of too low efficiency in purest MCTrack "
349 <<
"(mcTrackId=" << mcTrackId <<
").");
354 prTrack->addRelationTo(mcTrack, purity);
356 B2DEBUG(100,
"Classified PRTrack " << prTrackId <<
" as match -> mcTrackId "
357 << mcTrackId <<
" : " << purity);
364 for (TrackId mcTrackId = 0; mcTrackId < nMCTracks; ++mcTrackId) {
367 const pair<TrackId, Efficiency>& mostEfficiencyPRTrackId =
368 mostEfficientPRTrackId_by_mcTrackId[mcTrackId];
369 const TrackId& prTrackId = mostEfficiencyPRTrackId.first;
370 const Efficiency& efficiency = mostEfficiencyPRTrackId.second;
372 if (!(efficiency >= m_minEfficiency)) {
374 B2DEBUG(100,
"Classified MCTrack " << mcTrackId <<
" as missing.");
380 const pair<TrackId, Purity>& purestMCTrackId = purestMCTrackId_by_prTrackId[prTrackId];
381 const TrackId& mcTrackIdCompare = purestMCTrackId.first;
383 if (mcTrackId != mcTrackIdCompare) {
385 if (m_relateClonesAndMerged) {
390 B2DEBUG(100,
"Classifid MCTrack " << mcTrackId <<
" as merge -> prTrackId "
391 << prTrackId <<
" : " << -efficiency);
398 B2DEBUG(100,
"Classified MCTrack " << mcTrackId <<
" as match -> prTrackId "
399 << prTrackId <<
" : " << efficiency);
404 B2DEBUG(100,
"########## end MC matching ############");