CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/PhysicsTools/PatExamples/plugins/PatTauAnalyzer.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/PatExamples/plugins/PatTauAnalyzer.h"
00002 
00003 #include "FWCore/Framework/interface/Frameworkfwd.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "FWCore/ServiceRegistry/interface/Service.h"
00006 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00007 #include "DataFormats/PatCandidates/interface/Tau.h"
00008 
00009 #include <TMath.h>
00010 
00011 const reco::GenParticle* getGenTau(const pat::Tau& patTau)
00012 {
00013   std::vector<reco::GenParticleRef> associatedGenParticles = patTau.genParticleRefs();
00014   for ( std::vector<reco::GenParticleRef>::const_iterator it = associatedGenParticles.begin(); 
00015         it != associatedGenParticles.end(); ++it ) {
00016     if ( it->isAvailable() ) {
00017       const reco::GenParticleRef& genParticle = (*it);
00018       if ( genParticle->pdgId() == -15 || genParticle->pdgId() == +15 ) return genParticle.get();
00019     }
00020   }
00021 
00022   return 0;
00023 }
00024 
00025 PatTauAnalyzer::PatTauAnalyzer(const edm::ParameterSet& cfg)
00026 {
00027   //std::cout << "<PatTauAnalyzer::PatTauAnalyzer>:" << std::endl;
00028 
00029 //--- read name of pat::Tau collection
00030   src_ = cfg.getParameter<edm::InputTag>("src");
00031   //std::cout << " src = " << src_ << std::endl;
00032 
00033 //--- fill histograms for all tau-jet candidates or for "real" taus only ?
00034   requireGenTauMatch_ = cfg.getParameter<bool>("requireGenTauMatch");
00035   //std::cout << " requireGenTauMatch = " << requireGenTauMatch_ << std::endl;
00036 
00037 //--- read names of tau id. discriminators
00038   discrByLeadTrack_ = cfg.getParameter<std::string>("discrByLeadTrack");
00039   //std::cout << " discrByLeadTrack = " << discrByLeadTrack_ << std::endl;
00040 
00041   discrByIso_ = cfg.getParameter<std::string>("discrByIso");
00042   //std::cout << " discrByIso = " << discrByIso_ << std::endl;
00043 
00044   discrByTaNC_ = cfg.getParameter<std::string>("discrByTaNC");
00045   //std::cout << " discrByTaNC = " << discrByTaNC_ << std::endl;
00046 }
00047 
00048 PatTauAnalyzer::~PatTauAnalyzer()
00049 {
00050   //std::cout << "<PatTauAnalyzer::~PatTauAnalyzer>:" << std::endl;
00051 
00052 //--- clean-up memory;
00053 //    delete all histograms
00054 /*
00055   deletion of histograms taken care of by TFileService;
00056   do not delete them here (if the histograms are deleted here,
00057   they will not appear in the ROOT file written by TFileService)
00058 
00059   delete hGenTauEnergy_;
00060   delete hGenTauPt_;
00061   delete hGenTauEta_;
00062   delete hGenTauPhi_;
00063   delete hTauJetEnergy_;
00064   delete hTauJetPt_;
00065   delete hTauJetEta_;
00066   delete hTauJetPhi_;
00067   delete hNumTauJets_;
00068   delete hTauLeadTrackPt_;
00069   delete hTauNumSigConeTracks_;
00070   delete hTauNumIsoConeTracks_;
00071   delete hTauDiscrByIso_;
00072   delete hTauDiscrByTaNC_;
00073   delete hTauDiscrAgainstElectrons_;
00074   delete hTauDiscrAgainstMuons_;
00075   delete hTauJetEnergyIsoPassed_;
00076   delete hTauJetPtIsoPassed_;
00077   delete hTauJetEtaIsoPassed_;
00078   delete hTauJetPhiIsoPassed_;
00079   delete hTauJetEnergyTaNCpassed_;
00080   delete hTauJetPtTaNCpassed_;
00081   delete hTauJetEtaTaNCpassed_;
00082   delete hTauJetPhiTaNCpassed_;
00083  */
00084 }
00085 
00086 void PatTauAnalyzer::beginJob()
00087 {
00088 //--- retrieve handle to auxiliary service
00089 //    used for storing histograms into ROOT file
00090   edm::Service<TFileService> fs;
00091 
00092 //--- book generator level histograms
00093   hGenTauEnergy_ = fs->make<TH1F>("GenTauEnergy", "GenTauEnergy", 30, 0., 150.);
00094   hGenTauPt_ = fs->make<TH1F>("GenTauPt", "GenTauPt", 30, 0., 150.);
00095   hGenTauEta_ = fs->make<TH1F>("GenTauEta", "GenTauEta", 24, -3., +3.);
00096   hGenTauPhi_ = fs->make<TH1F>("GenTauPhi", "GenTauPhi", 18, -TMath::Pi(), +TMath::Pi());
00097   
00098 //--- book reconstruction level histograms
00099 //    for tau-jet Energy, Pt, Eta, Phi
00100   hTauJetEnergy_ = fs->make<TH1F>("TauJetEnergy", "TauJetEnergy", 30, 0., 150.);
00101   hTauJetPt_ = fs->make<TH1F>("TauJetPt", "TauJetPt", 30, 0., 150.);
00102 // 
00103 // TO-DO: add histograms for eta and phi of the tau-jet candidate
00104 //      
00105 // NOTE:
00106 //  1.) please use 
00107 //       "TauJetEta" and "TauJetPhi" 
00108 //      for the names of the histograms and choose the exact same binning 
00109 //      as is used for the histograms 
00110 //       "TauJetEtaIsoPassed" and "TauJetPhiIsoPassed" 
00111 //      below
00112 //
00113 //  2.) please check the histograms
00114 //       hTauJetEta_ and hTauJetPt_
00115 //      have already been defined in PatTauAnalyzer.h
00116 // 
00117 //hTauJetEta_ =...
00118 //hTauJetPt_ =...
00119 
00120 //... for number of tau-jet candidates
00121   hNumTauJets_ = fs->make<TH1F>("NumTauJets", "NumTauJets", 10, -0.5, 9.5);
00122 
00123 //... for Pt of highest Pt track within signal cone tau-jet...
00124   hTauLeadTrackPt_ = fs->make<TH1F>("TauLeadTrackPt", "TauLeadTrackPt", 40, 0., 100.);
00125   
00126 //... for total number of tracks within signal/isolation cones
00127   hTauNumSigConeTracks_ = fs->make<TH1F>("TauNumSigConeTracks", "TauNumSigConeTracks", 10, -0.5,  9.5);
00128   hTauNumIsoConeTracks_ = fs->make<TH1F>("TauNumIsoConeTracks", "TauNumIsoConeTracks", 20, -0.5, 19.5);
00129 
00130 //... for values of tau id. discriminators based on track isolation cut/
00131 //    neural network-based tau id.
00132   hTauDiscrByIso_ = fs->make<TH1F>("TauDiscrByIso", "TauDiscrByIso", 103, -0.015, 1.015);
00133   hTauDiscrByTaNC_ = fs->make<TH1F>("TauDiscrByTaNC", "TauDiscrByTaNC", 103, -0.015, 1.015);
00134   
00135 //... for values of tau id. discriminators against (unidentified) electrons and muons
00136   hTauDiscrAgainstElectrons_ = fs->make<TH1F>("TauDiscrAgainstElectrons", "TauDiscrAgainstElectrons", 103, -0.015, 1.015);
00137   hTauDiscrAgainstMuons_ = fs->make<TH1F>("TauDiscrAgainstMuons", "TauDiscrAgainstMuons", 103, -0.015, 1.015);
00138 
00139 //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByIsolation selection
00140   hTauJetEnergyIsoPassed_ = fs->make<TH1F>("TauJetEnergyIsoPassed", "TauJetEnergyIsoPassed", 30, 0., 150.);
00141   hTauJetPtIsoPassed_ = fs->make<TH1F>("TauJetPtIsoPassed", "TauJetPtIsoPassed", 30, 0., 150.);
00142   hTauJetEtaIsoPassed_ = fs->make<TH1F>("TauJetEtaIsoPassed", "TauJetEtaIsoPassed", 24, -3., +3.);
00143   hTauJetPhiIsoPassed_ = fs->make<TH1F>("TauJetPhiIsoPassed", "TauJetPhiIsoPassed", 18, -TMath::Pi(), +TMath::Pi());
00144 
00145 //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByTaNC ("Tau Neural Classifier") selection
00146   hTauJetEnergyTaNCpassed_ = fs->make<TH1F>("TauJetEnergyTaNCpassed", "TauJetEnergyTaNCpassed", 30, 0., 150.);
00147   hTauJetPtTaNCpassed_ = fs->make<TH1F>("TauJetPtTaNCpassed", "TauJetPtTaNCpassed", 30, 0., 150.);
00148   hTauJetEtaTaNCpassed_ = fs->make<TH1F>("TauJetEtaTaNCpassed", "TauJetEtaTaNCpassed", 24, -3., +3.);
00149   hTauJetPhiTaNCpassed_ = fs->make<TH1F>("TauJetPhiTaNCpassed", "TauJetPhiTaNCpassed", 18, -TMath::Pi(), +TMath::Pi());
00150 }
00151 
00152 void PatTauAnalyzer::analyze(const edm::Event& evt, const edm::EventSetup& es)
00153 {  
00154   //std::cout << "<PatTauAnalyzer::analyze>:" << std::endl; 
00155 
00156   edm::Handle<pat::TauCollection> patTaus;
00157   evt.getByLabel(src_, patTaus);
00158 
00159   hNumTauJets_->Fill(patTaus->size());
00160 
00161   for ( pat::TauCollection::const_iterator patTau = patTaus->begin(); 
00162         patTau != patTaus->end(); ++patTau ) {
00163 
00164 //--- skip fake taus in case configuration parameters set to do so...
00165     const reco::GenParticle* genTau = getGenTau(*patTau);
00166     if ( requireGenTauMatch_ && !genTau ) continue;
00167 
00168 //--- fill generator level histograms    
00169     if ( genTau ) {
00170       hGenTauEnergy_->Fill(genTau->energy());
00171       hGenTauPt_->Fill(genTau->pt());
00172       hGenTauEta_->Fill(genTau->eta());
00173       hGenTauPhi_->Fill(genTau->phi());
00174     }
00175 
00176 //--- fill reconstruction level histograms
00177 //    for Pt of highest Pt track within signal cone tau-jet...
00178     hTauJetEnergy_->Fill(patTau->energy());
00179     hTauJetPt_->Fill(patTau->pt());
00180 // 
00181 // TO-DO: 
00182 //  1.) fill histograms 
00183 //       hTauJetEta_ and hTauJetPhi_
00184 //      with the pseudo-rapidity and azimuthal angle
00185 //      of the tau-jet candidate respectively
00186 //  hTauJetEta_->...
00187 //  hTauJetPhi_->...
00188 //
00189 //  2.) fill histogram
00190 //       hTauLeadTrackPt_
00191 //      with the transverse momentum of the highest Pt ("leading") track within the tau-jet
00192 //       
00193 // NOTE: 
00194 //  1.) please have a look at
00195 //       https://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/Candidate/interface/Particle.h?revision=1.28&view=markup
00196 //      to find the methods for accessing eta and phi of the tau-jet
00197 //
00198 //  2.) please have a look at
00199 //       https://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/PatCandidates/interface/Tau.h?revision=1.25&view=markup
00200 //      to find the method for accessing the leading track
00201 //
00202 //  3.) the method pat::Tau::leadTrack returns a reference (reco::TrackRef) to a reco::Track object
00203 //      this reference can be null (in case no high Pt track has been reconstructed within the tau-jet),
00204 //      so a check if the leadTrack exists is needed before dereferencing the reco::TrackRef via operator->
00205 //
00206 //  if ( patTau->leadTrack().isAvailable() ) hTauLeadTrackPt_->Fill(patTau->leadTrack()->pt());
00207   
00208 //... for total number of tracks within signal/isolation cones
00209     hTauNumSigConeTracks_->Fill(patTau->signalTracks().size());
00210     hTauNumIsoConeTracks_->Fill(patTau->isolationTracks().size());
00211 
00212 //... for values of tau id. discriminators based on track isolation cut/
00213 //    neural network-based tau id.
00214 //    (combine with requirement of at least one "leading" track of Pt > 5. GeV 
00215 //     within the signal cone of the tau-jet)
00216     float discrByIso = ( patTau->tauID(discrByLeadTrack_.data()) > 0.5 ) ? patTau->tauID(discrByIso_.data()) : 0.;
00217     hTauDiscrByIso_->Fill(discrByIso);
00218     float discrByTaNC = ( patTau->tauID(discrByLeadTrack_.data()) > 0.5 ) ? patTau->tauID(discrByTaNC_.data()) : 0.;
00219     hTauDiscrByTaNC_->Fill(discrByTaNC);
00220 
00221 //... for values of tau id. discriminators against (unidentified) electrons and muons
00222 //
00223 // TO-DO: fill histogram
00224 //         hTauDiscrAgainstElectrons_
00225 //        with the value of the discriminatorAgainstElectrons
00226 //
00227 // NOTE:
00228 //  1.) please have a look at
00229 //       https://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/PatCandidates/interface/Tau.h?revision=1.25&view=markup
00230 //      to find the method for accessing the tau id. information
00231 //  
00232 //  2.) please have a look at
00233 //       https://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/PhysicsTools/PatAlgos/python/tools/tauTools.py?revision=1.18&view=markup
00234 //      and convince yourself that the string "againstElectron" needs to be passed as argument 
00235 //      of the pat::Tau::tauID method
00236 //
00237 //  hTauDiscrAgainstElectrons_->Fill...
00238     hTauDiscrAgainstMuons_->Fill(patTau->tauID("againstMuon"));
00239 
00240 //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByIsolation selection
00241     if ( discrByIso > 0.5 ) {
00242       hTauJetEnergyIsoPassed_->Fill(patTau->energy());
00243       hTauJetPtIsoPassed_->Fill(patTau->pt());
00244       hTauJetEtaIsoPassed_->Fill(patTau->eta());
00245       hTauJetPhiIsoPassed_->Fill(patTau->phi());
00246     }
00247 
00248 //... for Energy, Pt, Eta, Phi of tau-jets passing the discriminatorByTaNC ("Tau Neural Classifier") selection
00249     if ( discrByTaNC > 0.5 ) {
00250       hTauJetEnergyTaNCpassed_->Fill(patTau->energy());
00251       hTauJetPtTaNCpassed_->Fill(patTau->pt());
00252       hTauJetEtaTaNCpassed_->Fill(patTau->eta());
00253       hTauJetPhiTaNCpassed_->Fill(patTau->phi());
00254     }
00255   }
00256 }
00257 
00258 void PatTauAnalyzer::endJob()
00259 {
00260 //--- nothing to be done yet...
00261 }
00262 
00263 #include "FWCore/Framework/interface/MakerMacros.h"
00264 
00265 DEFINE_FWK_MODULE(PatTauAnalyzer);
00266