CMS 3D CMS Logo

MCSingleParticleFilter.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 betaBoost(iConfig.getUntrackedParameter("BetaBoost",0.))
15 {
16  //here do whatever other initialization is needed
17  vector<int> defpid ;
18  defpid.push_back(0) ;
19  particleID = iConfig.getUntrackedParameter< vector<int> >("ParticleID",defpid);
20  vector<double> defptmin ;
21  defptmin.push_back(0.);
22  ptMin = iConfig.getUntrackedParameter< vector<double> >("MinPt", defptmin);
23 
24  vector<double> defetamin ;
25  defetamin.push_back(-10.);
26  etaMin = iConfig.getUntrackedParameter< vector<double> >("MinEta", defetamin);
27  vector<double> defetamax ;
28  defetamax.push_back(10.);
29  etaMax = iConfig.getUntrackedParameter< vector<double> >("MaxEta", defetamax);
30  vector<int> defstat ;
31  defstat.push_back(0);
32  status = iConfig.getUntrackedParameter< vector<int> >("Status", defstat);
33 
34 
35  // check for same size
36  if ( (ptMin.size() > 1 && particleID.size() != ptMin.size())
37  || (etaMin.size() > 1 && particleID.size() != etaMin.size())
38  || (etaMax.size() > 1 && particleID.size() != etaMax.size())
39  || (status.size() > 1 && particleID.size() != status.size()) ) {
40  edm::LogInfo("MCSingleParticleFilter") << "WARNING: size of MinPthat and/or MaxPthat not matching with ProcessID size!!" << endl;
41  }
42 
43  // if ptMin size smaller than particleID , fill up further with defaults
44  if (particleID.size() > ptMin.size() ){
45  vector<double> defptmin2 ;
46  for (unsigned int i = 0; i < particleID.size(); i++){ defptmin2.push_back(0.);}
47  ptMin = defptmin2;
48  }
49  // if etaMin size smaller than particleID , fill up further with defaults
50  if (particleID.size() > etaMin.size() ){
51  vector<double> defetamin2 ;
52  for (unsigned int i = 0; i < particleID.size(); i++){ defetamin2.push_back(-10.);}
53  etaMin = defetamin2;
54  }
55  // if etaMax size smaller than particleID , fill up further with defaults
56  if (particleID.size() > etaMax.size() ){
57  vector<double> defetamax2 ;
58  for (unsigned int i = 0; i < particleID.size(); i++){ defetamax2.push_back(10.);}
59  etaMax = defetamax2;
60  }
61  // if status size smaller than particleID , fill up further with defaults
62  if (particleID.size() > status.size() ){
63  vector<int> defstat2 ;
64  for (unsigned int i = 0; i < particleID.size(); i++){ defstat2.push_back(0);}
65  status = defstat2;
66  }
67 
68  // check if beta is smaller than 1
69  if (std::abs(betaBoost) >= 1 ){
70  edm::LogError("MCSingleParticleFilter") << "Input beta boost is >= 1 !";
71  }
72 
73 }
74 
75 
77 {
78 
79  // do anything here that needs to be done at desctruction time
80  // (e.g. close files, deallocate resources etc.)
81 
82 }
83 
84 
85 // ------------ method called to skim the data ------------
87 {
88  using namespace edm;
89  bool accepted = false;
91  iEvent.getByToken(token_, evt);
92 
93  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
94 
95 
96  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
97  p != myGenEvent->particles_end(); ++p ) {
98 
99 
100  for (unsigned int i = 0; i < particleID.size(); i++){
101  if (particleID[i] == (*p)->pdg_id() || particleID[i] == 0) {
102 
103  if ( (*p)->momentum().perp() > ptMin[i]
104  && ((*p)->status() == status[i] || status[i] == 0)) {
105 
106  HepMC::FourVector mom = MCFilterZboostHelper::zboost((*p)->momentum(),betaBoost);
107  if ( mom.eta() > etaMin[i] && mom.eta() < etaMax[i] ) {
108  accepted = true;
109  }
110 
111  }
112 
113  }
114  }
115 
116 
117  }
118 
119  return accepted;
120 
121 }
HepMC::FourVector zboost(const HepMC::FourVector &, double)
T getUntrackedParameter(std::string const &, T const &) const
std::vector< double > ptMin
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const edm::EDGetTokenT< edm::HepMCProduct > token_
std::vector< double > etaMin
int iEvent
Definition: GenABIO.cc:224
std::vector< double > etaMax
std::vector< int > particleID
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
bool accepted(std::vector< std::string_view > const &, std::string_view)
HLT enums.
MCSingleParticleFilter(const edm::ParameterSet &)