Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

MTD simulation: fix association tracking particle - MtdSimLayerCluster #47018

Merged
merged 2 commits into from
Jan 7, 2025
Merged
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Original file line number Diff line number Diff line change
Expand Up @@ -29,11 +29,13 @@ reco::SimToTPCollectionMtd MtdSimLayerClusterToTPAssociatorByTrackIdImpl::associ
std::map<std::pair<unsigned int, uint32_t>, TrackingParticleRef> tpIdMap;
for (auto tpIt = trackingParticles.begin(); tpIt != trackingParticles.end(); tpIt++) {
const auto& tp = *tpIt;
unsigned int tpTrackId = tp.g4Tracks()[0].trackId();
EncodedEventId tpEventId = tp.eventId();
TrackingParticleRef tpRef =
edm::Ref<TrackingParticleCollection>(trackingParticleH, tpIt - trackingParticles.begin());
tpIdMap[std::make_pair(tpTrackId, tpEventId.rawId())] = tpRef;
for (unsigned int igt = 0; igt < tp.g4Tracks().size(); igt++) {
unsigned int tpTrackId = tp.g4Tracks()[igt].trackId();
TrackingParticleRef tpRef =
edm::Ref<TrackingParticleCollection>(trackingParticleH, tpIt - trackingParticles.begin());
tpIdMap[std::make_pair(tpTrackId, tpEventId.rawId())] = tpRef;
}
}

// -- loop over sim clusters and get the trackId, eventId
Expand All @@ -45,31 +47,26 @@ reco::SimToTPCollectionMtd MtdSimLayerClusterToTPAssociatorByTrackIdImpl::associ
const auto& simClus = *simClusIt;
size_t simClusIndex = simClusIt - simClusters.begin();
MtdSimLayerClusterRef simClusterRef = edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIndex);
unsigned int simClusTrackId = simClus.g4Tracks()[0].trackId();
EncodedEventId simClusEventId = simClus.eventId();

// -- Check the trackId offset of the sim hits and keep only clusters with "direct" hits (offset == 0)
/*
if (simClus.trackIdOffset() != 0 ) continue;
*/

std::pair uniqueId = std::make_pair(simClusTrackId, simClusEventId.rawId());
auto it = tpIdMap.find(uniqueId);

if (it != tpIdMap.end()) {
TrackingParticleRef tpRef = tpIdMap[uniqueId];
outputCollection.insert(simClusterRef, tpRef);

LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
<< "MtdSimLayerCluster: index = " << simClusIndex << " simClus TrackId = " << simClusTrackId
<< " simClus EventId = " << simClusEventId.rawId() << " simClus Eta = " << simClus.eta()
<< " simClus Phi = " << simClus.phi() << " simClus Time = " << simClus.simLCTime()
<< " simClus Energy = " << simClus.simLCEnergy() << std::endl;
LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
<< " --> Found associated tracking particle: tp TrackId = " << (*tpRef).g4Tracks()[0].trackId()
<< " tp EventId = " << (*tpRef).eventId().rawId() << std::endl;
for (unsigned int igt = 0; igt < simClus.g4Tracks().size(); igt++) {
unsigned int simClusTrackId = simClus.g4Tracks()[igt].trackId();
std::pair uniqueId = std::make_pair(simClusTrackId, simClusEventId.rawId());
auto it = tpIdMap.find(uniqueId);

if (it != tpIdMap.end()) {
TrackingParticleRef tpRef = tpIdMap[uniqueId];
outputCollection.insert(simClusterRef, tpRef);

LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
<< "MtdSimLayerCluster: index = " << simClusIndex << " simClus TrackId = " << simClusTrackId
<< " simClus EventId = " << simClusEventId.rawId() << " simClus Eta = " << simClus.eta()
<< " simClus Phi = " << simClus.phi() << " simClus Time = " << simClus.simLCTime()
<< " simClus Energy = " << simClus.simLCEnergy() << std::endl;
LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
<< " --> Found associated tracking particle: tp TrackId = " << (*tpRef).g4Tracks()[0].trackId()
<< " tp EventId = " << (*tpRef).eventId().rawId() << std::endl;
}
}

} // -- end loop over sim clus

