CMS 3D CMS Logo

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