CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/PhysicsTools/PatExamples/plugins/PatElectronAnalyzer.cc

Go to the documentation of this file.
00001 #include "TH1.h"
00002 #include "TH2.h"
00003 #include "TMath.h"
00004 
00005 #include "FWCore/Framework/interface/Event.h"
00006 #include "FWCore/Framework/interface/EDAnalyzer.h"
00007 #include "FWCore/Utilities/interface/InputTag.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "FWCore/ServiceRegistry/interface/Service.h"
00010 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00011 
00012 class PatElectronAnalyzer : public edm::EDAnalyzer {
00013 
00014  public:
00015   
00016   explicit PatElectronAnalyzer(const edm::ParameterSet&);
00017   ~PatElectronAnalyzer();
00018   
00019  private:
00020   
00021   virtual void beginJob() ;
00022   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00023   virtual void endJob() ;
00024 
00025   // restrictions for the electron to be 
00026   // considered
00027   double minPt_;
00028   double maxEta_;
00029 
00030   // decide in which mode to run the analyzer
00031   // 0 : genMatch, 1 : tagAndProbe are 
00032   // supported depending on the line comments
00033   unsigned int mode_;
00034 
00035   // choose a given electronID for the electron
00036   // in consideration; the following types are 
00037   // available: 
00038   // * eidRobustLoose
00039   // * eidRobustTight
00040   // * eidLoose
00041   // * eidTight
00042   // * eidRobustHighEnergy
00043   std::string electronID_;
00044 
00045   // source of electrons
00046   edm::InputTag electronSrc_;
00047   // source of generator particles
00048   edm::InputTag particleSrc_;
00049 
00050   edm::ParameterSet genMatchMode_;
00051   edm::ParameterSet tagAndProbeMode_;
00052 
00053   // internal variables for genMatchMode and
00054   // tagAndProbeMode 
00055   double maxDeltaR_;
00056   double maxDeltaM_;
00057   double maxTagIso_;
00058   
00059   // book histograms of interest
00060   TH1I *nr_;
00061   TH1F *pt_;
00062   TH1F *eta_;
00063   TH1F *phi_;
00064   TH1F *genPt_;
00065   TH1F *genEta_;
00066   TH1F *genPhi_;
00067   TH1F *deltaR_;
00068   TH1F *isoTag_;
00069   TH1F *invMass_;
00070   TH1F *deltaPhi_;
00071 };  
00072 
00073 #include "Math/VectorUtil.h"
00074 #include "DataFormats/PatCandidates/interface/Electron.h"
00075 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00076 
00077 PatElectronAnalyzer::PatElectronAnalyzer(const edm::ParameterSet& cfg):
00078   minPt_ (cfg.getParameter<double>("minPt")),
00079   maxEta_ (cfg.getParameter<double>("maxEta")),
00080   mode_ (cfg.getParameter<unsigned int>("mode")),
00081   electronID_  (cfg.getParameter<std::string>("electronID")), 
00082   electronSrc_ (cfg.getParameter<edm::InputTag>("electronSrc")),
00083   particleSrc_ (cfg.getParameter<edm::InputTag>("particleSrc")),
00084   genMatchMode_(cfg.getParameter<edm::ParameterSet>("genMatchMode")),
00085   tagAndProbeMode_(cfg.getParameter<edm::ParameterSet>("tagAndProbeMode"))
00086 {
00087   // complete the configuration of the analyzer
00088   maxDeltaR_ = genMatchMode_   .getParameter<double>("maxDeltaR"); 
00089   maxDeltaM_ = tagAndProbeMode_.getParameter<double>("maxDeltaM"); 
00090   maxTagIso_ = tagAndProbeMode_.getParameter<double>("maxTagIso"); 
00091 
00092 
00093   // register histograms to the TFileService
00094   edm::Service<TFileService> fs;
00095   nr_      = fs->make<TH1I>("nr",       "nr",        10,   0 ,  10 );
00096   pt_      = fs->make<TH1F>("pt",       "pt",        20,   0., 100.);
00097   eta_     = fs->make<TH1F>("eta",      "eta",       30,  -3.,   3.);
00098   phi_     = fs->make<TH1F>("phi",      "phi",       35, -3.5,  3.5);
00099   genPt_   = fs->make<TH1F>("genPt",    "pt",        20,   0., 100.);
00100   genEta_  = fs->make<TH1F>("genEta",   "eta",       30,  -3.,   3.);
00101   genPhi_  = fs->make<TH1F>("genPhi",   "phi",       35, -3.5,  3.5);
00102   deltaR_  = fs->make<TH1F>("deltaR",   "log(dR)",   50,  -5.,   0.);
00103   isoTag_  = fs->make<TH1F>("isoTag",   "iso",       50,   0.,  10.);
00104   invMass_ = fs->make<TH1F>("invMass",  "m",        100,  50., 150.);
00105   deltaPhi_= fs->make<TH1F>("deltaPhi", "deltaPhi", 100, -3.5,  3.5);
00106 }
00107 
00108 PatElectronAnalyzer::~PatElectronAnalyzer()
00109 {
00110 }
00111 
00112 void
00113 PatElectronAnalyzer::analyze(const edm::Event& evt, const edm::EventSetup& setup)
00114 {       
00115   // get electron collection
00116   edm::Handle<std::vector<pat::Electron> > electrons;
00117   evt.getByLabel(electronSrc_, electrons); 
00118   // get generator particle collection
00119   edm::Handle<reco::GenParticleCollection> particles;
00120   evt.getByLabel(particleSrc_, particles); 
00121 
00122   nr_->Fill( electrons->size() );
00123 
00124   // ----------------------------------------------------------------------
00125   //
00126   // First Part Mode 0: genMatch
00127   //
00128   // ----------------------------------------------------------------------
00129   if( mode_==0 ){
00130     // loop generator particles 
00131     for(reco::GenParticleCollection::const_iterator part=particles->begin(); 
00132         part!=particles->end(); ++part){
00133       // only loop stable electrons
00134       if( part->status()==1  && abs(part->pdgId())==11 ){
00135         if( part->pt()>minPt_ && fabs(part->eta())<maxEta_ ){
00136           genPt_ ->Fill( part->pt()  );
00137           genEta_->Fill( part->eta() );
00138           genPhi_->Fill( part->phi() );      
00139         }
00140       }
00141     }
00142     
00143     // loop electrons
00144     for( std::vector<pat::Electron>::const_iterator elec=electrons->begin(); elec!=electrons->end(); ++elec ){
00145       if( elec->genLepton() ){
00146         float deltaR = ROOT::Math::VectorUtil::DeltaR(elec->genLepton()->p4(), elec->p4());
00147         deltaR_->Fill(TMath::Log10(deltaR));    
00148         if( deltaR<maxDeltaR_ ){
00149           if( electronID_.compare("none")!=0 ){
00150             if( elec->electronID(electronID_)<0.5 ) 
00151               continue;
00152           }
00153           if( elec->pt()>minPt_ && fabs(elec->eta())<maxEta_ ){
00154             pt_ ->Fill( elec->pt()  );
00155             eta_->Fill( elec->eta() );
00156             phi_->Fill( elec->phi() );
00157           }
00158         }
00159       }
00160     }
00161   }
00162 
00163   // ----------------------------------------------------------------------
00164   //
00165   // Second Part Mode 1: tagAndProbe
00166   //
00167   // ----------------------------------------------------------------------
00168   if( mode_==1 ){
00169     // loop tag electron
00170     for( std::vector<pat::Electron>::const_iterator elec=electrons->begin(); elec!=electrons->end(); ++elec ){
00171       isoTag_->Fill(elec->trackIso());  
00172       if( elec->trackIso()<maxTagIso_  && elec->electronID("eidTight")>0.5 ){
00173         // loop probe electron
00174         for( std::vector<pat::Electron>::const_iterator probe=electrons->begin(); probe!=electrons->end(); ++probe ){
00175           // skip the tag electron itself
00176           if( probe==elec ) continue;
00177           
00178           float zMass = (probe->p4()+elec->p4()).mass();
00179           invMass_ ->Fill(zMass);       
00180           float deltaPhi = ROOT::Math::VectorUtil::DeltaPhi(elec->p4(), probe->p4());
00181           deltaPhi_->Fill(deltaPhi);    
00182           
00183           // check for the Z mass
00184           if( fabs( zMass-90. )<maxDeltaM_ ){
00185             if( electronID_.compare("none")!=0 ){
00186               if( probe->electronID(electronID_)<0.5 ) 
00187                 continue;
00188             }
00189             if( probe->pt()>minPt_ && fabs(probe->eta())<maxEta_ ){
00190               pt_ ->Fill( probe->pt()  );
00191               eta_->Fill( probe->eta() );
00192               phi_->Fill( probe->phi() );
00193             }
00194           }
00195         }
00196       }
00197     }
00198   }
00199 }
00200 
00201 void PatElectronAnalyzer::beginJob()
00202 {
00203 }
00204 
00205 void PatElectronAnalyzer::endJob()
00206 {
00207 }
00208 
00209 #include "FWCore/Framework/interface/MakerMacros.h"
00210 DEFINE_FWK_MODULE(PatElectronAnalyzer);