Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:12

0001 // Z->ee Filter
0002 
0003 // user includes
0004 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0007 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0008 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0009 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0010 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0011 #include "DataFormats/TrackReco/interface/HitPattern.h"
0012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/Frameworkfwd.h"
0015 #include "FWCore/Framework/interface/stream/EDFilter.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0020 
0021 // ROOT includes
0022 #include "TLorentzVector.h"
0023 
0024 class ZtoEEEventSelector : public edm::stream::EDFilter<> {
0025 public:
0026   explicit ZtoEEEventSelector(const edm::ParameterSet&);
0027   bool filter(edm::Event&, edm::EventSetup const&) override;
0028   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0029 
0030 private:
0031   // module config parameters
0032   const edm::InputTag electronTag_;
0033   const edm::InputTag bsTag_;
0034   const edm::EDGetTokenT<reco::GsfElectronCollection> electronToken_;
0035   const edm::EDGetTokenT<reco::BeamSpot> bsToken_;
0036 
0037   const double maxEta_;
0038   const double minPt_;
0039   const double maxDeltaPhiInEB_;
0040   const double maxDeltaEtaInEB_;
0041   const double maxHOEEB_;
0042   const double maxSigmaiEiEEB_;
0043   const double maxDeltaPhiInEE_;
0044   const double maxDeltaEtaInEE_;
0045   const double maxHOEEE_;
0046   const double maxSigmaiEiEEE_;
0047   const double maxNormChi2_;
0048   const double maxD0_;
0049   const double maxDz_;
0050   const int minPixelHits_;
0051   const int minStripHits_;
0052   const double maxIso_;
0053   const double minPtHighest_;
0054   const double minInvMass_;
0055   const double maxInvMass_;
0056 };
0057 
0058 using namespace std;
0059 using namespace edm;
0060 
0061 void ZtoEEEventSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0062   edm::ParameterSetDescription desc;
0063   desc.addUntracked<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"));
0064   desc.addUntracked<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"));
0065   desc.addUntracked<double>("maxEta", 2.4);
0066   desc.addUntracked<double>("minPt", 5);
0067   desc.addUntracked<double>("maxDeltaPhiInEB", .15);
0068   desc.addUntracked<double>("maxDeltaEtaInEB", .007);
0069   desc.addUntracked<double>("maxHOEEB", .12);
0070   desc.addUntracked<double>("maxSigmaiEiEEB", .01);
0071   desc.addUntracked<double>("maxDeltaPhiInEE", .1);
0072   desc.addUntracked<double>("maxDeltaEtaInEE", .009);
0073   desc.addUntracked<double>("maxHOEEB_", .10);
0074   desc.addUntracked<double>("maxSigmaiEiEEE", .03);
0075   desc.addUntracked<double>("maxNormChi2", 1000);
0076   desc.addUntracked<double>("maxD0", 0.02);
0077   desc.addUntracked<double>("maxDz", 20.);
0078   desc.addUntracked<uint32_t>("minPixelHits", 1);
0079   desc.addUntracked<uint32_t>("minStripHits", 8);
0080   desc.addUntracked<double>("maxIso", 0.3);
0081   desc.addUntracked<double>("minPtHighest", 24);
0082   desc.addUntracked<double>("minInvMass", 75);
0083   desc.addUntracked<double>("maxInvMass", 105);
0084   descriptions.addWithDefaultLabel(desc);
0085 }
0086 
0087 ZtoEEEventSelector::ZtoEEEventSelector(const edm::ParameterSet& ps)
0088     : electronTag_(ps.getUntrackedParameter<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"))),
0089       bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
0090       electronToken_(consumes<reco::GsfElectronCollection>(electronTag_)),
0091       bsToken_(consumes<reco::BeamSpot>(bsTag_)),
0092       maxEta_(ps.getUntrackedParameter<double>("maxEta", 2.4)),
0093       minPt_(ps.getUntrackedParameter<double>("minPt", 5)),
0094       maxDeltaPhiInEB_(ps.getUntrackedParameter<double>("maxDeltaPhiInEB", .15)),
0095       maxDeltaEtaInEB_(ps.getUntrackedParameter<double>("maxDeltaEtaInEB", .007)),
0096       maxHOEEB_(ps.getUntrackedParameter<double>("maxHOEEB", .12)),
0097       maxSigmaiEiEEB_(ps.getUntrackedParameter<double>("maxSigmaiEiEEB", .01)),
0098       maxDeltaPhiInEE_(ps.getUntrackedParameter<double>("maxDeltaPhiInEE", .1)),
0099       maxDeltaEtaInEE_(ps.getUntrackedParameter<double>("maxDeltaEtaInEE", .009)),
0100       maxHOEEE_(ps.getUntrackedParameter<double>("maxHOEEB_", .10)),
0101       maxSigmaiEiEEE_(ps.getUntrackedParameter<double>("maxSigmaiEiEEE", .03)),
0102       maxNormChi2_(ps.getUntrackedParameter<double>("maxNormChi2", 1000)),
0103       maxD0_(ps.getUntrackedParameter<double>("maxD0", 0.02)),
0104       maxDz_(ps.getUntrackedParameter<double>("maxDz", 20.)),
0105       minPixelHits_(ps.getUntrackedParameter<uint32_t>("minPixelHits", 1)),
0106       minStripHits_(ps.getUntrackedParameter<uint32_t>("minStripHits", 8)),
0107       maxIso_(ps.getUntrackedParameter<double>("maxIso", 0.3)),
0108       minPtHighest_(ps.getUntrackedParameter<double>("minPtHighest", 24)),
0109       minInvMass_(ps.getUntrackedParameter<double>("minInvMass", 75)),
0110       maxInvMass_(ps.getUntrackedParameter<double>("maxInvMass", 105)) {}
0111 
0112 bool ZtoEEEventSelector::filter(edm::Event& iEvent, edm::EventSetup const& iSetup) {
0113   // Read Electron Collection
0114   edm::Handle<reco::GsfElectronCollection> electronColl;
0115   iEvent.getByToken(electronToken_, electronColl);
0116 
0117   edm::Handle<reco::BeamSpot> beamSpot;
0118   iEvent.getByToken(bsToken_, beamSpot);
0119 
0120   std::vector<TLorentzVector> list;
0121   std::vector<int> chrgeList;
0122 
0123   if (electronColl.isValid()) {
0124     for (auto const& ele : *electronColl) {
0125       if (!ele.ecalDriven())
0126         continue;
0127       if (ele.pt() < minPt_)
0128         continue;
0129       // set a max Eta cut
0130       if (!(ele.isEB() || ele.isEE()))
0131         continue;
0132 
0133       double hOverE = ele.hadronicOverEm();
0134       double sigmaee = ele.sigmaIetaIeta();
0135       double deltaPhiIn = ele.deltaPhiSuperClusterTrackAtVtx();
0136       double deltaEtaIn = ele.deltaEtaSuperClusterTrackAtVtx();
0137 
0138       // separate cut for barrel and endcap
0139       if (ele.isEB()) {
0140         if (fabs(deltaPhiIn) >= maxDeltaPhiInEB_ && fabs(deltaEtaIn) >= maxDeltaEtaInEB_ && hOverE >= maxHOEEB_ &&
0141             sigmaee >= maxSigmaiEiEEB_)
0142           continue;
0143       } else if (ele.isEE()) {
0144         if (fabs(deltaPhiIn) >= maxDeltaPhiInEE_ && fabs(deltaEtaIn) >= maxDeltaEtaInEE_ && hOverE >= maxHOEEE_ &&
0145             sigmaee >= maxSigmaiEiEEE_)
0146           continue;
0147       }
0148 
0149       reco::GsfTrackRef trk = ele.gsfTrack();
0150       if (!trk.isNonnull())
0151         continue;  // only electrons with tracks
0152       double chi2 = trk->chi2();
0153       double ndof = trk->ndof();
0154       double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
0155       if (chbyndof >= maxNormChi2_)
0156         continue;
0157 
0158       double trkd0 = trk->d0();
0159       if (beamSpot.isValid()) {
0160         trkd0 = -(trk->dxy(beamSpot->position()));
0161       } else {
0162         edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get BeamSpot for label: " << bsTag_;
0163       }
0164       if (std::fabs(trkd0) >= maxD0_)
0165         continue;
0166 
0167       const reco::HitPattern& hitp = trk->hitPattern();
0168       int nPixelHits = hitp.numberOfValidPixelHits();
0169       if (nPixelHits < minPixelHits_)
0170         continue;
0171 
0172       int nStripHits = hitp.numberOfValidStripHits();
0173       if (nStripHits < minStripHits_)
0174         continue;
0175 
0176       // PF Isolation
0177       reco::GsfElectron::PflowIsolationVariables pfIso = ele.pfIsolationVariables();
0178       float absiso =
0179           pfIso.sumChargedHadronPt + std::max(0.0, pfIso.sumNeutralHadronEt + pfIso.sumPhotonEt - 0.5 * pfIso.sumPUPt);
0180       float eiso = absiso / (ele.pt());
0181       if (eiso > maxIso_)
0182         continue;
0183 
0184       TLorentzVector lv;
0185       lv.SetPtEtaPhiE(ele.pt(), ele.eta(), ele.phi(), ele.energy());
0186       list.push_back(lv);
0187       chrgeList.push_back(ele.charge());
0188     }
0189   } else {
0190     edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get ElectronCollection for label: " << electronTag_;
0191   }
0192   if (list.size() < 2)
0193     return false;
0194   if (chrgeList[0] + chrgeList[1] != 0)
0195     return false;
0196 
0197   if (list[0].Pt() < minPtHighest_)
0198     return false;
0199   TLorentzVector zv = list[0] + list[1];
0200   if (zv.M() < minInvMass_ || zv.M() > maxInvMass_)
0201     return false;
0202 
0203   return true;
0204 }
0205 // Define this as a plug-in
0206 #include "FWCore/Framework/interface/MakerMacros.h"
0207 DEFINE_FWK_MODULE(ZtoEEEventSelector);