CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TwoVBGenFilter.cc
Go to the documentation of this file.
1 
3 #include <iostream>
4 using namespace std;
5 
7 {
8  //now do what ever initialization is needed
9  src_ = iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag("generator", "unsmeared"));
10 
11  eejj_ = iConfig.getParameter<bool>("eejj");
12  enujj_ = iConfig.getParameter<bool>("enujj");
13 
14  mumujj_ = iConfig.getParameter<bool>("mumujj");
15  munujj_ = iConfig.getParameter<bool>("munujj");
16 
17  tautaujj_ = iConfig.getParameter<bool>("tautaujj");
18  taunujj_ = iConfig.getParameter<bool>("taunujj");
19 
20  nunujj_ = iConfig.getParameter<bool>("nunujj");
21 
22  //cout << eejj_ << endl;
23 
24 }
25 
26 
28 {
29 
30  // do anything here that needs to be done at desctruction time
31  // (e.g. close files, deallocate resources etc.)
32 
33 }
34 
35 
36 //
37 // member functions
38 //
39 
40 // ------------ method called on each new Event ------------
41 bool
43 {
44  int nj = 0;
45  int ne = 0;
46  int nnu = 0;
47  int nmu = 0;
48  int ntau = 0;
49 
51  iEvent.getByLabel(src_, evt);
52 
53  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
54 
55  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
56  p != myGenEvent->particles_end(); ++p ) {
57 
58  if(abs((*p)->pdg_id()) !=23 && abs((*p)->pdg_id()) !=24) continue; // If it is not Z or W, go to the next particle.
59 
60  if ( (*p)->end_vertex() ) {
61  for ( HepMC::GenVertex::particle_iterator des=(*p)->end_vertex()->particles_begin(HepMC::children);
62  des != (*p)->end_vertex()->particles_end(HepMC::children); ++des ) {
63 
64  const HepMC::GenParticle* theDaughter = *des;
65  if(isQuark(theDaughter)) ++nj;
66  else if(isNeutrino(theDaughter)) ++nnu;
67  else if(isElectron(theDaughter)) ++ne;
68  else if(isMuon(theDaughter)) ++nmu;
69  else if(isTau(theDaughter)) ++ntau;
70 
71  }
72  }
73  }
74 
75  if (ne==2 && nj == 2 && eejj_) return true;
76  else if (ne==1 && nj == 2 && nnu==1 && enujj_) return true;
77  else if (nmu==2 && nj == 2 && mumujj_) return true;
78  else if (nmu==1 && nj == 2 && nnu==1 && munujj_) return true;
79  else if (ntau==2 && nj == 2 && tautaujj_) return true;
80  else if (ntau==1 && nj == 2 && nnu==1 && taunujj_) return true;
81  else if (nnu==2 && nj == 2 && nunujj_) return true;
82  else return false;
83 }
84 
85 // ------------ method called once each job just before starting event loop ------------
86 void
88 {
89 }
90 
91 // ------------ method called once each job just after ending the event loop ------------
92 void
94 }
95 
97  bool result;
98  int pdgid = std::abs(p->pdg_id());
99  if(pdgid == 1 || pdgid == 2 || pdgid == 3 ||
100  pdgid == 4 || pdgid == 5 || pdgid == 6)
101  result = true;
102  else
103  result = false;
104  return result;
105 }
106 
108  bool result;
109  int pdgid = std::abs(p->pdg_id());
110  if(pdgid == 12 || pdgid == 14 || pdgid == 16)
111  result = true;
112  else
113  result = false;
114  return result;
115 }
116 
118  bool result;
119  int pdgid = std::abs(p->pdg_id());
120  if(pdgid == 11)
121  result = true;
122  else
123  result = false;
124  return result;
125 }
126 
128  bool result;
129  int pdgid = std::abs(p->pdg_id());
130  if(pdgid == 13)
131  result = true;
132  else
133  result = false;
134  return result;
135 }
136 
138  bool result;
139  int pdgid = std::abs(p->pdg_id());
140  if(pdgid == 15)
141  result = true;
142  else
143  result = false;
144  return result;
145 }
146 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool isMuon(const Candidate &part)
Definition: pdgIdUtils.h:11
bool isQuark(const HepMC::GenParticle *)
bool isNeutrino(const HepMC::GenParticle *)
int iEvent
Definition: GenABIO.cc:230
bool isNeutrino(const Candidate &part)
Definition: pdgIdUtils.h:25
bool isElectron(const Candidate &part)
Definition: pdgIdUtils.h:7
virtual void beginJob()
tuple result
Definition: query.py:137
bool isTau(const HepMC::GenParticle *)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual bool filter(edm::Event &, const edm::EventSetup &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
bool isElectron(const HepMC::GenParticle *)
virtual void endJob()
bool isMuon(const HepMC::GenParticle *)
TwoVBGenFilter(const edm::ParameterSet &)
bool isTau(const Candidate &part)
Definition: pdgIdUtils.h:15