CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTTauDQMPathPlotter.cc
Go to the documentation of this file.
5 
6 namespace {
7  std::string stripVersion(const std::string& pathName) {
8  size_t versionStart = pathName.rfind("_v");
9  if(versionStart == std::string::npos)
10  return pathName;
11  return pathName.substr(0, versionStart);
12  }
13 }
14 
16  bool doRefAnalysis, const std::string& dqmBaseFolder,
17  const std::string& hltProcess, int ptbins, int etabins, int phibins,
18  double ptmax, double highptmax,
19  double l1MatchDr, double hltMatchDr):
20  HLTTauDQMPlotter(stripVersion(pathName), dqmBaseFolder),
21  ptbins_(ptbins),
22  etabins_(etabins),
23  phibins_(phibins),
24  ptmax_(ptmax),
25  highptmax_(highptmax),
26  l1MatchDr_(l1MatchDr),
27  hltMatchDr_(hltMatchDr),
28  doRefAnalysis_(doRefAnalysis),
29  hltPath_(pathName, hltProcess, doRefAnalysis_, HLTCP)
30 {
32 }
33 
34 #include <algorithm>
36  if(!isValid())
37  return;
38 
39  // Book histograms
40  iBooker.setCurrentFolder(triggerTag());
41 
42  hAcceptedEvents_ = iBooker.book1D("EventsPerFilter", "Accepted Events per filter;;entries", hltPath_.filtersSize(), 0, hltPath_.filtersSize());
43  for(size_t i=0; i<hltPath_.filtersSize(); ++i) {
45  }
46 
47  // Efficiency helpers
48  if(doRefAnalysis_) {
49  iBooker.setCurrentFolder(triggerTag()+"/helpers");
50  if(hltPath_.hasL2Taus()) {
51  hL2TrigTauEtEffNum_ = iBooker.book1D("L2TrigTauEtEffNum", "L2 #tau p_{T} efficiency;Ref #tau p_{T};entries", ptbins_, 0, ptmax_);
52  hL2TrigTauEtEffDenom_ = iBooker.book1D("L2TrigTauEtEffDenom", "L2 #tau p_{T} denominator;Ref #tau p_{T};Efficiency", ptbins_, 0, ptmax_);
53  hL2TrigTauEtaEffNum_ = iBooker.book1D("L2TrigTauEtaEffNum", "L2 #tau #eta efficiency;Ref #tau #eta;entries", etabins_, -2.5, 2.5);
54  hL2TrigTauEtaEffDenom_ = iBooker.book1D("L2TrigTauEtaEffDenom", "L2 #tau #eta denominator;Ref #tau #eta;Efficiency", etabins_, -2.5, 2.5);
55  hL2TrigTauPhiEffNum_ = iBooker.book1D("L2TrigTauPhiEffNum", "L2 #tau #phi efficiency;Ref #tau #phi;entries", phibins_, -3.2, 3.2);
56  hL2TrigTauPhiEffDenom_ = iBooker.book1D("L2TrigTauPhiEffDenom", "L2 #tau #phi denominator;Ref #tau #phi;Efficiency", phibins_, -3.2, 3.2);
57  hL2TrigTauHighEtEffNum_ = iBooker.book1D("L2TrigTauHighEtEffNum", "L2 #tau p_{T} efficiency (high p_{T});Ref #tau p_{T};entries", ptbins_, 0, highptmax_);
58  hL2TrigTauHighEtEffDenom_ = iBooker.book1D("L2TrigTauHighEtEffDenom", "L2 #tau p_{T} denominator (high p_{T});Ref #tau p_{T};Efficiency", ptbins_, 0, highptmax_);
59  }
60 
61  if(hltPath_.hasL3Taus()) {
62  hL3TrigTauEtEffNum_ = iBooker.book1D("L3TrigTauEtEffNum", "L3 #tau p_{T} efficiency;Ref #tau p_{T};entries", ptbins_, 0, ptmax_);
63  hL3TrigTauEtEffDenom_ = iBooker.book1D("L3TrigTauEtEffDenom", "L3 #tau p_{T} denominator;Ref #tau p_{T};Efficiency", ptbins_, 0, ptmax_);
64  hL3TrigTauEtaEffNum_ = iBooker.book1D("L3TrigTauEtaEffNum", "L3 #tau #eta efficiency;Ref #tau #eta;entries", etabins_, -2.5, 2.5);
65  hL3TrigTauEtaEffDenom_ = iBooker.book1D("L3TrigTauEtaEffDenom", "L3 #tau #eta denominator;Ref #tau #eta;Efficiency", etabins_, -2.5, 2.5);
66  hL3TrigTauPhiEffNum_ = iBooker.book1D("L3TrigTauPhiEffNum", "L3 #tau #phi efficiency;Ref #tau #phi;entries", phibins_, -3.2, 3.2);
67  hL3TrigTauPhiEffDenom_ = iBooker.book1D("L3TrigTauPhiEffDenom", "L3 #tau #phi denominator;Ref #tau #phi;Efficiency", phibins_, -3.2, 3.2);
68  hL3TrigTauHighEtEffNum_ = iBooker.book1D("L3TrigTauHighEtEffNum", "L3 #tau p_{T} efficiency (high p_{T});Ref #tau p_{T};entries", ptbins_, 0, highptmax_);
69  hL3TrigTauHighEtEffDenom_ = iBooker.book1D("L3TrigTauHighEtEffDenom", "L3 #tau p_{T} denominator (high p_{T});Ref #tau p_{T};Efficiency", ptbins_, 0, highptmax_);
70  }
71 
72  if(hltPath_.hasL2Electrons()) {
73  hL2TrigElectronEtEffNum_ = iBooker.book1D("L2TrigElectronEtEffNum", "L2 electron p_{T} efficiency;Ref electron p_{T};entries", ptbins_, 0, ptmax_);
74  hL2TrigElectronEtEffDenom_ = iBooker.book1D("L2TrigElectronEtEffDenom", "L2 electron p_{T} denominator;Ref electron p_{T};Efficiency", ptbins_, 0, ptmax_);
75  hL2TrigElectronEtaEffNum_ = iBooker.book1D("L2TrigElectronEtaEffNum", "L2 electron #eta efficiency;Ref electron #eta;entries", etabins_, -2.5, 2.5);
76  hL2TrigElectronEtaEffDenom_ = iBooker.book1D("L2TrigElectronEtaEffDenom", "L2 electron #eta denominator;Ref electron #eta;Efficiency", etabins_, -2.5, 2.5);
77  hL2TrigElectronPhiEffNum_ = iBooker.book1D("L2TrigElectronPhiEffNum", "L2 electron #phi efficiency;Ref electron #phi;entries", phibins_, -3.2, 3.2);
78  hL2TrigElectronPhiEffDenom_ = iBooker.book1D("L2TrigElectronPhiEffDenom", "L2 electron #phi denominator;Ref electron #phi;Efficiency", phibins_, -3.2, 3.2);
79  }
80 
81  if(hltPath_.hasL3Electrons()) {
82  hL3TrigElectronEtEffNum_ = iBooker.book1D("L3TrigElectronEtEffNum", "L3 electron p_{T} efficiency;Ref electron p_{T};entries", ptbins_, 0, ptmax_);
83  hL3TrigElectronEtEffDenom_ = iBooker.book1D("L3TrigElectronEtEffDenom", "L3 electron p_{T} denominator;Ref electron p_{T};Efficiency", ptbins_, 0, ptmax_);
84  hL3TrigElectronEtaEffNum_ = iBooker.book1D("L3TrigElectronEtaEffNum", "L3 electron #eta efficiency;Ref electron #eta;entries", etabins_, -2.5, 2.5);
85  hL3TrigElectronEtaEffDenom_ = iBooker.book1D("L3TrigElectronEtaEffDenom", "L3 electron #eta denominator;Ref electron #eta;Efficiency", etabins_, -2.5, 2.5);
86  hL3TrigElectronPhiEffNum_ = iBooker.book1D("L3TrigElectronPhiEffNum", "L3 electron #phi efficiency;Ref electron #phi;entries", phibins_, -3.2, 3.2);
87  hL3TrigElectronPhiEffDenom_ = iBooker.book1D("L3TrigElectronPhiEffDenom", "L3 electron #phi denominator;Ref electron #phi;Efficiency", phibins_, -3.2, 3.2);
88  }
89 
90  if(hltPath_.hasL2Muons()) {
91  hL2TrigMuonEtEffNum_ = iBooker.book1D("L2TrigMuonEtEffNum", "L2 muon p_{T} efficiency;Ref muon p_{T};entries", ptbins_, 0, ptmax_);
92  hL2TrigMuonEtEffDenom_ = iBooker.book1D("L2TrigMuonEtEffDenom", "L2 muon p_{T} denominator;Ref muon p_{T};Efficiency", ptbins_, 0, ptmax_);
93  hL2TrigMuonEtaEffNum_ = iBooker.book1D("L2TrigMuonEtaEffNum", "L2 muon #eta efficiency;Ref muon #eta;entries", etabins_, -2.5, 2.5);
94  hL2TrigMuonEtaEffDenom_ = iBooker.book1D("L2TrigMuonEtaEffDenom", "L2 muon #eta denominator;Ref muon #eta;Efficiency", etabins_, -2.5, 2.5);
95  hL2TrigMuonPhiEffNum_ = iBooker.book1D("L2TrigMuonPhiEffNum", "L2 muon #phi efficiency;Ref muon #phi;entries", phibins_, -3.2, 3.2);
96  hL2TrigMuonPhiEffDenom_ = iBooker.book1D("L2TrigMuonPhiEffDenom", "L2 muon #phi denominator;Ref muon #phi;Efficiency", phibins_, -3.2, 3.2);
97  }
98 
99  if(hltPath_.hasL3Muons()) {
100  hL3TrigMuonEtEffNum_ = iBooker.book1D("L3TrigMuonEtEffNum", "L3 muon p_{T} efficiency;Ref muon p_{T};entries", ptbins_, 0, ptmax_);
101  hL3TrigMuonEtEffDenom_ = iBooker.book1D("L3TrigMuonEtEffDenom", "L3 muon p_{T} denominator;Ref muon p_{T};Efficiency", ptbins_, 0, ptmax_);
102  hL3TrigMuonEtaEffNum_ = iBooker.book1D("L3TrigMuonEtaEffNum", "L3 muon #eta efficiency;Ref muon #eta;entries", etabins_, -2.5, 2.5);
103  hL3TrigMuonEtaEffDenom_ = iBooker.book1D("L3TrigMuonEtaEffDenom", "L3 muon #eta denominator;Ref muon #eta;Efficiency", etabins_, -2.5, 2.5);
104  hL3TrigMuonPhiEffNum_ = iBooker.book1D("L3TrigMuonPhiEffNum", "L3 muon #phi efficiency;Ref muon #phi;entries", phibins_, -3.2, 3.2);
105  hL3TrigMuonPhiEffDenom_ = iBooker.book1D("L3TrigMuonPhiEffDenom", "L3 muon #phi denominator;Ref muon #phi;Efficiency", phibins_, -3.2, 3.2);
106  }
107 
108  if(hltPath_.hasL2CaloMET()) {
109  hL2TrigMETEtEffNum_ = iBooker.book1D("L2TrigMETEtEffNum", "L2 MET efficiency;Ref MET;entries", ptbins_, 0, ptmax_);
110  hL2TrigMETEtEffDenom_ = iBooker.book1D("L2TrigMETEtEffDenom", "L2 MET denominator;Ref MET;Efficiency", ptbins_, 0, ptmax_);
111  }
112 
113  iBooker.setCurrentFolder(triggerTag());
114  }
115 
116  // Book di-object invariant mass histogram only for mu+tau, ele+tau, and di-tau paths
117  hMass_ = nullptr;
118  if(doRefAnalysis_) {
122 
123  int nmet = 0;
124  int lastMatchedMETFilter = -1;
125  for(size_t i = 0; i < hltPath_.filtersSize(); ++i){
126  if(hltPath_.getFilterName(i).find("hltMET") < hltPath_.getFilterName(i).length()) lastMatchedMETFilter = i;
127  }
128  if(lastMatchedMETFilter >= 0) nmet = hltPath_.getFilterMET(lastMatchedMETFilter);
129  auto create = [&](const std::string& name) {
130  if(name == "tau-met"){
131  this->hMass_ = iBooker.book1D("ReferenceMass", "Transverse mass of reference "+name+";Reference transverse mass;entries", 100, 0, 500);
132  }else{
133  this->hMass_ = iBooker.book1D("ReferenceMass", "Invariant mass of reference "+name+";Reference invariant mass;entries", 100, 0, 500);
134  }
135  };
136  LogDebug("HLTTauDQMOffline") << "Path " << hltPath_.getPathName() << " number of taus " << ntaus << " electrons " << neles << " muons " << nmus;
137  if(ntaus > 0) {
138  hTrigTauEt_ = iBooker.book1D("TrigTauEt", "Triggered #tau p_{T};#tau p_{T};entries", ptbins_, 0, ptmax_);
139  hTrigTauEta_ = iBooker.book1D("TrigTauEta", "Triggered #tau #eta;#tau #eta;entries", etabins_, -2.5, 2.5);
140  hTrigTauPhi_ = iBooker.book1D("TrigTauPhi", "Triggered #tau #phi;#tau #phi;entries", phibins_, -3.2, 3.2);
141  }
142  if(neles > 0) {
143  hTrigElectronEt_ = iBooker.book1D("TrigElectronEt", "Triggered electron p_{T};electron p_{T};entries", ptbins_, 0, ptmax_);
144  hTrigElectronEta_ = iBooker.book1D("TrigElectronEta", "Triggered electron #eta;electron #eta;entries", etabins_, -2.5, 2.5);
145  hTrigElectronPhi_ = iBooker.book1D("TrigElectronPhi", "Triggered electron #phi;electron #phi;entries", phibins_, -3.2, 3.2);
146  }
147  if(nmus > 0) {
148  hTrigMuonEt_ = iBooker.book1D("TrigMuonEt", "Triggered muon p_{T};muon p_{T};entries", ptbins_, 0, ptmax_);
149  hTrigMuonEta_ = iBooker.book1D("TrigMuonEta", "Triggered muon #eta;muon #eta;entries", etabins_, -2.5, 2.5);
150  hTrigMuonPhi_ = iBooker.book1D("TrigMuonPhi", "Triggered muon #phi;muon #phi;entries", phibins_, -3.2, 3.2);
151  }
152  if(nmet > 0) {
153  hTrigMETEt_ = iBooker.book1D("TrigMETEt", "Triggered MET E_{T};MET E_{T};entries", ptbins_, 0, ptmax_);
154  hTrigMETPhi_ = iBooker.book1D("TrigMETPhi", "Triggered MET #phi;MET #phi;entries", phibins_, -3.2, 3.2);
155  }
156  if(ntaus == 2 && neles == 0 && nmus == 0 && nmet == 0) create("di-tau");
157  if(ntaus == 1 && neles == 1 && nmus == 0 && nmet == 0) create("electron-tau");
158  if(ntaus == 1 && neles == 0 && nmus == 1 && nmet == 0) create("muon-tau");
159  if(ntaus == 1 && neles == 0 && nmus == 0 && nmet == 1) create("tau-met");
160  }
161 }
162 
163 
165 
167  std::vector<HLTTauDQMPath::Object> triggerObjs;
168  std::vector<HLTTauDQMPath::Object> matchedTriggerObjs;
169  HLTTauDQMOfflineObjects matchedOfflineObjs;
170 
171  // Events per filter
172  const int lastPassedFilter = hltPath_.lastPassedFilter(triggerResults);
173  int lastMatchedFilter = -1;
174  int lastMatchedMETFilter = -1;
175  int lastMatchedElectronFilter = -1;
176  int lastMatchedMuonFilter = -1;
177  int lastMatchedTauFilter = -1;
178  int firstMatchedMETFilter = -1;
179 
180  if(doRefAnalysis_) {
181  double matchDr = hltPath_.isFirstFilterL1Seed() ? l1MatchDr_ : hltMatchDr_;
182  for(int i=0; i<=lastPassedFilter; ++i) {
183  triggerObjs.clear();
184  matchedTriggerObjs.clear();
185  matchedOfflineObjs.clear();
186  hltPath_.getFilterObjects(triggerEvent, i, triggerObjs);
187  //std::cout << "Filter name " << hltPath_.getFilterName(i) << " nobjs " << triggerObjs.size() << " " << "ref size " << refCollection.taus.size() << std::endl;
188  bool matched = hltPath_.offlineMatching(i, triggerObjs, refCollection, matchDr, matchedTriggerObjs, matchedOfflineObjs);
189  //std::cout << " offline matching: " << matched << " " << matchedTriggerObjs.size() << std::endl;
190  matchDr = hltMatchDr_;
191  if(!matched)
192  break;
193 
194  hAcceptedEvents_->Fill(i+0.5);
195  lastMatchedFilter = i;
196  if(hltPath_.getFilterName(i).find("hltMET") < hltPath_.getFilterName(i).length()) lastMatchedMETFilter = i;
197  if(hltPath_.getFilterType(i) == "HLTMuonL3PreFilter" || hltPath_.getFilterType(i) == "HLTMuonIsoFilter") lastMatchedMuonFilter = i;
198  if(hltPath_.getFilterName(i).find("hltEle") < hltPath_.getFilterName(i).length()) lastMatchedElectronFilter = i;
199  if(hltPath_.getFilterName(i).find("hltPFTau") < hltPath_.getFilterName(i).length() ||
200  hltPath_.getFilterName(i).find("hltDoublePFTau") < hltPath_.getFilterName(i).length()) lastMatchedTauFilter = i;
201  if(firstMatchedMETFilter < 0 && hltPath_.getFilterName(i).find("hltMET") < hltPath_.getFilterName(i).length()) firstMatchedMETFilter = i;
202  }
203  }
204  else {
205  for(int i=0; i<=lastPassedFilter; ++i) {
206  hAcceptedEvents_->Fill(i+0.5);
207  }
208  }
209 
210  // Efficiency plots
211  if(doRefAnalysis_ && lastMatchedFilter >= 0) {
212  // L2 taus
213  if(hltPath_.hasL2Taus()) {
214  // Denominators
215  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL2TauIndex()) {
216  for(const LV& tau: refCollection.taus) {
217  hL2TrigTauEtEffDenom_->Fill(tau.pt());
218  hL2TrigTauHighEtEffDenom_->Fill(tau.pt());
219  hL2TrigTauEtaEffDenom_->Fill(tau.eta());
220  hL2TrigTauPhiEffDenom_->Fill(tau.phi());
221  }
222  }
223 
224  // Numerators
225  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL2TauFilterIndex()) {
226  triggerObjs.clear();
227  matchedTriggerObjs.clear();
228  matchedOfflineObjs.clear();
229  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL2TauFilterIndex(), triggerObjs);
230  bool matched = hltPath_.offlineMatching(hltPath_.getLastL2TauFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
231  if(matched) {
232  for(const LV& tau: matchedOfflineObjs.taus) {
233  hL2TrigTauEtEffNum_->Fill(tau.pt());
234  hL2TrigTauHighEtEffNum_->Fill(tau.pt());
235  hL2TrigTauEtaEffNum_->Fill(tau.eta());
236  hL2TrigTauPhiEffNum_->Fill(tau.phi());
237  }
238  }
239  }
240  }
241 
242  // L3 taus
243  if(hltPath_.hasL3Taus()) {
244  // Denominators
245  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL3TauIndex()) {
246  for(const LV& tau: refCollection.taus) {
247  hL3TrigTauEtEffDenom_->Fill(tau.pt());
248  hL3TrigTauHighEtEffDenom_->Fill(tau.pt());
249  hL3TrigTauEtaEffDenom_->Fill(tau.eta());
250  hL3TrigTauPhiEffDenom_->Fill(tau.phi());
251  }
252  }
253 
254  // Numerators
255  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL3TauFilterIndex()) {
256  triggerObjs.clear();
257  matchedTriggerObjs.clear();
258  matchedOfflineObjs.clear();
259  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL3TauFilterIndex(), triggerObjs);
260  bool matched = hltPath_.offlineMatching(hltPath_.getLastL3TauFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
261  if(matched) {
262  for(const LV& tau: matchedOfflineObjs.taus) {
263  hL3TrigTauEtEffNum_->Fill(tau.pt());
264  hL3TrigTauHighEtEffNum_->Fill(tau.pt());
265  hL3TrigTauEtaEffNum_->Fill(tau.eta());
266  hL3TrigTauPhiEffNum_->Fill(tau.phi());
267  }
268  }
269  }
270  }
271 
272  // L2 Electrons
273  if(hltPath_.hasL2Electrons()) {
274  // Denominators
275  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL2ElectronIndex()) {
276  for(const LV& electron: refCollection.electrons) {
277  hL2TrigElectronEtEffDenom_->Fill(electron.pt());
278  hL2TrigElectronEtaEffDenom_->Fill(electron.eta());
279  hL2TrigElectronPhiEffDenom_->Fill(electron.phi());
280  }
281  }
282 
283  // Numerators
284  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL2ElectronFilterIndex()) {
285  triggerObjs.clear();
286  matchedTriggerObjs.clear();
287  matchedOfflineObjs.clear();
288  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL2ElectronFilterIndex(), triggerObjs);
289  bool matched = hltPath_.offlineMatching(hltPath_.getLastL2ElectronFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
290  if(matched) {
291  for(const LV& electron: matchedOfflineObjs.electrons) {
292  hL2TrigElectronEtEffNum_->Fill(electron.pt());
293  hL2TrigElectronEtaEffNum_->Fill(electron.eta());
294  hL2TrigElectronPhiEffNum_->Fill(electron.phi());
295  }
296  }
297  }
298  }
299 
300  // L3 electron
301  if(hltPath_.hasL3Electrons()) {
302  // Denominators
303  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL3ElectronIndex()) {
304  for(const LV& electron: refCollection.electrons) {
305  hL3TrigElectronEtEffDenom_->Fill(electron.pt());
306  hL3TrigElectronEtaEffDenom_->Fill(electron.eta());
307  hL3TrigElectronPhiEffDenom_->Fill(electron.phi());
308  }
309  }
310 
311  // Numerators
312  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL3ElectronFilterIndex()) {
313  triggerObjs.clear();
314  matchedTriggerObjs.clear();
315  matchedOfflineObjs.clear();
316  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL3ElectronFilterIndex(), triggerObjs);
317  bool matched = hltPath_.offlineMatching(hltPath_.getLastL3ElectronFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
318  if(matched) {
319  for(const LV& electron: matchedOfflineObjs.electrons) {
320  hL3TrigElectronEtEffNum_->Fill(electron.pt());
321  hL3TrigElectronEtaEffNum_->Fill(electron.eta());
322  hL3TrigElectronPhiEffNum_->Fill(electron.phi());
323  }
324  }
325  }
326  }
327 
328  // L2 Muons
329  if(hltPath_.hasL2Muons()) {
330  // Denominators
331  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL2MuonIndex()) {
332  for(const LV& muon: refCollection.muons) {
333  hL2TrigMuonEtEffDenom_->Fill(muon.pt());
334  hL2TrigMuonEtaEffDenom_->Fill(muon.eta());
335  hL2TrigMuonPhiEffDenom_->Fill(muon.phi());
336  }
337  }
338 
339  // Numerators
340  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL2MuonFilterIndex()) {
341  triggerObjs.clear();
342  matchedTriggerObjs.clear();
343  matchedOfflineObjs.clear();
344  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL2MuonFilterIndex(), triggerObjs);
345  bool matched = hltPath_.offlineMatching(hltPath_.getLastL2MuonFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
346  if(matched) {
347  for(const LV& muon: matchedOfflineObjs.muons) {
348  hL2TrigMuonEtEffNum_->Fill(muon.pt());
349  hL2TrigMuonEtaEffNum_->Fill(muon.eta());
350  hL2TrigMuonPhiEffNum_->Fill(muon.phi());
351  }
352  }
353  }
354  }
355 
356  // L3 muon
357  if(hltPath_.hasL3Muons()) {
358  // Denominators
359  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL3MuonIndex()) {
360  for(const LV& muon: refCollection.muons) {
361  hL3TrigMuonEtEffDenom_->Fill(muon.pt());
362  hL3TrigMuonEtaEffDenom_->Fill(muon.eta());
363  hL3TrigMuonPhiEffDenom_->Fill(muon.phi());
364  }
365  }
366 
367  // Numerators
368  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL3MuonFilterIndex()) {
369  triggerObjs.clear();
370  matchedTriggerObjs.clear();
371  matchedOfflineObjs.clear();
372  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL3MuonFilterIndex(), triggerObjs);
373  bool matched = hltPath_.offlineMatching(hltPath_.getLastL3MuonFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
374  if(matched) {
375  for(const LV& muon: matchedOfflineObjs.muons) {
376  hL3TrigMuonEtEffNum_->Fill(muon.pt());
377  hL3TrigMuonEtaEffNum_->Fill(muon.eta());
378  hL3TrigMuonPhiEffNum_->Fill(muon.phi());
379  }
380  }
381  }
382  }
383 
384  // L2 CaloMET
385  if(hltPath_.hasL2CaloMET()) {
386  // Denominators
387  if(static_cast<size_t>(firstMatchedMETFilter) >= hltPath_.getFirstFilterBeforeL2CaloMETIndex()) {
388  hL2TrigMETEtEffDenom_->Fill(refCollection.met[0].pt());
389  }
390 
391  // Numerators
392  if(static_cast<size_t>(lastMatchedMETFilter) >= hltPath_.getLastL2CaloMETFilterIndex()) {
393  triggerObjs.clear();
394  matchedTriggerObjs.clear();
395  matchedOfflineObjs.clear();
396  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL2CaloMETFilterIndex(), triggerObjs);
397  bool matched = hltPath_.offlineMatching(hltPath_.getLastL2CaloMETFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
398  if(matched) {
399  hL2TrigMETEtEffNum_->Fill(matchedOfflineObjs.met[0].pt());
400  }
401  }
402  }
403  }
404 
405  if(hltPath_.fired(triggerResults)) {
406  triggerObjs.clear();
407  matchedTriggerObjs.clear();
408  matchedOfflineObjs.clear();
409 
410  if(lastMatchedMETFilter >= 0) hltPath_.getFilterObjects(triggerEvent, lastMatchedMETFilter, triggerObjs);
411  if(lastMatchedMuonFilter >= 0) hltPath_.getFilterObjects(triggerEvent, lastMatchedMuonFilter, triggerObjs);
412  if(lastMatchedElectronFilter >= 0) hltPath_.getFilterObjects(triggerEvent, lastMatchedElectronFilter, triggerObjs);
413 
414  if(lastMatchedTauFilter >= 0) hltPath_.getFilterObjects(triggerEvent, lastMatchedTauFilter, triggerObjs);
415 
416  if(doRefAnalysis_) {
417  bool matched = hltPath_.offlineMatching(lastPassedFilter, triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
418  if(matched) {
419  // Di-object invariant mass
420  if(hMass_) {
421  const int ntaus = hltPath_.getFilterNTaus(lastPassedFilter);
422  if(ntaus == 2 && hltPath_.getFilterNElectrons(lastMatchedElectronFilter) == 0 && hltPath_.getFilterNMuons(lastMatchedMuonFilter) == 0) {
423  // Di-tau (matchedOfflineObjs are already sorted)
424  hMass_->Fill( (matchedOfflineObjs.taus[0]+matchedOfflineObjs.taus[1]).M() );
425  }
426  // Electron+tau
427  else if(ntaus == 1 && hltPath_.getFilterNElectrons(lastPassedFilter) == 1) {
428  hMass_->Fill( (matchedOfflineObjs.taus[0]+matchedOfflineObjs.electrons[0]).M() );
429  }
430  // Muon+tau
431  else if(ntaus == 1 && hltPath_.getFilterNMuons(lastPassedFilter) == 1) {
432  hMass_->Fill( (matchedOfflineObjs.taus[0]+matchedOfflineObjs.muons[0]).M() );
433  }
434  // Tau+MET
435  if(hltPath_.getFilterNTaus(lastPassedFilter) == 1 && hltPath_.getFilterMET(lastMatchedMETFilter) == 1) {
436  double taupt = matchedOfflineObjs.taus[0].Pt();
437  double tauphi = matchedOfflineObjs.taus[0].Phi();
438  double met = matchedOfflineObjs.met[0].Pt();
439  double metphi = matchedOfflineObjs.met[0].Phi();
440  double mT = sqrt( 2 * taupt * met * (1-cos(tauphi-metphi)) );
441 
442  hMass_->Fill( mT );
443  }
444  }
445  }
446 
447  // Triggered object kinematics
448  for(const HLTTauDQMPath::Object& obj: triggerObjs) {
449  if(obj.id == trigger::TriggerTau){
450  hTrigTauEt_->Fill(obj.object.pt());
451  hTrigTauEta_->Fill(obj.object.eta());
452  hTrigTauPhi_->Fill(obj.object.phi());
453  }
455  hTrigElectronEt_->Fill(obj.object.pt());
456  hTrigElectronEta_->Fill(obj.object.eta());
457  hTrigElectronPhi_->Fill(obj.object.phi());
458  }
459  if(obj.id == trigger::TriggerMuon){
460  hTrigMuonEt_->Fill(obj.object.pt());
461  hTrigMuonEta_->Fill(obj.object.eta());
462  hTrigMuonPhi_->Fill(obj.object.phi());
463  }
464  if(obj.id == trigger::TriggerMET){
465  hTrigMETEt_->Fill(obj.object.pt());
466  hTrigMETPhi_->Fill(obj.object.phi());
467  }
468  }
469  }
470  }
471 }
#define LogDebug(id)
bool isValid() const
size_t getLastL2TauFilterIndex() const
Definition: HLTTauDQMPath.h:72
int i
Definition: DBlmapReader.cc:9
MonitorElement * hL3TrigTauEtEffDenom_
MonitorElement * hTrigElectronEta_
bool hasL3Electrons() const
Definition: HLTTauDQMPath.h:67
MonitorElement * hL3TrigTauPhiEffDenom_
MonitorElement * hL2TrigTauPhiEffDenom_
MonitorElement * hL3TrigMuonPhiEffNum_
MonitorElement * hL3TrigMuonEtaEffDenom_
MonitorElement * hL2TrigTauEtEffDenom_
MonitorElement * hTrigTauEt_
std::vector< LV > electrons
The single EDProduct to be saved for each event (AOD case)
Definition: TriggerEvent.h:25
size_t getLastL3ElectronFilterIndex() const
Definition: HLTTauDQMPath.h:79
std::vector< LV > taus
HLTTauDQMPathPlotter(const std::string &pathName, const HLTConfigProvider &HLTCP, bool doRefAnalysis, const std::string &dqmBaseFolder, const std::string &hltProcess, int ptbins, int etabins, int phibins, double ptmax, double highptmax, double l1MatchDr, double hltMatchDr)
size_t getLastL3MuonFilterIndex() const
Definition: HLTTauDQMPath.h:84
MonitorElement * hTrigElectronPhi_
bool isValid() const
Definition: HLTTauDQMPath.h:41
bool hasL2Muons() const
Definition: HLTTauDQMPath.h:68
MonitorElement * hL3TrigTauEtEffNum_
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
int lastPassedFilter(const edm::TriggerResults &triggerResults) const
MonitorElement * hL2TrigMuonEtaEffNum_
const std::string & triggerTag() const
MonitorElement * hL3TrigMuonEtaEffNum_
MonitorElement * hL2TrigTauEtEffNum_
MonitorElement * hAcceptedEvents_
MonitorElement * hL2TrigElectronEtEffNum_
MonitorElement * hL2TrigMuonPhiEffNum_
int getFilterNElectrons(size_t i) const
Definition: HLTTauDQMPath.h:56
int getFilterNTaus(size_t i) const
Definition: HLTTauDQMPath.h:55
bool isFirstFilterL1Seed() const
Definition: HLTTauDQMPath.h:61
MonitorElement * hL2TrigTauPhiEffNum_
void Fill(long long x)
size_t getLastFilterBeforeL3TauIndex() const
Definition: HLTTauDQMPath.h:73
MonitorElement * hL3TrigMuonEtEffNum_
math::XYZTLorentzVectorD LV
MonitorElement * hL2TrigMuonEtEffDenom_
void getFilterObjects(const trigger::TriggerEvent &triggerEvent, size_t i, std::vector< Object > &retval) const
const std::string & getFilterType(size_t i) const
Definition: HLTTauDQMPath.h:54
bool hasL2CaloMET() const
Definition: HLTTauDQMPath.h:70
MonitorElement * hL3TrigTauHighEtEffNum_
MonitorElement * hL3TrigElectronEtaEffNum_
MonitorElement * hL2TrigElectronPhiEffDenom_
MonitorElement * hL3TrigTauEtaEffDenom_
int getFilterNMuons(size_t i) const
Definition: HLTTauDQMPath.h:57
MonitorElement * hL2TrigMuonEtaEffDenom_
MonitorElement * hL2TrigMETEtEffDenom_
T sqrt(T t)
Definition: SSEVec.h:18
size_t getLastFilterBeforeL3ElectronIndex() const
Definition: HLTTauDQMPath.h:78
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
size_t filtersSize() const
Definition: HLTTauDQMPath.h:52
size_t getLastFilterBeforeL2MuonIndex() const
Definition: HLTTauDQMPath.h:81
MonitorElement * hL3TrigElectronEtEffDenom_
MonitorElement * hL3TrigTauPhiEffNum_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * hTrigMETEt_
MonitorElement * hL3TrigTauEtaEffNum_
size_t getLastL2CaloMETFilterIndex() const
Definition: HLTTauDQMPath.h:87
void analyze(const edm::TriggerResults &triggerResults, const trigger::TriggerEvent &triggerEvent, const HLTTauDQMOfflineObjects &refCollection)
MonitorElement * hTrigMETPhi_
size_t getLastL2MuonFilterIndex() const
Definition: HLTTauDQMPath.h:82
bool hasL3Taus() const
Definition: HLTTauDQMPath.h:65
static std::string const triggerResults
Definition: EdmProvDump.cc:41
MonitorElement * hTrigMuonEt_
MonitorElement * hL2TrigMETEtEffNum_
size_t getLastFilterBeforeL3MuonIndex() const
Definition: HLTTauDQMPath.h:83
size_t getLastL3TauFilterIndex() const
Definition: HLTTauDQMPath.h:74
std::vector< LV > met
MonitorElement * hL2TrigTauHighEtEffNum_
bool offlineMatching(size_t i, const std::vector< Object > &triggerObjects, const HLTTauDQMOfflineObjects &offlineObjects, double dR, std::vector< Object > &matchedTriggerObjects, HLTTauDQMOfflineObjects &matchedOfflineObjects) const
MonitorElement * hL3TrigElectronEtEffNum_
MonitorElement * hL3TrigElectronEtaEffDenom_
void bookHistograms(DQMStore::IBooker &iBooker)
MonitorElement * hL2TrigElectronPhiEffNum_
MonitorElement * hTrigTauPhi_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
MonitorElement * hL2TrigElectronEtaEffDenom_
bool hasL3Muons() const
Definition: HLTTauDQMPath.h:69
MonitorElement * hTrigElectronEt_
MonitorElement * hL3TrigMuonPhiEffDenom_
const std::string & getFilterName(size_t i) const
Definition: HLTTauDQMPath.h:53
bool fired(const edm::TriggerResults &triggerResults) const
size_t getLastL2ElectronFilterIndex() const
Definition: HLTTauDQMPath.h:77
MonitorElement * hTrigMuonEta_
size_t getLastFilterBeforeL2TauIndex() const
Definition: HLTTauDQMPath.h:71
int getFilterMET(size_t i) const
Definition: HLTTauDQMPath.h:58
size_t getLastFilterBeforeL2ElectronIndex() const
Definition: HLTTauDQMPath.h:76
const std::string & getPathName() const
Definition: HLTTauDQMPath.h:49
MonitorElement * hL2TrigTauEtaEffNum_
MonitorElement * hL2TrigTauHighEtEffDenom_
bool hasL2Taus() const
Definition: HLTTauDQMPath.h:64
size_t getFirstFilterBeforeL2CaloMETIndex() const
Definition: HLTTauDQMPath.h:88
MonitorElement * hL3TrigElectronPhiEffDenom_
bool hasL2Electrons() const
Definition: HLTTauDQMPath.h:66
MonitorElement * hTrigTauEta_
MonitorElement * hL2TrigTauEtaEffDenom_
MonitorElement * hL3TrigMuonEtEffDenom_
MonitorElement * hL2TrigElectronEtEffDenom_
MonitorElement * hL2TrigElectronEtaEffNum_
MonitorElement * hTrigMuonPhi_
MonitorElement * hL3TrigElectronPhiEffNum_
MonitorElement * hL2TrigMuonPhiEffDenom_
std::vector< LV > muons
MonitorElement * hL3TrigTauHighEtEffDenom_
MonitorElement * hL2TrigMuonEtEffNum_