CMS 3D CMS Logo

PythiaDauFilter.cc
Go to the documentation of this file.
1 
3 
4 
6 #include <iostream>
7 
8 using namespace edm;
9 using namespace std;
10 using namespace Pythia8;
11 
12 
14 token_(consumes<edm::HepMCProduct>(edm::InputTag(iConfig.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared"))),
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  // create pythia8 instance to access particle data
28  edm::LogInfo("PythiaDauFilter::PythiaDauFilter") << "Creating pythia8 instance for particle properties" << endl;
29  if(!fLookupGen.get()) fLookupGen.reset(new Pythia());
30 }
31 
32 
34 {
35 
36  // do anything here that needs to be done at desctruction time
37  // (e.g. close files, deallocate resources etc.)
38 
39 }
40 
41 
42 //
43 // member functions
44 //
45 
46 // ------------ method called to produce the data ------------
48 {
49  using namespace edm;
50  bool accepted = false;
52  iEvent.getByToken(token_, evt);
53 
54  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
55 
56 
57  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
58  p != myGenEvent->particles_end(); ++p ) {
59 
60  if( (*p)->pdg_id() != particleID ) continue ;
61  int ndauac = 0;
62  int ndau = 0;
63  if ( (*p)->end_vertex() ) {
64  for ( HepMC::GenVertex::particle_iterator
65  des=(*p)->end_vertex()->particles_begin(HepMC::children);
66  des != (*p)->end_vertex()->particles_end(HepMC::children);
67  ++des ) {
68  ++ndau;
69  for( unsigned int i=0; i<dauIDs.size(); ++i) {
70  if( (*des)->pdg_id() != dauIDs[i] ) continue ;
71  if( (*des)->momentum().perp() > minptcut &&
72  (*des)->momentum().perp() < maxptcut &&
73  (*des)->momentum().eta() > minetacut &&
74  (*des)->momentum().eta() < maxetacut ) {
75  ++ndauac;
76  break;
77  }
78  }
79  }
80  }
81  if( ndau == ndaughters && ndauac == ndaughters ) {
82  accepted = true;
83  break;
84  }
85 
86  }
87 
88 
89  if( !accepted && chargeconju ) {
90 
91  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
92  p != myGenEvent->particles_end(); ++p ) {
93 
94  if( (*p)->pdg_id() != -particleID ) continue ;
95  int ndauac = 0;
96  int ndau = 0;
97  if ( (*p)->end_vertex() ) {
98  for ( HepMC::GenVertex::particle_iterator
99  des=(*p)->end_vertex()->particles_begin(HepMC::children);
100  des != (*p)->end_vertex()->particles_end(HepMC::children);
101  ++des ) {
102  ++ndau;
103  for( unsigned int i=0; i<dauIDs.size(); ++i) {
104  int IDanti = -dauIDs[i];
105  if ( !(fLookupGen->particleData.isParticle( IDanti )) ) IDanti = dauIDs[i];
106  if( (*des)->pdg_id() != IDanti ) continue ;
107  if( (*des)->momentum().perp() > minptcut &&
108  (*des)->momentum().perp() < maxptcut &&
109  (*des)->momentum().eta() > minetacut &&
110  (*des)->momentum().eta() < maxetacut ) {
111  ++ndauac;
112  break;
113  }
114  }
115  }
116  }
117  if( ndau == ndaughters && ndauac == ndaughters ) {
118  accepted = true;
119  break;
120  }
121  }
122 
123  }
124 
125  if (accepted){
126  return true; } else {return false;}
127 
128 }
std::unique_ptr< Pythia8::Pythia > fLookupGen
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > dauIDs
const double maxetacut
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
const double minptcut
PythiaDauFilter(const edm::ParameterSet &)
const double maxptcut
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< edm::HepMCProduct > token_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
~PythiaDauFilter() override
bool accepted(std::vector< std::string_view > const &, std::string_view)
const bool chargeconju
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
HLT enums.
const int particleID
const double minetacut
const int ndaughters