Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 #include "HLTMuonL1RegionalFilter.h"
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0005 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "DataFormats/L1Trigger/interface/L1MuonParticle.h"
0008 #include "DataFormats/L1Trigger/interface/L1MuonParticleFwd.h"
0009 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTCand.h"
0010 #include "FWCore/Utilities/interface/EDMException.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 
0013 HLTMuonL1RegionalFilter::HLTMuonL1RegionalFilter(const edm::ParameterSet& iConfig)
0014     : HLTFilter(iConfig),
0015       candTag_(iConfig.getParameter<edm::InputTag>("CandTag")),
0016       candToken_(consumes<l1extra::L1MuonParticleCollection>(candTag_)),
0017       previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
0018       previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
0019       minN_(iConfig.getParameter<int>("MinN")) {
0020   using namespace std;
0021   using namespace edm;
0022 
0023   // read in the eta-range dependent parameters
0024   const vector<ParameterSet> cuts = iConfig.getParameter<vector<ParameterSet> >("Cuts");
0025   size_t ranges = cuts.size();
0026   if (ranges == 0) {
0027     throw edm::Exception(errors::Configuration) << "Please provide at least one PSet in the Cuts VPSet!";
0028   }
0029   etaBoundaries_.reserve(ranges + 1);
0030   minPts_.reserve(ranges);
0031   qualityBitMasks_.reserve(ranges);
0032   for (size_t i = 0; i < ranges; i++) {
0033     //set the eta range
0034     vector<double> etaRange = cuts[i].getParameter<vector<double> >("EtaRange");
0035     if (etaRange.size() != 2 || etaRange[0] >= etaRange[1]) {
0036       throw edm::Exception(errors::Configuration) << "EtaRange must have two non-equal values in increasing order!";
0037     }
0038     if (i == 0) {
0039       etaBoundaries_.push_back(etaRange[0]);
0040     } else if (etaBoundaries_[i] != etaRange[0]) {
0041       throw edm::Exception(errors::Configuration)
0042           << "EtaRanges must be disjoint without gaps and listed in increasing eta order!";
0043     }
0044     etaBoundaries_.push_back(etaRange[1]);
0045 
0046     //set the minPt
0047     minPts_.push_back(cuts[i].getParameter<double>("MinPt"));
0048 
0049     //set the quality bit masks
0050     qualityBitMasks_.push_back(0);
0051     vector<unsigned int> qualities = cuts[i].getParameter<vector<unsigned int> >("QualityBits");
0052     for (unsigned int qualitie : qualities) {
0053       if (7U < qualitie) {  // qualities[j] >= 0, since qualities[j] is unsigned
0054         throw edm::Exception(errors::Configuration) << "QualityBits must be between 0 and 7 !";
0055       }
0056       qualityBitMasks_[i] |= 1 << qualitie;
0057     }
0058   }
0059 
0060   // dump parameters for debugging
0061   if (edm::isDebugEnabled()) {
0062     ostringstream ss;
0063     ss << "Constructed with parameters:" << endl;
0064     ss << "    CandTag = " << candTag_.encode() << endl;
0065     ss << "    PreviousCandTag = " << previousCandTag_.encode() << endl;
0066     ss << "    EtaBoundaries = \t" << etaBoundaries_[0];
0067     for (size_t i = 1; i < etaBoundaries_.size(); i++) {
0068       ss << '\t' << etaBoundaries_[i];
0069     }
0070     ss << endl;
0071     ss << "    MinPts =        \t    " << minPts_[0];
0072     for (size_t i = 1; i < minPts_.size(); i++) {
0073       ss << "\t    " << minPts_[i];
0074     }
0075     ss << endl;
0076     ss << "    QualityBitMasks =  \t    " << qualityBitMasks_[0];
0077     for (size_t i = 1; i < qualityBitMasks_.size(); i++) {
0078       ss << "\t    " << qualityBitMasks_[i];
0079     }
0080     ss << endl;
0081     ss << "    MinN = " << minN_ << endl;
0082     ss << "    saveTags= " << saveTags();
0083     LogDebug("HLTMuonL1RegionalFilter") << ss.str();
0084   }
0085 }
0086 
0087 HLTMuonL1RegionalFilter::~HLTMuonL1RegionalFilter() = default;
0088 
0089 void HLTMuonL1RegionalFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0090   edm::ParameterSetDescription desc;
0091   makeHLTFilterDescription(desc);
0092   desc.add<edm::InputTag>("CandTag", edm::InputTag("hltL1extraParticles"));
0093   desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag("hltL1sL1SingleMu20"));
0094   desc.add<int>("MinN", 1);
0095 
0096   edm::ParameterSetDescription validator;
0097   std::vector<edm::ParameterSet> defaults(3);
0098 
0099   std::vector<double> etaRange;
0100   double minPt;
0101   std::vector<unsigned int> qualityBits;
0102 
0103   etaRange.clear();
0104   etaRange.push_back(-2.5);
0105   etaRange.push_back(+2.5);
0106   minPt = 20.0;
0107   qualityBits.clear();
0108   qualityBits.push_back(6);
0109   qualityBits.push_back(7);
0110   validator.add<std::vector<double> >("EtaRange", etaRange);
0111   validator.add<double>("MinPt", minPt);
0112   validator.add<std::vector<unsigned int> >("QualityBits", qualityBits);
0113 
0114   etaRange.clear();
0115   etaRange.push_back(-2.5);
0116   etaRange.push_back(-1.6);
0117   minPt = 20.0;
0118   qualityBits.clear();
0119   qualityBits.push_back(6);
0120   qualityBits.push_back(7);
0121   defaults[0].addParameter<std::vector<double> >("EtaRange", etaRange);
0122   defaults[0].addParameter<double>("MinPt", minPt);
0123   defaults[0].addParameter<std::vector<unsigned int> >("QualityBits", qualityBits);
0124 
0125   etaRange.clear();
0126   etaRange.push_back(-1.6);
0127   etaRange.push_back(+1.6);
0128   minPt = 20.0;
0129   qualityBits.clear();
0130   qualityBits.push_back(7);
0131   defaults[1].addParameter<std::vector<double> >("EtaRange", etaRange);
0132   defaults[1].addParameter<double>("MinPt", minPt);
0133   defaults[1].addParameter<std::vector<unsigned int> >("QualityBits", qualityBits);
0134 
0135   etaRange.clear();
0136   etaRange.push_back(+1.6);
0137   etaRange.push_back(+2.5);
0138   minPt = 20.0;
0139   qualityBits.clear();
0140   qualityBits.push_back(6);
0141   qualityBits.push_back(7);
0142   edm::ParameterSetDescription element2;
0143   defaults[2].addParameter<std::vector<double> >("EtaRange", etaRange);
0144   defaults[2].addParameter<double>("MinPt", minPt);
0145   defaults[2].addParameter<std::vector<unsigned int> >("QualityBits", qualityBits);
0146 
0147   desc.addVPSet("Cuts", validator, defaults);
0148 
0149   descriptions.add("hltMuonL1RegionalFilter", desc);
0150 }
0151 
0152 bool HLTMuonL1RegionalFilter::hltFilter(edm::Event& iEvent,
0153                                         const edm::EventSetup& iSetup,
0154                                         trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0155   using namespace std;
0156   using namespace edm;
0157   using namespace trigger;
0158   using namespace l1extra;
0159 
0160   // All HLT filters must create and fill an HLT filter object,
0161   // recording any reconstructed physics objects satisfying (or not)
0162   // this HLT filter, and place it in the Event.
0163 
0164   // get hold of all muons
0165   Handle<L1MuonParticleCollection> allMuons;
0166   iEvent.getByToken(candToken_, allMuons);
0167 
0168   // get hold of muons that fired the previous level
0169   Handle<TriggerFilterObjectWithRefs> previousLevelCands;
0170   iEvent.getByToken(previousCandToken_, previousLevelCands);
0171   vector<L1MuonParticleRef> prevMuons;
0172   previousLevelCands->getObjects(TriggerL1Mu, prevMuons);
0173 
0174   // look at all mucands,  check cuts and add to filter object
0175   int n = 0;
0176   for (size_t i = 0; i < allMuons->size(); i++) {
0177     L1MuonParticleRef muon(allMuons, i);
0178 
0179     //check if triggered by the previous level
0180     if (find(prevMuons.begin(), prevMuons.end(), muon) == prevMuons.end())
0181       continue;
0182 
0183     //find eta region, continue otherwise
0184     float eta = muon->eta();
0185     int region = -1;
0186     for (size_t r = 0; r < etaBoundaries_.size() - 1; r++) {
0187       if (etaBoundaries_[r] <= eta && eta <= etaBoundaries_[r + 1]) {
0188         region = r;
0189         break;
0190       }
0191     }
0192     if (region == -1)
0193       continue;
0194 
0195     //check pT cut
0196     if (muon->pt() < minPts_[region])
0197       continue;
0198 
0199     //check quality cut
0200     if (qualityBitMasks_[region]) {
0201       int quality = muon->gmtMuonCand().empty() ? 0 : (1 << muon->gmtMuonCand().quality());
0202       if ((quality & qualityBitMasks_[region]) == 0)
0203         continue;
0204     }
0205 
0206     //we have a good candidate
0207     filterproduct.addObject(TriggerL1Mu, muon);
0208 
0209     n++;
0210   }
0211 
0212   if (saveTags())
0213     filterproduct.addCollectionTag(candTag_);
0214 
0215   // filter decision
0216   const bool accept(n >= minN_);
0217 
0218   // dump event for debugging
0219   if (edm::isDebugEnabled()) {
0220     ostringstream ss;
0221     ss.precision(2);
0222     ss << "L1mu#" << '\t' << "q*pt" << '\t' << '\t' << "eta" << '\t' << "phi" << '\t' << "quality" << '\t' << "isPrev"
0223        << '\t' << "isFired" << endl;
0224     ss << "---------------------------------------------------------------" << endl;
0225 
0226     vector<L1MuonParticleRef> firedMuons;
0227     filterproduct.getObjects(TriggerL1Mu, firedMuons);
0228     for (size_t i = 0; i < allMuons->size(); i++) {
0229       L1MuonParticleRef mu(allMuons, i);
0230       int quality = mu->gmtMuonCand().empty() ? 0 : mu->gmtMuonCand().quality();
0231       bool isPrev = find(prevMuons.begin(), prevMuons.end(), mu) != prevMuons.end();
0232       bool isFired = find(firedMuons.begin(), firedMuons.end(), mu) != firedMuons.end();
0233       ss << i << '\t' << scientific << mu->charge() * mu->pt() << fixed << '\t' << mu->eta() << '\t' << mu->phi()
0234          << '\t' << quality << '\t' << isPrev << '\t' << isFired << endl;
0235     }
0236     ss << "---------------------------------------------------------------" << endl;
0237     LogDebug("HLTMuonL1RegionalFilter") << ss.str() << "Decision of filter is " << accept
0238                                         << ", number of muons passing = " << filterproduct.l1muonSize();
0239   }
0240 
0241   return accept;
0242 }
0243 
0244 // declare this class as a framework plugin
0245 #include "FWCore/Framework/interface/MakerMacros.h"
0246 DEFINE_FWK_MODULE(HLTMuonL1RegionalFilter);