CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PythiaMomDauFilter.cc
Go to the documentation of this file.
1 
3 
4 
6 #include "HepMC/PythiaWrapper6_4.h"
7 #include <iostream>
8 
9 using namespace edm;
10 using namespace std;
11 
12 
14 label_(iConfig.getUntrackedParameter("moduleLabel",std::string("generator"))),
15 particleID(iConfig.getUntrackedParameter("ParticleID", 0)),
16 daughterID(iConfig.getUntrackedParameter("DaughterID", 0)),
17 chargeconju(iConfig.getUntrackedParameter("ChargeConjugation", true)),
18 ndaughters(iConfig.getUntrackedParameter("NumberDaughters", 0)),
19 ndescendants(iConfig.getUntrackedParameter("NumberDescendants", 0)),
20 minptcut(iConfig.getUntrackedParameter("MinPt", 0.)),
21 maxptcut(iConfig.getUntrackedParameter("MaxPt", 14000.)),
22 minetacut(iConfig.getUntrackedParameter("MinEta", -10.)),
23 maxetacut(iConfig.getUntrackedParameter("MaxEta", 10.)),
24 mom_minptcut(iConfig.getUntrackedParameter("MomMinPt", 0.)),
25 mom_maxptcut(iConfig.getUntrackedParameter("MomMaxPt", 14000.)),
26 mom_minetacut(iConfig.getUntrackedParameter("MomMinEta", -10.)),
27 mom_maxetacut(iConfig.getUntrackedParameter("MomMaxEta", 10.))
28 {
29  //now do what ever initialization is needed
30  vector<int> defdauID;
31  defdauID.push_back(0);
32  vector<int> defdesID;
33  defdesID.push_back(0);
34  dauIDs = iConfig.getUntrackedParameter< vector<int> >("DaughterIDs", defdauID);
35  desIDs = iConfig.getUntrackedParameter< vector<int> >("DescendantsIDs", defdesID);
36 }
37 
38 
40 {
41 
42  // do anything here that needs to be done at desctruction time
43  // (e.g. close files, deallocate resources etc.)
44 
45 }
46 
47 
48 //
49 // member functions
50 //
51 
52 // ------------ method called to produce the data ------------
54 {
55  using namespace edm;
56  bool accepted = false;
57  bool mom_accepted = false;
59  iEvent.getByLabel(label_, evt);
60 
61  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
62  int ndauac = 0;
63  int ndau = 0;
64  int ndesac = 0;
65  int ndes = 0;
66 
67  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
68 
69  if( (*p)->pdg_id() != particleID ) continue ;
70  if( (*p)->momentum().perp() > mom_minptcut && (*p)->momentum().perp() < mom_maxptcut && (*p)->momentum().eta() > mom_minetacut && (*p)->momentum().eta() < mom_maxetacut ){
71  mom_accepted = true;
72  ndauac = 0;
73  ndau = 0;
74  ndesac = 0;
75  ndes = 0;
76  if ( (*p)->end_vertex() ) {
77  for ( HepMC::GenVertex::particle_iterator dau =(*p)->end_vertex()->particles_begin(HepMC::children); dau != (*p)->end_vertex()->particles_end(HepMC::children); ++dau ) {
78  ++ndau;
79  for( unsigned int i=0; i<dauIDs.size(); ++i) {
80  if( (*dau)->pdg_id() != dauIDs[i] ) continue ;
81  ++ndauac;
82  break;
83  }
84  if((*dau)->pdg_id() == daughterID){
85  for(HepMC::GenVertex::particle_iterator des = (*dau)->end_vertex()->particles_begin(HepMC::children); des != (*des)->end_vertex()->particles_end(HepMC::children); ++des ){
86  ++ndes;
87  for( unsigned int i=0; i<desIDs.size(); ++i) {
88  if( (*des)->pdg_id() != desIDs[i] ) continue ;
89  if( (*des)->momentum().perp() > minptcut && (*des)->momentum().perp() < maxptcut && (*des)->momentum().eta() > minetacut && (*des)->momentum().eta() < maxetacut ) {
90  ++ndesac;
91  break;
92  }
93  }
94  }
95  }
96  }
97  }
98  }
99  if( ndau == ndaughters && ndauac == ndaughters && mom_accepted && ndes == ndescendants && ndesac == ndescendants ) {
100  accepted = true;
101  break;
102  }
103 
104  }
105 
106 
107  if( !accepted && chargeconju ) {
108 
109  mom_accepted = false;
110  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
111 
112  if( (*p)->pdg_id() != -particleID ) continue ;
113  if( (*p)->momentum().perp() > mom_minptcut && (*p)->momentum().perp() < mom_maxptcut && (*p)->momentum().eta() > mom_minetacut && (*p)->momentum().eta() < mom_maxetacut ){
114  mom_accepted = true;
115  ndauac = 0;
116  ndau = 0;
117  ndesac = 0;
118  ndes = 0;
119  if ( (*p)->end_vertex() ) {
120  for ( HepMC::GenVertex::particle_iterator dau =(*p)->end_vertex()->particles_begin(HepMC::children); dau != (*p)->end_vertex()->particles_end(HepMC::children); ++dau ) {
121  ++ndau;
122  for( unsigned int i=0; i<dauIDs.size(); ++i) {
123  int IDanti = -dauIDs[i];
124  int pythiaCode = PYCOMP(dauIDs[i]);
125  int has_antipart = pydat2.kchg[3-1][pythiaCode-1];
126  if( has_antipart == 0 ) IDanti = dauIDs[i];
127  if( (*dau)->pdg_id() != IDanti ) continue ;
128  ++ndauac;
129  break;
130  }
131  if((*dau)->pdg_id() == -daughterID) {
132  for( HepMC::GenVertex::particle_iterator des = (*dau)->end_vertex()->particles_begin(HepMC::children); des != (*des)->end_vertex()->particles_end(HepMC::children); ++des ){
133  ++ndes;
134  for( unsigned int i=0; i<desIDs.size(); ++i) {
135  int IDanti = -desIDs[i];
136  int pythiaCode = PYCOMP(desIDs[i]);
137  int has_antipart = pydat2.kchg[3-1][pythiaCode-1];
138  if( has_antipart == 0 ) IDanti = desIDs[i];
139  if( (*des)->pdg_id() != IDanti ) continue ;
140  if( (*des)->momentum().perp() > minptcut && (*des)->momentum().perp() < maxptcut && (*des)->momentum().eta() > minetacut && (*des)->momentum().eta() < maxetacut ) {
141  ++ndesac;
142  break;
143  }
144  }
145  }
146  }
147  }
148  }
149  }
150  if( ndau == ndaughters && ndauac == ndaughters && mom_accepted && ndes == ndescendants && ndesac == ndescendants ) {
151  accepted = true;
152  break;
153  }
154  }
155 
156  }
157 
158  if (accepted){
159  return true; } else {return false;}
160 
161 }
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< int > desIDs
PythiaMomDauFilter(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:230
#define PYCOMP
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
virtual bool filter(edm::Event &, const edm::EventSetup &)
std::vector< int > dauIDs