Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:38:11

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/TTStubAssociator.h"
0011 
0012 /// Implement the producer
0013 template <>
0014 void TTStubAssociator<Ref_Phase2TrackerDigi_>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0015   /// Exit if real data
0016   if (iEvent.isRealData())
0017     return;
0018 
0019   /// Exit if the vectors are uncorrectly dimensioned
0020   if (ttClusterTruthInputTags_.size() != ttStubsInputTags_.size()) {
0021     edm::LogError("TTStubAsso ") << "E R R O R! the InputTag vectors have different size!";
0022     return;
0023   }
0024 
0025   int ncont1 = 0;
0026 
0027   const TrackerGeometry* const theTrackerGeom = &iSetup.getData(theTrackerGeometryToken_);
0028   const TrackerTopology* const tTopo = &iSetup.getData(theTrackerTopologyToken_);
0029 
0030   /// Loop over the InputTags to handle multiple collections
0031 
0032   for (const auto& iTag : ttStubsTokens_) {
0033     /// Prepare output
0034     auto associationMapForOutput = std::make_unique<TTStubAssociationMap<Ref_Phase2TrackerDigi_>>();
0035 
0036     /// Get the Stubs already stored away
0037     edm::Handle<TTStubDetSetVec> ttStubHandle;
0038     iEvent.getByToken(iTag, ttStubHandle);
0039 
0040     /// Get the Cluster MC truth
0041     edm::Handle<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> ttClusterAssociationMapHandle;
0042     iEvent.getByToken(ttClusterTruthTokens_.at(ncont1), ttClusterAssociationMapHandle);
0043 
0044     /// Prepare the necessary maps
0045     std::map<TTStubRef, TrackingParticlePtr> stubToTrackingParticleMap;
0046     std::map<TrackingParticlePtr, std::vector<TTStubRef>> trackingParticleToStubVectorMap;
0047 
0048     /// Loop over the input Stubs
0049 
0050     if (not ttStubHandle->empty()) {
0051       for (const auto& gd : theTrackerGeom->dets()) {
0052         DetId detid = gd->geographicalId();
0053         if (detid.subdetId() != StripSubdetector::TOB && detid.subdetId() != StripSubdetector::TID)
0054           continue;  // only run on OT
0055 
0056         if (!tTopo->isLower(detid))
0057           continue;  // loop on the stacks: choose the lower arbitrarily
0058 
0059         DetId stackDetid = tTopo->stack(detid);  // Stub module detid
0060 
0061         if (ttStubHandle->find(stackDetid) == ttStubHandle->end())
0062           continue;
0063 
0064         /// Get the DetSets of the Clusters
0065         edmNew::DetSet<TTStub<Ref_Phase2TrackerDigi_>> stubs = (*ttStubHandle)[stackDetid];
0066 
0067         for (auto contentIter = stubs.begin(); contentIter != stubs.end(); ++contentIter) {
0068           /// Make the reference to be put in the map
0069           TTStubRef tempStubRef = edmNew::makeRefTo(ttStubHandle, contentIter);
0070 
0071           /// Fill the inclusive map which is careless of the stub classification
0072           for (unsigned int ic = 0; ic < 2; ic++) {
0073             const std::vector<TrackingParticlePtr>& tempTPs =
0074                 ttClusterAssociationMapHandle->findTrackingParticlePtrs(tempStubRef->clusterRef(ic));
0075 
0076             for (const TrackingParticlePtr& testTP : tempTPs) {
0077               if (testTP.isNull())
0078                 continue;
0079 
0080               /// Prepare the maps wrt TrackingParticle
0081               if (trackingParticleToStubVectorMap.find(testTP) == trackingParticleToStubVectorMap.end()) {
0082                 std::vector<TTStubRef> stubVector;
0083                 trackingParticleToStubVectorMap.emplace(testTP, stubVector);
0084               }
0085               trackingParticleToStubVectorMap.find(testTP)->second.push_back(tempStubRef);  /// Fill the auxiliary map
0086             }
0087           }
0088 
0089           /// GENUINE for clusters means not combinatoric and
0090           /// not unknown: same MC truth content MUST be found
0091           /// in both clusters composing the stub
0092           if (ttClusterAssociationMapHandle->isUnknown(tempStubRef->clusterRef(0)) ||
0093               ttClusterAssociationMapHandle->isUnknown(tempStubRef->clusterRef(1))) {
0094             /// If at least one cluster is unknown, it means
0095             /// either unknown, either combinatoric
0096             /// Do nothing, and go to the next Stub
0097             continue;
0098           } else {
0099             /// Here both are clusters are genuine/combinatoric
0100             /// If both clusters have some known SimTrack content
0101             /// they must be compared to each other
0102             if (ttClusterAssociationMapHandle->isGenuine(tempStubRef->clusterRef(0)) &&
0103                 ttClusterAssociationMapHandle->isGenuine(tempStubRef->clusterRef(1))) {
0104               /// If both clusters are genuine, they must be associated to the same TrackingParticle
0105               /// in order to return a genuine stub. Period. Note we can perform safely
0106               /// this comparison because, if both clusters are genuine, their TrackingParticle shall NEVER be NULL
0107               if (ttClusterAssociationMapHandle->findTrackingParticlePtr(tempStubRef->clusterRef(0)).get() ==
0108                   ttClusterAssociationMapHandle->findTrackingParticlePtr(tempStubRef->clusterRef(1)).get()) {
0109                 /// Two genuine clusters with same SimTrack content mean genuine
0110                 const TrackingParticlePtr& testTP =
0111                     ttClusterAssociationMapHandle->findTrackingParticlePtr(tempStubRef->clusterRef(0));
0112 
0113                 /// Fill the map: by construction, this will always be the first time the
0114                 /// stub is inserted into the map: no need for "find"
0115                 stubToTrackingParticleMap.emplace(tempStubRef, testTP);
0116 
0117                 /// At this point, go to the next Stub
0118                 continue;
0119               } else {
0120                 /// It means combinatoric
0121                 continue;
0122               }
0123             }  /// End of two genuine clusters
0124             else {
0125               /// Here, at least one cluster is combinatoric
0126               TrackingParticle* prevTPAddress = nullptr;
0127               unsigned int whichTP = 0;
0128 
0129               const std::vector<TrackingParticlePtr>& trackingParticles0 =
0130                   ttClusterAssociationMapHandle->findTrackingParticlePtrs(tempStubRef->clusterRef(0));
0131               const std::vector<TrackingParticlePtr>& trackingParticles1 =
0132                   ttClusterAssociationMapHandle->findTrackingParticlePtrs(tempStubRef->clusterRef(1));
0133 
0134               bool escape = false;
0135 
0136               for (unsigned int i = 0; i < trackingParticles0.size() && !escape; i++) {
0137                 /// Skip NULL pointers
0138                 if (trackingParticles0.at(i).isNull())
0139                   continue;
0140 
0141                 for (unsigned int k = 0; k < trackingParticles1.size() && !escape; k++) {
0142                   /// Skip NULL pointers
0143                   if (trackingParticles1.at(k).isNull())
0144                     continue;
0145 
0146                   if (trackingParticles0.at(i).get() == trackingParticles1.at(k).get()) {
0147                     /// Same SimTrack is present in both clusters
0148                     if (prevTPAddress == nullptr) {
0149                       prevTPAddress = const_cast<TrackingParticle*>(trackingParticles1.at(k).get());
0150                       whichTP = k;
0151                     }
0152 
0153                     /// If two different SimTracks are found in both clusters,
0154                     /// then the stub is for sure combinatoric
0155                     if (prevTPAddress != const_cast<TrackingParticle*>(trackingParticles1.at(k).get())) {
0156                       escape = true;
0157                       continue;
0158                     }
0159                   }
0160                 }
0161               }  /// End of double loop over SimTracks of both clusters
0162 
0163               /// If two different SimTracks are found in both clusters,
0164               /// go to the next stub
0165               if (escape)
0166                 continue;
0167 
0168               if (prevTPAddress == nullptr) {
0169                 /// No SimTracks were found to be in both clusters
0170                 continue;
0171               } else {
0172                 /// Only one SimTrack was found to be present in both clusters
0173                 /// even if one of the clusters (or both) are combinatoric:
0174                 /// this means there is only one track that participates in
0175                 /// both clusters, hence the stub is genuine
0176                 TrackingParticlePtr testTP = trackingParticles1.at(whichTP);
0177 
0178                 /// Fill the map: by construction, this will always be the first time the
0179                 /// stub is inserted into the map: no need for "find"
0180                 stubToTrackingParticleMap.emplace(tempStubRef, testTP);
0181 
0182                 /// At this point, go to the next Stub
0183                 continue;
0184               }  /// End of one single SimTrack in both clusters
0185             }  /// End of "at least one cluster is combinatoric"
0186           }  /// End of "both clusters are known, somehow..."
0187         }
0188       }  /// End of loop over Stubs
0189     }
0190 
0191     /// Clean the only map that needs cleaning
0192     /// Prepare the output map wrt TrackingParticle
0193     for (auto& p : trackingParticleToStubVectorMap) {
0194       /// Get the vector of references to TTStub
0195       std::vector<TTStubRef>& tempVector = p.second;
0196 
0197       /// Sort and remove duplicates
0198       std::sort(tempVector.begin(), tempVector.end());
0199       tempVector.erase(std::unique(tempVector.begin(), tempVector.end()), tempVector.end());
0200     }
0201 
0202     /// Also, create the pointer to the TTClusterAssociationMap
0203     edm::RefProd<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> theCluAssoMap(ttClusterAssociationMapHandle);
0204 
0205     /// Put the maps in the association object
0206     associationMapForOutput->setTTStubToTrackingParticleMap(stubToTrackingParticleMap);
0207     associationMapForOutput->setTrackingParticleToTTStubsMap(trackingParticleToStubVectorMap);
0208     associationMapForOutput->setTTClusterAssociationMap(theCluAssoMap);
0209 
0210     /// Put output in the event
0211     iEvent.put(std::move(associationMapForOutput), ttStubsInputTags_.at(ncont1).instance());
0212 
0213     ++ncont1;
0214   }  /// End of loop over InputTags
0215 }