Back to home page

Project CMSSW displayed by LXR

 
 

    


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/Utilities/interface/EDGetToken.h"
0008 #include "FWCore/Utilities/interface/ESGetToken.h"
0009 #include "FWCore/Utilities/interface/InputTag.h"
0010 #include "FWCore/Utilities/interface/Exception.h"
0011 #include "DataFormats/Common/interface/Handle.h"
0012 
0013 #include "L1Trigger/TrackerTFP/interface/Demonstrator.h"
0014 #include "L1Trigger/TrackFindingTracklet/interface/ChannelAssignment.h"
0015 
0016 #include <sstream>
0017 #include <utility>
0018 #include <numeric>
0019 
0020 using namespace std;
0021 using namespace edm;
0022 using namespace tt;
0023 using namespace trackerTFP;
0024 
0025 namespace trklet {
0026 
0027   /*! \class  trklet::AnalyzerDemonstrator
0028    *  \brief  Class to demontrate correctness of track trigger emulators
0029    *          by comparing FW with SW
0030    *  \author Thomas Schuh
0031    *  \date   2022, March
0032    */
0033   class AnalyzerDemonstrator : public one::EDAnalyzer<one::WatchRuns> {
0034   public:
0035     AnalyzerDemonstrator(const ParameterSet& iConfig);
0036     void beginJob() override {}
0037     void beginRun(const Run& iEvent, const EventSetup& iSetup) override;
0038     void analyze(const Event& iEvent, const EventSetup& iSetup) override;
0039     void endRun(const Run& iEvent, const EventSetup& iSetup) override {}
0040     void endJob() override;
0041 
0042   private:
0043     //
0044     void convert(const Event& iEvent,
0045                  const EDGetTokenT<StreamsTrack>& tokenTracks,
0046                  const EDGetTokenT<StreamsStub>& tokenStubs,
0047                  vector<vector<Frame>>& bits) const;
0048     //
0049     template <typename T>
0050     void convert(const T& collection, vector<vector<Frame>>& bits) const;
0051     // ED input token of Tracks
0052     EDGetTokenT<StreamsStub> edGetTokenStubsIn_;
0053     EDGetTokenT<StreamsStub> edGetTokenStubsOut_;
0054     // ED input token of Stubs
0055     EDGetTokenT<StreamsTrack> edGetTokenTracksIn_;
0056     EDGetTokenT<StreamsTrack> edGetTokenTracksOut_;
0057     // Setup token
0058     ESGetToken<Setup, SetupRcd> esGetTokenSetup_;
0059     // ChannelAssignment token
0060     ESGetToken<ChannelAssignment, ChannelAssignmentRcd> esGetTokenChannelAssignment_;
0061     // Demonstrator token
0062     ESGetToken<Demonstrator, DemonstratorRcd> esGetTokenDemonstrator_;
0063     //
0064     const Setup* setup_ = nullptr;
0065     //
0066     const ChannelAssignment* channelAssignment_ = nullptr;
0067     //
0068     const Demonstrator* demonstrator_ = nullptr;
0069     //
0070     int nEvents_ = 0;
0071     //
0072     int nEventsSuccessful_ = 0;
0073   };
0074 
0075   AnalyzerDemonstrator::AnalyzerDemonstrator(const ParameterSet& iConfig) {
0076     // book in- and output ED products
0077     const string& labelIn = iConfig.getParameter<string>("LabelIn");
0078     const string& labelOut = iConfig.getParameter<string>("LabelOut");
0079     const string& branchStubs = iConfig.getParameter<string>("BranchAcceptedStubs");
0080     const string& branchTracks = iConfig.getParameter<string>("BranchAcceptedTracks");
0081     edGetTokenStubsIn_ = consumes<StreamsStub>(InputTag(labelIn, branchStubs));
0082     edGetTokenStubsOut_ = consumes<StreamsStub>(InputTag(labelOut, branchStubs));
0083     if (labelOut != "TrackFindingTrackletProducerKFout")
0084       edGetTokenStubsOut_ = consumes<StreamsStub>(InputTag(labelOut, branchStubs));
0085     if (labelIn != "TrackFindingTrackletProducerIRin")
0086       edGetTokenTracksIn_ = consumes<StreamsTrack>(InputTag(labelIn, branchTracks));
0087     if (labelOut != "TrackFindingTrackletProducerIRin")
0088       edGetTokenTracksOut_ = consumes<StreamsTrack>(InputTag(labelOut, branchTracks));
0089     // book ES products
0090     esGetTokenSetup_ = esConsumes<Setup, SetupRcd, Transition::BeginRun>();
0091     esGetTokenChannelAssignment_ = esConsumes<ChannelAssignment, ChannelAssignmentRcd, Transition::BeginRun>();
0092     esGetTokenDemonstrator_ = esConsumes<Demonstrator, DemonstratorRcd, Transition::BeginRun>();
0093   }
0094 
0095   void AnalyzerDemonstrator::beginRun(const Run& iEvent, const EventSetup& iSetup) {
0096     //
0097     setup_ = &iSetup.getData(esGetTokenSetup_);
0098     //
0099     channelAssignment_ = &iSetup.getData(esGetTokenChannelAssignment_);
0100     //
0101     demonstrator_ = &iSetup.getData(esGetTokenDemonstrator_);
0102   }
0103 
0104   void AnalyzerDemonstrator::analyze(const Event& iEvent, const EventSetup& iSetup) {
0105     nEvents_++;
0106     vector<vector<Frame>> input;
0107     vector<vector<Frame>> output;
0108     convert(iEvent, edGetTokenTracksIn_, edGetTokenStubsIn_, input);
0109     convert(iEvent, edGetTokenTracksOut_, edGetTokenStubsOut_, output);
0110     if (demonstrator_->analyze(input, output))
0111       nEventsSuccessful_++;
0112   }
0113 
0114   //
0115   void AnalyzerDemonstrator::convert(const Event& iEvent,
0116                                      const EDGetTokenT<StreamsTrack>& tokenTracks,
0117                                      const EDGetTokenT<StreamsStub>& tokenStubs,
0118                                      vector<vector<Frame>>& bits) const {
0119     const bool tracks = !tokenTracks.isUninitialized();
0120     const bool stubs = !tokenStubs.isUninitialized();
0121     Handle<StreamsStub> handleStubs;
0122     Handle<StreamsTrack> handleTracks;
0123     int numChannelTracks(0);
0124     if (tracks) {
0125       iEvent.getByToken<StreamsTrack>(tokenTracks, handleTracks);
0126       numChannelTracks = handleTracks->size();
0127     }
0128     numChannelTracks /= setup_->numRegions();
0129     int numChannelStubs(0);
0130     vector<int> numChannelsStubs(numChannelTracks, 0);
0131     if (stubs) {
0132       iEvent.getByToken<StreamsStub>(tokenStubs, handleStubs);
0133       numChannelStubs = handleStubs->size() / setup_->numRegions();
0134       const int numChannel = tracks ? numChannelStubs / numChannelTracks : numChannelStubs;
0135       numChannelsStubs = vector<int>(numChannelTracks, numChannel);
0136     }
0137     bits.reserve(numChannelTracks + numChannelStubs);
0138     for (int region = 0; region < setup_->numRegions(); region++) {
0139       if (tracks) {
0140         const int offsetTracks = region * numChannelTracks;
0141         for (int channelTracks = 0; channelTracks < numChannelTracks; channelTracks++) {
0142           const int offsetStubs =
0143               region * numChannelStubs +
0144               accumulate(numChannelsStubs.begin(), next(numChannelsStubs.begin(), channelTracks), 0);
0145           convert(handleTracks->at(offsetTracks + channelTracks), bits);
0146           if (stubs) {
0147             for (int channelStubs = 0; channelStubs < numChannelsStubs[channelTracks]; channelStubs++)
0148               convert(handleStubs->at(offsetStubs + channelStubs), bits);
0149           }
0150         }
0151       } else {
0152         const int offsetStubs = region * numChannelStubs;
0153         for (int channelStubs = 0; channelStubs < numChannelStubs; channelStubs++)
0154           convert(handleStubs->at(offsetStubs + channelStubs), bits);
0155       }
0156     }
0157   }
0158 
0159   //
0160   template <typename T>
0161   void AnalyzerDemonstrator::convert(const T& collection, vector<vector<Frame>>& bits) const {
0162     bits.emplace_back();
0163     vector<Frame>& bvs = bits.back();
0164     bvs.reserve(collection.size());
0165     transform(collection.begin(), collection.end(), back_inserter(bvs), [](const auto& frame) { return frame.second; });
0166   }
0167 
0168   void AnalyzerDemonstrator::endJob() {
0169     stringstream log;
0170     log << "Successrate: " << nEventsSuccessful_ << " / " << nEvents_ << " = " << nEventsSuccessful_ / (double)nEvents_;
0171     LogPrint("L1Trigger/TrackerTFP") << log.str();
0172   }
0173 
0174 }  // namespace trklet
0175 
0176 DEFINE_FWK_MODULE(trklet::AnalyzerDemonstrator);