Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-12 02:41:52

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