CMS 3D CMS Logo

HLTMuonL1TRegionalFilter.cc
Go to the documentation of this file.
1 
9 
11  candTag_ (iConfig.getParameter<edm::InputTag>("CandTag") ),
12  candToken_(consumes<l1t::MuonBxCollection>(candTag_)),
13  previousCandTag_ (iConfig.getParameter<edm::InputTag>("PreviousCandTag") ),
14  previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
15  minN_( iConfig.getParameter<int>("MinN") ),
16  centralBxOnly_( iConfig.getParameter<bool>("CentralBxOnly") )
17 {
18  using namespace std;
19  using namespace edm;
20 
21  // read in the eta-range dependent parameters
22  const vector<ParameterSet> cuts = iConfig.getParameter<vector<ParameterSet> >("Cuts");
23  size_t ranges = cuts.size();
24  if(ranges==0){
25  throw edm::Exception(errors::Configuration) << "Please provide at least one PSet in the Cuts VPSet!";
26  }
27  etaBoundaries_.reserve(ranges+1);
28  minPts_.reserve(ranges);
29  qualityBitMasks_.reserve(ranges);
30  for(size_t i=0; i<ranges; i++){
31  //set the eta range
32  vector<double> etaRange = cuts[i].getParameter<vector<double> >("EtaRange");
33  if(etaRange.size() != 2 || etaRange[0] >= etaRange[1]){
34  throw edm::Exception(errors::Configuration) << "EtaRange must have two non-equal values in increasing order!";
35  }
36  if(i==0){
37  etaBoundaries_.push_back( etaRange[0] );
38  }else if(etaBoundaries_[i] != etaRange[0]){
39  throw edm::Exception(errors::Configuration) << "EtaRanges must be disjoint without gaps and listed in increasing eta order!";
40  }
41  etaBoundaries_.push_back( etaRange[1] );
42 
43  //set the minPt
44  minPts_.push_back( cuts[i].getParameter<double>("MinPt") );
45 
46  //set the quality bit masks
47  qualityBitMasks_.push_back( 0 );
48  vector<unsigned int> qualities = cuts[i].getParameter<vector<unsigned int> >("QualityBits");
49  for(unsigned int qualitie : qualities){
50 // if(7U < qualities[j]){ // qualities[j] >= 0, since qualities[j] is unsigned //FIXME: this will be updated once we have info from L1
51 // throw edm::Exception(errors::Configuration) << "QualityBits must be between 0 and 7 !";
52 // }
53  qualityBitMasks_[i] |= 1 << qualitie;
54  }
55  }
56 
57  // dump parameters for debugging
58  if(edm::isDebugEnabled()){
59  ostringstream ss;
60  ss<<"Constructed with parameters:"<<endl;
61  ss<<" CandTag = "<<candTag_.encode()<<endl;
62  ss<<" PreviousCandTag = "<<previousCandTag_.encode()<<endl;
63  ss<<" EtaBoundaries = \t"<<etaBoundaries_[0];
64  for(size_t i=1; i<etaBoundaries_.size(); i++){
65  ss<<'\t'<<etaBoundaries_[i];
66  }
67  ss<<endl;
68  ss<<" MinPts = \t "<<minPts_[0];
69  for(size_t i=1; i<minPts_.size(); i++){
70  ss<<"\t "<<minPts_[i];
71  }
72  ss<<endl;
73  ss<<" QualityBitMasks = \t "<<qualityBitMasks_[0];
74  for(size_t i=1; i<qualityBitMasks_.size(); i++){
75  ss<<"\t "<<qualityBitMasks_[i];
76  }
77  ss<<endl;
78  ss<<" MinN = "<<minN_<<endl;
79  ss<<" saveTags= "<<saveTags();
80  LogDebug("HLTMuonL1TRegionalFilter")<<ss.str();
81  }
82 }
83 
85 
86 void
90  desc.add<edm::InputTag>("CandTag",edm::InputTag("hltGmtStage2Digis"));
91  desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag("hltL1sL1SingleMu20"));
92  desc.add<int>("MinN",1);
93  desc.add<bool>("CentralBxOnly", true);
94 
96  std::vector<edm::ParameterSet> defaults(3);
97 
98  std::vector<double> etaRange;
99  double minPt;
100  std::vector<unsigned int> qualityBits;
101 
102  etaRange.clear();
103  etaRange.push_back(-2.5);
104  etaRange.push_back(+2.5);
105  minPt=20.0;
106  qualityBits.clear();
107  qualityBits.push_back(6);
108  qualityBits.push_back(7);
109  validator.add<std::vector<double> >("EtaRange",etaRange);
110  validator.add<double>("MinPt",minPt);
111  validator.add<std::vector<unsigned int> >("QualityBits",qualityBits);
112 
113 
114  etaRange.clear();
115  etaRange.push_back(-2.5);
116  etaRange.push_back(-1.6);
117  minPt=20.0;
118  qualityBits.clear();
119  qualityBits.push_back(6);
120  qualityBits.push_back(7);
121  defaults[0].addParameter<std::vector<double> >("EtaRange",etaRange);
122  defaults[0].addParameter<double>("MinPt",minPt);
123  defaults[0].addParameter<std::vector<unsigned int> >("QualityBits",qualityBits);
124 
125 
126  etaRange.clear();
127  etaRange.push_back(-1.6);
128  etaRange.push_back(+1.6);
129  minPt=20.0;
130  qualityBits.clear();
131  qualityBits.push_back(7);
132  defaults[1].addParameter<std::vector<double> >("EtaRange",etaRange);
133  defaults[1].addParameter<double>("MinPt",minPt);
134  defaults[1].addParameter<std::vector<unsigned int> >("QualityBits",qualityBits);
135 
136 
137  etaRange.clear();
138  etaRange.push_back(+1.6);
139  etaRange.push_back(+2.5);
140  minPt=20.0;
141  qualityBits.clear();
142  qualityBits.push_back(6);
143  qualityBits.push_back(7);
145  defaults[2].addParameter<std::vector<double> >("EtaRange",etaRange);
146  defaults[2].addParameter<double>("MinPt",minPt);
147  defaults[2].addParameter<std::vector<unsigned int> >("QualityBits",qualityBits);
148 
149  desc.addVPSet("Cuts",validator,defaults);
150 
151  descriptions.add("hltMuonL1TRegionalFilter", desc);
152 }
153 
155  using namespace std;
156  using namespace edm;
157  using namespace trigger;
158  using namespace l1t;
159 
160  // All HLT filters must create and fill an HLT filter object,
161  // recording any reconstructed physics objects satisfying (or not)
162  // this HLT filter, and place it in the Event.
163 
164  // get hold of all muons
166  iEvent.getByToken(candToken_, allMuons);
167 
168  // get hold of muons that fired the previous level
169  Handle<TriggerFilterObjectWithRefs> previousLevelCands;
170  iEvent.getByToken(previousCandToken_, previousLevelCands);
171 
172  vector<MuonRef> prevMuons;
173  previousLevelCands->getObjects(TriggerL1Mu, prevMuons);
174 
175  // look at all mucands, check cuts and add to filter object
176  int n = 0;
177 
178 
179  for (int ibx = allMuons->getFirstBX(); ibx <= allMuons->getLastBX(); ++ibx) {
180  if (centralBxOnly_ && (ibx != 0)) continue;
181  for (auto it = allMuons->begin(ibx); it != allMuons->end(ibx); it++){
182 
183  MuonRef muon(allMuons, distance(allMuons->begin(allMuons->getFirstBX()),it) );
184 
185  // Only select muons that were selected in the previous level
186  if(find(prevMuons.begin(), prevMuons.end(), muon) == prevMuons.end()) continue;
187 
188  //check maxEta cut
189  float eta = muon->eta();
190  int region = -1;
191  for(size_t r=0; r<etaBoundaries_.size()-1; r++){
192  if(etaBoundaries_[r]<=eta && eta<=etaBoundaries_[r+1]){
193  region = r;
194  break;
195  }
196  }
197  if(region == -1) continue;
198 
199  //check pT cut
200  if(muon->pt() < minPts_[region]) continue;
201 
202  //check quality cut
203  if(qualityBitMasks_[region]){
204  int quality = (it->hwQual() == 0 ? 0 : (1 << it->hwQual()));
205  if((quality & qualityBitMasks_[region]) == 0) continue;
206  }
207 
208  //we have a good candidate
209  n++;
210  filterproduct.addObject(TriggerL1Mu,muon);
211  }
212  }
213 
214 
215  if (saveTags()) filterproduct.addCollectionTag(candTag_);
216 
217  // filter decision
218  const bool accept (n >= minN_);
219 
220  // dump event for debugging
221  if(edm::isDebugEnabled()){
222 
223  LogTrace("HLTMuonL1TRegionalFilter")<< "\nHLTMuonL1TRegionalFilter -----------------------------------------------" << endl;
224  LogTrace("HLTMuonL1TRegionalFilter")<< "L1mu#" << '\t'
225  << "q*pt" << '\t' << '\t'
226  << "eta" << '\t' << "phi" << '\t'
227  << "quality" << '\t' << "isPrev\t " << endl;
228  LogTrace("HLTMuonL1TRegionalFilter")<< "--------------------------------------------------------------------------" << endl;
229 
230  vector<MuonRef> firedMuons;
231  filterproduct.getObjects(TriggerL1Mu, firedMuons);
232  for(size_t i=0; i<firedMuons.size(); i++){
233  l1t::MuonRef mu = firedMuons[i];
234  bool isPrev = find(prevMuons.begin(), prevMuons.end(), mu) != prevMuons.end();
235  LogTrace("HLTMuonL1TRegionalFilter")<< i << '\t' << setprecision(2) << scientific
236  << mu->charge()* mu->pt() << '\t' << fixed
237  << mu->eta() << '\t' << mu->phi() << '\t'
238  << mu->hwQual() << '\t' << isPrev << endl;
239  }
240  LogTrace("HLTMuonL1TRegionalFilter")<< "--------------------------------------------------------------------------" << endl;
241  LogTrace("HLTMuonL1TRegionalFilter")<< "Decision of this filter is " << accept << ", number of muons passing = " << filterproduct.l1tmuonSize();
242 
243  }
244 
245  return accept;
246 }
247 
248 
249 // declare this class as a framework plugin
#define LogDebug(id)
const_iterator end(int bx) const
T getParameter(std::string const &) const
bool isDebugEnabled()
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
~HLTMuonL1TRegionalFilter() override
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
HLTMuonL1TRegionalFilter(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
enum start value shifted to 81 so as to avoid clashes with PDG codes
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
delete x;
Definition: CaloConfig.h:22
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
std::string encode() const
Definition: InputTag.cc:159
bool centralBxOnly_
use central bx only muons
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
std::vector< double > minPts_
the vector of MinPt values, one for eta each region
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
def defaults(locpath, dataType, var)
int minN_
required number of passing candidates to pass the filter
const int mu
Definition: Constants.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
#define LogTrace(id)
string ranges
Definition: diffTwoXMLs.py:79
BXVector< Muon > MuonBxCollection
Definition: Muon.h:11
edm::EDGetTokenT< l1t::MuonBxCollection > candToken_
int getFirstBX() const
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< double > etaBoundaries_
the vector of eta region boundaries; note: # of boundaries = # of regions + 1
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::InputTag candTag_
input tag identifying the product containing muons
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
edm::InputTag previousCandTag_
input tag identifying the product containing refs to muons passing the previous level ...
int getLastBX() const
const_iterator begin(int bx) const