return outputCollection;
Expand All @@ -88,43 +85,42 @@ reco::TPToSimCollectionMtd MtdSimLayerClusterToTPAssociatorByTrackIdImpl::associ
std::map<std::pair<unsigned int, uint32_t>, std::vector<MtdSimLayerClusterRef>> simClusIdMap;
for (auto simClusIt = simClusters.begin(); simClusIt != simClusters.end(); simClusIt++) {
const auto& simClus = *simClusIt;
unsigned int simClusTrackId = simClus.g4Tracks()[0].trackId();
EncodedEventId simClusEventId = simClus.eventId();
MtdSimLayerClusterRef simClusterRef =
edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIt - simClusters.begin());
simClusIdMap[std::make_pair(simClusTrackId, simClusEventId.rawId())].push_back(simClusterRef);
for (unsigned int igt = 0; igt < simClus.g4Tracks().size(); igt++) {
unsigned int simClusTrackId = simClus.g4Tracks()[igt].trackId();
MtdSimLayerClusterRef simClusterRef =
edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIt - simClusters.begin());
simClusIdMap[std::make_pair(simClusTrackId, simClusEventId.rawId())].push_back(simClusterRef);
}
}

// -- Loop over the tracking particles
for (auto tpIt = trackingParticles.begin(); tpIt != trackingParticles.end(); tpIt++) {
const auto& tp = *tpIt;
size_t tpIndex = tpIt - trackingParticles.begin();
TrackingParticleRef tpRef = edm::Ref<TrackingParticleCollection>(trackingParticleH, tpIndex);
unsigned int tpTrackId = tp.g4Tracks()[0].trackId();
EncodedEventId tpEventId = tp.eventId();

std::pair uniqueId = std::make_pair(tpTrackId, tpEventId.rawId());
auto it = simClusIdMap.find(uniqueId);

if (it != simClusIdMap.end()) {
for (unsigned int i = 0; i < simClusIdMap[uniqueId].size(); i++) {
MtdSimLayerClusterRef simClusterRef = simClusIdMap[uniqueId][i];
// -- Check the trackId offset of the sim hits and keep only clusters with "direct" hits (offset == 0)
/*
if (simClus.trackIdOffset() != 0 ) continue;
*/

outputCollection.insert(tpRef, simClusterRef);

LogDebug("MtdSimLayerClusterToTPAssociator")
<< "Tracking particle: index = " << tpIndex << " tp TrackId = " << tpTrackId
<< " tp EventId = " << tpEventId.rawId();
LogDebug("MtdSimLayerClusterToTPAssociator")
<< " --> Found associated MtdSimLayerCluster: simClus TrackId = "
<< (*simClusterRef).g4Tracks()[0].trackId() << " simClus EventId = " << (*simClusterRef).eventId().rawId()
<< " simClus Eta = " << (*simClusterRef).eta() << " simClus Phi = " << (*simClusterRef).phi()
<< " simClus Time = " << (*simClusterRef).simLCTime()
<< " simClus Energy = " << (*simClusterRef).simLCEnergy() << std::endl;
for (unsigned int igt = 0; igt < tp.g4Tracks().size(); igt++) {
unsigned int tpTrackId = tp.g4Tracks()[igt].trackId();
std::pair uniqueId = std::make_pair(tpTrackId, tpEventId.rawId());
auto it = simClusIdMap.find(uniqueId);

if (it != simClusIdMap.end()) {
for (unsigned int i = 0; i < simClusIdMap[uniqueId].size(); i++) {
MtdSimLayerClusterRef simClusterRef = simClusIdMap[uniqueId][i];

outputCollection.insert(tpRef, simClusterRef);

LogDebug("MtdSimLayerClusterToTPAssociator")
<< "Tracking particle: index = " << tpIndex << " tp TrackId = " << tpTrackId
<< " tp EventId = " << tpEventId.rawId();
LogDebug("MtdSimLayerClusterToTPAssociator")
<< " --> Found associated MtdSimLayerCluster: simClus TrackId = "
<< (*simClusterRef).g4Tracks()[0].trackId() << " simClus EventId = " << (*simClusterRef).eventId().rawId()
<< " simClus Eta = " << (*simClusterRef).eta() << " simClus Phi = " << (*simClusterRef).phi()
<< " simClus Time = " << (*simClusterRef).simLCTime()
<< " simClus Energy = " << (*simClusterRef).simLCEnergy() << std::endl;
}
}
}
}
Expand Down