Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:13:40

0001 /*
0002 Package:    GeneralInterface/GenFilters/ModelpMSSMFilter
0003 Class:      ModelpMSSMFilter
0004 
0005 class ModelpMSSMFilter ModelpMSSMFilter.cc GeneratorInterface/GenFilters/plugins/ModelpMSSMFilter.cc
0006 
0007 Description: EDFilter which checks the event passes a baseline selection for the run-II pMSSM effort.
0008 
0009 Implementation: 
0010 
0011 The following input parameters are used 
0012   gpssrc = cms.InputTag("X") : gen particle collection label as input
0013   jetsrc  = cms.InputTag("X") : genjet collection collection label as input
0014   jetPtCut = cms.double(#) : GenJet pT cut for HT
0015   jetEtaCut = cms.double(#) : GenJet eta cut for HT
0016   genHTcut = cms.double(#) : GenHT cut
0017   muPtCut = cms.double(#) : muon pT cut
0018   muEtaCut = cms.double(#) : muon eta cut
0019   elPtCut = cms.double(#) : electron pT cut
0020   elEtaCut = cms.double(#) : electron eta cut
0021   gammaPtCut = cms.double(#) : photon pT cut
0022   gammaEtaCut = cms.double(#) : photon eta cut
0023   loosemuPtCut = cms.double(#) : loose muon pT cut
0024   looseelPtCut = cms.double(#) : loose electron pT cut
0025   loosegammaPtCut = cms.double(#) : loose photon pT cut
0026   veryloosegammaPtCut = cms.double(#) : even looser photon pT cut
0027 Original Author:  Malte Mrowietz
0028          Created:  Jun 2019
0029 */
0030 
0031 //System include files
0032 #include <cmath>
0033 #include <memory>
0034 #include <vector>
0035 //User include files
0036 #include "FWCore/Framework/interface/global/EDFilter.h"
0037 #include "FWCore/Framework/interface/Frameworkfwd.h"
0038 #include "FWCore/Framework/interface/Event.h"
0039 #include "FWCore/Framework/interface/MakerMacros.h"
0040 
0041 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0042 #include "FWCore/Utilities/interface/InputTag.h"
0043 
0044 #include "DataFormats/Common/interface/Handle.h"
0045 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0046 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0047 
0048 //Class declaration
0049 class ModelpMSSMFilter : public edm::global::EDFilter<> {
0050 public:
0051   explicit ModelpMSSMFilter(const edm::ParameterSet&);
0052   ~ModelpMSSMFilter() override;
0053 
0054 private:
0055   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0056 
0057   //Member data
0058   edm::EDGetTokenT<reco::GenParticleCollection> token_;
0059   edm::EDGetTokenT<reco::GenJetCollection> token2_;
0060   double muPtCut_, muEtaCut_, tauPtCut_, tauEtaCut_, elPtCut_, elEtaCut_, gammaPtCut_, gammaEtaCut_, loosemuPtCut_,
0061       looseelPtCut_, loosegammaPtCut_, veryloosegammaPtCut_, jetPtCut_, jetEtaCut_, genHTcut_;
0062 };
0063 
0064 //Constructor
0065 ModelpMSSMFilter::ModelpMSSMFilter(const edm::ParameterSet& params)
0066     : token_(consumes<reco::GenParticleCollection>(params.getParameter<edm::InputTag>("gpssrc"))),
0067       token2_(consumes<reco::GenJetCollection>(params.getParameter<edm::InputTag>("jetsrc"))),
0068       muPtCut_(params.getParameter<double>("muPtCut")),
0069       muEtaCut_(params.getParameter<double>("muEtaCut")),
0070       tauPtCut_(params.getParameter<double>("tauPtCut")),
0071       tauEtaCut_(params.getParameter<double>("tauEtaCut")),
0072       elPtCut_(params.getParameter<double>("elPtCut")),
0073       elEtaCut_(params.getParameter<double>("elEtaCut")),
0074       gammaPtCut_(params.getParameter<double>("gammaPtCut")),
0075       gammaEtaCut_(params.getParameter<double>("gammaEtaCut")),
0076       loosemuPtCut_(params.getParameter<double>("loosemuPtCut")),
0077       looseelPtCut_(params.getParameter<double>("looseelPtCut")),
0078       loosegammaPtCut_(params.getParameter<double>("loosegammaPtCut")),
0079       veryloosegammaPtCut_(params.getParameter<double>("veryloosegammaPtCut")),
0080       jetPtCut_(params.getParameter<double>("jetPtCut")),
0081       jetEtaCut_(params.getParameter<double>("jetEtaCut")),
0082       genHTcut_(params.getParameter<double>("genHTcut")) {}
0083 
0084 //Destructor
0085 ModelpMSSMFilter::~ModelpMSSMFilter() {}
0086 
0087 bool ModelpMSSMFilter::filter(edm::StreamID, edm::Event& evt, const edm::EventSetup& params) const {
0088   using namespace std;
0089   using namespace edm;
0090   using namespace reco;
0091   edm::Handle<reco::GenParticleCollection> gps;
0092   evt.getByToken(token_, gps);
0093   edm::Handle<reco::GenJetCollection> generatedJets;
0094   evt.getByToken(token2_, generatedJets);
0095   int looseel = 0;
0096   int loosemu = 0;
0097   int loosegamma = 0;
0098   int veryloosegamma = 0;
0099   float decaylength;
0100   for (std::vector<reco::GenParticle>::const_iterator it = gps->begin(); it != gps->end(); ++it) {
0101     const reco::GenParticle& gp = *it;
0102     if (gp.isLastCopy()) {
0103       if (fabs(gp.pdgId()) == 15) {
0104         if (gp.pt() > tauPtCut_ && fabs(gp.eta()) < tauEtaCut_) {
0105           return true;
0106         }
0107       }
0108       if (fabs(gp.pdgId()) == 13) {
0109         if (gp.pt() > muPtCut_ && fabs(gp.eta()) < muEtaCut_) {
0110           return true;
0111         }
0112         if (gp.pt() > loosemuPtCut_ && fabs(gp.eta()) < muEtaCut_) {
0113           loosemu += 1;
0114         }
0115       }
0116       if (fabs(gp.pdgId()) == 11) {
0117         if (gp.pt() > elPtCut_ && fabs(gp.eta()) < elEtaCut_) {
0118           return true;
0119         }
0120         if (gp.pt() > looseelPtCut_ && fabs(gp.eta()) < elEtaCut_) {
0121           looseel += 1;
0122         }
0123       }
0124       if (fabs(gp.pdgId()) == 22) {
0125         if (gp.pt() > gammaPtCut_ && fabs(gp.eta()) < gammaEtaCut_) {
0126           return true;
0127         }
0128         if (gp.pt() > loosegammaPtCut_ && fabs(gp.eta()) < gammaEtaCut_) {
0129           loosegamma += 1;
0130         } else {
0131           if (gp.pt() > veryloosegammaPtCut_ && fabs(gp.eta()) < gammaEtaCut_) {
0132             veryloosegamma += 1;
0133           }
0134         }
0135       }
0136       if (fabs(gp.pdgId()) == 1000024) {
0137         if (gp.numberOfDaughters() > 0) {
0138           decaylength = sqrt(pow(gp.vx() - gp.daughter(0)->vx(), 2) + pow(gp.vy() - gp.daughter(0)->vy(), 2));
0139           if (decaylength > 300) {
0140             return true;
0141           }
0142         } else {
0143           return true;
0144         }
0145       }
0146     }
0147   }
0148   if (looseel + loosemu + loosegamma > 1) {
0149     return true;
0150   }
0151   if (loosegamma > 0 && veryloosegamma > 0) {
0152     return true;
0153   }
0154   double genHT = 0.0;
0155   for (std::vector<reco::GenJet>::const_iterator it = generatedJets->begin(); it != generatedJets->end(); ++it) {
0156     const reco::GenJet& gjet = *it;
0157     //Add GenJet pt to genHT if GenJet complies with given HT definition
0158     if (gjet.pt() > jetPtCut_ && fabs(gjet.eta()) < jetEtaCut_) {
0159       genHT += gjet.pt();
0160     }
0161   }
0162   return (genHT > genHTcut_);
0163 }
0164 
0165 // Define module as a plug-in
0166 DEFINE_FWK_MODULE(ModelpMSSMFilter);