Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:47

0001 //
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "FWCore/Utilities/interface/Exception.h"
0004 
0005 #include "MtdSimLayerClusterToTPAssociatorByTrackIdImpl.h"
0006 
0007 using namespace reco;
0008 using namespace std;
0009 
0010 /* Constructor */
0011 MtdSimLayerClusterToTPAssociatorByTrackIdImpl::MtdSimLayerClusterToTPAssociatorByTrackIdImpl(
0012     edm::EDProductGetter const& productGetter)
0013     : productGetter_(&productGetter) {}
0014 
0015 //
0016 //---member functions
0017 //
0018 
0019 reco::SimToTPCollectionMtd MtdSimLayerClusterToTPAssociatorByTrackIdImpl::associateSimToTP(
0020     const edm::Handle<MtdSimLayerClusterCollection>& simClusH,
0021     const edm::Handle<TrackingParticleCollection>& trackingParticleH) const {
0022   SimToTPCollectionMtd outputCollection(productGetter_);
0023 
0024   // -- get the collections
0025   const auto& simClusters = *simClusH.product();
0026   const auto& trackingParticles = *trackingParticleH.product();
0027 
0028   // -- Loop over tracking particles and build a temporary map of trackId, eventId  --> tpRef
0029   std::map<std::pair<unsigned int, uint32_t>, TrackingParticleRef> tpIdMap;
0030   for (auto tpIt = trackingParticles.begin(); tpIt != trackingParticles.end(); tpIt++) {
0031     const auto& tp = *tpIt;
0032     unsigned int tpTrackId = tp.g4Tracks()[0].trackId();
0033     EncodedEventId tpEventId = tp.eventId();
0034     TrackingParticleRef tpRef =
0035         edm::Ref<TrackingParticleCollection>(trackingParticleH, tpIt - trackingParticles.begin());
0036     tpIdMap[std::make_pair(tpTrackId, tpEventId.rawId())] = tpRef;
0037   }
0038 
0039   // -- loop over sim clusters and get the trackId, eventId
0040 
0041   LogDebug("MtdSimLayerClusterToTPAssociator")
0042       << " Found " << simClusters.size() << " MtdSimLayerClusters in the event";
0043 
0044   for (auto simClusIt = simClusters.begin(); simClusIt != simClusters.end(); simClusIt++) {
0045     const auto& simClus = *simClusIt;
0046     size_t simClusIndex = simClusIt - simClusters.begin();
0047     MtdSimLayerClusterRef simClusterRef = edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIndex);
0048     unsigned int simClusTrackId = simClus.g4Tracks()[0].trackId();
0049     EncodedEventId simClusEventId = simClus.eventId();
0050 
0051     // -- Check the trackId offset of the sim hits and keep only clusters with "direct" hits (offset == 0)
0052     /* 
0053        if (simClus.trackIdOffset() != 0 ) continue; 
0054     */
0055 
0056     std::pair uniqueId = std::make_pair(simClusTrackId, simClusEventId.rawId());
0057     auto it = tpIdMap.find(uniqueId);
0058 
0059     if (it != tpIdMap.end()) {
0060       TrackingParticleRef tpRef = tpIdMap[uniqueId];
0061       outputCollection.insert(simClusterRef, tpRef);
0062 
0063       LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
0064           << "MtdSimLayerCluster: index = " << simClusIndex << "   simClus TrackId = " << simClusTrackId
0065           << " simClus EventId = " << simClusEventId.rawId() << " simClus Eta = " << simClus.eta()
0066           << " simClus Phi = " << simClus.phi() << "  simClus Time = " << simClus.simLCTime()
0067           << "  simClus Energy = " << simClus.simLCEnergy() << std::endl;
0068       LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
0069           << "  --> Found associated tracking particle:  tp TrackId = " << (*tpRef).g4Tracks()[0].trackId()
0070           << " tp EventId = " << (*tpRef).eventId().rawId() << std::endl;
0071     }
0072 
0073   }  // -- end loop over sim clus
0074 
0075   return outputCollection;
0076 }
0077 
0078 reco::TPToSimCollectionMtd MtdSimLayerClusterToTPAssociatorByTrackIdImpl::associateTPToSim(
0079     const edm::Handle<MtdSimLayerClusterCollection>& simClusH,
0080     const edm::Handle<TrackingParticleCollection>& trackingParticleH) const {
0081   TPToSimCollectionMtd outputCollection(productGetter_);
0082 
0083   // -- get the collections
0084   const auto& simClusters = *simClusH.product();
0085   const auto& trackingParticles = *trackingParticleH.product();
0086 
0087   // -- Loop over MtdSimLayerClusters and build a temporary map of trackId, eventId --> simClusterRef
0088   std::map<std::pair<unsigned int, uint32_t>, std::vector<MtdSimLayerClusterRef>> simClusIdMap;
0089   for (auto simClusIt = simClusters.begin(); simClusIt != simClusters.end(); simClusIt++) {
0090     const auto& simClus = *simClusIt;
0091     unsigned int simClusTrackId = simClus.g4Tracks()[0].trackId();
0092     EncodedEventId simClusEventId = simClus.eventId();
0093     MtdSimLayerClusterRef simClusterRef =
0094         edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIt - simClusters.begin());
0095     simClusIdMap[std::make_pair(simClusTrackId, simClusEventId.rawId())].push_back(simClusterRef);
0096   }
0097 
0098   // -- Loop over the tracking particles
0099   for (auto tpIt = trackingParticles.begin(); tpIt != trackingParticles.end(); tpIt++) {
0100     const auto& tp = *tpIt;
0101     size_t tpIndex = tpIt - trackingParticles.begin();
0102     TrackingParticleRef tpRef = edm::Ref<TrackingParticleCollection>(trackingParticleH, tpIndex);
0103     unsigned int tpTrackId = tp.g4Tracks()[0].trackId();
0104     EncodedEventId tpEventId = tp.eventId();
0105 
0106     std::pair uniqueId = std::make_pair(tpTrackId, tpEventId.rawId());
0107     auto it = simClusIdMap.find(uniqueId);
0108 
0109     if (it != simClusIdMap.end()) {
0110       for (unsigned int i = 0; i < simClusIdMap[uniqueId].size(); i++) {
0111         MtdSimLayerClusterRef simClusterRef = simClusIdMap[uniqueId][i];
0112         // -- Check the trackId offset of the sim hits and keep only clusters with "direct" hits (offset == 0)
0113         /*
0114       if (simClus.trackIdOffset() != 0 ) continue; 
0115     */
0116 
0117         outputCollection.insert(tpRef, simClusterRef);
0118 
0119         LogDebug("MtdSimLayerClusterToTPAssociator")
0120             << "Tracking particle:  index = " << tpIndex << "  tp TrackId = " << tpTrackId
0121             << "  tp EventId = " << tpEventId.rawId();
0122         LogDebug("MtdSimLayerClusterToTPAssociator")
0123             << " --> Found associated MtdSimLayerCluster:  simClus TrackId = "
0124             << (*simClusterRef).g4Tracks()[0].trackId() << " simClus EventId = " << (*simClusterRef).eventId().rawId()
0125             << " simClus Eta = " << (*simClusterRef).eta() << " simClus Phi = " << (*simClusterRef).phi()
0126             << "  simClus Time = " << (*simClusterRef).simLCTime()
0127             << "  simClus Energy = " << (*simClusterRef).simLCEnergy() << std::endl;
0128       }
0129     }
0130   }
0131 
0132   return outputCollection;
0133 }