CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTTauDQML1Plotter.cc
Go to the documentation of this file.
2 
4 
5 #include<cstring>
6 
7 namespace {
8  double getMaxEta(int binsEta, double widthEta) {
9  if(widthEta <= 0.0) {
10  edm::LogWarning("HLTTauDQMOffline") << "HLTTauDQML1Plotter::HLTTauDQML1Plotter: EtaHistoBinWidth = " << widthEta << " <= 0, using default value 0.348 instead";
11  widthEta = 0.348;
12  }
13  return binsEta/2*widthEta;
14  }
15 }
16 
17 HLTTauDQML1Plotter::HLTTauDQML1Plotter(const edm::ParameterSet& ps, edm::ConsumesCollector&& cc, int phibins, double maxpt, double maxhighpt, bool ref, double dr, const std::string& dqmBaseFolder):
18  HLTTauDQMPlotter(ps, dqmBaseFolder),
19  doRefAnalysis_(ref),
20  matchDeltaR_(dr),
21  maxPt_(maxpt),
22  maxHighPt_(maxhighpt),
23  binsEt_(ps.getUntrackedParameter<int>("EtHistoBins", 25)),
24  binsEta_(ps.getUntrackedParameter<int>("EtaHistoBins", 14)),
25  binsPhi_(phibins),
26  maxEta_(getMaxEta(binsEta_, ps.getUntrackedParameter<double>("EtaHistoBinWidth", 0.348)))
27 {
28  if(!configValid_)
29  return;
30 
31  //Process PSet
36  l1JetMinEt_ = ps.getUntrackedParameter<double>("L1JetMinEt");
37  configValid_ = true;
38 }
39 
41  if(!configValid_)
42  return;
43 
44  // The L1 phi plot is asymmetric around 0 because of the discrete nature of L1 phi
45  constexpr float pi = 3.1416f;
46  constexpr float phiShift = pi/18; // half of 2pi/18 bin
47  constexpr float minPhi = -pi+phiShift;
48  constexpr float maxPhi = pi+phiShift;
49 
50  constexpr int BUFMAX = 256;
51  char buffer[BUFMAX] = "";
52 
53  //Create the histograms
54  iBooker.setCurrentFolder(triggerTag());
55 
56  l1tauEt_ = iBooker.book1D("L1TauEt","L1 #tau E_{T};L1 #tau E_{T};entries",binsEt_,0,maxPt_);
57  l1tauEta_ = iBooker.book1D("L1TauEta","L1 #tau #eta;L1 #tau #eta;entries",binsEta_,-maxEta_,maxEta_);
58  l1tauPhi_ = iBooker.book1D("L1TauPhi","L1 #tau #phi;L1 #tau #phi;entries",binsPhi_,minPhi,maxPhi);
59 
60  l1jetEt_ = iBooker.book1D("L1JetEt","L1 central jet E_{T};L1 jet E_{T};entries",binsEt_,0,maxPt_);
61  snprintf(buffer, BUFMAX, "L1 central jet #eta (E_{T} > %.1f);L1 jet #eta;entries", l1JetMinEt_);
62  l1jetEta_ = iBooker.book1D("L1JetEta", buffer, binsEta_, -maxEta_, maxEta_);
63  snprintf(buffer, BUFMAX, "L1 central jet #phi (E_{T} > %.1f);L1 jet #phi;entries", l1JetMinEt_);
64  l1jetPhi_ = iBooker.book1D("L1JetPhi", buffer, binsPhi_, minPhi, maxPhi);
65 
66  snprintf(buffer, BUFMAX, "L1 leading (#tau OR central jet E_{T} > %.1f) E_{T};L1 (#tau or central jet) E_{T};entries", l1JetMinEt_);
67  firstTauEt_ = iBooker.book1D("L1LeadTauEt", buffer, binsEt_, 0, maxPt_);
68  snprintf(buffer, BUFMAX, "L1 leading (#tau OR central jet E_{T} > %.1f) #eta;L1 (#tau or central jet) #eta;entries", l1JetMinEt_);
69  firstTauEta_ = iBooker.book1D("L1LeadTauEta", buffer, binsEta_, -maxEta_, maxEta_);
70  snprintf(buffer, BUFMAX, "L1 leading (#tau OR central jet E_{T} > %.1f) #phi;L1 (#tau or central jet) #phi;entries", l1JetMinEt_);
71  firstTauPhi_ = iBooker.book1D("L1LeadTauPhi", buffer, binsPhi_, minPhi, maxPhi);
72 
73  snprintf(buffer, BUFMAX, "L1 second-leading (#tau OR central jet E_{T} > %.1f) E_{T};L1 (#tau or central jet) E_{T};entries", l1JetMinEt_);
74  secondTauEt_ = iBooker.book1D("L1SecondTauEt", buffer, binsEt_, 0, maxPt_);
75  snprintf(buffer, BUFMAX, "L1 second-leading (#tau OR central jet E_{T} > %.1f) #eta;L1 (#tau or central jet) #eta;entries", l1JetMinEt_);
76  secondTauEta_ = iBooker.book1D("L1SecondTauEta", buffer, binsEta_, -maxEta_, maxEta_);
77  snprintf(buffer, BUFMAX, "L1 second-leading (#tau OR central jet E_{T} > %.1f) #phi;L1 (#tau or central jet) #phi;entries", l1JetMinEt_);
78  secondTauPhi_ = iBooker.book1D("L1SecondTauPhi", buffer, binsPhi_, minPhi, maxPhi);
79 
80  if (doRefAnalysis_) {
81  l1tauEtRes_ = iBooker.book1D("L1TauEtResol","L1 #tau E_{T} resolution;[L1 #tau E_{T}-Ref #tau E_{T}]/Ref #tau E_{T};entries",60,-1,4);
82  snprintf(buffer, BUFMAX, "L1 central jet E_{T} resolution (E_{T} > %.1f);[L1 jet E_{T}-Ref #tau E_{T}]/Ref #tau E_{T};entries", l1JetMinEt_);
83  l1jetEtRes_ = iBooker.book1D("L1JetEtResol", buffer, 60, -1, 4);
84 
85  iBooker.setCurrentFolder(triggerTag()+"/helpers");
86 
87  l1tauEtEffNum_ = iBooker.book1D("L1TauEtEffNum","L1 #tau E_{T} Efficiency;Ref #tau E_{T};entries",binsEt_,0,maxPt_);
88  l1tauHighEtEffNum_ = iBooker.book1D("L1TauHighEtEffNum","L1 #tau E_{T} Efficiency (high E_{T});Ref #tau E_{T};entries",binsEt_,0,maxHighPt_);
89 
90  l1tauEtEffDenom_ = iBooker.book1D("L1TauEtEffDenom","L1 #tau E_{T} Denominator;Ref #tau E_{T};entries",binsEt_,0,maxPt_);
91  l1tauHighEtEffDenom_ = iBooker.book1D("L1TauHighEtEffDenom","L1 #tau E_{T} Denominator (high E_{T});Ref #tau E_{T};Efficiency",binsEt_,0,maxHighPt_);
92 
93  l1tauEtaEffNum_ = iBooker.book1D("L1TauEtaEffNum","L1 #tau #eta Efficiency;Ref #tau #eta;entries",binsEta_,-maxEta_,maxEta_);
94  l1tauEtaEffDenom_ = iBooker.book1D("L1TauEtaEffDenom","L1 #tau #eta Denominator;Ref #tau #eta;entries",binsEta_,-maxEta_,maxEta_);
95 
96  l1tauPhiEffNum_ = iBooker.book1D("L1TauPhiEffNum","L1 #tau #phi Efficiency;Ref #tau #phi;entries",binsPhi_,minPhi,maxPhi);
97  l1tauPhiEffDenom_ = iBooker.book1D("L1TauPhiEffDenom","L1 #tau #phi Denominator;Ref #tau #phi;Efficiency",binsPhi_,minPhi,maxPhi);
98 
99  l1jetEtEffNum_ = iBooker.book1D("L1JetEtEffNum","L1 central jet E_{T} Efficiency;Ref #tau E_{T};entries",binsEt_,0,maxPt_);
100  l1jetHighEtEffNum_ = iBooker.book1D("L1JetHighEtEffNum","L1 central jet E_{T} Efficiency (high E_{T});Ref #tau E_{T};entries",binsEt_,0,maxHighPt_);
101 
102  l1jetEtEffDenom_ = iBooker.book1D("L1JetEtEffDenom","L1 central jet E_{T} Denominator;Ref #tau E_{T};entries",binsEt_,0,maxPt_);
103  l1jetHighEtEffDenom_ = iBooker.book1D("L1JetHighEtEffDenom","L1 central jet E_{T} Denominator (high E_{T});Ref #tau E_{T};Efficiency",binsEt_,0,maxHighPt_);
104 
105  snprintf(buffer, BUFMAX, "L1 central jet #eta Efficiency (E_{T} > %.1f);Ref #tau #eta;entries", l1JetMinEt_);
106  l1jetEtaEffNum_ = iBooker.book1D("L1JetEtaEffNum", buffer, binsEta_, -maxEta_, maxEta_);
107 
108  snprintf(buffer, BUFMAX, "L1 central jet #eta Denominator (E_{T} > %.1f);Ref #tau #eta;Efficiency", l1JetMinEt_);
109  l1jetEtaEffDenom_ = iBooker.book1D("L1JetEtaEffDenom", buffer, binsEta_, -maxEta_, maxEta_);
110 
111  snprintf(buffer, BUFMAX, "L1 central jet #phi Efficiency (E_{T} > %.1f);Ref #tau #phi;entries", l1JetMinEt_);
112  l1jetPhiEffNum_ = iBooker.book1D("L1JetPhiEffNum", buffer, binsPhi_, minPhi, maxPhi);
113 
114  snprintf(buffer, BUFMAX, "L1 central jet #phi Efficiency (E_{T} > %.1f);Ref #tau #phi;Efficiency", l1JetMinEt_);
115  l1jetPhiEffDenom_ = iBooker.book1D("L1JetPhiEffDenom", buffer, binsPhi_, minPhi, maxPhi);
116  }
117 }
118 
119 
121 }
122 
123 //
124 // member functions
125 //
126 
128  if ( doRefAnalysis_ ) {
129  //Tau reference
130  for ( LVColl::const_iterator iter = refC.taus.begin(); iter != refC.taus.end(); ++iter ) {
131  l1tauEtEffDenom_->Fill(iter->pt());
132  l1jetEtEffDenom_->Fill(iter->pt());
135 
136  l1tauEtaEffDenom_->Fill(iter->eta());
137  l1jetEtaEffDenom_->Fill(iter->eta());
138 
139  l1tauPhiEffDenom_->Fill(iter->phi());
140  l1jetPhiEffDenom_->Fill(iter->phi());
141  }
142  }
143 
144  //Analyze L1 Objects (Tau+Jets)
147  iEvent.getByToken(l1ExtraTausToken_, taus);
148  iEvent.getByToken(l1ExtraJetsToken_, jets);
149 
150  LVColl pathTaus;
151 
152  //Set Variables for the threshold plot
153  LVColl l1taus;
154  LVColl l1jets;
155 
156  if(taus.isValid()) {
157  for(l1extra::L1JetParticleCollection::const_iterator i = taus->begin(); i != taus->end(); ++i) {
158  l1taus.push_back(i->p4());
159  if(!doRefAnalysis_) {
160  l1tauEt_->Fill(i->et());
161  l1tauEta_->Fill(i->eta());
162  l1tauPhi_->Fill(i->phi());
163  pathTaus.push_back(i->p4());
164  }
165  }
166  }
167  else {
168  edm::LogWarning("HLTTauDQMOffline") << "HLTTauDQML1Plotter::analyze: unable to read L1 tau collection " << l1ExtraTaus_.encode();
169  }
170 
171  if(jets.isValid()) {
172  for(l1extra::L1JetParticleCollection::const_iterator i = jets->begin(); i != jets->end(); ++i) {
173  l1jets.push_back(i->p4());
174  if(!doRefAnalysis_) {
175  l1jetEt_->Fill(i->et());
176  if(i->et() >= l1JetMinEt_) {
177  l1jetEta_->Fill(i->eta());
178  l1jetPhi_->Fill(i->phi());
179  pathTaus.push_back(i->p4());
180  }
181  }
182  }
183  }
184  else {
185  edm::LogWarning("HLTTauDQMOffline") << "HLTTauDQML1Plotter::analyze: unable to read L1 jet collection " << l1ExtraJets_.encode();
186  }
187 
188  //Now do the efficiency matching
189  if ( doRefAnalysis_ ) {
190  for ( LVColl::const_iterator i = refC.taus.begin(); i != refC.taus.end(); ++i ) {
191  std::pair<bool,LV> m = match(*i,l1taus,matchDeltaR_);
192  if ( m.first ) {
193  l1tauEt_->Fill(m.second.pt());
194  l1tauEta_->Fill(m.second.eta());
195  l1tauPhi_->Fill(m.second.phi());
196 
197  l1tauEtEffNum_->Fill(i->pt());
198  l1tauHighEtEffNum_->Fill(i->pt());
199  l1tauEtaEffNum_->Fill(i->eta());
200  l1tauPhiEffNum_->Fill(i->phi());
201 
202  l1tauEtRes_->Fill((m.second.pt()-i->pt())/i->pt());
203 
204  pathTaus.push_back(m.second);
205  }
206  }
207 
208  for ( LVColl::const_iterator i = refC.taus.begin(); i != refC.taus.end(); ++i ) {
209  std::pair<bool,LV> m = match(*i,l1jets,matchDeltaR_);
210  if ( m.first ) {
211  l1jetEt_->Fill(m.second.pt());
212  if(m.second.pt() >= l1JetMinEt_) {
213  l1jetEta_->Fill(m.second.eta());
214  l1jetPhi_->Fill(m.second.phi());
215 
216  l1jetEtEffNum_->Fill(i->pt());
217  l1jetHighEtEffNum_->Fill(i->pt());
218  l1jetEtaEffNum_->Fill(i->eta());
219  l1jetPhiEffNum_->Fill(i->phi());
220 
221  l1jetEtRes_->Fill((m.second.pt()-i->pt())/i->pt());
222 
223  pathTaus.push_back(m.second);
224  }
225  }
226  }
227  }
228 
229 
230  //Fill the Threshold Monitoring
231  if(pathTaus.size() > 1) std::sort(pathTaus.begin(), pathTaus.end(), [](const LV& a, const LV& b) { return a.pt() > b.pt(); });
232 
233  if ( pathTaus.size() > 0 ) {
234  firstTauEt_->Fill(pathTaus[0].pt());
235  firstTauEta_->Fill(pathTaus[0].eta());
236  firstTauPhi_->Fill(pathTaus[0].phi());
237  }
238  if ( pathTaus.size() > 1 ) {
239  secondTauEt_->Fill(pathTaus[1].pt());
240  secondTauEta_->Fill(pathTaus[1].eta());
241  secondTauPhi_->Fill(pathTaus[1].phi());
242  }
243 }
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
const double matchDeltaR_
MonitorElement * l1jetEta_
edm::EDGetTokenT< l1extra::L1JetParticleCollection > l1ExtraJetsToken_
MonitorElement * secondTauEta_
std::vector< LV > taus
void bookHistograms(DQMStore::IBooker &iBooker)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
MonitorElement * firstTauPhi_
MonitorElement * l1jetEtEffNum_
MonitorElement * l1jetPhiEffDenom_
MonitorElement * l1jetHighEtEffDenom_
std::vector< L1JetParticle > L1JetParticleCollection
MonitorElement * l1jetEtaEffNum_
MonitorElement * l1jetHighEtEffNum_
const std::string & triggerTag() const
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup, const HLTTauDQMOfflineObjects &refC)
T eta() const
edm::EDGetTokenT< l1extra::L1JetParticleCollection > l1ExtraTausToken_
MonitorElement * firstTauEta_
std::string encode() const
Definition: InputTag.cc:164
#define constexpr
MonitorElement * l1tauEt_
MonitorElement * secondTauEt_
MonitorElement * l1tauEtEffDenom_
void Fill(long long x)
const Double_t pi
MonitorElement * l1jetEtRes_
MonitorElement * l1jetPhiEffNum_
math::XYZTLorentzVectorD LV
int iEvent
Definition: GenABIO.cc:230
MonitorElement * l1tauEtaEffNum_
MonitorElement * firstTauEt_
MonitorElement * l1tauHighEtEffDenom_
std::vector< LV > LVColl
vector< PseudoJet > jets
MonitorElement * l1jetEtEffDenom_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:113
MonitorElement * l1tauEtEffNum_
MonitorElement * l1tauPhi_
edm::InputTag l1ExtraTaus_
bool isValid() const
Definition: HandleBase.h:76
MonitorElement * l1tauEta_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * secondTauPhi_
MonitorElement * l1tauPhiEffDenom_
HLTTauDQML1Plotter(const edm::ParameterSet &, edm::ConsumesCollector &&cc, int phibins, double maxpt, double maxhighpt, bool ref, double dr, const std::string &dqmBaseFolder)
double b
Definition: hdecay.h:120
edm::InputTag l1ExtraJets_
MonitorElement * l1tauEtaEffDenom_
MonitorElement * l1tauEtRes_
MonitorElement * l1jetPhi_
double a
Definition: hdecay.h:121
MonitorElement * l1jetEt_
MonitorElement * l1jetEtaEffDenom_
std::pair< bool, LV > match(const LV &, const LVColl &, double)
MonitorElement * l1tauPhiEffNum_
MonitorElement * l1tauHighEtEffNum_
Definition: DDAxes.h:10