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 
35  if(!isValid())
36  return;
37 
38  // Book histograms
39  iBooker.setCurrentFolder(triggerTag());
40 
41  hAcceptedEvents_ = iBooker.book1D("EventsPerFilter", "Accepted Events per filter;;entries", hltPath_.filtersSize(), 0, hltPath_.filtersSize());
42  for(size_t i=0; i<hltPath_.filtersSize(); ++i) {
44  }
45 
46  hTrigTauEt_ = iBooker.book1D("TrigTauEt", "Triggered #tau p_{T};#tau p_{T};entries", ptbins_, 0, ptmax_);
47  hTrigTauEta_ = iBooker.book1D("TrigTauEta", "Triggered #tau #eta;#tau #eta;entries", etabins_, -2.5, 2.5);
48  hTrigTauPhi_ = iBooker.book1D("TrigTauPhi", "Triggered #tau #phi;#tau #phi;entries", phibins_, -3.2, 3.2);
49 
50  // Efficiency helpers
51  if(doRefAnalysis_) {
52  iBooker.setCurrentFolder(triggerTag()+"/helpers");
53  if(hltPath_.hasL2Taus()) {
54  hL2TrigTauEtEffNum_ = iBooker.book1D("L2TrigTauEtEffNum", "L2 #tau p_{T} efficiency;Ref #tau p_{T};entries", ptbins_, 0, ptmax_);
55  hL2TrigTauEtEffDenom_ = iBooker.book1D("L2TrigTauEtEffDenom", "L2 #tau p_{T} denominator;Ref #tau p_{T};Efficiency", ptbins_, 0, ptmax_);
56  hL2TrigTauEtaEffNum_ = iBooker.book1D("L2TrigTauEtaEffNum", "L2 #tau #eta efficiency;Ref #tau #eta;entries", etabins_, -2.5, 2.5);
57  hL2TrigTauEtaEffDenom_ = iBooker.book1D("L2TrigTauEtaEffDenom", "L2 #tau #eta denominator;Ref #tau #eta;Efficiency", etabins_, -2.5, 2.5);
58  hL2TrigTauPhiEffNum_ = iBooker.book1D("L2TrigTauPhiEffNum", "L2 #tau #phi efficiency;Ref #tau #phi;entries", phibins_, -3.2, 3.2);
59  hL2TrigTauPhiEffDenom_ = iBooker.book1D("L2TrigTauPhiEffDenom", "L2 #tau #phi denominator;Ref #tau #phi;Efficiency", phibins_, -3.2, 3.2);
60  hL2TrigTauHighEtEffNum_ = iBooker.book1D("L2TrigTauHighEtEffNum", "L2 #tau p_{T} efficiency (high p_{T});Ref #tau p_{T};entries", ptbins_, 0, highptmax_);
61  hL2TrigTauHighEtEffDenom_ = iBooker.book1D("L2TrigTauHighEtEffDenom", "L2 #tau p_{T} denominator (high p_{T});Ref #tau p_{T};Efficiency", ptbins_, 0, highptmax_);
62  }
63 
64  if(hltPath_.hasL3Taus()) {
65  hL3TrigTauEtEffNum_ = iBooker.book1D("L3TrigTauEtEffNum", "L3 #tau p_{T} efficiency;Ref #tau p_{T};entries", ptbins_, 0, ptmax_);
66  hL3TrigTauEtEffDenom_ = iBooker.book1D("L3TrigTauEtEffDenom", "L3 #tau p_{T} denominator;Ref #tau p_{T};Efficiency", ptbins_, 0, ptmax_);
67  hL3TrigTauEtaEffNum_ = iBooker.book1D("L3TrigTauEtaEffNum", "L3 #tau #eta efficiency;Ref #tau #eta;entries", etabins_, -2.5, 2.5);
68  hL3TrigTauEtaEffDenom_ = iBooker.book1D("L3TrigTauEtaEffDenom", "L3 #tau #eta denominator;Ref #tau #eta;Efficiency", etabins_, -2.5, 2.5);
69  hL3TrigTauPhiEffNum_ = iBooker.book1D("L3TrigTauPhiEffNum", "L3 #tau #phi efficiency;Ref #tau #phi;entries", phibins_, -3.2, 3.2);
70  hL3TrigTauPhiEffDenom_ = iBooker.book1D("L3TrigTauPhiEffDenom", "L3 #tau #phi denominator;Ref #tau #phi;Efficiency", phibins_, -3.2, 3.2);
71  hL3TrigTauHighEtEffNum_ = iBooker.book1D("L3TrigTauHighEtEffNum", "L3 #tau p_{T} efficiency (high p_{T});Ref #tau p_{T};entries", ptbins_, 0, highptmax_);
72  hL3TrigTauHighEtEffDenom_ = iBooker.book1D("L3TrigTauHighEtEffDenom", "L3 #tau p_{T} denominator (high p_{T});Ref #tau p_{T};Efficiency", ptbins_, 0, highptmax_);
73  }
74  iBooker.setCurrentFolder(triggerTag());
75  }
76 
77  // Book di-object invariant mass histogram only for mu+tau, ele+tau, and di-tau paths
78  hMass_ = nullptr;
79  if(doRefAnalysis_) {
80  const int lastFilter = hltPath_.filtersSize()-1;
81  const int ntaus = hltPath_.getFilterNTaus(lastFilter);
82  const int neles = hltPath_.getFilterNElectrons(lastFilter);
83  const int nmus = hltPath_.getFilterNMuons(lastFilter);
84  auto create = [&](const std::string& name) {
85  this->hMass_ = iBooker.book1D("ReferenceMass", "Invariant mass of reference "+name+";Reference invariant mass;entries", 100, 0, 500);
86  };
87 
88  LogDebug("HLTTauDQMOffline") << "Path " << hltPath_.getPathName() << " number of taus " << ntaus << " electrons " << neles << " muons " << nmus;
89 
90  if(ntaus == 2)
91  create("di-tau");
92  else if(ntaus == 1) {
93  if(neles == 1)
94  create("electron-tau");
95  else if(nmus == 1)
96  create("muon-tau");
97  }
98  }
99 }
100 
101 
103 
105 
106  std::vector<HLTTauDQMPath::Object> triggerObjs;
107  std::vector<HLTTauDQMPath::Object> matchedTriggerObjs;
108  HLTTauDQMOfflineObjects matchedOfflineObjs;
109 
110  // Events per filter
111  const int lastPassedFilter = hltPath_.lastPassedFilter(triggerResults);
112  int lastMatchedFilter = -1;
113  //std::cout << "Last passed filter " << lastPassedFilter << " " << (lastPassedFilter >= 0 ? hltPath_.getFilterName(lastPassedFilter) : "") << std::endl;
114  if(doRefAnalysis_) {
115  double matchDr = hltPath_.isFirstFilterL1Seed() ? l1MatchDr_ : hltMatchDr_;
116  for(int i=0; i<=lastPassedFilter; ++i) {
117  triggerObjs.clear();
118  matchedTriggerObjs.clear();
119  matchedOfflineObjs.clear();
120  hltPath_.getFilterObjects(triggerEvent, i, triggerObjs);
121  //std::cout << "Filter name " << hltPath_.getFilterName(i) << " nobjs " << triggerObjs.size() << std::endl;
122  bool matched = hltPath_.offlineMatching(i, triggerObjs, refCollection, matchDr, matchedTriggerObjs, matchedOfflineObjs);
123  //std::cout << " offline matching: " << matched << std::endl;
124  matchDr = hltMatchDr_;
125  if(!matched)
126  break;
127 
128  hAcceptedEvents_->Fill(i+0.5);
129  lastMatchedFilter = i;
130  }
131  }
132  else {
133  for(int i=0; i<=lastPassedFilter; ++i) {
134  hAcceptedEvents_->Fill(i+0.5);
135  }
136  }
137 
138  // Efficiency plots
139  if(doRefAnalysis_ && lastMatchedFilter >= 0) {
140  // L2 taus
141  if(hltPath_.hasL2Taus()) {
142  // Denominators
143  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL2TauIndex()) {
144  for(const LV& tau: refCollection.taus) {
145  hL2TrigTauEtEffDenom_->Fill(tau.pt());
146  hL2TrigTauHighEtEffDenom_->Fill(tau.pt());
147  hL2TrigTauEtaEffDenom_->Fill(tau.eta());
148  hL2TrigTauPhiEffDenom_->Fill(tau.phi());
149  }
150  }
151 
152  // Numerators
153  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL2TauFilterIndex()) {
154  triggerObjs.clear();
155  matchedTriggerObjs.clear();
156  matchedOfflineObjs.clear();
157  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL2TauFilterIndex(), triggerObjs);
158  bool matched = hltPath_.offlineMatching(hltPath_.getLastL2TauFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
159  if(matched) {
160  for(const LV& tau: matchedOfflineObjs.taus) {
161  hL2TrigTauEtEffNum_->Fill(tau.pt());
162  hL2TrigTauHighEtEffNum_->Fill(tau.pt());
163  hL2TrigTauEtaEffNum_->Fill(tau.eta());
164  hL2TrigTauPhiEffNum_->Fill(tau.phi());
165  }
166  }
167  }
168  }
169 
170  // L3 taus
171  if(hltPath_.hasL3Taus()) {
172  // Denominators
173  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastFilterBeforeL3TauIndex()) {
174  for(const LV& tau: refCollection.taus) {
175  hL3TrigTauEtEffDenom_->Fill(tau.pt());
176  hL3TrigTauHighEtEffDenom_->Fill(tau.pt());
177  hL3TrigTauEtaEffDenom_->Fill(tau.eta());
178  hL3TrigTauPhiEffDenom_->Fill(tau.phi());
179  }
180  }
181 
182  // Numerators
183  if(static_cast<size_t>(lastMatchedFilter) >= hltPath_.getLastL3TauFilterIndex()) {
184  triggerObjs.clear();
185  matchedTriggerObjs.clear();
186  matchedOfflineObjs.clear();
187  hltPath_.getFilterObjects(triggerEvent, hltPath_.getLastL3TauFilterIndex(), triggerObjs);
188  bool matched = hltPath_.offlineMatching(hltPath_.getLastL3TauFilterIndex(), triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
189  if(matched) {
190  for(const LV& tau: matchedOfflineObjs.taus) {
191  hL3TrigTauEtEffNum_->Fill(tau.pt());
192  hL3TrigTauHighEtEffNum_->Fill(tau.pt());
193  hL3TrigTauEtaEffNum_->Fill(tau.eta());
194  hL3TrigTauPhiEffNum_->Fill(tau.phi());
195  }
196  }
197  }
198  }
199  }
200 
201  if(hltPath_.fired(triggerResults)) {
202  triggerObjs.clear();
203  matchedTriggerObjs.clear();
204  matchedOfflineObjs.clear();
205  hltPath_.getFilterObjects(triggerEvent, lastPassedFilter, triggerObjs);
206  if(doRefAnalysis_) {
207  bool matched = hltPath_.offlineMatching(lastPassedFilter, triggerObjs, refCollection, hltMatchDr_, matchedTriggerObjs, matchedOfflineObjs);
208  if(matched) {
209  // Di-object invariant mass
210  if(hMass_) {
211  const int ntaus = hltPath_.getFilterNTaus(lastPassedFilter);
212  if(ntaus == 2) {
213  // Di-tau (matchedOfflineObjs are already sorted)
214  hMass_->Fill( (matchedOfflineObjs.taus[0]+matchedOfflineObjs.taus[1]).M() );
215  }
216  // Electron+tau
217  else if(ntaus == 1 && hltPath_.getFilterNElectrons(lastPassedFilter) == 1) {
218  hMass_->Fill( (matchedOfflineObjs.taus[0]+matchedOfflineObjs.electrons[0]).M() );
219  }
220  // Muon+tau
221  else if(ntaus == 1 && hltPath_.getFilterNMuons(lastPassedFilter) == 1) {
222  hMass_->Fill( (matchedOfflineObjs.taus[0]+matchedOfflineObjs.muons[0]).M() );
223  }
224  }
225  }
226 
227  if(matched)
228  triggerObjs.swap(matchedTriggerObjs);
229  else
230  triggerObjs.clear();
231  }
232 
233  // Triggered tau kinematics
234  for(const HLTTauDQMPath::Object& obj: triggerObjs) {
235  if(obj.id != trigger::TriggerTau)
236  continue;
237  hTrigTauEt_->Fill(obj.object.pt());
238  hTrigTauEta_->Fill(obj.object.eta());
239  hTrigTauPhi_->Fill(obj.object.phi());
240  }
241  }
242 }
#define LogDebug(id)
bool isValid() const
size_t getLastL2TauFilterIndex() const
Definition: HLTTauDQMPath.h:63
int i
Definition: DBlmapReader.cc:9
MonitorElement * hL3TrigTauEtEffDenom_
MonitorElement * hL3TrigTauPhiEffDenom_
MonitorElement * hL2TrigTauPhiEffDenom_
MonitorElement * hL2TrigTauEtEffDenom_
MonitorElement * hTrigTauEt_
std::vector< LV > electrons
The single EDProduct to be saved for each event (AOD case)
Definition: TriggerEvent.h:25
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)
bool isValid() const
Definition: HLTTauDQMPath.h:40
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
const std::string & triggerTag() const
MonitorElement * hL2TrigTauEtEffNum_
MonitorElement * hAcceptedEvents_
int getFilterNElectrons(size_t i) const
Definition: HLTTauDQMPath.h:54
int getFilterNTaus(size_t i) const
Definition: HLTTauDQMPath.h:53
bool isFirstFilterL1Seed() const
Definition: HLTTauDQMPath.h:57
MonitorElement * hL2TrigTauPhiEffNum_
void Fill(long long x)
size_t getLastFilterBeforeL3TauIndex() const
Definition: HLTTauDQMPath.h:64
math::XYZTLorentzVectorD LV
void getFilterObjects(const trigger::TriggerEvent &triggerEvent, size_t i, std::vector< Object > &retval) const
MonitorElement * hL3TrigTauHighEtEffNum_
MonitorElement * hL3TrigTauEtaEffDenom_
int getFilterNMuons(size_t i) const
Definition: HLTTauDQMPath.h:55
size_t filtersSize() const
Definition: HLTTauDQMPath.h:51
MonitorElement * hL3TrigTauPhiEffNum_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * hL3TrigTauEtaEffNum_
void analyze(const edm::TriggerResults &triggerResults, const trigger::TriggerEvent &triggerEvent, const HLTTauDQMOfflineObjects &refCollection)
bool hasL3Taus() const
Definition: HLTTauDQMPath.h:61
static std::string const triggerResults
Definition: EdmProvDump.cc:41
size_t getLastL3TauFilterIndex() const
Definition: HLTTauDQMPath.h:65
MonitorElement * hL2TrigTauHighEtEffNum_
bool offlineMatching(size_t i, const std::vector< Object > &triggerObjects, const HLTTauDQMOfflineObjects &offlineObjects, double dR, std::vector< Object > &matchedTriggerObjects, HLTTauDQMOfflineObjects &matchedOfflineObjects) const
void bookHistograms(DQMStore::IBooker &iBooker)
MonitorElement * hTrigTauPhi_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
const std::string & getFilterName(size_t i) const
Definition: HLTTauDQMPath.h:52
bool fired(const edm::TriggerResults &triggerResults) const
size_t getLastFilterBeforeL2TauIndex() const
Definition: HLTTauDQMPath.h:62
const std::string & getPathName() const
Definition: HLTTauDQMPath.h:48
MonitorElement * hL2TrigTauEtaEffNum_
MonitorElement * hL2TrigTauHighEtEffDenom_
bool hasL2Taus() const
Definition: HLTTauDQMPath.h:60
MonitorElement * hTrigTauEta_
MonitorElement * hL2TrigTauEtaEffDenom_
SurfaceDeformation * create(int type, const std::vector< double > &params)
std::vector< LV > muons
MonitorElement * hL3TrigTauHighEtEffDenom_