Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:01

0001 #include "TH1.h"
0002 #include "TH2.h"
0003 #include "TMath.h"
0004 
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0007 #include "FWCore/Utilities/interface/InputTag.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ServiceRegistry/interface/Service.h"
0010 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0011 
0012 #include "DataFormats/PatCandidates/interface/Electron.h"
0013 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0014 
0015 class PatElectronAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0016 public:
0017   explicit PatElectronAnalyzer(const edm::ParameterSet &);
0018   ~PatElectronAnalyzer() override;
0019 
0020 private:
0021   void beginJob() override;
0022   void analyze(const edm::Event &, const edm::EventSetup &) override;
0023   void endJob() override;
0024 
0025   // restrictions for the electron to be
0026   // considered
0027   double minPt_;
0028   double maxEta_;
0029 
0030   // decide in which mode to run the analyzer
0031   // 0 : genMatch, 1 : tagAndProbe are
0032   // supported depending on the line comments
0033   unsigned int mode_;
0034 
0035   // choose a given electronID for the electron
0036   // in consideration; the following types are
0037   // available:
0038   // * eidRobustLoose
0039   // * eidRobustTight
0040   // * eidLoose
0041   // * eidTight
0042   // * eidRobustHighEnergy
0043   std::string electronID_;
0044 
0045   // source of electrons
0046   edm::EDGetTokenT<std::vector<pat::Electron> > electronSrcToken_;
0047   // source of generator particles
0048   edm::EDGetTokenT<reco::GenParticleCollection> particleSrcToken_;
0049 
0050   edm::ParameterSet genMatchMode_;
0051   edm::ParameterSet tagAndProbeMode_;
0052 
0053   // internal variables for genMatchMode and
0054   // tagAndProbeMode
0055   double maxDeltaR_;
0056   double maxDeltaM_;
0057   double maxTagIso_;
0058 
0059   // book histograms of interest
0060   TH1I *nr_;
0061   TH1F *pt_;
0062   TH1F *eta_;
0063   TH1F *phi_;
0064   TH1F *genPt_;
0065   TH1F *genEta_;
0066   TH1F *genPhi_;
0067   TH1F *deltaR_;
0068   TH1F *isoTag_;
0069   TH1F *invMass_;
0070   TH1F *deltaPhi_;
0071 };
0072 
0073 #include "Math/VectorUtil.h"
0074 
0075 PatElectronAnalyzer::PatElectronAnalyzer(const edm::ParameterSet &cfg)
0076     : minPt_(cfg.getParameter<double>("minPt")),
0077       maxEta_(cfg.getParameter<double>("maxEta")),
0078       mode_(cfg.getParameter<unsigned int>("mode")),
0079       electronID_(cfg.getParameter<std::string>("electronID")),
0080       electronSrcToken_(consumes<std::vector<pat::Electron> >(cfg.getParameter<edm::InputTag>("electronSrc"))),
0081       particleSrcToken_(consumes<reco::GenParticleCollection>(cfg.getParameter<edm::InputTag>("particleSrc"))),
0082       genMatchMode_(cfg.getParameter<edm::ParameterSet>("genMatchMode")),
0083       tagAndProbeMode_(cfg.getParameter<edm::ParameterSet>("tagAndProbeMode")) {
0084   usesResource(TFileService::kSharedResource);
0085 
0086   // complete the configuration of the analyzer
0087   maxDeltaR_ = genMatchMode_.getParameter<double>("maxDeltaR");
0088   maxDeltaM_ = tagAndProbeMode_.getParameter<double>("maxDeltaM");
0089   maxTagIso_ = tagAndProbeMode_.getParameter<double>("maxTagIso");
0090 
0091   // register histograms to the TFileService
0092   edm::Service<TFileService> fs;
0093   nr_ = fs->make<TH1I>("nr", "nr", 10, 0, 10);
0094   pt_ = fs->make<TH1F>("pt", "pt", 20, 0., 100.);
0095   eta_ = fs->make<TH1F>("eta", "eta", 30, -3., 3.);
0096   phi_ = fs->make<TH1F>("phi", "phi", 35, -3.5, 3.5);
0097   genPt_ = fs->make<TH1F>("genPt", "pt", 20, 0., 100.);
0098   genEta_ = fs->make<TH1F>("genEta", "eta", 30, -3., 3.);
0099   genPhi_ = fs->make<TH1F>("genPhi", "phi", 35, -3.5, 3.5);
0100   deltaR_ = fs->make<TH1F>("deltaR", "log(dR)", 50, -5., 0.);
0101   isoTag_ = fs->make<TH1F>("isoTag", "iso", 50, 0., 10.);
0102   invMass_ = fs->make<TH1F>("invMass", "m", 100, 50., 150.);
0103   deltaPhi_ = fs->make<TH1F>("deltaPhi", "deltaPhi", 100, -3.5, 3.5);
0104 }
0105 
0106 PatElectronAnalyzer::~PatElectronAnalyzer() {}
0107 
0108 void PatElectronAnalyzer::analyze(const edm::Event &evt, const edm::EventSetup &setup) {
0109   // get electron collection
0110   edm::Handle<std::vector<pat::Electron> > electrons;
0111   evt.getByToken(electronSrcToken_, electrons);
0112   // get generator particle collection
0113   edm::Handle<reco::GenParticleCollection> particles;
0114   evt.getByToken(particleSrcToken_, particles);
0115 
0116   nr_->Fill(electrons->size());
0117 
0118   // ----------------------------------------------------------------------
0119   //
0120   // First Part Mode 0: genMatch
0121   //
0122   // ----------------------------------------------------------------------
0123   if (mode_ == 0) {
0124     // loop generator particles
0125     for (reco::GenParticleCollection::const_iterator part = particles->begin(); part != particles->end(); ++part) {
0126       // only loop stable electrons
0127       if (part->status() == 1 && abs(part->pdgId()) == 11) {
0128         if (part->pt() > minPt_ && fabs(part->eta()) < maxEta_) {
0129           genPt_->Fill(part->pt());
0130           genEta_->Fill(part->eta());
0131           genPhi_->Fill(part->phi());
0132         }
0133       }
0134     }
0135 
0136     // loop electrons
0137     for (std::vector<pat::Electron>::const_iterator elec = electrons->begin(); elec != electrons->end(); ++elec) {
0138       if (elec->genLepton()) {
0139         float deltaR = ROOT::Math::VectorUtil::DeltaR(elec->genLepton()->p4(), elec->p4());
0140         deltaR_->Fill(TMath::Log10(deltaR));
0141         if (deltaR < maxDeltaR_) {
0142           if (electronID_ != "none") {
0143             if (elec->electronID(electronID_) < 0.5)
0144               continue;
0145           }
0146           if (elec->pt() > minPt_ && fabs(elec->eta()) < maxEta_) {
0147             pt_->Fill(elec->pt());
0148             eta_->Fill(elec->eta());
0149             phi_->Fill(elec->phi());
0150           }
0151         }
0152       }
0153     }
0154   }
0155 
0156   // ----------------------------------------------------------------------
0157   //
0158   // Second Part Mode 1: tagAndProbe
0159   //
0160   // ----------------------------------------------------------------------
0161   if (mode_ == 1) {
0162     // loop tag electron
0163     for (std::vector<pat::Electron>::const_iterator elec = electrons->begin(); elec != electrons->end(); ++elec) {
0164       isoTag_->Fill(elec->trackIso());
0165       if (elec->trackIso() < maxTagIso_ && elec->electronID("eidTight") > 0.5) {
0166         // loop probe electron
0167         for (std::vector<pat::Electron>::const_iterator probe = electrons->begin(); probe != electrons->end();
0168              ++probe) {
0169           // skip the tag electron itself
0170           if (probe == elec)
0171             continue;
0172 
0173           float zMass = (probe->p4() + elec->p4()).mass();
0174           invMass_->Fill(zMass);
0175           float deltaPhi = ROOT::Math::VectorUtil::DeltaPhi(elec->p4(), probe->p4());
0176           deltaPhi_->Fill(deltaPhi);
0177 
0178           // check for the Z mass
0179           if (fabs(zMass - 90.) < maxDeltaM_) {
0180             if (electronID_ != "none") {
0181               if (probe->electronID(electronID_) < 0.5)
0182                 continue;
0183             }
0184             if (probe->pt() > minPt_ && fabs(probe->eta()) < maxEta_) {
0185               pt_->Fill(probe->pt());
0186               eta_->Fill(probe->eta());
0187               phi_->Fill(probe->phi());
0188             }
0189           }
0190         }
0191       }
0192     }
0193   }
0194 }
0195 
0196 void PatElectronAnalyzer::beginJob() {}
0197 
0198 void PatElectronAnalyzer::endJob() {}
0199 
0200 #include "FWCore/Framework/interface/MakerMacros.h"
0201 DEFINE_FWK_MODULE(PatElectronAnalyzer);