CMS 3D CMS Logo

PythiaFilter.cc
Go to the documentation of this file.
1 
4 
6 #include <iostream>
7 
8 using namespace edm;
9 using namespace std;
10 
11 
13 token_(consumes<edm::HepMCProduct>(edm::InputTag(iConfig.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared"))),
14 particleID(iConfig.getUntrackedParameter("ParticleID", 0)),
15 minpcut(iConfig.getUntrackedParameter("MinP", 0.)),
16 maxpcut(iConfig.getUntrackedParameter("MaxP", 10000.)),
17 minptcut(iConfig.getUntrackedParameter("MinPt", 0.)),
18 maxptcut(iConfig.getUntrackedParameter("MaxPt", 10000.)),
19 minetacut(iConfig.getUntrackedParameter("MinEta", -10.)),
20 maxetacut(iConfig.getUntrackedParameter("MaxEta", 10.)),
21 minrapcut(iConfig.getUntrackedParameter("MinRapidity", -20.)),
22 maxrapcut(iConfig.getUntrackedParameter("MaxRapidity", 20.)),
23 minphicut(iConfig.getUntrackedParameter("MinPhi", -3.5)),
24 maxphicut(iConfig.getUntrackedParameter("MaxPhi", 3.5)),
25 status(iConfig.getUntrackedParameter("Status", 0)),
26 motherID(iConfig.getUntrackedParameter("MotherID", 0)),
27 processID(iConfig.getUntrackedParameter("ProcessID", 0)),
28 betaBoost(iConfig.getUntrackedParameter("BetaBoost",0.))
29 {
30  //now do what ever initialization is needed
31 
32 }
33 
34 
36 {
37 
38  // do anything here that needs to be done at desctruction time
39  // (e.g. close files, deallocate resources etc.)
40 
41 }
42 
43 
44 //
45 // member functions
46 //
47 
48 // ------------ method called to produce the data ------------
50 {
51  using namespace edm;
52  bool accepted = false;
54  iEvent.getByToken(token_, evt);
55 
56  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
57 
58  if(processID == 0 || processID == myGenEvent->signal_process_id()) {
59 
60  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
61  p != myGenEvent->particles_end(); ++p ) {
62  HepMC::FourVector mom = MCFilterZboostHelper::zboost((*p)->momentum(),betaBoost);
63  double rapidity = 0.5*log( (mom.e()+mom.pz()) / (mom.e()-mom.pz()) );
64 
65  if ( abs((*p)->pdg_id()) == particleID
66  && mom.rho() > minpcut
67  && mom.rho() < maxpcut
68  && (*p)->momentum().perp() > minptcut
69  && (*p)->momentum().perp() < maxptcut
70  && mom.eta() > minetacut
71  && mom.eta() < maxetacut
72  && rapidity > minrapcut
73  && rapidity < maxrapcut
74  && (*p)->momentum().phi() > minphicut
75  && (*p)->momentum().phi() < maxphicut ) {
76 
77 
78 
79  if (status == 0 && motherID == 0){
80  accepted = true;
81  }
82  if (status != 0 && motherID == 0){
83  if ((*p)->status() == status)
84  accepted = true;
85  }
86 
87  HepMC::GenParticle* mother = (*((*p)->production_vertex()->particles_in_const_begin()));
88 
89  if (status == 0 && motherID != 0){
90  if (abs(mother->pdg_id()) == abs(motherID)) {
91  accepted = true;
92  }
93  }
94  if (status != 0 && motherID != 0){
95 
96  if ((*p)->status() == status && abs(mother->pdg_id()) == abs(motherID)){
97  accepted = true;
98 
99  }
100  }
101 
102  /*
103  if (status == 0 && motherID != 0){
104  if (abs(((*p)->mother())->pdg_id()) == abs(motherID)) {
105  accepted = true;
106  }
107  }
108  if (status != 0 && motherID != 0){
109 
110  if ((*p)->status() == status && abs(((*p)->mother())->pdg_id()) == abs(motherID)){
111  accepted = true;
112 
113  }
114  }
115  */
116 
117  }
118  }
119 
120  } else { accepted = true; }
121 
122  if (accepted){
123  return true; } else {return false;}
124 
125 }
HepMC::FourVector zboost(const HepMC::FourVector &, double)
~PythiaFilter() override
Definition: PythiaFilter.cc:35
const double maxetacut
Definition: PythiaFilter.h:58
PythiaFilter(const edm::ParameterSet &)
Definition: PythiaFilter.cc:12
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
const double minptcut
Definition: PythiaFilter.h:55
const int status
Definition: PythiaFilter.h:64
const double maxpcut
Definition: PythiaFilter.h:54
const double maxrapcut
Definition: PythiaFilter.h:60
const double minrapcut
Definition: PythiaFilter.h:59
const double maxptcut
Definition: PythiaFilter.h:56
int iEvent
Definition: GenABIO.cc:224
const double minpcut
Definition: PythiaFilter.h:53
const int motherID
Definition: PythiaFilter.h:65
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const int processID
Definition: PythiaFilter.h:66
const double maxphicut
Definition: PythiaFilter.h:62
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
bool accepted(std::vector< std::string_view > const &, std::string_view)
const double minphicut
Definition: PythiaFilter.h:61
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: PythiaFilter.cc:49
HLT enums.
const double betaBoost
Definition: PythiaFilter.h:68
const double minetacut
Definition: PythiaFilter.h:57
const edm::EDGetTokenT< edm::HepMCProduct > token_
Definition: PythiaFilter.h:51
const int particleID
Definition: PythiaFilter.h:52