8#include <tracking/modules/mcMatcher/MCRecoTracksMatcherModule.h>
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>
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)
47 return begin() == end();
52 inline iter_pair_range<Iter> as_range(std::pair<Iter, Iter>
const& x)
54 return iter_pair_range<Iter>(x);
59 using RecoTrackId = int;
60 struct WeightedRecoTrackId {
61 operator int()
const {
return id; }
66 using DetHitIdPair = std::pair<DetId, HitId>;
68 struct CompDetHitIdPair {
69 bool operator()(
const std::pair<DetId, HitId>& lhs,
70 const std::pair<std::pair<DetId, HitId>, WeightedRecoTrackId>& rhs)
72 return lhs < rhs.first;
75 bool operator()(
const std::pair<std::pair<DetId, HitId>, WeightedRecoTrackId>& lhs,
76 const std::pair<DetId, HitId>& rhs)
78 return lhs.first < rhs;
83 template <
class AMapOrSet>
84 void fillIDsFromStoreArray(AMapOrSet& recoTrackID_by_hitID,
87 RecoTrackId recoTrackId = -1;
88 for (
const RecoTrack& recoTrack : storedRecoTracks) {
90 std::vector<std::pair<DetHitIdPair, WeightedRecoTrackId> > hitIDsInTrack;
91 double totalWeight = 0;
93 const OriginTrackFinder c_MCTrackFinderAuxiliaryHit =
94 OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
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}});
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}});
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}});
116 if (totalWeight == 0) {
117 for (std::pair<DetHitIdPair, WeightedRecoTrackId>& recoTrack_for_hitID : hitIDsInTrack) {
118 recoTrack_for_hitID.second.weight = 1;
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);
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.");
142 addParam(
"prRecoTracksStoreArrayName",
144 "Name of the collection containing the tracks as generate a pattern recognition algorithm to be evaluated ",
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"));
154 "Name of the Tracks StoreArray to be used when checking fitted tracks.",
160 "Set true if PXDHits or PXDClusters should be used in the matching in case they are present",
165 "Set true if SVDHits or SVDClusters should be used in the matching in case they are present",
170 "Set true if CDCHits should be used in the matching in case they are present",
175 "Set true if only the axial CDCHits should be used",
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.",
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.",
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.",
244 B2DEBUG(23,
"Skipping MC Track Matcher as there are no MC Particles registered in the DataStore.");
248 B2DEBUG(23,
"########## MCRecoTracksMatcherModule ############");
253 B2DEBUG(23,
"Number pattern recognition tracks is " << nPRRecoTracks);
254 B2DEBUG(23,
"Number Monte-Carlo tracks is " << nMCRecoTracks);
256 if (not nMCRecoTracks or not nPRRecoTracks) {
264 std::multimap<DetHitIdPair, WeightedRecoTrackId > mcId_by_hitId;
272 std::set<std::pair<DetHitIdPair, WeightedRecoTrackId>> prId_by_hitId;
278 std::map<DetId, int> nHits_by_detId;
292 nHits_by_detId[Const::CDC] =
m_CDCHits.getEntries();
299 Eigen::MatrixXd confusionMatrix = Eigen::MatrixXd::Zero(nPRRecoTracks, nMCRecoTracks + 1);
300 Eigen::MatrixXd weightedConfusionMatrix = Eigen::MatrixXd::Zero(nPRRecoTracks, nMCRecoTracks + 1);
304 Eigen::RowVectorXd totalNDF_by_mcId = Eigen::RowVectorXd::Zero(nMCRecoTracks + 1);
305 Eigen::RowVectorXd totalWeight_by_mcId = Eigen::RowVectorXd::Zero(nMCRecoTracks + 1);
309 Eigen::VectorXd totalNDF_by_prId = Eigen::VectorXd::Zero(nPRRecoTracks);
312 const int mcBkgId = nMCRecoTracks;
316 for (
const std::pair<const DetId, NDF>& detId_nHits_pair : nHits_by_detId) {
318 DetId detId = detId_nHits_pair.first;
319 int nHits = detId_nHits_pair.second;
322 for (HitId hitId = 0; hitId < nHits; ++hitId) {
323 DetHitIdPair detId_hitId_pair(detId, hitId);
327 const CDCHit* cdcHit = cdcHits[hitId];
335 const auto mcIds_for_detId_hitId_pair =
336 as_range(mcId_by_hitId.equal_range(detId_hitId_pair));
339 const auto prIds_for_detId_hitId_pair =
340 as_range(std::equal_range(prId_by_hitId.begin(),
343 CompDetHitIdPair()));
347 if (mcIds_for_detId_hitId_pair.empty()) {
350 RecoTrackId mcId = mcBkgId;
352 totalNDF_by_mcId(mcId) += ndfForOneHit;
353 totalWeight_by_mcId(mcId) += ndfForOneHit * mcWeight;
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;
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;
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()) {
376 confusionMatrix(prId, mcId) += ndfForOneHit;
377 weightedConfusionMatrix(prId, mcId) += ndfForOneHit * mcWeight;
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;
390 B2DEBUG(24,
"Confusion matrix of the event : " << std::endl << confusionMatrix);
391 B2DEBUG(24,
"Weighted confusion matrix of the event : " << std::endl << weightedConfusionMatrix);
393 B2DEBUG(24,
"totalNDF_by_mcId : " << std::endl << totalNDF_by_mcId);
394 B2DEBUG(24,
"totalWeight_by_mcId : " << std::endl << totalWeight_by_mcId);
396 B2DEBUG(24,
"totalNDF_by_prId : " << std::endl << totalNDF_by_prId);
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();
402 B2DEBUG(23,
"Purities");
403 B2DEBUG(23, purityMatrix);
405 B2DEBUG(23,
"Efficiencies");
406 B2DEBUG(23, efficiencyMatrix);
408 B2DEBUG(23,
"Weighted efficiencies");
409 B2DEBUG(23, weightedEfficiencyMatrix);
413 using Efficiency = float;
414 using Purity = float;
416 struct MostWeightEfficientPRId {
418 Efficiency weightedEfficiency;
420 Efficiency efficiency;
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);
427 RecoTrackId bestPrId = 0;
428 Efficiency bestWeightedEfficiency = weightedEfficiencyCol(0);
429 Efficiency bestEfficiency = efficiencyCol(0);
430 Purity bestPurity = purityMatrix.row(0)(mcId);
434 bestWeightedEfficiency = 0;
438 for (RecoTrackId prId = 1; prId < nPRRecoTracks; ++prId) {
439 Eigen::RowVectorXd purityRow = purityMatrix.row(prId);
441 Efficiency currentWeightedEfficiency = weightedEfficiencyCol(prId);
442 Efficiency currentEfficiency = efficiencyCol(prId);
443 Purity currentPurity = purityRow(mcId);
447 currentWeightedEfficiency = 0;
450 if (std::tie(currentWeightedEfficiency, currentEfficiency, currentPurity) >
451 std::tie(bestWeightedEfficiency, bestEfficiency, bestPurity)) {
453 bestEfficiency = currentEfficiency;
454 bestWeightedEfficiency = currentWeightedEfficiency;
455 bestPurity = currentPurity;
459 bestWeightedEfficiency = weightedEfficiencyCol(bestPrId);
460 bestEfficiency = efficiencyCol(bestPrId);
461 mostWeightEfficientPRId_by_mcId[mcId] = {bestPrId, bestWeightedEfficiency, bestEfficiency};
466 struct MostPureMCId {
471 std::vector<MostPureMCId> mostPureMCId_by_prId(nPRRecoTracks);
472 for (
int prId = 0; prId < nPRRecoTracks; ++prId) {
473 Eigen::RowVectorXd purityRow = purityMatrix.row(prId);
476 Purity highestPurity = purityRow.maxCoeff(&mcId);
478 mostPureMCId_by_prId[prId] = {mcId, highestPurity};
484 RecoTrackId mcId = -1;
485 B2DEBUG(24,
"MCTrack to highest weighted efficiency PRTrack relation");
486 for (
const auto& mostWeightEfficientPRId_for_mcId : mostWeightEfficientPRId_by_mcId) {
488 const Efficiency& weightedEfficiency = mostWeightEfficientPRId_for_mcId.weightedEfficiency;
489 const RecoTrackId& prId = mostWeightEfficientPRId_for_mcId.id;
491 "mcId : " << mcId <<
" -> prId : " << prId <<
" with weighted efficiency "
492 << weightedEfficiency);
499 RecoTrackId prId = -1;
500 B2DEBUG(24,
"PRTrack to highest purity MCTrack relation");
501 for (
const auto& mostPureMCId_for_prId : mostPureMCId_by_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);
510 int nMatched{}, nWrongCharge{}, nBackground{}, nClones{}, nClonesWrongCharge{}, nGhost{};
516 for (RecoTrackId prId = 0; prId < nPRRecoTracks; ++prId) {
519 const MostPureMCId& mostPureMCId = mostPureMCId_by_prId[prId];
521 const RecoTrackId& mcId = mostPureMCId.id;
522 const Purity& purity = mostPureMCId.purity;
528 B2DEBUG(23,
"Stored PRTrack " << prId <<
" as ghost because of too low purity");
534 if (mcId == mcBkgId) {
537 B2DEBUG(23,
"Stored PRTrack " << prId <<
" as background because of too low purity.");
550 B2ASSERT(
"No relation from MCRecoTrack to MCParticle.", mcParticle);
552 const MostWeightEfficientPRId& mostWeightEfficientPRId_for_mcId =
553 mostWeightEfficientPRId_by_mcId[mcId];
555 const RecoTrackId& mostWeightEfficientPRId = mostWeightEfficientPRId_for_mcId.id;
556 const Efficiency& weightedEfficiency = mostWeightEfficientPRId_for_mcId.weightedEfficiency;
560 const short MCParticleTrackCharge = mcParticle->
getCharge() > 0 ? 1 : -1;
564 short nPositiveCharges = 0;
565 short nNegativeCharges = 0;
569 if (fittedTracks.
size() > 0) {
571 for (
const auto& fittedTrack : fittedTracks) {
577 trackFitResult->
getChargeSign() > 0 ? nPositiveCharges++ : nNegativeCharges++;
580 trackFitResult->
getChargeSign() > 0 ? nPositiveCharges++ : nNegativeCharges++;
584 if (nPositiveCharges > 0 and nNegativeCharges > 0) {
586 "There are different charges attributed to the same track, this shouldn't happen. Continue with the majority of positive or negative charges");
590 if (fittedTracks.
size() > 0 and (nPositiveCharges > 0 or nNegativeCharges > 0)) {
591 foundTrackCharge = nPositiveCharges > nNegativeCharges ? 1 : -1;
599 if (prId == mostWeightEfficientPRId) {
600 if (foundTrackCharge != MCParticleTrackCharge) {
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);
626 B2DEBUG(23,
"Stored PRTrack " << prId <<
" as ghost because of too low efficiency.");
633 if (foundTrackCharge != MCParticleTrackCharge) {
635 ++nClonesWrongCharge;
646 B2DEBUG(23,
"Stored PRTrack " << prId <<
" as clone.");
647 B2DEBUG(23,
"MC Match prId " << prId <<
" to mcPartId " << mcParticle->
getArrayIndex());
648 B2DEBUG(23,
"Purity rel: prId " << prId <<
" -> mcId " << mcId <<
" : " << -purity);
653 B2DEBUG(23,
"Number of matches " << nMatched);
654 B2DEBUG(23,
"Number of wrongCharge " << nWrongCharge);
655 B2DEBUG(23,
"Number of clones " << nClones);
656 B2DEBUG(23,
"Number of clones wrongCharge " << nClonesWrongCharge);
657 B2DEBUG(23,
"Number of bkg " << nBackground);
658 B2DEBUG(23,
"Number of ghost " << nGhost);
662 for (RecoTrackId mcId = 0; mcId < nMCRecoTracks; ++mcId) {
666 const MostWeightEfficientPRId& mostWeightEfficiencyPRId = mostWeightEfficientPRId_by_mcId[mcId];
668 const RecoTrackId& prId = mostWeightEfficiencyPRId.id;
669 const Efficiency& weightedEfficiency = mostWeightEfficiencyPRId.weightedEfficiency;
672 B2ASSERT(
"Index of pattern recognition tracks out of range.", prId < nPRRecoTracks and prId >= 0);
676 const MostPureMCId& mostPureMCId_for_prId = mostPureMCId_by_prId[prId];
677 const RecoTrackId& mostPureMCId = mostPureMCId_for_prId.id;
680 if (mcId == mostPureMCId and
682 prRecoTrack->
getMatchingStatus() == RecoTrack::MatchingStatus::c_matchedWrongCharge or
684 prRecoTrack->
getMatchingStatus() == RecoTrack::MatchingStatus::c_cloneWrongCharge)) {
688 B2DEBUG(23,
"Efficiency rel: mcId " << mcId <<
" -> prId " << prId <<
" : " << weightedEfficiency);
696 bool isMergedMCRecoTrack =
699 prRecoTrack->
getMatchingStatus() == RecoTrack::MatchingStatus::c_matchedWrongCharge or
701 prRecoTrack->
getMatchingStatus() == RecoTrack::MatchingStatus::c_cloneWrongCharge);
703 if (isMergedMCRecoTrack) {
704 mcRecoTrack->
addRelationTo(prRecoTrack, -weightedEfficiency);
706 B2DEBUG(23,
"Efficiency rel: mcId " << mcId <<
" -> prId " << prId <<
" : " << -weightedEfficiency);
713 B2DEBUG(23,
"mcId " << mcId <<
" is missing. No relation created.");
714 B2DEBUG(23,
"is Primary " <<
m_MCRecoTracks[mcId]->getRelatedTo<MCParticle>()->isPrimaryParticle());
715 B2DEBUG(23,
"best prId " << prId <<
" with purity " << mostPureMCId_for_prId.purity <<
" -> " << mostPureMCId);
716 B2DEBUG(23,
"MC Total ndf " << totalNDF_by_mcId[mcId]);
717 B2DEBUG(23,
"MC Total weight " << totalWeight_by_mcId[mcId]);
718 B2DEBUG(23,
"MC Overlap ndf\n " << confusionMatrix.col(mcId).transpose());
719 B2DEBUG(23,
"MC Overlap weight\n " << weightedConfusionMatrix.col(mcId).transpose());
720 B2DEBUG(23,
"MC Efficiencies for the track\n" << efficiencyMatrix.col(mcId).transpose());
721 B2DEBUG(23,
"MC Weighted efficiencies for the track\n" << weightedEfficiencyMatrix.col(mcId).transpose());
724 B2DEBUG(23,
"########## End MCRecoTracksMatcherModule ############");
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
unsigned short getISuperLayer() const
Getter for iSuperLayer.
Provides a type-safe way to pass members of the chargedStableSet set.
EDetector
Enum for identifying the detector components (detector and subdetector).
static const ChargedStable pion
charged pion particle
A Class to store the Monte Carlo particle information.
int getArrayIndex() const
Get 0-based index of the particle in the corresponding MCParticle list.
float getCharge() const
Return the particle charge defined in TDatabasePDG.
int getPDG() const
Return PDG code of particle.
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 definition 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.
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
This is the Reconstruction Event-Data Model Track.
void setMatchingStatus(MatchingStatus matchingStatus)
Set the matching status (used by the TrackMatcher module)
MatchingStatus getMatchingStatus() const
Return the matching status set by the TrackMatcher module.
short int getChargeSeed() const
Return the charge seed stored in the reco track. ATTENTION: This is not the fitted charge.
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.
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.
Accessor to arrays stored in the data store.
int getEntries() const
Get the number of objects in the array.
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.
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.
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.