Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-07-10 23:57:00

0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/stream/EDProducer.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0008 #include "FWCore/Utilities/interface/EDGetToken.h"
0009 #include "FWCore/Utilities/interface/EDPutToken.h"
0010 #include "FWCore/Utilities/interface/Exception.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 #include "FWCore/Utilities/interface/StreamID.h"
0013 
0014 // L1 scouting
0015 #include "DataFormats/L1Scouting/interface/L1ScoutingBMTFStub.h"
0016 #include "DataFormats/L1Scouting/interface/OrbitCollection.h"
0017 
0018 #include <memory>
0019 #include <utility>
0020 #include <vector>
0021 #include <string>
0022 #include <numeric>
0023 #include <set>
0024 
0025 namespace l1ScoutingRun3 {
0026   enum BMTFSelectorConditionType { Simple, Wheel };
0027 }
0028 
0029 class BMTFStubMultiBxSelector : public edm::stream::EDProducer<> {
0030 public:
0031   explicit BMTFStubMultiBxSelector(const edm::ParameterSet&);
0032   ~BMTFStubMultiBxSelector() override {}
0033   static void fillDescriptions(edm::ConfigurationDescriptions&);
0034 
0035   unsigned int makeWheelPattern(const edm::Handle<OrbitCollection<l1ScoutingRun3::BMTFStub>>& stubs, unsigned int bx);
0036   bool windowHasCloseWheels(const std::vector<unsigned int>& bxs);
0037 
0038 private:
0039   void produce(edm::Event&, const edm::EventSetup&) override;
0040 
0041   // Tokens for scouting data
0042   edm::EDGetTokenT<OrbitCollection<l1ScoutingRun3::BMTFStub>> stubsTokenData_;
0043 
0044   // Condition
0045   l1ScoutingRun3::BMTFSelectorConditionType condition_;
0046 
0047   // Selection thresholds
0048   unsigned int bxWindowLength_;
0049   unsigned int minNBMTFStub_;
0050 };
0051 
0052 BMTFStubMultiBxSelector::BMTFStubMultiBxSelector(const edm::ParameterSet& iPSet)
0053     : stubsTokenData_(consumes(iPSet.getParameter<edm::InputTag>("stubsTag"))),
0054       bxWindowLength_(iPSet.getParameter<unsigned int>("bxWindowLength")),
0055       minNBMTFStub_(iPSet.getParameter<unsigned int>("minNBMTFStub")) {
0056   std::string conditionStr = iPSet.getParameter<std::string>("condition");
0057   if (conditionStr == "simple")
0058     condition_ = l1ScoutingRun3::Simple;
0059   else if (conditionStr == "wheel")
0060     condition_ = l1ScoutingRun3::Wheel;
0061   else
0062     throw cms::Exception("BMTFStubMultiBxSelector::BMTFStubMultiBxSelector")
0063         << "Condition '" << conditionStr << "' not supported or not found";
0064 
0065   produces<std::vector<unsigned int>>("SelBx").setBranchAlias("MultiBxStubsSelectedBx");
0066 }
0067 
0068 // ------------ method called for each ORBIT  ------------
0069 void BMTFStubMultiBxSelector::produce(edm::Event& iEvent, const edm::EventSetup&) {
0070   edm::Handle<OrbitCollection<l1ScoutingRun3::BMTFStub>> stubsCollection;
0071 
0072   iEvent.getByToken(stubsTokenData_, stubsCollection);
0073 
0074   // This is necessary to have a clean 'MultiBxStubsSelectedBx' collection
0075   // Indeed, if there is overlap betweem two valid windows, there will be duplicated BXs
0076   std::unique_ptr<std::set<unsigned int>> uniqueStubSelectedBxs = std::make_unique<std::set<unsigned int>>();
0077 
0078   // Loop over valid bunch crossings
0079   std::vector<unsigned int> vNumStubBx(bxWindowLength_, 0);
0080   std::vector<unsigned int> vWheelPatternBx(bxWindowLength_, 0);
0081   for (const unsigned int& bx : stubsCollection->getFilledBxs()) {
0082     // Get number of stubs in current window of size 'ws': [BX-ws+1, BX]
0083     for (unsigned int i = 0; i < std::min(bxWindowLength_, bx); ++i)
0084       vNumStubBx[i] = stubsCollection->getBxSize(bx - i);
0085 
0086     // Get number of stubs in current BX (last of the window) and enque
0087     unsigned int numStubsWindow = std::reduce(vNumStubBx.begin(), vNumStubBx.end());
0088 
0089     // Simple condition: just number of stubs (already checked)
0090     if (condition_ == l1ScoutingRun3::Simple) {
0091       // If there are enough stubs in window...
0092       if (numStubsWindow >= minNBMTFStub_) {
0093         // ...loop over window to add BXs (std::min to include edge case of first BXs)
0094         for (unsigned int i = 0; i < std::min(bxWindowLength_, bx); ++i)
0095           uniqueStubSelectedBxs->insert(bx - i);
0096       }
0097     }
0098     // Wheel condition: enough longitudinally "neighbouring" stubs
0099     else if (condition_ == l1ScoutingRun3::Wheel) {
0100       for (unsigned int i = 0; i < std::min(bxWindowLength_, bx); ++i) {
0101         // Prepare pattern and add to window vector
0102         vWheelPatternBx[i] = makeWheelPattern(stubsCollection, bx - i);
0103       }
0104 
0105       // Check if there are stubs in different BXs with neighbouring wheels
0106       // For example (with window of 3 BXs, neighbouring condition abs(wheel_pair) <= 1)
0107       // BX-2         BX-1        BX0
0108       // s0(wh = 2)   s0(wh=0)    s1(wh=-2)
0109       // s1(wh = 1)
0110       //
0111       // => valid window! (neighbouring stubs s1 in BX-2 and s0 in BX-1, assuming nStub threshold is satisfied)
0112       bool validWindow = windowHasCloseWheels(vWheelPatternBx) && (numStubsWindow >= minNBMTFStub_);
0113 
0114       // If window is valid....
0115       if (validWindow) {
0116         // ...loop over window to add BXs (std::min to include edge case of first BXs)
0117         for (unsigned int i = 0; i < std::min(bxWindowLength_, bx); ++i)
0118           uniqueStubSelectedBxs->insert(bx - i);
0119       }
0120     }
0121   }  // end orbit loop
0122 
0123   // Convert set of selected BXs to a vector and put collection in event content
0124   std::unique_ptr<std::vector<unsigned int>> stubSelectedBx =
0125       std::make_unique<std::vector<unsigned int>>(uniqueStubSelectedBxs->begin(), uniqueStubSelectedBxs->end());
0126   iEvent.put(std::move(stubSelectedBx), "SelBx");
0127 }
0128 
0129 unsigned int BMTFStubMultiBxSelector::makeWheelPattern(
0130     const edm::Handle<OrbitCollection<l1ScoutingRun3::BMTFStub>>& stubs, unsigned int bx) {
0131   // 5 wheel numbers + 2 to handle boundaries in a more comfortable way:
0132   // 0b         0   X   X   X   X   X   0
0133   // wheels   bnd  +2  +1   0  -1  -2 bnd
0134   // position   6   5   4   3   2   1   0
0135   unsigned int wheelPatternBx = 0;
0136   for (const auto& s : stubs->bxIterator(bx))
0137     wheelPatternBx |= (1 << (s.wheel() + 3));
0138   return wheelPatternBx;
0139 }
0140 
0141 bool BMTFStubMultiBxSelector::windowHasCloseWheels(const std::vector<unsigned int>& bxs) {
0142   for (size_t i = 0; i < std::size(bxs); ++i) {
0143     for (size_t j = i + 1; j < std::size(bxs); ++j) {
0144       // Assume for example that we have the following patterns
0145       //          bwwwwwb (b = boundary, w = wheel fields)
0146       // BX-2 : 0b0010000 (stub with wheel=+2)
0147       // BX-1 : 0b0001000 (stub with wheel=+1)
0148       // Bitwise AND comparison between:
0149       // BX-2 : 0b0010000
0150       // BX-1 : 0b0011100
0151       //            bsb   (s = stub wheel, b = boundary wheels)
0152       unsigned int compare = bxs[i] & ((bxs[j] << 1) | bxs[j] | (bxs[j] >> 1));
0153       bool checkWindow = compare & 0b0111110;
0154       if (checkWindow)
0155         return true;
0156     }
0157   }
0158   return false;
0159 }
0160 
0161 void BMTFStubMultiBxSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0162   edm::ParameterSetDescription desc;
0163   desc.add<edm::InputTag>("stubsTag")->setComment("input BMTF stubs orbit collection");
0164   desc.add<unsigned int>("bxWindowLength")->setComment("number of consecutive bunch crossings of window to analyse");
0165   desc.add<unsigned int>("minNBMTFStub")
0166       ->setComment("minimum number of stubs that should be in bunch crossings window");
0167   desc.add<std::string>("condition")->setComment("condition to apply for bunch crossing window selection");
0168   descriptions.add("BMTFStubMultiBxSelector", desc);
0169 }
0170 
0171 DEFINE_FWK_MODULE(BMTFStubMultiBxSelector);