CMS 3D CMS Logo

DYToMuMuGenFilter.cc
Go to the documentation of this file.
7 
10 
14 
15 
16 
18  public:
19  explicit DYToMuMuGenFilter(const edm::ParameterSet&);
20  ~DYToMuMuGenFilter() override;
21 
22  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
23 
24  private:
25  void beginStream(edm::StreamID) override;
26  bool filter(edm::Event&, const edm::EventSetup&) override;
27  void endStream() override;
28 
31 
33 
34  //virtual void beginRun(edm::Run const&, edm::EventSetup const&) override;
35  //virtual void endRun(edm::Run const&, edm::EventSetup const&) override;
36  //virtual void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
37  //virtual void endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
38 
39  // ----------member data ---------------------------
40 };
41 
42 
43 
44 
46 {
47  inputTag_= iConfig.getParameter<edm::InputTag>("inputTag");
48  genParticleCollection_ = consumes<reco::GenParticleCollection>(inputTag_);
49 }
50 
51 
53 }
54 
55 
57 
58 
60 
61  for(unsigned int i = 0; i < gen_handle->size(); i++)
62  {
63  const reco::GenParticle gen_particle = (*gen_handle)[i];
64  // Check if Z Boson decayed into two leptons
65  if (gen_particle.pdgId() == 23 && gen_particle.numberOfDaughters() == 2)
66  {
67  //Debug output
68  //std::cout << "pdgId" << gen_particle.pdgId() << std::endl;
69  //std::cout << "nDau" << gen_particle.numberOfDaughters() << std::endl;
70  //std::cout << "Dau1" << gen_particle->daughters->at(0).pdgId() << std::endl;
71  //std::cout << "Dau2" << gen_particle.numberOfDaughters() << std::endl;
72  //std::cout << "Dau1 " << gen_particle.daughter(0)->pdgId() << std::endl;
73  //std::cout << "Dau2 " << gen_particle.daughter(1)->pdgId() << std::endl;
74  //std::cout << gen_particle.daughter(1)->pdgId()+gen_particle.daughter(0)->pdgId() << std::endl;
75 
76  // Check if daugther particles are muons
77  if (std::abs(gen_particle.daughter(0)->pdgId()) == 13
78  && std::abs(gen_particle.daughter(0)->eta())<2.6
79  && std::abs(gen_particle.daughter(1)->eta())<2.6
80  && gen_particle.daughter(0)->pt()>7
81  && gen_particle.daughter(1)->pt()>7)
82  {
83  //std::cout << "pdgId" << gen_particle.pdgId() << std::endl;
84  //std::cout << "nDau" << gen_particle.numberOfDaughters() << std::endl;
85  //std::cout << "Dau1 " << gen_particle.daughter(0)->pdgId() << std::endl;
86  //std::cout << "Dau1 pt " << gen_particle.daughter(0)->pt() << std::endl;
87  //std::cout << "Dau1 pt " << gen_particle.daughter(0)->eta() << std::endl;
88  //std::cout << "Dau2 " << gen_particle.daughter(1)->pdgId() << std::endl;
89  //std::cout << "Dau2 pt " << gen_particle.daughter(1)->pt() << std::endl;
90  //std::cout << "Dau2 pt " << gen_particle.daughter(1)->eta() << std::endl;
91  //std::cout << gen_particle.daughter(1)->pdgId()+gen_particle.daughter(0)->pdgId() << std::endl;
92  return true;
93  }
94  else
95  {
96  return false;
97  }
98  }
99  }
100  return false;
101 }
102 // ------------ method called once each stream before processing any runs, lumis or events ------------
103 void
105 {
106 }
107 
108 // ------------ method called once each stream after processing all runs, lumis and events ------------
109 void
111 }
112 
113 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
114 void
116  //The following says we do not know what parameters are allowed so do no validation
117  // Please change this to state exactly what you do use, even if it is no parameters
119  desc.setUnknown();
120  descriptions.addDefault(desc);
121 }
122 
123 
124 
T getParameter(std::string const &) const
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) ...
void endStream() override
bool filter(edm::Event &, const edm::EventSetup &) override
edm::Handle< reco::GenParticleCollection > gen_handle
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
virtual int pdgId() const =0
PDG identifier.
~DYToMuMuGenFilter() override
edm::InputTag inputTag_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual double eta() const =0
momentum pseudorapidity
void beginStream(edm::StreamID) override
edm::EDGetTokenT< reco::GenParticleCollection > genParticleCollection_
virtual double pt() const =0
transverse momentum
size_t numberOfDaughters() const override
number of daughters
DYToMuMuGenFilter(const edm::ParameterSet &)