CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LQGenFilter.cc
Go to the documentation of this file.
1 
3 
5 {
6  //now do what ever initialization is needed
7  src_ = iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag("source"));
8 
9  eejj_ = iConfig.getParameter<bool>("eejj");
10  enuejj_ = iConfig.getParameter<bool>("enuejj");
11  nuenuejj_ = iConfig.getParameter<bool>("nuenuejj");
12 
13  mumujj_ = iConfig.getParameter<bool>("mumujj");
14  munumujj_ = iConfig.getParameter<bool>("munumujj");
15  numunumujj_ = iConfig.getParameter<bool>("numunumujj");
16 }
17 
18 
20 {
21 
22  // do anything here that needs to be done at desctruction time
23  // (e.g. close files, deallocate resources etc.)
24 
25 }
26 
27 
28 //
29 // member functions
30 //
31 
32 // ------------ method called on each new Event ------------
33 bool
35 {
36  int ne = 0;
37  int nnue = 0;
38  int nmu = 0;
39  int nnumu = 0;
40 
42  iEvent.getByLabel(src_, evt);
43 
44  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
45 
46  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
47  p != myGenEvent->particles_end(); ++p ) {
48 
49  if( abs((*p)->pdg_id()) != 42 ) continue; // skip if not a leptoquark
50 
51  if ( (*p)->end_vertex() ) {
52  for ( HepMC::GenVertex::particle_iterator des=(*p)->end_vertex()->particles_begin(HepMC::children);
53  des != (*p)->end_vertex()->particles_end(HepMC::children); ++des ) {
54 
55  if ( abs((*des)->pdg_id()) == 11 ) ++ne;
56  else if ( abs((*des)->pdg_id()) == 12 ) ++nnue;
57  else if ( abs((*des)->pdg_id()) == 13 ) ++nmu;
58  else if ( abs((*des)->pdg_id()) == 14 ) ++nnumu;
59  }
60  }
61  }
62 
63  if (ne==2 && eejj_) return true;
64  else if (ne==1 && nnue==1 && enuejj_) return true;
65  else if (nnue==2 && nuenuejj_) return true;
66  else if (nmu==2 && mumujj_) return true;
67  else if (nmu==1 && nnumu==1 && munumujj_) return true;
68  else if (nnumu==2 && numunumujj_) return true;
69  else return false;
70 }
71 
72 // ------------ method called once each job just before starting event loop ------------
73 void
75 {
76 }
77 
78 // ------------ method called once each job just after ending the event loop ------------
79 void
81 }
82 
T getParameter(std::string const &) const
bool numunumujj_
Definition: LQGenFilter.h:50
T getUntrackedParameter(std::string const &, T const &) const
bool enuejj_
Definition: LQGenFilter.h:50
virtual void endJob()
Definition: LQGenFilter.cc:80
LQGenFilter(const edm::ParameterSet &)
Definition: LQGenFilter.cc:4
int iEvent
Definition: GenABIO.cc:230
bool mumujj_
Definition: LQGenFilter.h:50
edm::InputTag src_
Definition: LQGenFilter.h:49
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
virtual bool filter(edm::Event &, const edm::EventSetup &)
Definition: LQGenFilter.cc:34
bool nuenuejj_
Definition: LQGenFilter.h:50
virtual void beginJob()
Definition: LQGenFilter.cc:74
bool munumujj_
Definition: LQGenFilter.h:50