Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-13 02:58:23

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/L1ScoutingMuon.h"
0016 #include "DataFormats/L1Scouting/interface/L1ScoutingCalo.h"
0017 #include "DataFormats/L1Scouting/interface/OrbitCollection.h"
0018 #include "L1TriggerScouting/Utilities/interface/conversion.h"
0019 
0020 // root libraries
0021 #include "TLorentzVector.h"
0022 #include "Math/VectorUtil.h"
0023 
0024 #include <memory>
0025 #include <utility>
0026 #include <vector>
0027 
0028 using namespace l1ScoutingRun3;
0029 
0030 class MuTagJetBxSelector : public edm::stream::EDProducer<> {
0031 public:
0032   explicit MuTagJetBxSelector(const edm::ParameterSet&);
0033   ~MuTagJetBxSelector() override {}
0034   static void fillDescriptions(edm::ConfigurationDescriptions&);
0035 
0036 private:
0037   void produce(edm::Event&, const edm::EventSetup&) override;
0038 
0039   // tokens for scouting data
0040   edm::EDGetTokenT<OrbitCollection<l1ScoutingRun3::Muon>> muonsTokenData_;
0041   edm::EDGetTokenT<OrbitCollection<l1ScoutingRun3::Jet>> jetsTokenData_;
0042 
0043   // SELECTION THRESHOLDS
0044   int minNJet_;
0045   std::vector<double> minJetEt_;
0046   std::vector<double> maxJetEta_;
0047   std::vector<double> minMuPt_;
0048   std::vector<double> maxMuEta_;
0049   std::vector<int> minMuTfIndex_;
0050   std::vector<int> maxMuTfIndex_;
0051   std::vector<int> minMuHwQual_;
0052   std::vector<double> maxDR_;
0053 };
0054 
0055 MuTagJetBxSelector::MuTagJetBxSelector(const edm::ParameterSet& iPSet)
0056     : muonsTokenData_(consumes(iPSet.getParameter<edm::InputTag>("muonsTag"))),
0057       jetsTokenData_(consumes(iPSet.getParameter<edm::InputTag>("jetsTag"))),
0058       minNJet_(iPSet.getParameter<int>("minNJet")),
0059       minJetEt_(iPSet.getParameter<std::vector<double>>("minJetEt")),
0060       maxJetEta_(iPSet.getParameter<std::vector<double>>("maxJetEta")),
0061       minMuPt_(iPSet.getParameter<std::vector<double>>("minMuPt")),
0062       maxMuEta_(iPSet.getParameter<std::vector<double>>("maxMuEta")),
0063       minMuTfIndex_(iPSet.getParameter<std::vector<int>>("minMuTfIndex")),
0064       maxMuTfIndex_(iPSet.getParameter<std::vector<int>>("maxMuTfIndex")),
0065       minMuHwQual_(iPSet.getParameter<std::vector<int>>("minMuHwQual")),
0066       maxDR_(iPSet.getParameter<std::vector<double>>("maxDR")) {
0067   if ((minJetEt_.size() != (size_t)minNJet_) || (maxJetEta_.size() != (size_t)minNJet_) ||
0068       (minMuPt_.size() != (size_t)minNJet_) || (maxMuEta_.size() != (size_t)minNJet_) ||
0069       (minMuTfIndex_.size() != (size_t)minNJet_) || (maxMuTfIndex_.size() != (size_t)minNJet_) ||
0070       (minMuHwQual_.size() != (size_t)minNJet_) || (maxDR_.size() != (size_t)minNJet_))
0071     throw cms::Exception("MuTagJetBxSelector::MuTagJetBxSelector")
0072         << "size mismatch: size of minJetEt or maxJetEta or  minMuPt or maxMuEta or minMuTfIndex or maxMuTfIndex or "
0073            "minMuHwQual or maxDR != minNMu.";
0074 
0075   produces<std::vector<unsigned>>("SelBx").setBranchAlias("MuTagJetSelectedBx");
0076 }
0077 
0078 // ------------ method called for each ORBIT  ------------
0079 void MuTagJetBxSelector::produce(edm::Event& iEvent, const edm::EventSetup&) {
0080   edm::Handle<OrbitCollection<l1ScoutingRun3::Muon>> muonsCollection;
0081   edm::Handle<OrbitCollection<l1ScoutingRun3::Jet>> jetsCollection;
0082 
0083   iEvent.getByToken(muonsTokenData_, muonsCollection);
0084   iEvent.getByToken(jetsTokenData_, jetsCollection);
0085 
0086   std::unique_ptr<std::vector<unsigned>> muTagJetBx(new std::vector<unsigned>);
0087 
0088   // loop over valid bunch crossings
0089   for (const unsigned& bx : jetsCollection->getFilledBxs()) {
0090     const auto& jets = jetsCollection->bxIterator(bx);
0091     const auto& muons = muonsCollection->bxIterator(bx);
0092 
0093     // we have at least N jets and N muons
0094     if (std::ssize(jets) < minNJet_ && std::ssize(muons) < minNJet_)
0095       continue;
0096 
0097     // it must satisfy certain requirements
0098     bool jetCond = false;
0099     bool muCond = false;
0100     int nAccJets = 0;
0101     for (const auto& jet : jets) {
0102       jetCond = (std::abs(demux::fEta(jet.hwEta())) < maxJetEta_[nAccJets]) &&
0103                 (demux::fEt(jet.hwEt()) >= minJetEt_[nAccJets]);
0104       if (!jetCond)
0105         continue;  // jet does not satisfy requirements, next one
0106       ROOT::Math::PtEtaPhiMVector jetLV(demux::fEt(jet.hwEt()), demux::fEta(jet.hwEta()), demux::fPhi(jet.hwPhi()), 0);
0107 
0108       for (const auto& muon : muons) {
0109         muCond = (std::abs(ugmt::fEta(muon.hwEta())) < maxMuEta_[nAccJets]) &&
0110                  (muon.tfMuonIndex() <= maxMuTfIndex_[nAccJets]) && (muon.tfMuonIndex() >= minMuTfIndex_[nAccJets]) &&
0111                  (ugmt::fPt(muon.hwPt()) >= minMuPt_[nAccJets]) && (muon.hwQual() >= minMuHwQual_[nAccJets]);
0112         if (!muCond)
0113           continue;  // muon does not satisfy requirements, next one
0114         ROOT::Math::PtEtaPhiMVector muLV(
0115             ugmt::fPt(muon.hwPt()), ugmt::fEta(muon.hwEta()), ugmt::fPhi(muon.hwPhi()), 0.1057);
0116 
0117         float dr = ROOT::Math::VectorUtil::DeltaR(jetLV, muLV);
0118         if (dr < maxDR_[nAccJets]) {
0119           nAccJets++;  // found mu-tag for current jet, end muon loop
0120           break;
0121         }
0122       }
0123 
0124       if (nAccJets == minNJet_)
0125         break;  // found all requested mu-tagged jets
0126     }
0127 
0128     if (nAccJets < minNJet_)
0129       continue;
0130 
0131     muTagJetBx->push_back(bx);
0132 
0133   }  // end orbit loop
0134 
0135   iEvent.put(std::move(muTagJetBx), "SelBx");
0136 }
0137 
0138 void MuTagJetBxSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0139   edm::ParameterSetDescription desc;
0140   desc.setUnknown();
0141   descriptions.addDefault(desc);
0142 }
0143 
0144 DEFINE_FWK_MODULE(MuTagJetBxSelector);