Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:36

0001 /* \class HLTMuonL1TFilter
0002  *
0003  * This is an HLTFilter implementing filtering on L1T Stage2 GMT objects
0004  * 
0005  * \author:  V. Rekovic
0006  *
0007  */
0008 
0009 #include "HLTMuonL1TFilter.h"
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0012 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/Utilities/interface/EDMException.h"
0015 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0017 #include "FWCore/Utilities/interface/InputTag.h"
0018 #include "DataFormats/Math/interface/deltaR.h"
0019 #include "TMath.h"
0020 
0021 #include <vector>
0022 
0023 HLTMuonL1TFilter::HLTMuonL1TFilter(const edm::ParameterSet& iConfig)
0024     : HLTFilter(iConfig),
0025       candTag_(iConfig.getParameter<edm::InputTag>("CandTag")),
0026       candToken_(consumes<l1t::MuonBxCollection>(candTag_)),
0027       previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
0028       previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
0029       maxEta_(iConfig.getParameter<double>("MaxEta")),
0030       minPt_(iConfig.getParameter<double>("MinPt")),
0031       maxDR_(iConfig.getParameter<double>("MaxDeltaR")),
0032       maxDR2_(maxDR_ * maxDR_),
0033       minN_(iConfig.getParameter<int>("MinN")),
0034       centralBxOnly_(iConfig.getParameter<bool>("CentralBxOnly")) {
0035   //set the quality bit mask
0036   qualityBitMask_ = 0;
0037   vector<int> selectQualities = iConfig.getParameter<vector<int> >("SelectQualities");
0038   for (int selectQualitie : selectQualities) {
0039     //     if(selectQualities[i] > 7){  // FIXME: this will be updated once we have info from L1
0040     //       throw edm::Exception(edm::errors::Configuration) << "QualityBits must be smaller than 8!";
0041     //     }
0042     qualityBitMask_ |= 1 << selectQualitie;
0043   }
0044   //make sure cut parameter for candidate matching is strictly positive
0045   if (maxDR_ <= 0.) {
0046     throw cms::Exception("HLTMuonL1TFilterConfiguration")
0047         << "invalid value for parameter \"MaxDeltaR\" (must be > 0): " << maxDR_;
0048   }
0049   // dump parameters for debugging
0050   if (edm::isDebugEnabled()) {
0051     ostringstream ss;
0052     ss << "Constructed with parameters:" << endl;
0053     ss << "    CandTag = " << candTag_.encode() << endl;
0054     ss << "    PreviousCandTag = " << previousCandTag_.encode() << endl;
0055     ss << "    MaxEta = " << maxEta_ << endl;
0056     ss << "    MinPt = " << minPt_ << endl;
0057     ss << "    SelectQualities =";
0058     //     for(size_t i = 0; i < 8; i++){
0059     //       if((qualityBitMask_>>i) % 2) ss << " " << i;
0060     //     }
0061     ss << endl;
0062     ss << "    MinN = " << minN_ << endl;
0063     ss << "    saveTags= " << saveTags();
0064     LogDebug("HLTMuonL1TFilter") << ss.str();
0065   }
0066 }
0067 
0068 HLTMuonL1TFilter::~HLTMuonL1TFilter() = default;
0069 
0070 void HLTMuonL1TFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0071   edm::ParameterSetDescription desc;
0072   makeHLTFilterDescription(desc);
0073   desc.add<edm::InputTag>("CandTag", edm::InputTag("hltGmtStage2Digis"));
0074   desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
0075   desc.add<double>("MaxEta", 2.5);
0076   desc.add<double>("MinPt", 0.0);
0077   desc.add<double>("MaxDeltaR", 0.3);
0078   desc.add<int>("MinN", 1);
0079   desc.add<bool>("CentralBxOnly", true);
0080   {
0081     std::vector<int> temp1;
0082     temp1.reserve(0);
0083     desc.add<std::vector<int> >("SelectQualities", temp1);
0084   }
0085   descriptions.add("hltMuonL1TFilter", desc);
0086 }
0087 
0088 bool HLTMuonL1TFilter::hltFilter(edm::Event& iEvent,
0089                                  const edm::EventSetup& iSetup,
0090                                  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0091   using namespace std;
0092   using namespace edm;
0093   using namespace trigger;
0094   using namespace l1t;
0095 
0096   // All HLT filters must create and fill an HLT filter object,
0097   // recording any reconstructed physics objects satisfying (or not)
0098   // this HLT filter, and place it in the Event.
0099 
0100   // get hold of all muons
0101   Handle<MuonBxCollection> allMuons;
0102   iEvent.getByToken(candToken_, allMuons);
0103 
0104   // get hold of TFOWRs that fired the previous level
0105   Handle<TriggerFilterObjectWithRefs> previousLevelTFOWR;
0106   iEvent.getByToken(previousCandToken_, previousLevelTFOWR);
0107 
0108   vector<MuonRef> prevMuons;
0109   previousLevelTFOWR->getObjects(TriggerL1Mu, prevMuons);
0110 
0111   // look at all muon candidates, check cuts and add to filter object
0112   int n = 0;
0113 
0114   for (int ibx = allMuons->getFirstBX(); ibx <= allMuons->getLastBX(); ++ibx) {
0115     if (centralBxOnly_ && (ibx != 0))
0116       continue;
0117     for (auto it = allMuons->begin(ibx); it != allMuons->end(ibx); it++) {
0118       MuonRef muon(allMuons, distance(allMuons->begin(allMuons->getFirstBX()), it));
0119 
0120       //check maxEta cut
0121       if (fabs(muon->eta()) > maxEta_)
0122         continue;
0123 
0124       //check pT cut
0125       if (muon->pt() < minPt_)
0126         continue;
0127 
0128       //check quality cut
0129       if (qualityBitMask_) {
0130         int quality = (it->hwQual() == 0 ? 0 : (1 << it->hwQual()));
0131         if ((quality & qualityBitMask_) == 0)
0132           continue;
0133       }
0134 
0135       // Only select muons that were selected in the previous level
0136       bool matchPrevL1 = false;
0137       int prevSize = prevMuons.size();
0138       for (int it2 = 0; it2 < prevSize; it2++) {
0139         if (deltaR2(muon->eta(), muon->phi(), prevMuons[it2]->eta(), prevMuons[it2]->phi()) < maxDR2_) {
0140           matchPrevL1 = true;
0141           break;
0142         }
0143       }
0144       if (!matchPrevL1)
0145         continue;
0146 
0147       //we have a good candidate
0148       n++;
0149       filterproduct.addObject(TriggerL1Mu, muon);
0150     }
0151   }
0152 
0153   if (saveTags())
0154     filterproduct.addCollectionTag(candTag_);
0155 
0156   // filter decision
0157   const bool accept(n >= minN_);
0158 
0159   // dump event for debugging
0160   if (edm::isDebugEnabled()) {
0161     LogTrace("HLTMuonL1TFilter") << "\nHLTMuonL1TFilter -----------------------------------------------" << endl;
0162     LogTrace("HLTMuonL1TFilter") << "L1mu#" << '\t' << "q" << '\t' << "pt" << '\t' << '\t' << "eta" << '\t' << "phi"
0163                                  << '\t' << "quality" << '\t' << "isPrev\t (|maxEta| = " << maxEta_ << ")" << endl;
0164     LogTrace("HLTMuonL1TFilter") << "--------------------------------------------------------------------------"
0165                                  << endl;
0166 
0167     vector<MuonRef> firedMuons;
0168     filterproduct.getObjects(TriggerL1Mu, firedMuons);
0169     for (size_t i = 0; i < firedMuons.size(); i++) {
0170       l1t::MuonRef mu = firedMuons[i];
0171       bool isPrev = find(prevMuons.begin(), prevMuons.end(), mu) != prevMuons.end();
0172       LogTrace("HLTMuonL1TFilter") << i << '\t' << setprecision(2) << scientific << mu->charge() << '\t' << mu->pt()
0173                                    << '\t' << fixed << mu->eta() << '\t' << mu->phi() << '\t' << mu->hwQual() << '\t'
0174                                    << isPrev << endl;
0175     }
0176     LogTrace("HLTMuonL1TFilter") << "--------------------------------------------------------------------------"
0177                                  << endl;
0178     LogTrace("HLTMuonL1TFilter") << "Decision of this filter is " << accept
0179                                  << ", number of muons passing = " << filterproduct.l1tmuonSize();
0180   }
0181 
0182   return accept;
0183 }
0184 
0185 // declare this class as a framework plugin
0186 #include "FWCore/Framework/interface/MakerMacros.h"
0187 DEFINE_FWK_MODULE(HLTMuonL1TFilter);