Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:07

0001 /*! \brief   
0002  *  \details Here, in the source file, the methods which do depend
0003  *           on the specific type <T> that can fit the template.
0004  *
0005  *  \author Nicola Pozzobon
0006  *  \date   2013, Jul 19
0007  *
0008  */
0009 
0010 #include "SimTracker/TrackTriggerAssociation/plugins/TTClusterAssociator.h"
0011 
0012 /// Implement the producer
0013 template <>
0014 void TTClusterAssociator<Ref_Phase2TrackerDigi_>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0015   /// Exit if real data
0016   if (iEvent.isRealData())
0017     return;
0018 
0019   /// Get the PixelDigiSimLink
0020   iEvent.getByToken(digisimLinkToken_, thePixelDigiSimLinkHandle_);
0021 
0022   /// Get the TrackingParticles
0023 
0024   iEvent.getByToken(tpToken_, trackingParticleHandle_);
0025 
0026   //  const TrackerTopology* const tTopo = theTrackerTopology_.product();
0027   const TrackerGeometry* const theTrackerGeom = &iSetup.getData(theTrackerGeometryToken_);
0028 
0029   /// Preliminary task: map SimTracks by TrackingParticle
0030   /// Prepare the map
0031   std::map<std::pair<unsigned int, EncodedEventId>, TrackingParticlePtr> simTrackUniqueToTPMap;
0032 
0033   if (not trackingParticleHandle_->empty()) {
0034     /// Loop over TrackingParticles
0035     for (unsigned int tpCnt = 0; tpCnt < trackingParticleHandle_->size(); tpCnt++) {
0036       /// Make the pointer to the TrackingParticle
0037       TrackingParticlePtr tempTPPtr(trackingParticleHandle_, tpCnt);
0038 
0039       /// Get the EncodedEventId
0040       EncodedEventId eventId = EncodedEventId(tempTPPtr->eventId());
0041 
0042       /// Loop over SimTracks inside TrackingParticle
0043       for (const auto& simTrack : tempTPPtr->g4Tracks()) {
0044         /// Use the unique SimTrack Id (which is SimTrack ID + EncodedEventId)
0045         simTrackUniqueToTPMap.emplace(std::make_pair(simTrack.trackId(), eventId), tempTPPtr);
0046       }
0047     }  /// End of loop over TrackingParticles
0048   }
0049 
0050   /// Loop over InputTags to handle multiple collections
0051 
0052   int ncont1 = 0;
0053 
0054   for (const auto& iTag : ttClustersTokens_) {
0055     /// Prepare output
0056     auto associationMapForOutput = std::make_unique<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>>();
0057 
0058     /// Get the Clusters already stored away
0059     edm::Handle<TTClusterDetSetVec> TTClusterHandle;
0060 
0061     iEvent.getByToken(iTag, TTClusterHandle);
0062 
0063     /// Prepare the necessary maps
0064     std::map<TTClusterRef, std::vector<TrackingParticlePtr>> clusterToTrackingParticleVectorMap;
0065     std::map<TrackingParticlePtr, std::vector<TTClusterRef>> trackingParticleToClusterVectorMap;
0066 
0067     /// Loop over the input Clusters
0068     for (const auto& gd : theTrackerGeom->dets()) {
0069       DetId detid = gd->geographicalId();
0070       if (detid.subdetId() != StripSubdetector::TOB && detid.subdetId() != StripSubdetector::TID)
0071         continue;  // only run on OT
0072 
0073       if (TTClusterHandle->find(detid) == TTClusterHandle->end())
0074         continue;
0075 
0076       /// Get the DetSets of the Clusters
0077       edmNew::DetSet<TTCluster<Ref_Phase2TrackerDigi_>> clusters = (*TTClusterHandle)[detid];
0078 
0079       for (auto contentIter = clusters.begin(); contentIter != clusters.end(); ++contentIter) {
0080         /// Make the reference to be put in the map
0081         TTClusterRef tempCluRef = edmNew::makeRefTo(TTClusterHandle, contentIter);
0082 
0083         /// Prepare the maps wrt TTCluster
0084         if (clusterToTrackingParticleVectorMap.find(tempCluRef) == clusterToTrackingParticleVectorMap.end()) {
0085           std::vector<TrackingParticlePtr> tpVector;
0086           clusterToTrackingParticleVectorMap.emplace(tempCluRef, tpVector);
0087         }
0088 
0089         /// Get the PixelDigiSimLink
0090         /// Safety check added after new digitizer (Oct 2014)
0091         if (thePixelDigiSimLinkHandle_->find(detid) == thePixelDigiSimLinkHandle_->end()) {
0092           /// Sensor is not found in DigiSimLink.
0093           /// Set MC truth to NULL for all hits in this sensor. Period.
0094 
0095           /// Get the Digis and loop over them
0096           std::vector<Ref_Phase2TrackerDigi_> theseHits = tempCluRef->getHits();
0097           for (unsigned int i = 0; i < theseHits.size(); i++) {
0098             /// No SimLink is found by definition
0099             /// Then store NULL MC truth for all the digis
0100             TrackingParticlePtr tempTPPtr;
0101             clusterToTrackingParticleVectorMap.find(tempCluRef)->second.push_back(tempTPPtr);
0102           }
0103 
0104           /// Go to the next sensor
0105           continue;
0106         }
0107 
0108         edm::DetSet<PixelDigiSimLink> thisDigiSimLink = (*(thePixelDigiSimLinkHandle_))[detid];
0109         edm::DetSet<PixelDigiSimLink>::const_iterator iterSimLink;
0110 
0111         /// Get the Digis and loop over them
0112         std::vector<Ref_Phase2TrackerDigi_> theseHits = tempCluRef->getHits();
0113         for (unsigned int i = 0; i < theseHits.size(); i++) {
0114           /// Loop over PixelDigiSimLink
0115           for (iterSimLink = thisDigiSimLink.data.begin(); iterSimLink != thisDigiSimLink.data.end(); iterSimLink++) {
0116             /// Find the link and, if there's not, skip
0117             if (static_cast<int>(iterSimLink->channel()) != static_cast<int>(theseHits.at(i)->channel()))
0118               continue;
0119 
0120             /// Get SimTrack Id and type
0121             unsigned int curSimTrkId = iterSimLink->SimTrackId();
0122             EncodedEventId curSimEvId = iterSimLink->eventId();
0123 
0124             /// Prepare the SimTrack Unique ID
0125             std::pair<unsigned int, EncodedEventId> thisUniqueId = std::make_pair(curSimTrkId, curSimEvId);
0126 
0127             /// Get the corresponding TrackingParticle
0128             if (simTrackUniqueToTPMap.find(thisUniqueId) != simTrackUniqueToTPMap.end()) {
0129               TrackingParticlePtr thisTrackingParticle = simTrackUniqueToTPMap.find(thisUniqueId)->second;
0130 
0131               /// Store the TrackingParticle
0132               clusterToTrackingParticleVectorMap.find(tempCluRef)->second.push_back(thisTrackingParticle);
0133 
0134               /// Prepare the maps wrt TrackingParticle
0135               if (trackingParticleToClusterVectorMap.find(thisTrackingParticle) ==
0136                   trackingParticleToClusterVectorMap.end()) {
0137                 std::vector<TTClusterRef> clusterVector;
0138                 trackingParticleToClusterVectorMap.emplace(thisTrackingParticle, clusterVector);
0139               }
0140               trackingParticleToClusterVectorMap.find(thisTrackingParticle)
0141                   ->second.push_back(tempCluRef);  /// Fill the auxiliary map
0142             } else {
0143               /// In case no TrackingParticle is found, store a NULL pointer
0144 
0145               TrackingParticlePtr tempTPPtr;
0146               clusterToTrackingParticleVectorMap.find(tempCluRef)->second.push_back(tempTPPtr);
0147             }
0148           }  /// End of loop over PixelDigiSimLink
0149         }    /// End of loop over all the hits composing the Cluster
0150 
0151         /// Check that the cluster has a non-NULL TP pointer
0152         const std::vector<TrackingParticlePtr>& theseClusterTrackingParticlePtrs =
0153             clusterToTrackingParticleVectorMap.find(tempCluRef)->second;
0154         bool allOfThemAreNull = true;
0155         for (unsigned int tpi = 0; tpi < theseClusterTrackingParticlePtrs.size() && allOfThemAreNull; tpi++) {
0156           if (theseClusterTrackingParticlePtrs.at(tpi).isNull() == false)
0157             allOfThemAreNull = false;
0158         }
0159 
0160         if (allOfThemAreNull) {
0161           /// In case no TrackingParticle is found at all, drop the map element
0162           clusterToTrackingParticleVectorMap.erase(tempCluRef);  /// Use "erase by key"
0163         }
0164       }
0165     }  /// End of loop over all the TTClusters of the event
0166 
0167     /// Clean the maps that need cleaning
0168     /// Prepare the output map wrt TrackingParticle
0169     for (auto& p : trackingParticleToClusterVectorMap) {
0170       /// Get the vector of references to TTCluster
0171       std::vector<TTClusterRef>& tempVector = p.second;
0172 
0173       /// Sort and remove duplicates
0174       std::sort(tempVector.begin(), tempVector.end());
0175       tempVector.erase(std::unique(tempVector.begin(), tempVector.end()), tempVector.end());
0176     }
0177 
0178     /// Put the maps in the association object
0179     associationMapForOutput->setTTClusterToTrackingParticlesMap(clusterToTrackingParticleVectorMap);
0180     associationMapForOutput->setTrackingParticleToTTClustersMap(trackingParticleToClusterVectorMap);
0181 
0182     /// Put output in the event
0183     iEvent.put(std::move(associationMapForOutput), ttClustersInputTags_.at(ncont1).instance());
0184 
0185     ++ncont1;
0186 
0187   }  /// End of loop over input tags
0188 }