CMS 3D CMS Logo

ModelpMSSMFilter.cc
Go to the documentation of this file.
1 /*
2 Package: GeneralInterface/GenFilters/ModelpMSSMFilter
3 Class: ModelpMSSMFilter
4 
5 class ModelpMSSMFilter ModelpMSSMFilter.cc GeneratorInterface/GenFilters/plugins/ModelpMSSMFilter.cc
6 
7 Description: EDFilter which checks the event passes a baseline selection for the run-II pMSSM effort.
8 
9 Implementation:
10 
11 The following input parameters are used
12  gpssrc = cms.InputTag("X") : gen particle collection label as input
13  jetsrc = cms.InputTag("X") : genjet collection collection label as input
14  jetPtCut = cms.double(#) : GenJet pT cut for HT
15  jetEtaCut = cms.double(#) : GenJet eta cut for HT
16  genHTcut = cms.double(#) : GenHT cut
17  muPtCut = cms.double(#) : muon pT cut
18  muEtaCut = cms.double(#) : muon eta cut
19  elPtCut = cms.double(#) : electron pT cut
20  elEtaCut = cms.double(#) : electron eta cut
21  gammaPtCut = cms.double(#) : photon pT cut
22  gammaEtaCut = cms.double(#) : photon eta cut
23  loosemuPtCut = cms.double(#) : loose muon pT cut
24  looseelPtCut = cms.double(#) : loose electron pT cut
25  loosegammaPtCut = cms.double(#) : loose photon pT cut
26  veryloosegammaPtCut = cms.double(#) : even looser photon pT cut
27 Original Author: Malte Mrowietz
28  Created: Jun 2019
29 */
30 
31 //System include files
32 #include <cmath>
33 #include <memory>
34 #include <vector>
35 //User include files
40 
43 
47 
48 //Class declaration
50 public:
51  explicit ModelpMSSMFilter(const edm::ParameterSet&);
52  ~ModelpMSSMFilter() override;
53 
54 private:
55  bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
56 
57  //Member data
62 };
63 
64 //Constructor
66  : token_(consumes<reco::GenParticleCollection>(params.getParameter<edm::InputTag>("gpssrc"))),
67  token2_(consumes<reco::GenJetCollection>(params.getParameter<edm::InputTag>("jetsrc"))),
68  muPtCut_(params.getParameter<double>("muPtCut")),
69  muEtaCut_(params.getParameter<double>("muEtaCut")),
70  tauPtCut_(params.getParameter<double>("tauPtCut")),
71  tauEtaCut_(params.getParameter<double>("tauEtaCut")),
72  elPtCut_(params.getParameter<double>("elPtCut")),
73  elEtaCut_(params.getParameter<double>("elEtaCut")),
74  gammaPtCut_(params.getParameter<double>("gammaPtCut")),
75  gammaEtaCut_(params.getParameter<double>("gammaEtaCut")),
76  loosemuPtCut_(params.getParameter<double>("loosemuPtCut")),
77  looseelPtCut_(params.getParameter<double>("looseelPtCut")),
78  loosegammaPtCut_(params.getParameter<double>("loosegammaPtCut")),
79  veryloosegammaPtCut_(params.getParameter<double>("veryloosegammaPtCut")),
80  jetPtCut_(params.getParameter<double>("jetPtCut")),
81  jetEtaCut_(params.getParameter<double>("jetEtaCut")),
82  genHTcut_(params.getParameter<double>("genHTcut")) {}
83 
84 //Destructor
86 
88  using namespace std;
89  using namespace edm;
90  using namespace reco;
92  evt.getByToken(token_, gps);
94  evt.getByToken(token2_, generatedJets);
95  int looseel = 0;
96  int loosemu = 0;
97  int loosegamma = 0;
98  int veryloosegamma = 0;
99  float decaylength;
100  for (std::vector<reco::GenParticle>::const_iterator it = gps->begin(); it != gps->end(); ++it) {
101  const reco::GenParticle& gp = *it;
102  if (gp.isLastCopy()) {
103  if (fabs(gp.pdgId()) == 15) {
104  if (gp.pt() > tauPtCut_ && fabs(gp.eta()) < tauEtaCut_) {
105  return true;
106  }
107  }
108  if (fabs(gp.pdgId()) == 13) {
109  if (gp.pt() > muPtCut_ && fabs(gp.eta()) < muEtaCut_) {
110  return true;
111  }
112  if (gp.pt() > loosemuPtCut_ && fabs(gp.eta()) < muEtaCut_) {
113  loosemu += 1;
114  }
115  }
116  if (fabs(gp.pdgId()) == 11) {
117  if (gp.pt() > elPtCut_ && fabs(gp.eta()) < elEtaCut_) {
118  return true;
119  }
120  if (gp.pt() > looseelPtCut_ && fabs(gp.eta()) < elEtaCut_) {
121  looseel += 1;
122  }
123  }
124  if (fabs(gp.pdgId()) == 22) {
125  if (gp.pt() > gammaPtCut_ && fabs(gp.eta()) < gammaEtaCut_) {
126  return true;
127  }
128  if (gp.pt() > loosegammaPtCut_ && fabs(gp.eta()) < gammaEtaCut_) {
129  loosegamma += 1;
130  } else {
131  if (gp.pt() > veryloosegammaPtCut_ && fabs(gp.eta()) < gammaEtaCut_) {
132  veryloosegamma += 1;
133  }
134  }
135  }
136  if (fabs(gp.pdgId()) == 1000024) {
137  if (gp.numberOfDaughters() > 0) {
138  decaylength = sqrt(pow(gp.vx() - gp.daughter(0)->vx(), 2) + pow(gp.vy() - gp.daughter(0)->vy(), 2));
139  if (decaylength > 300) {
140  return true;
141  }
142  } else {
143  return true;
144  }
145  }
146  }
147  }
148  if (looseel + loosemu + loosegamma > 1) {
149  return true;
150  }
151  if (loosegamma > 0 && veryloosegamma > 0) {
152  return true;
153  }
154  double genHT = 0.0;
155  for (std::vector<reco::GenJet>::const_iterator it = generatedJets->begin(); it != generatedJets->end(); ++it) {
156  const reco::GenJet& gjet = *it;
157  //Add GenJet pt to genHT if GenJet complies with given HT definition
158  if (gjet.pt() > jetPtCut_ && fabs(gjet.eta()) < jetEtaCut_) {
159  genHT += gjet.pt();
160  }
161  }
162  return (genHT > genHTcut_);
163 }
164 
165 // Define module as a plug-in
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int pdgId() const final
PDG identifier.
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
double eta() const final
momentum pseudorapidity
virtual double vx() const =0
x coordinate of vertex position
double vy() const override
y coordinate of vertex position
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool isLastCopy() const
Definition: GenParticle.h:98
std::vector< GenJet > GenJetCollection
collection of GenJet objects
double pt() const final
transverse momentum
virtual double vy() const =0
y coordinate of vertex position
~ModelpMSSMFilter() override
edm::EDGetTokenT< reco::GenJetCollection > token2_
edm::EDGetTokenT< reco::GenParticleCollection > token_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
Jets made from MC generator particles.
Definition: GenJet.h:25
ModelpMSSMFilter(const edm::ParameterSet &)
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
size_t numberOfDaughters() const override
number of daughters
fixed size matrix
HLT enums.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double vx() const override
x coordinate of vertex position