CMS 3D CMS Logo

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