CMS 3D CMS Logo

HLTMuonL1TFilter.cc
Go to the documentation of this file.
1 /* \class HLTMuonL1TFilter
2  *
3  * This is an HLTFilter implementing filtering on L1T Stage2 GMT objects
4  *
5  * \author: V. Rekovic
6  *
7  */
8 
9 #include "HLTMuonL1TFilter.h"
10 
14 
16 
18 
20 #include "TMath.h"
21 
25 
26 #include <vector>
27 
29  candTag_( iConfig.getParameter<edm::InputTag>("CandTag") ),
30  candToken_( consumes<l1t::MuonBxCollection>(candTag_)),
31  previousCandTag_( iConfig.getParameter<edm::InputTag>("PreviousCandTag") ),
32  previousCandToken_( consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
33  maxEta_( iConfig.getParameter<double>("MaxEta") ),
34  minPt_( iConfig.getParameter<double>("MinPt") ),
35  minN_( iConfig.getParameter<int>("MinN") ),
36  centralBxOnly_( iConfig.getParameter<bool>("CentralBxOnly") )
37 {
38 
39  //set the quality bit mask
40  qualityBitMask_ = 0;
41  vector<int> selectQualities = iConfig.getParameter<vector<int> >("SelectQualities");
42  for(int selectQualitie : selectQualities){
43 // if(selectQualities[i] > 7){ // FIXME: this will be updated once we have info from L1
44 // throw edm::Exception(edm::errors::Configuration) << "QualityBits must be smaller than 8!";
45 // }
46  qualityBitMask_ |= 1<<selectQualitie;
47  }
48 
49  // dump parameters for debugging
50  if(edm::isDebugEnabled()){
51  ostringstream ss;
52  ss << "Constructed with parameters:" << endl;
53  ss << " CandTag = " << candTag_.encode() << endl;
54  ss << " PreviousCandTag = " << previousCandTag_.encode()<< endl;
55  ss << " MaxEta = " << maxEta_ << endl;
56  ss << " MinPt = " << minPt_ << endl;
57  ss << " SelectQualities =";
58 // for(size_t i = 0; i < 8; i++){
59 // if((qualityBitMask_>>i) % 2) ss << " " << i;
60 // }
61  ss << endl;
62  ss << " MinN = " << minN_ << endl;
63  ss << " saveTags= " << saveTags();
64  LogDebug("HLTMuonL1TFilter") << ss.str();
65  }
66 
67 
68 }
69 
70 
72 
73 
74 
75 void
79  desc.add<edm::InputTag>("CandTag",edm::InputTag("hltGmtStage2Digis"));
80  desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag(""));
81  desc.add<double>("MaxEta",2.5);
82  desc.add<double>("MinPt",0.0);
83  desc.add<int>("MinN",1);
84  desc.add<bool>("CentralBxOnly", true);
85  {
86  std::vector<int> temp1;
87  temp1.reserve(0);
88  desc.add<std::vector<int> >("SelectQualities",temp1);
89  }
90  descriptions.add("hltMuonL1TFilter",desc);
91 }
92 
94  using namespace std;
95  using namespace edm;
96  using namespace trigger;
97  using namespace l1t;
98 
99  // All HLT filters must create and fill an HLT filter object,
100  // recording any reconstructed physics objects satisfying (or not)
101  // this HLT filter, and place it in the Event.
102 
103  // get hold of all muons
105  iEvent.getByToken(candToken_, allMuons);
106 
107  // get hold of TFOWRs that fired the previous level
108  Handle<TriggerFilterObjectWithRefs> previousLevelTFOWR;
109  iEvent.getByToken(previousCandToken_, previousLevelTFOWR);
110 
111  vector<MuonRef> prevMuons;
112  previousLevelTFOWR->getObjects(TriggerL1Mu, prevMuons);
113 
114  // look at all muon candidates, check cuts and add to filter object
115  int n = 0;
116 
117  for (int ibx = allMuons->getFirstBX(); ibx <= allMuons->getLastBX(); ++ibx) {
118  if (centralBxOnly_ && (ibx != 0)) continue;
119  for (auto it = allMuons->begin(ibx); it != allMuons->end(ibx); it++){
120 
121  MuonRef muon(allMuons, distance(allMuons->begin(allMuons->getFirstBX()),it) );
122 
123  // Only select muons that were selected in the previous level
124  if(find(prevMuons.begin(), prevMuons.end(), muon) == prevMuons.end()) continue;
125 
126  //check maxEta cut
127  if(fabs(muon->eta()) > maxEta_) continue;
128 
129  //check pT cut
130  if(muon->pt() < minPt_) continue;
131 
132  //check quality cut
133  if(qualityBitMask_){
134  int quality = (it->hwQual() == 0 ? 0 : (1 << it->hwQual()));
135  if((quality & qualityBitMask_) == 0) continue;
136  }
137 
138  //we have a good candidate
139  n++;
140  filterproduct.addObject(TriggerL1Mu,muon);
141  }
142  }
143 
144 
145 
146  if (saveTags()) filterproduct.addCollectionTag(candTag_);
147 
148  // filter decision
149  const bool accept(n >= minN_);
150 
151  // dump event for debugging
152  if(edm::isDebugEnabled()){
153 
154  LogTrace("HLTMuonL1TFilter")<< "\nHLTMuonL1TFilter -----------------------------------------------" << endl;
155  LogTrace("HLTMuonL1TFilter")<< "L1mu#" << '\t'
156  << "q" << '\t' << "pt" << '\t' << '\t'
157  << "eta" << '\t' << "phi" << '\t'
158  << "quality" << '\t' << "isPrev\t (|maxEta| = " << maxEta_ << ")" << endl;
159  LogTrace("HLTMuonL1TFilter")<< "--------------------------------------------------------------------------" << endl;
160 
161  vector<MuonRef> firedMuons;
162  filterproduct.getObjects(TriggerL1Mu, firedMuons);
163  for(size_t i=0; i<firedMuons.size(); i++){
164  l1t::MuonRef mu = firedMuons[i];
165  bool isPrev = find(prevMuons.begin(), prevMuons.end(), mu) != prevMuons.end();
166  LogTrace("HLTMuonL1TFilter")<< i << '\t' << setprecision(2) << scientific
167  << mu->charge() << '\t' << mu->pt() << '\t' << fixed
168  << mu->eta() << '\t' << mu->phi() << '\t'
169  << mu->hwQual() << '\t' << isPrev << endl;
170  }
171  LogTrace("HLTMuonL1TFilter")<< "--------------------------------------------------------------------------" << endl;
172  LogTrace("HLTMuonL1TFilter")<< "Decision of this filter is " << accept << ", number of muons passing = " << filterproduct.l1tmuonSize();
173  }
174 
175  return accept;
176 }
177 
178 // declare this class as a framework plugin
#define LogDebug(id)
const_iterator end(int bx) const
T getParameter(std::string const &) const
bool isDebugEnabled()
HLTMuonL1TFilter(const edm::ParameterSet &)
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
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
double maxEta_
max Eta cut
edm::InputTag previousCandTag_
input tag identifying the product containing refs to muons passing the previous level ...
delete x;
Definition: CaloConfig.h:22
double minN_
min N objects
bool centralBxOnly_
use central bx only muons
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
~HLTMuonL1TFilter() override
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< l1t::MuonBxCollection > candToken_
const int mu
Definition: Constants.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define LogTrace(id)
BXVector< Muon > MuonBxCollection
Definition: Muon.h:11
int getFirstBX() const
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
double minPt_
pT threshold
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::InputTag candTag_
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
int getLastBX() const
const_iterator begin(int bx) const