CMS 3D CMS Logo

MCSingleParticleYPt.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MCSingleParticleYPt
4 // Class: MCSingleParticleYPt
5 //
6 /*
7 
8  Description: filter events based on the Pythia particleID, Pt, Y and status. It is based on MCSingleParticleFilter.
9  It will used to filter a b-hadron with the given kinematics, only one b-hadron is required to match.
10  Implementation: inherits from generic EDFilter
11 
12 */
13 //
14 // Author: Alberto Sanchez-Hernandez
15 // Adapted on: August 2016
16 //
17 //
18 // Filter based on MCSingleParticleFilter.cc, but using rapidity instead of eta
19 
30 
31 #include <cmath>
32 #include <ostream>
33 #include <string>
34 #include <vector>
35 
37 public:
38  explicit MCSingleParticleYPt(const edm::ParameterSet&);
39 
40  bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
41 
42 private:
44  const int fVerbose;
45  const bool fchekantiparticle;
46  std::vector<int> particleID;
47  std::vector<double> ptMin;
48  std::vector<double> rapMin;
49  std::vector<double> rapMax;
50  std::vector<int> status;
51 };
52 
53 using namespace edm;
54 using namespace std;
55 
57  : token_(consumes<edm::HepMCProduct>(
58  edm::InputTag(iConfig.getUntrackedParameter("moduleLabel", std::string("generator")), "unsmeared"))),
59  fVerbose(iConfig.getUntrackedParameter("verbose", 0)),
60  fchekantiparticle(iConfig.getUntrackedParameter("CheckAntiparticle", true)) {
61  vector<int> defpid;
62  defpid.push_back(0);
63  particleID = iConfig.getUntrackedParameter<vector<int> >("ParticleID", defpid);
64  vector<double> defptmin;
65  defptmin.push_back(0.);
66  ptMin = iConfig.getUntrackedParameter<vector<double> >("MinPt", defptmin);
67  vector<double> defrapmin;
68  defrapmin.push_back(-10.);
69  rapMin = iConfig.getUntrackedParameter<vector<double> >("MinY", defrapmin);
70  vector<double> defrapmax;
71  defrapmax.push_back(10.);
72  rapMax = iConfig.getUntrackedParameter<vector<double> >("MaxY", defrapmax);
73  vector<int> defstat;
74  defstat.push_back(0);
75  status = iConfig.getUntrackedParameter<vector<int> >("Status", defstat);
76 
77  // check for same size
78  if ((ptMin.size() > 1 && particleID.size() != ptMin.size()) ||
79  (rapMin.size() > 1 && particleID.size() != rapMin.size()) ||
80  (rapMax.size() > 1 && particleID.size() != rapMax.size()) ||
81  (status.size() > 1 && particleID.size() != status.size())) {
82  edm::LogWarning("MCSingleParticleYPt")
83  << "WARNING: MCSingleParticleYPt : size of vector cuts do not match!!" << endl;
84  }
85 
86  // if ptMin size smaller than particleID , fill up further with defaults
87  if (particleID.size() > ptMin.size()) {
88  vector<double> defptmin2;
89  for (unsigned int i = 0; i < particleID.size(); i++) {
90  defptmin2.push_back(0.);
91  }
92  ptMin = defptmin2;
93  }
94  // if etaMin size smaller than particleID , fill up further with defaults
95  if (particleID.size() > rapMin.size()) {
96  vector<double> defrapmin2;
97  for (unsigned int i = 0; i < particleID.size(); i++) {
98  defrapmin2.push_back(-10.);
99  }
100  rapMin = defrapmin2;
101  }
102  // if etaMax size smaller than particleID , fill up further with defaults
103  if (particleID.size() > rapMax.size()) {
104  vector<double> defrapmax2;
105  for (unsigned int i = 0; i < particleID.size(); i++) {
106  defrapmax2.push_back(10.);
107  }
108  rapMax = defrapmax2;
109  }
110  // if status size smaller than particleID , fill up further with defaults
111  if (particleID.size() > status.size()) {
112  vector<int> defstat2;
113  for (unsigned int i = 0; i < particleID.size(); i++) {
114  defstat2.push_back(0);
115  }
116  status = defstat2;
117  }
118 
119  if (fVerbose > 0) {
120  edm::LogInfo("MCSingleParticleYPt") << "----------------------------------------------------------------------"
121  << std::endl;
122  edm::LogInfo("MCSingleParticleYPt") << "----- MCSingleParticleYPt" << std::endl;
123  for (unsigned int i = 0; i < particleID.size(); ++i) {
124  edm::LogInfo("MCSingleParticleYPt") << " ID: " << particleID[i] << " pT > " << ptMin[i] << ", " << rapMin[i]
125  << " < y < " << rapMax[i] << ", status = " << status[i] << std::endl;
126  }
127  if (fchekantiparticle)
128  edm::LogInfo("MCSingleParticleYPt") << " anti-particles will be tested as well." << std::endl;
129  edm::LogInfo("MCSingleParticleYPt") << "----------------------------------------------------------------------"
130  << std::endl;
131  }
132 }
133 
135  bool accepted = false;
137  iEvent.getByToken(token_, evt);
138 
139  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
140  for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
141  ++p) {
142  if (fVerbose > 3)
143  edm::LogInfo("MCSingleParticleYPt")
144  << "Looking at particle : " << (*p)->pdg_id() << " status : " << (*p)->status() << std::endl;
145 
146  for (unsigned int i = 0; i < particleID.size(); i++) {
147  if (particleID[i] == (*p)->pdg_id() || (fchekantiparticle && (-particleID[i] == (*p)->pdg_id())) ||
148  particleID[i] == 0) {
149  // calculate rapidity just for the desired particle and make sure, this particles has enough energy
150  double rapidity = ((*p)->momentum().e() - (*p)->momentum().pz()) > 0.
151  ? 0.5 * log(((*p)->momentum().e() + (*p)->momentum().pz()) /
152  ((*p)->momentum().e() - (*p)->momentum().pz()))
153  : rapMax[i] + .1;
154  if (fVerbose > 2)
155  edm::LogInfo("MCSingleParticleYPt")
156  << "Testing particle : " << (*p)->pdg_id() << " pT: " << (*p)->momentum().perp() << " y: " << rapidity
157  << " status : " << (*p)->status() << endl;
158  if ((*p)->momentum().perp() > ptMin[i] && rapidity > rapMin[i] && rapidity < rapMax[i] &&
159  ((*p)->status() == status[i] || status[i] == 0)) {
160  accepted = true;
161  if (fVerbose > 1)
162  edm::LogInfo("MCSingleParticleYPt")
163  << "Accepted particle : " << (*p)->pdg_id() << " pT: " << (*p)->momentum().perp() << " y: " << rapidity
164  << " status : " << (*p)->status() << endl;
165  break;
166  }
167  }
168  }
169  if (accepted)
170  break;
171  }
172  return accepted;
173 }
174 
edm::StreamID
Definition: StreamID.h:30
MCSingleParticleYPt::fchekantiparticle
const bool fchekantiparticle
Definition: MCSingleParticleYPt.cc:45
Handle.h
mps_fire.i
i
Definition: mps_fire.py:428
MessageLogger.h
edm::EDGetTokenT< edm::HepMCProduct >
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
MCSingleParticleYPt::token_
const edm::EDGetTokenT< edm::HepMCProduct > token_
Definition: MCSingleParticleYPt.cc:43
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
edm::Handle
Definition: AssociativeIterator.h:50
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
MCSingleParticleYPt::particleID
std::vector< int > particleID
Definition: MCSingleParticleYPt.cc:46
MCSingleParticleYPt::filter
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: MCSingleParticleYPt.cc:134
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
MCSingleParticleYPt::fVerbose
const int fVerbose
Definition: MCSingleParticleYPt.cc:44
MakerMacros.h
MCSingleParticleYPt
Definition: MCSingleParticleYPt.cc:36
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
EDGetToken.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
MCSingleParticleYPt::rapMax
std::vector< double > rapMax
Definition: MCSingleParticleYPt.cc:49
funct::true
true
Definition: Factorize.h:173
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup
Definition: EventSetup.h:57
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
InputTag.h
MCSingleParticleYPt::MCSingleParticleYPt
MCSingleParticleYPt(const edm::ParameterSet &)
Definition: MCSingleParticleYPt.cc:56
std
Definition: JetResolutionObject.h:76
Frameworkfwd.h
MCSingleParticleYPt::status
std::vector< int > status
Definition: MCSingleParticleYPt.cc:50
MCSingleParticleYPt::rapMin
std::vector< double > rapMin
Definition: MCSingleParticleYPt.cc:48
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
edm::global::EDFilter
Definition: EDFilter.h:32
edm::HepMCProduct
Definition: HepMCProduct.h:18
ParameterSet.h
HepMCProduct.h
edm::Event
Definition: Event.h:73
EDFilter.h
edm::InputTag
Definition: InputTag.h:15
MCSingleParticleYPt::ptMin
std::vector< double > ptMin
Definition: MCSingleParticleYPt.cc:47