CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PythiaDauFilter.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 chargeconju(iConfig.getUntrackedParameter("ChargeConjugation", true)),
17 ndaughters(iConfig.getUntrackedParameter("NumberDaughters", 0)),
18 minptcut(iConfig.getUntrackedParameter("MinPt", 0.)),
19 maxptcut(iConfig.getUntrackedParameter("MaxPt", 14000.)),
20 minetacut(iConfig.getUntrackedParameter("MinEta", -10.)),
21 maxetacut(iConfig.getUntrackedParameter("MaxEta", 10.))
22 {
23  //now do what ever initialization is needed
24  vector<int> defdauID;
25  defdauID.push_back(0);
26  dauIDs = iConfig.getUntrackedParameter< vector<int> >("DaughterIDs",defdauID);
27 }
28 
29 
31 {
32 
33  // do anything here that needs to be done at desctruction time
34  // (e.g. close files, deallocate resources etc.)
35 
36 }
37 
38 
39 //
40 // member functions
41 //
42 
43 // ------------ method called to produce the data ------------
45 {
46  using namespace edm;
47  bool accepted = false;
49  iEvent.getByLabel(label_, evt);
50 
51  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
52 
53 
54  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
55  p != myGenEvent->particles_end(); ++p ) {
56 
57  if( (*p)->pdg_id() != particleID ) continue ;
58  int ndauac = 0;
59  int ndau = 0;
60  if ( (*p)->end_vertex() ) {
61  for ( HepMC::GenVertex::particle_iterator
62  des=(*p)->end_vertex()->particles_begin(HepMC::children);
63  des != (*p)->end_vertex()->particles_end(HepMC::children);
64  ++des ) {
65  ++ndau;
66  for( unsigned int i=0; i<dauIDs.size(); ++i) {
67  if( (*des)->pdg_id() != dauIDs[i] ) continue ;
68  if( (*des)->momentum().perp() > minptcut &&
69  (*des)->momentum().perp() < maxptcut &&
70  (*des)->momentum().eta() > minetacut &&
71  (*des)->momentum().eta() < maxetacut ) {
72  ++ndauac;
73  break;
74  }
75  }
76  }
77  }
78  if( ndau == ndaughters && ndauac == ndaughters ) {
79  accepted = true;
80  break;
81  }
82 
83  }
84 
85 
86  if( !accepted && chargeconju ) {
87 
88  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
89  p != myGenEvent->particles_end(); ++p ) {
90 
91  if( (*p)->pdg_id() != -particleID ) continue ;
92  int ndauac = 0;
93  int ndau = 0;
94  if ( (*p)->end_vertex() ) {
95  for ( HepMC::GenVertex::particle_iterator
96  des=(*p)->end_vertex()->particles_begin(HepMC::children);
97  des != (*p)->end_vertex()->particles_end(HepMC::children);
98  ++des ) {
99  ++ndau;
100  for( unsigned int i=0; i<dauIDs.size(); ++i) {
101  int IDanti = -dauIDs[i];
102  int pythiaCode = PYCOMP(dauIDs[i]);
103  int has_antipart = pydat2.kchg[3-1][pythiaCode-1];
104  if( has_antipart == 0 ) IDanti = dauIDs[i];
105  if( (*des)->pdg_id() != IDanti ) continue ;
106  if( (*des)->momentum().perp() > minptcut &&
107  (*des)->momentum().perp() < maxptcut &&
108  (*des)->momentum().eta() > minetacut &&
109  (*des)->momentum().eta() < maxetacut ) {
110  ++ndauac;
111  break;
112  }
113  }
114  }
115  }
116  if( ndau == ndaughters && ndauac == ndaughters ) {
117  accepted = true;
118  break;
119  }
120  }
121 
122  }
123 
124  if (accepted){
125  return true; } else {return false;}
126 
127 }
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< int > dauIDs
PythiaDauFilter(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:243
#define PYCOMP
virtual bool filter(edm::Event &, const edm::EventSetup &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
std::string label_