CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PatTauAnalyzer.cc
Go to the documentation of this file.
11 
12 #include <TH1.h>
13 #include <TMath.h>
14 
15 #include <string>
16 
18 public:
19  explicit PatTauAnalyzer(const edm::ParameterSet&);
20  ~PatTauAnalyzer() override;
21 
22  //--- methods inherited from EDAnalyzer base-class
23  void beginJob() override;
24  void analyze(const edm::Event&, const edm::EventSetup&) override;
25  void endJob() override;
26 
27 private:
28  //--- configuration parameters
31 
33 
37 
38  //--- generator level histograms
40  TH1* hGenTauPt_;
43 
44  //--- reconstruction level histograms
46  TH1* hTauJetPt_;
49 
51 
53 
56 
61 
66 
71 };
72 
73 const reco::GenParticle* getGenTau(const pat::Tau& patTau) {
74  std::vector<reco::GenParticleRef> associatedGenParticles = patTau.genParticleRefs();
75  for (std::vector<reco::GenParticleRef>::const_iterator it = associatedGenParticles.begin();
76  it != associatedGenParticles.end();
77  ++it) {
78  if (it->isAvailable()) {
79  const reco::GenParticleRef& genParticle = (*it);
80  if (genParticle->pdgId() == -15 || genParticle->pdgId() == +15)
81  return genParticle.get();
82  }
83  }
84 
85  return nullptr;
86 }
87 
89  //std::cout << "<PatTauAnalyzer::PatTauAnalyzer>:" << std::endl;
90 
91  //--- read name of pat::Tau collection
92  src_ = cfg.getParameter<edm::InputTag>("src");
93  //std::cout << " src = " << src_ << std::endl;
94  srcToken_ = consumes<pat::TauCollection>(src_);
95 
96  //--- fill histograms for all tau-jet candidates or for "real" taus only ?
97  requireGenTauMatch_ = cfg.getParameter<bool>("requireGenTauMatch");
98  //std::cout << " requireGenTauMatch = " << requireGenTauMatch_ << std::endl;
99 
100  //--- read names of tau id. discriminators
101  discrByLeadTrack_ = cfg.getParameter<std::string>("discrByLeadTrack");
102  //std::cout << " discrByLeadTrack = " << discrByLeadTrack_ << std::endl;
103 
104  discrByIso_ = cfg.getParameter<std::string>("discrByIso");
105  //std::cout << " discrByIso = " << discrByIso_ << std::endl;
106 
107  discrByTaNC_ = cfg.getParameter<std::string>("discrByTaNC");
108  //std::cout << " discrByTaNC = " << discrByTaNC_ << std::endl;
109 }
110 
112  //std::cout << "<PatTauAnalyzer::~PatTauAnalyzer>:" << std::endl;
113 
114  //--- clean-up memory;
115  // delete all histograms
116  /*
117  deletion of histograms taken care of by TFileService;
118  do not delete them here (if the histograms are deleted here,
119  they will not appear in the ROOT file written by TFileService)
120 
121  delete hGenTauEnergy_;
122  delete hGenTauPt_;
123  delete hGenTauEta_;
124  delete hGenTauPhi_;
125  delete hTauJetEnergy_;
126  delete hTauJetPt_;
127  delete hTauJetEta_;
128  delete hTauJetPhi_;
129  delete hNumTauJets_;
130  delete hTauLeadTrackPt_;
131  delete hTauNumSigConeTracks_;
132  delete hTauNumIsoConeTracks_;
133  delete hTauDiscrByIso_;
134  delete hTauDiscrByTaNC_;
135  delete hTauDiscrAgainstElectrons_;
136  delete hTauDiscrAgainstMuons_;
137  delete hTauJetEnergyIsoPassed_;
138  delete hTauJetPtIsoPassed_;
139  delete hTauJetEtaIsoPassed_;
140  delete hTauJetPhiIsoPassed_;
141  delete hTauJetEnergyTaNCpassed_;
142  delete hTauJetPtTaNCpassed_;
143  delete hTauJetEtaTaNCpassed_;
144  delete hTauJetPhiTaNCpassed_;
145  */
146 }
147 
149  //--- retrieve handle to auxiliary service
150  // used for storing histograms into ROOT file
152 
153  //--- book generator level histograms
154  hGenTauEnergy_ = fs->make<TH1F>("GenTauEnergy", "GenTauEnergy", 30, 0., 150.);
155  hGenTauPt_ = fs->make<TH1F>("GenTauPt", "GenTauPt", 30, 0., 150.);
156  hGenTauEta_ = fs->make<TH1F>("GenTauEta", "GenTauEta", 24, -3., +3.);
157  hGenTauPhi_ = fs->make<TH1F>("GenTauPhi", "GenTauPhi", 18, -TMath::Pi(), +TMath::Pi());
158 
159  //--- book reconstruction level histograms
160  // for tau-jet Energy, Pt, Eta, Phi
161  hTauJetEnergy_ = fs->make<TH1F>("TauJetEnergy", "TauJetEnergy", 30, 0., 150.);
162  hTauJetPt_ = fs->make<TH1F>("TauJetPt", "TauJetPt", 30, 0., 150.);
163  //
164  // TO-DO: add histograms for eta and phi of the tau-jet candidate
165  //
166  // NOTE:
167  // 1.) please use
168  // "TauJetEta" and "TauJetPhi"
169  // for the names of the histograms and choose the exact same binning
170  // as is used for the histograms
171  // "TauJetEtaIsoPassed" and "TauJetPhiIsoPassed"
172  // below
173  //
174  // 2.) please check the histograms
175  // hTauJetEta_ and hTauJetPt_
176  // have already been defined in PatTauAnalyzer.h
177  //
178  //hTauJetEta_ =...
179  //hTauJetPt_ =...
180 
181  //... for number of tau-jet candidates
182  hNumTauJets_ = fs->make<TH1F>("NumTauJets", "NumTauJets", 10, -0.5, 9.5);
183 
184  //... for Pt of highest Pt track within signal cone tau-jet...
185  hTauLeadTrackPt_ = fs->make<TH1F>("TauLeadTrackPt", "TauLeadTrackPt", 40, 0., 100.);
186 
187  //... for total number of tracks within signal/isolation cones
188  hTauNumSigConeTracks_ = fs->make<TH1F>("TauNumSigConeTracks", "TauNumSigConeTracks", 10, -0.5, 9.5);
189  hTauNumIsoConeTracks_ = fs->make<TH1F>("TauNumIsoConeTracks", "TauNumIsoConeTracks", 20, -0.5, 19.5);
190 
191  //... for values of tau id. discriminators based on track isolation cut/
192  // neural network-based tau id.
193  hTauDiscrByIso_ = fs->make<TH1F>("TauDiscrByIso", "TauDiscrByIso", 103, -0.015, 1.015);
194  hTauDiscrByTaNC_ = fs->make<TH1F>("TauDiscrByTaNC", "TauDiscrByTaNC", 103, -0.015, 1.015);
195 
196  //... for values of tau id. discriminators against (unidentified) electrons and muons
198  fs->make<TH1F>("TauDiscrAgainstElectrons", "TauDiscrAgainstElectrons", 103, -0.015, 1.015);
199  hTauDiscrAgainstMuons_ = fs->make<TH1F>("TauDiscrAgainstMuons", "TauDiscrAgainstMuons", 103, -0.015, 1.015);
200 
201  //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByIsolation selection
202  hTauJetEnergyIsoPassed_ = fs->make<TH1F>("TauJetEnergyIsoPassed", "TauJetEnergyIsoPassed", 30, 0., 150.);
203  hTauJetPtIsoPassed_ = fs->make<TH1F>("TauJetPtIsoPassed", "TauJetPtIsoPassed", 30, 0., 150.);
204  hTauJetEtaIsoPassed_ = fs->make<TH1F>("TauJetEtaIsoPassed", "TauJetEtaIsoPassed", 24, -3., +3.);
205  hTauJetPhiIsoPassed_ = fs->make<TH1F>("TauJetPhiIsoPassed", "TauJetPhiIsoPassed", 18, -TMath::Pi(), +TMath::Pi());
206 
207  //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByTaNC ("Tau Neural Classifier") selection
208  hTauJetEnergyTaNCpassed_ = fs->make<TH1F>("TauJetEnergyTaNCpassed", "TauJetEnergyTaNCpassed", 30, 0., 150.);
209  hTauJetPtTaNCpassed_ = fs->make<TH1F>("TauJetPtTaNCpassed", "TauJetPtTaNCpassed", 30, 0., 150.);
210  hTauJetEtaTaNCpassed_ = fs->make<TH1F>("TauJetEtaTaNCpassed", "TauJetEtaTaNCpassed", 24, -3., +3.);
211  hTauJetPhiTaNCpassed_ = fs->make<TH1F>("TauJetPhiTaNCpassed", "TauJetPhiTaNCpassed", 18, -TMath::Pi(), +TMath::Pi());
212 }
213 
215  //std::cout << "<PatTauAnalyzer::analyze>:" << std::endl;
216 
218  evt.getByToken(srcToken_, patTaus);
219 
220  hNumTauJets_->Fill(patTaus->size());
221 
222  for (pat::TauCollection::const_iterator patTau = patTaus->begin(); patTau != patTaus->end(); ++patTau) {
223  //--- skip fake taus in case configuration parameters set to do so...
224  const reco::GenParticle* genTau = getGenTau(*patTau);
225  if (requireGenTauMatch_ && !genTau)
226  continue;
227 
228  //--- fill generator level histograms
229  if (genTau) {
230  hGenTauEnergy_->Fill(genTau->energy());
231  hGenTauPt_->Fill(genTau->pt());
232  hGenTauEta_->Fill(genTau->eta());
233  hGenTauPhi_->Fill(genTau->phi());
234  }
235 
236  //--- fill reconstruction level histograms
237  // for Pt of highest Pt track within signal cone tau-jet...
238  hTauJetEnergy_->Fill(patTau->energy());
239  hTauJetPt_->Fill(patTau->pt());
240  //
241  // TO-DO:
242  // 1.) fill histograms
243  // hTauJetEta_ and hTauJetPhi_
244  // with the pseudo-rapidity and azimuthal angle
245  // of the tau-jet candidate respectively
246  // hTauJetEta_->...
247  // hTauJetPhi_->...
248  //
249  // 2.) fill histogram
250  // hTauLeadTrackPt_
251  // with the transverse momentum of the highest Pt ("leading") track within the tau-jet
252  //
253  // NOTE:
254  // 1.) please have a look at
255  // http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/Candidate/interface/Particle.h?revision=1.28&view=markup
256  // to find the methods for accessing eta and phi of the tau-jet
257  //
258  // 2.) please have a look at
259  // http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/PatCandidates/interface/Tau.h?revision=1.25&view=markup
260  // to find the method for accessing the leading track
261  //
262  // 3.) the method pat::Tau::leadTrack returns a reference (reco::TrackRef) to a reco::Track object
263  // this reference can be null (in case no high Pt track has been reconstructed within the tau-jet),
264  // so a check if the leadTrack exists is needed before dereferencing the reco::TrackRef via operator->
265  //
266  // if ( patTau->leadTrack().isAvailable() ) hTauLeadTrackPt_->Fill(patTau->leadTrack()->pt());
267 
268  //... for total number of tracks within signal/isolation cones
269  hTauNumSigConeTracks_->Fill(patTau->signalTracks().size());
270  hTauNumIsoConeTracks_->Fill(patTau->isolationTracks().size());
271 
272  //... for values of tau id. discriminators based on track isolation cut/
273  // neural network-based tau id.
274  // (combine with requirement of at least one "leading" track of Pt > 5. GeV
275  // within the signal cone of the tau-jet)
276  float discrByIso = (patTau->tauID(discrByLeadTrack_.data()) > 0.5) ? patTau->tauID(discrByIso_.data()) : 0.;
277  hTauDiscrByIso_->Fill(discrByIso);
278  float discrByTaNC = (patTau->tauID(discrByLeadTrack_.data()) > 0.5) ? patTau->tauID(discrByTaNC_.data()) : 0.;
279  hTauDiscrByTaNC_->Fill(discrByTaNC);
280 
281  //... for values of tau id. discriminators against (unidentified) electrons and muons
282  //
283  // TO-DO: fill histogram
284  // hTauDiscrAgainstElectrons_
285  // with the value of the discriminatorAgainstElectronsLoose
286  //
287  // NOTE:
288  // 1.) please have a look at
289  // http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/PatCandidates/interface/Tau.h?revision=1.25&view=markup
290  // to find the method for accessing the tau id. information
291  //
292  // 2.) please have a look at
293  // http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/PhysicsTools/PatAlgos/python/tools/tauTools.py?revision=1.43&view=markup
294  // and convince yourself that the string "againstElectronLoose" needs to be passed as argument
295  // of the pat::Tau::tauID method
296  //
297  // hTauDiscrAgainstElectrons_->Fill...
298  hTauDiscrAgainstMuons_->Fill(patTau->tauID("againstMuonLoose"));
299 
300  //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByIsolation selection
301  if (discrByIso > 0.5) {
302  hTauJetEnergyIsoPassed_->Fill(patTau->energy());
303  hTauJetPtIsoPassed_->Fill(patTau->pt());
304  hTauJetEtaIsoPassed_->Fill(patTau->eta());
305  hTauJetPhiIsoPassed_->Fill(patTau->phi());
306  }
307 
308  //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByTaNC ("Tau Neural Classifier") selection
309  if (discrByTaNC > 0.5) {
310  hTauJetEnergyTaNCpassed_->Fill(patTau->energy());
311  hTauJetPtTaNCpassed_->Fill(patTau->pt());
312  hTauJetEtaTaNCpassed_->Fill(patTau->eta());
313  hTauJetPhiTaNCpassed_->Fill(patTau->phi());
314  }
315  }
316 }
317 
319  //--- nothing to be done yet...
320 }
321 
323 
const double Pi
TH1 * hTauDiscrAgainstMuons_
~PatTauAnalyzer() override
std::string discrByLeadTrack_
tuple cfg
Definition: looper.py:296
double pt() const final
transverse momentum
TH1 * hTauNumIsoConeTracks_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const reco::GenParticle * getGenTau(const pat::Tau &patTau)
PatTauAnalyzer(const edm::ParameterSet &)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
std::vector< reco::GenParticleRef > genParticleRefs() const
Definition: PATObject.h:779
TH1 * hTauDiscrAgainstElectrons_
TH1 * hTauJetEnergyIsoPassed_
void beginJob() override
std::string discrByIso_
TH1 * hTauJetPtTaNCpassed_
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
TH1 * hTauJetEnergyTaNCpassed_
Analysis-level tau class.
Definition: Tau.h:53
TH1 * hTauJetPtIsoPassed_
edm::InputTag src_
TH1 * hTauNumSigConeTracks_
std::string discrByTaNC_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void endJob() override
TH1 * hTauJetPhiTaNCpassed_
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< pat::TauCollection > srcToken_
TH1 * hTauJetEtaIsoPassed_
TH1 * hTauJetPhiIsoPassed_
double phi() const final
momentum azimuthal angle
TH1 * hTauJetEtaTaNCpassed_
double energy() const final
energy
double eta() const final
momentum pseudorapidity