CMS 3D CMS Logo

HighMultiplicityGenFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HighMultiplicityGenFilter
4 // Class: HighMultiplicityGenFilter
5 //
13 //
14 // Original Author: Wei Li
15 // Created: Tue Dec 8 23:51:37 EST 2009
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
33 //
34 // class declaration
35 //
36 
38  public:
40  ~HighMultiplicityGenFilter() override;
41 
42  private:
43  void beginJob() override;
44  bool filter(edm::Event&, const edm::EventSetup&) override;
45  void endJob() override;
46 
47  // ----------member data ---------------------------
50  double etaMax;
51  double ptMin;
52  int nMin;
53  int nAccepted;
54 };
55 
56 //
57 // constants, enums and typedefs
58 //
59 
60 //
61 // static data member definitions
62 //
63 
64 //
65 // constructors and destructor
66 //
68  hepmcSrc(consumes<edm::HepMCProduct>(iConfig.getParameter< edm::InputTag > ("generatorSmeared"))),
69  etaMax(iConfig.getUntrackedParameter<double>("etaMax")),
70  ptMin(iConfig.getUntrackedParameter<double>("ptMin")),
71  nMin(iConfig.getUntrackedParameter<int>("nMin"))
72 {
73  //now do what ever initialization is needed
74  nAccepted = 0;
75 }
76 
77 
79 {
80 
81  // do anything here that needs to be done at desctruction time
82  // (e.g. close files, deallocate resources etc.)
83 
84 }
85 
86 
87 //
88 // member functions
89 //
90 
91 // ------------ method called on each new Event ------------
92 bool
94 {
95 
96  bool accepted = false;
98  iEvent.getByToken(hepmcSrc,evt);
99 
100  iSetup.getData(pdt);
101 
102  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
103 
104  int nMult=0;
105  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
106 
107  if((*p)->status()!=1) continue;
108 
109  double charge = 0;
110  int pid = (*p)->pdg_id();
111  if(abs(pid) > 100000) { std::cout<<"pid="<<pid<<" status="<<(*p)->status()<<std::endl; continue; }
112  const ParticleData* part = pdt->particle(pid);
113  if(part) charge = part->charge();
114  if(charge == 0) continue;
115 
116  if (
117  (*p)->momentum().perp() > ptMin
118  && fabs((*p)->momentum().eta()) < etaMax ) nMult++;
119  }
120  if(nMult>=nMin) { nAccepted++; accepted = true; }
121  return accepted;
122 }
123 
124 // ------------ method called once each job just before starting event loop ------------
125 void
127 {}
128 
129 // ------------ method called once each job just after ending the event loop ------------
130 void
132  std::cout<<"There are "<<nAccepted<<" events with multiplicity greater than "<<nMin<<std::endl;
133 }
134 
135 //define this as a plug-in
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool getData(T &iHolder) const
Definition: EventSetup.h:111
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::ESHandle< ParticleDataTable > pdt
HepPDT::ParticleData ParticleData
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
bool filter(edm::Event &, const edm::EventSetup &) override
part
Definition: HCALResponse.h:20
bool accepted(std::vector< std::string_view > const &, std::string_view)
HLT enums.
edm::EDGetTokenT< edm::HepMCProduct > hepmcSrc
HighMultiplicityGenFilter(const edm::ParameterSet &)