Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-08 03:36:35

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     EncodedEventId tpEventId = tp.eventId();
0033     for (unsigned int igt = 0; igt < tp.g4Tracks().size(); igt++) {
0034       unsigned int tpTrackId = tp.g4Tracks()[igt].trackId();
0035       TrackingParticleRef tpRef =
0036           edm::Ref<TrackingParticleCollection>(trackingParticleH, tpIt - trackingParticles.begin());
0037       tpIdMap[std::make_pair(tpTrackId, tpEventId.rawId())] = tpRef;
0038     }
0039   }
0040 
0041   // -- loop over sim clusters and get the trackId, eventId
0042 
0043   LogDebug("MtdSimLayerClusterToTPAssociator")
0044       << " Found " << simClusters.size() << " MtdSimLayerClusters in the event";
0045 
0046   for (auto simClusIt = simClusters.begin(); simClusIt != simClusters.end(); simClusIt++) {
0047     const auto& simClus = *simClusIt;
0048     size_t simClusIndex = simClusIt - simClusters.begin();
0049     MtdSimLayerClusterRef simClusterRef = edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIndex);
0050     EncodedEventId simClusEventId = simClus.eventId();
0051     for (unsigned int igt = 0; igt < simClus.g4Tracks().size(); igt++) {
0052       unsigned int simClusTrackId = simClus.g4Tracks()[igt].trackId();
0053       std::pair uniqueId = std::make_pair(simClusTrackId, simClusEventId.rawId());
0054       auto it = tpIdMap.find(uniqueId);
0055 
0056       if (it != tpIdMap.end()) {
0057         TrackingParticleRef tpRef = tpIdMap[uniqueId];
0058         outputCollection.insert(simClusterRef, tpRef);
0059 
0060         LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
0061             << "MtdSimLayerCluster: index = " << simClusIndex << "   simClus TrackId = " << simClusTrackId
0062             << " simClus EventId = " << simClusEventId.rawId() << " simClus Eta = " << simClus.eta()
0063             << " simClus Phi = " << simClus.phi() << "  simClus Time = " << simClus.simLCTime()
0064             << "  simClus Energy = " << simClus.simLCEnergy() << std::endl;
0065         LogDebug("MtdSimLayerClusterToTPAssociator::associateSimToTP")
0066             << "  --> Found associated tracking particle:  tp TrackId = " << (*tpRef).g4Tracks()[0].trackId()
0067             << " tp EventId = " << (*tpRef).eventId().rawId() << std::endl;
0068       }
0069     }
0070   }  // -- end loop over sim clus
0071 
0072   return outputCollection;
0073 }
0074 
0075 reco::TPToSimCollectionMtd MtdSimLayerClusterToTPAssociatorByTrackIdImpl::associateTPToSim(
0076     const edm::Handle<MtdSimLayerClusterCollection>& simClusH,
0077     const edm::Handle<TrackingParticleCollection>& trackingParticleH) const {
0078   TPToSimCollectionMtd outputCollection(productGetter_);
0079 
0080   // -- get the collections
0081   const auto& simClusters = *simClusH.product();
0082   const auto& trackingParticles = *trackingParticleH.product();
0083 
0084   // -- Loop over MtdSimLayerClusters and build a temporary map of trackId, eventId --> simClusterRef
0085   std::map<std::pair<unsigned int, uint32_t>, std::vector<MtdSimLayerClusterRef>> simClusIdMap;
0086   for (auto simClusIt = simClusters.begin(); simClusIt != simClusters.end(); simClusIt++) {
0087     const auto& simClus = *simClusIt;
0088     EncodedEventId simClusEventId = simClus.eventId();
0089     for (unsigned int igt = 0; igt < simClus.g4Tracks().size(); igt++) {
0090       unsigned int simClusTrackId = simClus.g4Tracks()[igt].trackId();
0091       MtdSimLayerClusterRef simClusterRef =
0092           edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIt - simClusters.begin());
0093       simClusIdMap[std::make_pair(simClusTrackId, simClusEventId.rawId())].push_back(simClusterRef);
0094     }
0095   }
0096 
0097   // -- Loop over the tracking particles
0098   for (auto tpIt = trackingParticles.begin(); tpIt != trackingParticles.end(); tpIt++) {
0099     const auto& tp = *tpIt;
0100     size_t tpIndex = tpIt - trackingParticles.begin();
0101     TrackingParticleRef tpRef = edm::Ref<TrackingParticleCollection>(trackingParticleH, tpIndex);
0102     EncodedEventId tpEventId = tp.eventId();
0103     for (unsigned int igt = 0; igt < tp.g4Tracks().size(); igt++) {
0104       unsigned int tpTrackId = tp.g4Tracks()[igt].trackId();
0105       std::pair uniqueId = std::make_pair(tpTrackId, tpEventId.rawId());
0106       auto it = simClusIdMap.find(uniqueId);
0107 
0108       if (it != simClusIdMap.end()) {
0109         for (unsigned int i = 0; i < simClusIdMap[uniqueId].size(); i++) {
0110           MtdSimLayerClusterRef simClusterRef = simClusIdMap[uniqueId][i];
0111 
0112           outputCollection.insert(tpRef, simClusterRef);
0113 
0114           LogDebug("MtdSimLayerClusterToTPAssociator")
0115               << "Tracking particle:  index = " << tpIndex << "  tp TrackId = " << tpTrackId
0116               << "  tp EventId = " << tpEventId.rawId();
0117           LogDebug("MtdSimLayerClusterToTPAssociator")
0118               << " --> Found associated MtdSimLayerCluster:  simClus TrackId = "
0119               << (*simClusterRef).g4Tracks()[0].trackId() << " simClus EventId = " << (*simClusterRef).eventId().rawId()
0120               << " simClus Eta = " << (*simClusterRef).eta() << " simClus Phi = " << (*simClusterRef).phi()
0121               << "  simClus Time = " << (*simClusterRef).simLCTime()
0122               << "  simClus Energy = " << (*simClusterRef).simLCEnergy() << std::endl;
0123         }
0124       }
0125     }
0126   }
0127 
0128   return outputCollection;
0129 }