Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-03 00:12:25

0001 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0002 #include "FWCore/Framework/interface/Run.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/ServiceRegistry/interface/Service.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/Utilities/interface/EDGetToken.h"
0010 #include "FWCore/Utilities/interface/InputTag.h"
0011 #include "FWCore/Utilities/interface/Exception.h"
0012 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0013 #include "DataFormats/Common/interface/Handle.h"
0014 
0015 #include "SimTracker/TrackTriggerAssociation/interface/StubAssociation.h"
0016 #include "L1Trigger/TrackTrigger/interface/Setup.h"
0017 #include "L1Trigger/TrackFindingTracklet/interface/DataFormats.h"
0018 #include "L1Trigger/TrackFindingTracklet/interface/ChannelAssignment.h"
0019 
0020 #include <TProfile.h>
0021 #include <TH1F.h>
0022 
0023 #include <vector>
0024 #include <deque>
0025 #include <map>
0026 #include <set>
0027 #include <cmath>
0028 #include <numeric>
0029 #include <sstream>
0030 
0031 namespace trklet {
0032 
0033   /*! \class  trklet::AnalyzerTM
0034    *  \brief  Class to analyze hardware like structured TTStub Collection generated by TM module
0035    *  \author Thomas Schuh
0036    *  \date   2023, Jan
0037    */
0038   class AnalyzerTM : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0039   public:
0040     AnalyzerTM(const edm::ParameterSet& iConfig);
0041     void beginJob() override {}
0042     void beginRun(const edm::Run& iEvent, const edm::EventSetup& iSetup) override;
0043     void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0044     void endRun(const edm::Run& iEvent, const edm::EventSetup& iSetup) override {}
0045     void endJob() override;
0046 
0047   private:
0048     //
0049     void formTracks(const tt::StreamsTrack& streamsTrack,
0050                     const tt::StreamsStub& streamsStubs,
0051                     std::vector<std::vector<TTStubRef>>& tracks,
0052                     int channel) const;
0053     //
0054     void associate(const std::vector<std::vector<TTStubRef>>& tracks,
0055                    const tt::StubAssociation* ass,
0056                    std::set<TPPtr>& tps,
0057                    int& sum,
0058                    bool perfect = false) const;
0059 
0060     // ED input token of stubs
0061     edm::EDGetTokenT<tt::StreamsStub> edGetTokenStubs_;
0062     // ED input token of tracks
0063     edm::EDGetTokenT<tt::StreamsTrack> edGetTokenTracks_;
0064     // ED input token of TTStubRef to TPPtr association for tracking efficiency
0065     edm::EDGetTokenT<tt::StubAssociation> edGetTokenSelection_;
0066     // ED input token of TTStubRef to recontructable TPPtr association
0067     edm::EDGetTokenT<tt::StubAssociation> edGetTokenReconstructable_;
0068     // Setup token
0069     edm::ESGetToken<tt::Setup, tt::SetupRcd> esGetTokenSetup_;
0070     // DataFormats token
0071     edm::ESGetToken<DataFormats, ChannelAssignmentRcd> esGetTokenDataFormats_;
0072     // ChannelAssignment token
0073     edm::ESGetToken<ChannelAssignment, ChannelAssignmentRcd> esGetTokenChannelAssignment_;
0074     // stores, calculates and provides run-time constants
0075     const tt::Setup* setup_ = nullptr;
0076     // helper class to extract structured data from tt::Frames
0077     const DataFormats* dataFormats_ = nullptr;
0078     // helper class to assign tracklet track to channel
0079     const ChannelAssignment* channelAssignment_ = nullptr;
0080     // enables analyze of TPs
0081     bool useMCTruth_;
0082     //
0083     int nEvents_ = 0;
0084 
0085     // Histograms
0086 
0087     TProfile* prof_;
0088     TProfile* profChannel_;
0089     TH1F* hisChannel_;
0090 
0091     // printout
0092     std::stringstream log_;
0093   };
0094 
0095   AnalyzerTM::AnalyzerTM(const edm::ParameterSet& iConfig) : useMCTruth_(iConfig.getParameter<bool>("UseMCTruth")) {
0096     usesResource("TFileService");
0097     // book in- and output ED products
0098     const std::string& label = iConfig.getParameter<std::string>("OutputLabelTM");
0099     const std::string& branchStubs = iConfig.getParameter<std::string>("BranchStubs");
0100     const std::string& branchTracks = iConfig.getParameter<std::string>("BranchTracks");
0101     edGetTokenStubs_ = consumes<tt::StreamsStub>(edm::InputTag(label, branchStubs));
0102     edGetTokenTracks_ = consumes<tt::StreamsTrack>(edm::InputTag(label, branchTracks));
0103     if (useMCTruth_) {
0104       const auto& inputTagSelecttion = iConfig.getParameter<edm::InputTag>("InputTagSelection");
0105       const auto& inputTagReconstructable = iConfig.getParameter<edm::InputTag>("InputTagReconstructable");
0106       edGetTokenSelection_ = consumes<tt::StubAssociation>(inputTagSelecttion);
0107       edGetTokenReconstructable_ = consumes<tt::StubAssociation>(inputTagReconstructable);
0108     }
0109     // book ES products
0110     esGetTokenSetup_ = esConsumes<edm::Transition::BeginRun>();
0111     esGetTokenDataFormats_ = esConsumes<edm::Transition::BeginRun>();
0112     esGetTokenChannelAssignment_ = esConsumes<edm::Transition::BeginRun>();
0113     // log config
0114     log_.setf(std::ios::fixed, std::ios::floatfield);
0115     log_.precision(4);
0116   }
0117 
0118   void AnalyzerTM::beginRun(const edm::Run& iEvent, const edm::EventSetup& iSetup) {
0119     // helper class to store configurations
0120     setup_ = &iSetup.getData(esGetTokenSetup_);
0121     // helper class to extract structured data from tt::Frames
0122     dataFormats_ = &iSetup.getData(esGetTokenDataFormats_);
0123     // helper class to assign tracklet track to channel
0124     channelAssignment_ = &iSetup.getData(esGetTokenChannelAssignment_);
0125     // book histograms
0126     edm::Service<TFileService> fs;
0127     TFileDirectory dir;
0128     dir = fs->mkdir("TM");
0129     prof_ = dir.make<TProfile>("Counts", ";", 10, 0.5, 10.5);
0130     prof_->GetXaxis()->SetBinLabel(1, "Stubs");
0131     prof_->GetXaxis()->SetBinLabel(2, "Tracks");
0132     prof_->GetXaxis()->SetBinLabel(4, "Matched Tracks");
0133     prof_->GetXaxis()->SetBinLabel(5, "All Tracks");
0134     prof_->GetXaxis()->SetBinLabel(6, "Found TPs");
0135     prof_->GetXaxis()->SetBinLabel(7, "Found selected TPs");
0136     prof_->GetXaxis()->SetBinLabel(9, "All TPs");
0137     prof_->GetXaxis()->SetBinLabel(10, "Perfect TPs");
0138     // channel occupancy
0139     constexpr int maxOcc = 180;
0140     const int numChannels = 1;
0141     hisChannel_ = dir.make<TH1F>("His Channel Occupancy", ";", maxOcc, -.5, maxOcc - .5);
0142     profChannel_ = dir.make<TProfile>("Prof Channel Occupancy", ";", numChannels, -.5, numChannels - .5);
0143   }
0144 
0145   void AnalyzerTM::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0146     // read in ht products
0147     edm::Handle<tt::StreamsStub> handleStubs;
0148     iEvent.getByToken<tt::StreamsStub>(edGetTokenStubs_, handleStubs);
0149     const tt::StreamsStub& streamsStub = *handleStubs;
0150     edm::Handle<tt::StreamsTrack> handleTracks;
0151     iEvent.getByToken<tt::StreamsTrack>(edGetTokenTracks_, handleTracks);
0152     const tt::StreamsTrack& streamsTrack = *handleTracks;
0153     // read in MCTruth
0154     const tt::StubAssociation* selection = nullptr;
0155     const tt::StubAssociation* reconstructable = nullptr;
0156     if (useMCTruth_) {
0157       edm::Handle<tt::StubAssociation> handleSelection;
0158       iEvent.getByToken<tt::StubAssociation>(edGetTokenSelection_, handleSelection);
0159       selection = handleSelection.product();
0160       prof_->Fill(9, selection->numTPs());
0161       edm::Handle<tt::StubAssociation> handleReconstructable;
0162       iEvent.getByToken<tt::StubAssociation>(edGetTokenReconstructable_, handleReconstructable);
0163       reconstructable = handleReconstructable.product();
0164     }
0165     // analyze ht products and associate found tracks with reconstrucable TrackingParticles
0166     std::set<TPPtr> tpPtrs;
0167     std::set<TPPtr> tpPtrsSelection;
0168     std::set<TPPtr> tpPtrsPerfect;
0169     int allMatched(0);
0170     int allTracks(0);
0171     for (int region = 0; region < setup_->numRegions(); region++) {
0172       int nStubs(0);
0173       int nTracks(0);
0174       std::vector<std::vector<TTStubRef>> tracks;
0175       formTracks(streamsTrack, streamsStub, tracks, region);
0176       nTracks += tracks.size();
0177       nStubs += std::accumulate(tracks.begin(), tracks.end(), 0, [](int sum, const std::vector<TTStubRef>& track) {
0178         return sum + track.size();
0179       });
0180       allTracks += tracks.size();
0181       if (!useMCTruth_)
0182         continue;
0183       int tmp(0);
0184       associate(tracks, selection, tpPtrsSelection, tmp);
0185       associate(tracks, selection, tpPtrsPerfect, tmp, true);
0186       associate(tracks, reconstructable, tpPtrs, allMatched);
0187       const tt::StreamTrack& stream = streamsTrack[region];
0188       const auto end = std::find_if(
0189           stream.rbegin(), stream.rend(), [](const tt::FrameTrack& frame) { return frame.first.isNonnull(); });
0190       const int size = std::distance(stream.begin(), end.base()) - 1;
0191       hisChannel_->Fill(size);
0192       profChannel_->Fill(1, size);
0193       prof_->Fill(1, nStubs);
0194       prof_->Fill(2, nTracks);
0195     }
0196     prof_->Fill(4, allMatched);
0197     prof_->Fill(5, allTracks);
0198     prof_->Fill(6, tpPtrs.size());
0199     prof_->Fill(7, tpPtrsSelection.size());
0200     prof_->Fill(10, tpPtrsPerfect.size());
0201     nEvents_++;
0202   }
0203 
0204   void AnalyzerTM::endJob() {
0205     if (nEvents_ == 0)
0206       return;
0207     // printout summary
0208     const double totalTPs = prof_->GetBinContent(9);
0209     const double numStubs = prof_->GetBinContent(1);
0210     const double numTracks = prof_->GetBinContent(2);
0211     const double totalTracks = prof_->GetBinContent(5);
0212     const double numTracksMatched = prof_->GetBinContent(4);
0213     const double numTPsAll = prof_->GetBinContent(6);
0214     const double numTPsEff = prof_->GetBinContent(7);
0215     const double numTPsEffPerfect = prof_->GetBinContent(10);
0216     const double errStubs = prof_->GetBinError(1);
0217     const double errTracks = prof_->GetBinError(2);
0218     const double fracFake = (totalTracks - numTracksMatched) / totalTracks;
0219     const double fracDup = (numTracksMatched - numTPsAll) / totalTracks;
0220     const double eff = numTPsEff / totalTPs;
0221     const double errEff = std::sqrt(eff * (1. - eff) / totalTPs / nEvents_);
0222     const double effPerfect = numTPsEffPerfect / totalTPs;
0223     const double errEffPerfect = std::sqrt(effPerfect * (1. - effPerfect) / totalTPs / nEvents_);
0224     const std::vector<double> nums = {numStubs, numTracks};
0225     const std::vector<double> errs = {errStubs, errTracks};
0226     const int wNums = std::ceil(std::log10(*std::max_element(nums.begin(), nums.end()))) + 5;
0227     const int wErrs = std::ceil(std::log10(*std::max_element(errs.begin(), errs.end()))) + 5;
0228     log_ << "                          TM  SUMMARY                        " << std::endl;
0229     log_ << "number of stubs       per TFP = " << std::setw(wNums) << numStubs << " +- " << std::setw(wErrs) << errStubs
0230          << std::endl;
0231     log_ << "number of tracks      per TFP = " << std::setw(wNums) << numTracks << " +- " << std::setw(wErrs)
0232          << errTracks << std::endl;
0233     log_ << "  current tracking efficiency = " << std::setw(wNums) << effPerfect << " +- " << std::setw(wErrs)
0234          << errEffPerfect << std::endl;
0235     log_ << "  max     tracking efficiency = " << std::setw(wNums) << eff << " +- " << std::setw(wErrs) << errEff
0236          << std::endl;
0237     log_ << "                    fake rate = " << std::setw(wNums) << fracFake << std::endl;
0238     log_ << "               duplicate rate = " << std::setw(wNums) << fracDup << std::endl;
0239     log_ << "=============================================================";
0240     edm::LogPrint(moduleDescription().moduleName()) << log_.str();
0241   }
0242 
0243   //
0244   void AnalyzerTM::formTracks(const tt::StreamsTrack& streamsTrack,
0245                               const tt::StreamsStub& streamsStubs,
0246                               std::vector<std::vector<TTStubRef>>& tracks,
0247                               int channel) const {
0248     static const int numLayers = channelAssignment_->tmNumLayers();
0249     const int offset = channel * numLayers;
0250     const tt::StreamTrack& streamTrack = streamsTrack[channel];
0251     const int numTracks =
0252         std::accumulate(streamTrack.begin(), streamTrack.end(), 0, [](int sum, const tt::FrameTrack& frame) {
0253           return sum + (frame.first.isNonnull() ? 1 : 0);
0254         });
0255     tracks.reserve(numTracks);
0256     for (int frame = 0; frame < static_cast<int>(streamTrack.size()); frame++) {
0257       const tt::FrameTrack& frameTrack = streamTrack[frame];
0258       if (frameTrack.first.isNull())
0259         continue;
0260       std::vector<TTStubRef> ttStubRefs;
0261       ttStubRefs.reserve(numLayers);
0262       for (int layer = 0; layer < numLayers; layer++) {
0263         const tt::FrameStub& stub = streamsStubs[offset + layer][frame];
0264         if (stub.first.isNonnull())
0265           ttStubRefs.push_back(stub.first);
0266       }
0267       tracks.push_back(ttStubRefs);
0268     }
0269   }
0270 
0271   //
0272   void AnalyzerTM::associate(const std::vector<std::vector<TTStubRef>>& tracks,
0273                              const tt::StubAssociation* ass,
0274                              std::set<TPPtr>& tps,
0275                              int& sum,
0276                              bool perfect) const {
0277     for (const std::vector<TTStubRef>& ttStubRefs : tracks) {
0278       const std::vector<TPPtr>& tpPtrs = perfect ? ass->associateFinal(ttStubRefs) : ass->associate(ttStubRefs);
0279       if (tpPtrs.empty())
0280         continue;
0281       sum++;
0282       std::copy(tpPtrs.begin(), tpPtrs.end(), std::inserter(tps, tps.begin()));
0283     }
0284   }
0285 
0286 }  // namespace trklet
0287 
0288 DEFINE_FWK_MODULE(trklet::AnalyzerTM);