File indexing completed on 2024-04-06 12:22:06
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/TrackerTFP/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 <set>
0026 #include <cmath>
0027 #include <numeric>
0028 #include <sstream>
0029
0030 using namespace std;
0031 using namespace edm;
0032 using namespace trackerTFP;
0033 using namespace tt;
0034
0035 namespace trklet {
0036
0037
0038
0039
0040
0041
0042 class AnalyzerKFin : public one::EDAnalyzer<one::WatchRuns, one::SharedResources> {
0043 public:
0044 AnalyzerKFin(const ParameterSet& iConfig);
0045 void beginJob() override {}
0046 void beginRun(const Run& iEvent, const EventSetup& iSetup) override;
0047 void analyze(const Event& iEvent, const EventSetup& iSetup) override;
0048 void endRun(const Run& iEvent, const EventSetup& iSetup) override {}
0049 void endJob() override;
0050
0051 private:
0052
0053 void formTracks(const StreamsTrack& streamsTrack,
0054 const StreamsStub& streamsStubs,
0055 vector<vector<TTStubRef>>& tracks,
0056 int channel) const;
0057
0058 void associate(const vector<vector<TTStubRef>>& tracks,
0059 const StubAssociation* ass,
0060 set<TPPtr>& tps,
0061 int& sum,
0062 bool perfect = false) const;
0063
0064 EDGetTokenT<StreamsStub> edGetTokenAcceptedStubs_;
0065
0066 EDGetTokenT<StreamsTrack> edGetTokenAcceptedTracks_;
0067
0068 EDGetTokenT<StreamsStub> edGetTokenLostStubs_;
0069
0070 EDGetTokenT<StreamsTrack> edGetTokenLostTracks_;
0071
0072 EDGetTokenT<StubAssociation> edGetTokenSelection_;
0073
0074 EDGetTokenT<StubAssociation> edGetTokenReconstructable_;
0075
0076 ESGetToken<Setup, SetupRcd> esGetTokenSetup_;
0077
0078 ESGetToken<DataFormats, DataFormatsRcd> esGetTokenDataFormats_;
0079
0080 const Setup* setup_ = nullptr;
0081
0082 const DataFormats* dataFormats_ = nullptr;
0083
0084 bool useMCTruth_;
0085
0086 int nEvents_ = 0;
0087
0088
0089
0090 TProfile* prof_;
0091 TProfile* profChannel_;
0092 TH1F* hisChannel_;
0093
0094
0095 stringstream log_;
0096 };
0097
0098 AnalyzerKFin::AnalyzerKFin(const ParameterSet& iConfig) : useMCTruth_(iConfig.getParameter<bool>("UseMCTruth")) {
0099 usesResource("TFileService");
0100
0101 const string& label = iConfig.getParameter<string>("LabelKFin");
0102 const string& branchAcceptedStubs = iConfig.getParameter<string>("BranchAcceptedStubs");
0103 const string& branchAcceptedTracks = iConfig.getParameter<string>("BranchAcceptedTracks");
0104 const string& branchLostStubs = iConfig.getParameter<string>("BranchLostStubs");
0105 const string& branchLostTracks = iConfig.getParameter<string>("BranchLostTracks");
0106 edGetTokenAcceptedStubs_ = consumes<StreamsStub>(InputTag(label, branchAcceptedStubs));
0107 edGetTokenAcceptedTracks_ = consumes<StreamsTrack>(InputTag(label, branchAcceptedTracks));
0108 edGetTokenLostStubs_ = consumes<StreamsStub>(InputTag(label, branchLostStubs));
0109 edGetTokenLostTracks_ = consumes<StreamsTrack>(InputTag(label, branchLostTracks));
0110 if (useMCTruth_) {
0111 const auto& inputTagSelecttion = iConfig.getParameter<InputTag>("InputTagSelection");
0112 const auto& inputTagReconstructable = iConfig.getParameter<InputTag>("InputTagReconstructable");
0113 edGetTokenSelection_ = consumes<StubAssociation>(inputTagSelecttion);
0114 edGetTokenReconstructable_ = consumes<StubAssociation>(inputTagReconstructable);
0115 }
0116
0117 esGetTokenSetup_ = esConsumes<Setup, SetupRcd, Transition::BeginRun>();
0118 esGetTokenDataFormats_ = esConsumes<DataFormats, DataFormatsRcd, Transition::BeginRun>();
0119
0120 log_.setf(ios::fixed, ios::floatfield);
0121 log_.precision(4);
0122 }
0123
0124 void AnalyzerKFin::beginRun(const Run& iEvent, const EventSetup& iSetup) {
0125
0126 setup_ = &iSetup.getData(esGetTokenSetup_);
0127
0128 dataFormats_ = &iSetup.getData(esGetTokenDataFormats_);
0129
0130 Service<TFileService> fs;
0131 TFileDirectory dir;
0132 dir = fs->mkdir("KFin");
0133 prof_ = dir.make<TProfile>("Counts", ";", 10, 0.5, 10.5);
0134 prof_->GetXaxis()->SetBinLabel(1, "Stubs");
0135 prof_->GetXaxis()->SetBinLabel(2, "Tracks");
0136 prof_->GetXaxis()->SetBinLabel(3, "Lost Tracks");
0137 prof_->GetXaxis()->SetBinLabel(4, "Matched Tracks");
0138 prof_->GetXaxis()->SetBinLabel(5, "All Tracks");
0139 prof_->GetXaxis()->SetBinLabel(6, "Found TPs");
0140 prof_->GetXaxis()->SetBinLabel(7, "Found selected TPs");
0141 prof_->GetXaxis()->SetBinLabel(8, "Lost TPs");
0142 prof_->GetXaxis()->SetBinLabel(9, "All TPs");
0143 prof_->GetXaxis()->SetBinLabel(10, "Perfect TPs");
0144
0145 constexpr int maxOcc = 180;
0146 const int numChannels = setup_->kfNumWorker();
0147 hisChannel_ = dir.make<TH1F>("His Channel Occupancy", ";", maxOcc, -.5, maxOcc - .5);
0148 profChannel_ = dir.make<TProfile>("Prof Channel Occupancy", ";", numChannels, -.5, numChannels - .5);
0149 }
0150
0151 void AnalyzerKFin::analyze(const Event& iEvent, const EventSetup& iSetup) {
0152
0153 Handle<StreamsStub> handleAcceptedStubs;
0154 iEvent.getByToken<StreamsStub>(edGetTokenAcceptedStubs_, handleAcceptedStubs);
0155 const StreamsStub& acceptedStubs = *handleAcceptedStubs;
0156 Handle<StreamsTrack> handleAcceptedTracks;
0157 iEvent.getByToken<StreamsTrack>(edGetTokenAcceptedTracks_, handleAcceptedTracks);
0158 const StreamsTrack& acceptedTracks = *handleAcceptedTracks;
0159 Handle<StreamsStub> handleLostStubs;
0160 iEvent.getByToken<StreamsStub>(edGetTokenLostStubs_, handleLostStubs);
0161 const StreamsStub& lostStubs = *handleLostStubs;
0162 Handle<StreamsTrack> handleLostTracks;
0163 iEvent.getByToken<StreamsTrack>(edGetTokenLostTracks_, handleLostTracks);
0164 const StreamsTrack& lostTracks = *handleLostTracks;
0165
0166 const StubAssociation* selection = nullptr;
0167 const StubAssociation* reconstructable = nullptr;
0168 if (useMCTruth_) {
0169 Handle<StubAssociation> handleSelection;
0170 iEvent.getByToken<StubAssociation>(edGetTokenSelection_, handleSelection);
0171 selection = handleSelection.product();
0172 prof_->Fill(9, selection->numTPs());
0173 Handle<StubAssociation> handleReconstructable;
0174 iEvent.getByToken<StubAssociation>(edGetTokenReconstructable_, handleReconstructable);
0175 reconstructable = handleReconstructable.product();
0176 }
0177
0178 set<TPPtr> tpPtrs;
0179 set<TPPtr> tpPtrsSelection;
0180 set<TPPtr> tpPtrsPerfect;
0181 set<TPPtr> tpPtrsLost;
0182 int allMatched(0);
0183 int allTracks(0);
0184 for (int region = 0; region < setup_->numRegions(); region++) {
0185 const int offset = region * setup_->kfNumWorker();
0186 int nStubs(0);
0187 int nTracks(0);
0188 int nLost(0);
0189 for (int channel = 0; channel < setup_->kfNumWorker(); channel++) {
0190 vector<vector<TTStubRef>> tracks;
0191 formTracks(acceptedTracks, acceptedStubs, tracks, offset + channel);
0192 vector<vector<TTStubRef>> lost;
0193 formTracks(lostTracks, lostStubs, lost, offset + channel);
0194 nTracks += tracks.size();
0195 nStubs += accumulate(tracks.begin(), tracks.end(), 0, [](int sum, const vector<TTStubRef>& track) {
0196 return sum + (int)track.size();
0197 });
0198 nLost += lost.size();
0199 allTracks += tracks.size();
0200 if (!useMCTruth_)
0201 continue;
0202 int tmp(0);
0203 associate(tracks, selection, tpPtrsSelection, tmp);
0204 associate(tracks, selection, tpPtrsPerfect, tmp, true);
0205 associate(lost, selection, tpPtrsLost, tmp);
0206 associate(tracks, reconstructable, tpPtrs, allMatched);
0207 const StreamTrack& stream = acceptedTracks[offset + channel];
0208 const auto end =
0209 find_if(stream.rbegin(), stream.rend(), [](const FrameTrack& frame) { return frame.first.isNonnull(); });
0210 const int size = distance(stream.begin(), end.base()) - 1;
0211 hisChannel_->Fill(size);
0212 profChannel_->Fill(channel, size);
0213 }
0214 prof_->Fill(1, nStubs);
0215 prof_->Fill(2, nTracks);
0216 prof_->Fill(3, nLost);
0217 }
0218 vector<TPPtr> recovered;
0219 recovered.reserve(tpPtrsLost.size());
0220 set_intersection(tpPtrsLost.begin(), tpPtrsLost.end(), tpPtrs.begin(), tpPtrs.end(), back_inserter(recovered));
0221 for (const TPPtr& tpPtr : recovered)
0222 tpPtrsLost.erase(tpPtr);
0223 prof_->Fill(4, allMatched);
0224 prof_->Fill(5, allTracks);
0225 prof_->Fill(6, tpPtrs.size());
0226 prof_->Fill(7, tpPtrsSelection.size());
0227 prof_->Fill(8, tpPtrsLost.size());
0228 prof_->Fill(10, tpPtrsPerfect.size());
0229 nEvents_++;
0230 }
0231
0232 void AnalyzerKFin::endJob() {
0233 if (nEvents_ == 0)
0234 return;
0235
0236 const double totalTPs = prof_->GetBinContent(9);
0237 const double numStubs = prof_->GetBinContent(1);
0238 const double numTracks = prof_->GetBinContent(2);
0239 const double numTracksLost = prof_->GetBinContent(3);
0240 const double totalTracks = prof_->GetBinContent(5);
0241 const double numTracksMatched = prof_->GetBinContent(4);
0242 const double numTPsAll = prof_->GetBinContent(6);
0243 const double numTPsEff = prof_->GetBinContent(7);
0244 const double numTPsLost = prof_->GetBinContent(8);
0245 const double numTPsEffPerfect = prof_->GetBinContent(10);
0246 const double errStubs = prof_->GetBinError(1);
0247 const double errTracks = prof_->GetBinError(2);
0248 const double errTracksLost = prof_->GetBinError(3);
0249 const double fracFake = (totalTracks - numTracksMatched) / totalTracks;
0250 const double fracDup = (numTracksMatched - numTPsAll) / totalTracks;
0251 const double eff = numTPsEff / totalTPs;
0252 const double errEff = sqrt(eff * (1. - eff) / totalTPs / nEvents_);
0253 const double effLoss = numTPsLost / totalTPs;
0254 const double errEffLoss = sqrt(effLoss * (1. - effLoss) / totalTPs / nEvents_);
0255 const double effPerfect = numTPsEffPerfect / totalTPs;
0256 const double errEffPerfect = sqrt(effPerfect * (1. - effPerfect) / totalTPs / nEvents_);
0257 const vector<double> nums = {numStubs, numTracks, numTracksLost};
0258 const vector<double> errs = {errStubs, errTracks, errTracksLost};
0259 const int wNums = ceil(log10(*max_element(nums.begin(), nums.end()))) + 5;
0260 const int wErrs = ceil(log10(*max_element(errs.begin(), errs.end()))) + 5;
0261 log_ << " KFin SUMMARY " << endl;
0262 log_ << "number of stubs per TFP = " << setw(wNums) << numStubs << " +- " << setw(wErrs) << errStubs << endl;
0263 log_ << "number of tracks per TFP = " << setw(wNums) << numTracks << " +- " << setw(wErrs) << errTracks
0264 << endl;
0265 log_ << "number of lost tracks per TFP = " << setw(wNums) << numTracksLost << " +- " << setw(wErrs) << errTracksLost
0266 << endl;
0267 log_ << " current tracking efficiency = " << setw(wNums) << effPerfect << " +- " << setw(wErrs) << errEffPerfect
0268 << endl;
0269 log_ << " max tracking efficiency = " << setw(wNums) << eff << " +- " << setw(wErrs) << errEff << endl;
0270 log_ << " lost tracking efficiency = " << setw(wNums) << effLoss << " +- " << setw(wErrs) << errEffLoss << endl;
0271 log_ << " fake rate = " << setw(wNums) << fracFake << endl;
0272 log_ << " duplicate rate = " << setw(wNums) << fracDup << endl;
0273 log_ << "=============================================================";
0274 LogPrint("L1Trigger/TrackerTFP") << log_.str();
0275 }
0276
0277
0278 void AnalyzerKFin::formTracks(const StreamsTrack& streamsTrack,
0279 const StreamsStub& streamsStubs,
0280 vector<vector<TTStubRef>>& tracks,
0281 int channel) const {
0282 const int offset = channel * setup_->numLayers();
0283 const StreamTrack& streamTrack = streamsTrack[channel];
0284 const int numTracks = accumulate(streamTrack.begin(), streamTrack.end(), 0, [](int sum, const FrameTrack& frame) {
0285 return sum + (frame.first.isNonnull() ? 1 : 0);
0286 });
0287 tracks.reserve(numTracks);
0288 for (int frame = 0; frame < (int)streamTrack.size(); frame++) {
0289 const FrameTrack& frameTrack = streamTrack[frame];
0290 if (frameTrack.first.isNull())
0291 continue;
0292 vector<TTStubRef> ttStubRefs;
0293 ttStubRefs.reserve(setup_->numLayers());
0294 for (int layer = 0; layer < setup_->numLayers(); layer++) {
0295 const FrameStub& stub = streamsStubs[offset + layer][frame];
0296 if (stub.first.isNonnull())
0297 ttStubRefs.push_back(stub.first);
0298 }
0299 tracks.push_back(ttStubRefs);
0300 }
0301 }
0302
0303
0304 void AnalyzerKFin::associate(const vector<vector<TTStubRef>>& tracks,
0305 const StubAssociation* ass,
0306 set<TPPtr>& tps,
0307 int& sum,
0308 bool perfect) const {
0309 for (const vector<TTStubRef>& ttStubRefs : tracks) {
0310 const vector<TPPtr>& tpPtrs = perfect ? ass->associateFinal(ttStubRefs) : ass->associate(ttStubRefs);
0311 if (tpPtrs.empty())
0312 continue;
0313 sum++;
0314 copy(tpPtrs.begin(), tpPtrs.end(), inserter(tps, tps.begin()));
0315 }
0316 }
0317
0318 }
0319
0320 DEFINE_FWK_MODULE(trklet::AnalyzerKFin);