Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:51

0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/EventSetup.h"
0003 #include "SimDataFormats/Track/interface/SimTrack.h"
0004 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
0005 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0006 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0007 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0008 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0009 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0010 #include "DataFormats/JetReco/interface/GenJet.h"
0011 
0012 #include "L1Trigger/TrackFindingTMTT/interface/InputData.h"
0013 #include "L1Trigger/TrackFindingTMTT/interface/Settings.h"
0014 #include "L1Trigger/TrackFindingTMTT/interface/StubKiller.h"
0015 #include "L1Trigger/TrackFindingTMTT/interface/StubWindowSuggest.h"
0016 #include "L1Trigger/TrackFindingTMTT/interface/DegradeBend.h"
0017 
0018 #include <map>
0019 #include <memory>
0020 
0021 using namespace std;
0022 
0023 namespace tmtt {
0024 
0025   InputData::InputData(const edm::Event& iEvent,
0026                        const edm::EventSetup& iSetup,
0027                        const Settings* settings,
0028                        StubWindowSuggest* stubWindowSuggest,
0029                        const DegradeBend* degradeBend,
0030                        const TrackerGeometry* trackerGeometry,
0031                        const TrackerTopology* trackerTopology,
0032                        const list<TrackerModule>& listTrackerModule,
0033                        const edm::EDGetTokenT<TrackingParticleCollection> tpToken,
0034                        const edm::EDGetTokenT<TTStubDetSetVec> stubToken,
0035                        const edm::EDGetTokenT<TTStubAssMap> stubTruthToken,
0036                        const edm::EDGetTokenT<TTClusterAssMap> clusterTruthToken,
0037                        const edm::EDGetTokenT<reco::GenJetCollection> genJetToken)
0038       :  // Note if job will use MC truth info (or skip it to save CPU).
0039         enableMCtruth_(settings->enableMCtruth()) {
0040     edm::Handle<TrackingParticleCollection> tpHandle;
0041     edm::Handle<TTStubDetSetVec> ttStubHandle;
0042     edm::Handle<TTStubAssMap> mcTruthTTStubHandle;
0043     edm::Handle<TTClusterAssMap> mcTruthTTClusterHandle;
0044     edm::Handle<reco::GenJetCollection> genJetHandle;
0045     iEvent.getByToken(stubToken, ttStubHandle);
0046     if (enableMCtruth_) {
0047       iEvent.getByToken(tpToken, tpHandle);
0048       iEvent.getByToken(stubTruthToken, mcTruthTTStubHandle);
0049       iEvent.getByToken(clusterTruthToken, mcTruthTTClusterHandle);
0050       iEvent.getByToken(genJetToken, genJetHandle);
0051     }
0052 
0053     // Get TrackingParticle info
0054 
0055     if (enableMCtruth_) {
0056       unsigned int tpCount = 0;
0057       for (unsigned int i = 0; i < tpHandle->size(); i++) {
0058         const TrackingParticle& tPart = tpHandle->at(i);
0059         // Creating Ptr uses CPU, so apply Pt cut here, copied from TP::fillUse(), to avoid doing it too often.
0060         constexpr float ptMinScale = 0.7;
0061         const float ptMin = min(settings->genMinPt(), ptMinScale * settings->houghMinPt());
0062         if (tPart.pt() > ptMin) {
0063           TrackingParticlePtr tpPtr(tpHandle, i);
0064           // Store the TrackingParticle info, using class TP to provide easy access to the most useful info.
0065           TP tp(tpPtr, tpCount, settings);
0066           // Only bother storing tp if it could be useful for tracking efficiency or fake rate measurements.
0067           if (tp.use()) {
0068             if (genJetHandle.isValid()) {
0069               tp.fillNearestJetInfo(genJetHandle.product());
0070             }
0071 
0072             vTPs_.push_back(tp);
0073             tpCount++;
0074           }
0075         }
0076       }
0077     }
0078 
0079     // Also create map relating edm::Ptr<TrackingParticle> to TP.
0080 
0081     map<edm::Ptr<TrackingParticle>, const TP*> translateTP;
0082 
0083     if (enableMCtruth_) {
0084       for (const TP& tp : vTPs_) {
0085         const TrackingParticlePtr& tpPtr = tp.trackingParticlePtr();
0086         translateTP[tpPtr] = &tp;
0087       }
0088     }
0089 
0090     // Initialize code for killing some stubs to model detector problems.
0091     const StubKiller::KillOptions killOpt = static_cast<StubKiller::KillOptions>(settings->killScenario());
0092     std::unique_ptr<const StubKiller> stubKiller;
0093     if (killOpt != StubKiller::KillOptions::none) {
0094       stubKiller = std::make_unique<StubKiller>(killOpt, trackerTopology, trackerGeometry, iEvent);
0095     }
0096 
0097     // Loop over tracker modules to get module info & stubs.
0098 
0099     for (const TrackerModule& trackerModule : listTrackerModule) {
0100       const DetId& stackedDetId = trackerModule.stackedDetId();
0101       TTStubDetSetVec::const_iterator p_module = ttStubHandle->find(stackedDetId);
0102       if (p_module != ttStubHandle->end()) {
0103         for (TTStubDetSet::const_iterator p_ttstub = p_module->begin(); p_ttstub != p_module->end(); p_ttstub++) {
0104           TTStubRef ttStubRef = edmNew::makeRefTo(ttStubHandle, p_ttstub);
0105           const unsigned int stubIndex = vAllStubs_.size();
0106 
0107           // Store the Stub info, using class Stub to provide easy access to the most useful info.
0108           vAllStubs_.emplace_back(
0109               ttStubRef, stubIndex, settings, trackerTopology, &trackerModule, degradeBend, stubKiller.get());
0110 
0111           // Also fill truth associating stubs to tracking particles.
0112           if (enableMCtruth_) {
0113             Stub& stub = vAllStubs_.back();
0114             stub.fillTruth(translateTP, mcTruthTTStubHandle, mcTruthTTClusterHandle);
0115           }
0116         }
0117       }
0118     }
0119 
0120     // Produced reduced list containing only the subset of stubs that the user has declared will be
0121     // output by the front-end readout electronics.
0122     for (Stub& s : vAllStubs_) {
0123       if (s.frontendPass()) {
0124         vStubs_.push_back(&s);
0125         vStubsConst_.push_back(&s);
0126       }
0127     }
0128     // Optionally sort stubs according to bend, so highest Pt ones are sent from DTC to GP first.
0129     if (settings->orderStubsByBend()) {
0130       auto orderStubsByBend = [](const Stub* a, const Stub* b) { return (std::abs(a->bend()) < std::abs(b->bend())); };
0131       vStubs_.sort(orderStubsByBend);
0132     }
0133 
0134     // Note list of stubs produced by each tracking particle.
0135     // (By passing vAllStubs_ here instead of vStubs_, it means that any algorithmic efficiencies
0136     // measured will be reduced if the tightened frontend electronics cuts, specified in section StubCuts
0137     // of Analyze_Defaults_cfi.py, are not 100% efficient).
0138     if (enableMCtruth_) {
0139       for (TP& tp : vTPs_) {
0140         tp.fillTruth(vAllStubs_);
0141       }
0142     }
0143 
0144     // If requested, recommend better FE stub window cuts.
0145     if (settings->printStubWindows()) {
0146       for (const Stub& s : vAllStubs_) {
0147         stubWindowSuggest->process(trackerTopology, &s);
0148       }
0149     }
0150   }
0151 
0152 }  // namespace tmtt