CMS 3D CMS Logo

GenTtbarCategorizer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: TopQuarkAnalysis/TopTools
4 // Class: GenTtbarCategorizer
5 //
48 //
49 // Original Author: Johannes Hauk, Nazar Bartosik
50 // Created: Sun, 14 Jun 2015 19:42:58 GMT
51 //
52 //
53 
54 
55 // system include files
56 #include <memory>
57 #include <algorithm>
58 #include <functional>
59 
60 // user include files
69 
70 
71 //
72 // class declaration
73 //
74 
76  public:
77  explicit GenTtbarCategorizer(const edm::ParameterSet&);
78  ~GenTtbarCategorizer() override;
79 
80  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
81 
82  private:
83  void beginJob() override;
84  void produce(edm::Event&, const edm::EventSetup&) override;
85  void endJob() override;
86 
87  //virtual void beginRun(edm::Run const&, edm::EventSetup const&) override;
88  //virtual void endRun(edm::Run const&, edm::EventSetup const&) override;
89  //virtual void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
90  //virtual void endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
91 
92  std::vector<int> nHadronsOrderedJetIndices(const std::map<int, int>& m_jetIndex);
93 
94  // ----------member data ---------------------------
95 
96  // Jet configuration
97  const double genJetPtMin_;
98  const double genJetAbsEtaMax_;
99 
100  // Input tags
102 
111 
116 
117 
118 };
119 
120 //
121 // constants, enums and typedefs
122 //
123 
124 
125 //
126 // static data member definitions
127 //
128 
129 //
130 // constructors and destructor
131 //
133 genJetPtMin_(iConfig.getParameter<double>("genJetPtMin")),
134 genJetAbsEtaMax_(iConfig.getParameter<double>("genJetAbsEtaMax")),
135 genJetsToken_(consumes<reco::GenJetCollection>(iConfig.getParameter<edm::InputTag>("genJets"))),
136 genBHadJetIndexToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genBHadJetIndex"))),
137 genBHadFlavourToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genBHadFlavour"))),
138 genBHadFromTopWeakDecayToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genBHadFromTopWeakDecay"))),
139 genBHadPlusMothersToken_(consumes<std::vector<reco::GenParticle> >(iConfig.getParameter<edm::InputTag>("genBHadPlusMothers"))),
140 genBHadPlusMothersIndicesToken_(consumes<std::vector<std::vector<int> > >(iConfig.getParameter<edm::InputTag>("genBHadPlusMothersIndices"))),
141 genBHadIndexToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genBHadIndex"))),
142 genBHadLeptonHadronIndexToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genBHadLeptonHadronIndex"))),
143 genBHadLeptonViaTauToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genBHadLeptonViaTau"))),
144 genCHadJetIndexToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genCHadJetIndex"))),
145 genCHadFlavourToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genCHadFlavour"))),
146 genCHadFromTopWeakDecayToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genCHadFromTopWeakDecay"))),
147 genCHadBHadronIdToken_(consumes<std::vector<int> >(iConfig.getParameter<edm::InputTag>("genCHadBHadronId")))
148 {
149  produces<int>("genTtbarId");
150 }
151 
152 
154 {}
155 
156 
157 //
158 // member functions
159 //
160 
161 // ------------ method called to produce the data ------------
162 void
164 {
165  // Access gen jets
167  iEvent.getByToken(genJetsToken_, genJets);
168 
169 
170  // Access B hadrons information
172  iEvent.getByToken(genBHadFlavourToken_, genBHadFlavour);
173 
175  iEvent.getByToken(genBHadJetIndexToken_, genBHadJetIndex);
176 
178  iEvent.getByToken(genBHadFromTopWeakDecayToken_, genBHadFromTopWeakDecay);
179 
181  iEvent.getByToken(genBHadPlusMothersToken_, genBHadPlusMothers);
182 
184  iEvent.getByToken(genBHadPlusMothersIndicesToken_, genBHadPlusMothersIndices);
185 
187  iEvent.getByToken(genBHadIndexToken_, genBHadIndex);
188 
190  iEvent.getByToken(genBHadLeptonHadronIndexToken_, genBHadLeptonHadronIndex);
191 
193  iEvent.getByToken(genBHadLeptonViaTauToken_, genBHadLeptonViaTau);
194 
195 
196  // Access C hadrons information
198  iEvent.getByToken(genCHadFlavourToken_, genCHadFlavour);
199 
201  iEvent.getByToken(genCHadJetIndexToken_, genCHadJetIndex);
202 
204  iEvent.getByToken(genCHadFromTopWeakDecayToken_, genCHadFromTopWeakDecay);
205 
207  iEvent.getByToken(genCHadBHadronIdToken_, genCHadBHadronId);
208 
209 
210  // Map <jet index, number of specific hadrons in jet>
211  // B jets with b hadrons directly from t->b decay
212  std::map<int, int> bJetFromTopIds;
213  // B jets with b hadrons from W->b decay
214  std::map<int, int> bJetFromWIds;
215  // C jets with c hadrons from W->c decay
216  std::map<int, int> cJetFromWIds;
217  // B jets with b hadrons before top quark decay chain
218  std::map<int, int> bJetAdditionalIds;
219  // C jets with c hadrons before top quark decay chain
220  std::map<int, int> cJetAdditionalIds;
221 
222 
223  // Count number of specific b hadrons in each jet
224  for(size_t hadronId = 0; hadronId < genBHadIndex->size(); ++hadronId) {
225  // Index of jet associated to the hadron
226  const int jetIndex = genBHadJetIndex->at(hadronId);
227  // Skip hadrons which have no associated jet
228  if(jetIndex < 0) continue;
229  // Skip if jet is not in acceptance
230  if(genJets->at(jetIndex).pt() < genJetPtMin_) continue;
231  if(std::fabs(genJets->at(jetIndex).eta()) > genJetAbsEtaMax_) continue;
232  // Flavour of the hadron's origin
233  const int flavour = genBHadFlavour->at(hadronId);
234  // Jet from t->b decay [pdgId(top)=6]
235  if(std::abs(flavour) == 6) {
236  if(bJetFromTopIds.count(jetIndex) < 1) bJetFromTopIds[jetIndex] = 1;
237  else bJetFromTopIds[jetIndex]++;
238  continue;
239  }
240  // Jet from W->b decay [pdgId(W)=24]
241  if(std::abs(flavour) == 24) {
242  if(bJetFromWIds.count(jetIndex) < 1) bJetFromWIds[jetIndex] = 1;
243  else bJetFromWIds[jetIndex]++;
244  continue;
245  }
246  // Identify jets with b hadrons not from top-quark or W-boson decay
247  if(bJetAdditionalIds.count(jetIndex) < 1) bJetAdditionalIds[jetIndex] = 1;
248  else bJetAdditionalIds[jetIndex]++;
249  }
250 
251  // Cleaning up b jets from W->b decays
252  for(std::map<int, int>::iterator it = bJetFromWIds.begin(); it != bJetFromWIds.end(); ) {
253  // Cannot be a b jet from t->b decay
254  if(bJetFromTopIds.count(it->first) > 0) bJetFromWIds.erase(it++);
255  else ++it;
256  }
257 
258  // Cleaning up additional b jets
259  for(std::map<int, int>::iterator it = bJetAdditionalIds.begin(); it != bJetAdditionalIds.end(); ) {
260  // Cannot be a b jet from t->b decay
261  if(bJetFromTopIds.count(it->first) > 0) bJetAdditionalIds.erase(it++);
262  // Cannot be a b jet from W->b decay
263  else if(bJetFromWIds.count(it->first) > 0) bJetAdditionalIds.erase(it++);
264  else ++it;
265  }
266 
267  // Count number of specific c hadrons in each c jet
268  for(size_t hadronId = 0; hadronId < genCHadJetIndex->size(); ++hadronId) {
269  // Index of jet associated to the hadron
270  const int jetIndex = genCHadJetIndex->at(hadronId);
271  // Skip hadrons which have no associated jet
272  if(jetIndex < 0) continue;
273  // Skip c hadrons that are coming from b hadrons
274  if(genCHadBHadronId->at(hadronId) >= 0) continue;
275  // Skip if jet is not in acceptance
276  if(genJets->at(jetIndex).pt() < genJetPtMin_) continue;
277  if(std::fabs(genJets->at(jetIndex).eta()) > genJetAbsEtaMax_) continue;
278  // Skip if jet contains a b hadron
279  if(bJetFromTopIds.count(jetIndex) > 0) continue;
280  if(bJetFromWIds.count(jetIndex) > 0) continue;
281  if(bJetAdditionalIds.count(jetIndex) > 0) continue;
282  // Flavour of the hadron's origin
283  const int flavour = genCHadFlavour->at(hadronId);
284  // Jet from W->c decay [pdgId(W)=24]
285  if(std::abs(flavour) == 24) {
286  if(cJetFromWIds.count(jetIndex) < 1) cJetFromWIds[jetIndex] = 1;
287  else cJetFromWIds[jetIndex]++;
288  continue;
289  }
290  // Identify jets with c hadrons not from W-boson decay
291  if(cJetAdditionalIds.count(jetIndex) < 1) cJetAdditionalIds[jetIndex] = 1;
292  else cJetAdditionalIds[jetIndex]++;
293  }
294 
295  // Cleaning up additional c jets
296  for(std::map<int, int>::iterator it = cJetAdditionalIds.begin(); it != cJetAdditionalIds.end(); ) {
297  // Cannot be a c jet from W->c decay
298  if(cJetFromWIds.count(it->first) > 0) cJetAdditionalIds.erase(it++);
299  else ++it;
300  }
301 
302  // Categorize event based on number of additional b/c jets
303  // and number of corresponding hadrons in each of them
304  int additionalJetEventId = bJetFromTopIds.size()*100 + bJetFromWIds.size()*1000 + cJetFromWIds.size()*10000;
305  // tt + 1 additional b jet
306  if(bJetAdditionalIds.size() == 1){
307  const int nHadronsInJet = bJetAdditionalIds.begin()->second;
308  // tt + 1 additional b jet from 1 additional b hadron
309  if(nHadronsInJet == 1) additionalJetEventId += 51;
310  // tt + 1 additional b jet from >=2 additional b hadrons
311  else additionalJetEventId += 52;
312  }
313  // tt + >=2 additional b jets
314  else if(bJetAdditionalIds.size() > 1){
315  // Check first two additional b jets (rare cases could have more)
316  const std::vector<int> orderedJetIndices = nHadronsOrderedJetIndices(bJetAdditionalIds);
317  int nHadronsInJet1 = bJetAdditionalIds.at(orderedJetIndices.at(0));
318  int nHadronsInJet2 = bJetAdditionalIds.at(orderedJetIndices.at(1));
319  // tt + >=2 additional b jets each from 1 additional b hadron
320  if(std::max(nHadronsInJet1, nHadronsInJet2) == 1) additionalJetEventId += 53;
321  // tt + >=2 additional b jets one of which from >=2 additional b hadrons
322  else if(std::min(nHadronsInJet1, nHadronsInJet2) == 1 && std::max(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 54;
323  // tt + >=2 additional b jets each from >=2 additional b hadrons
324  else if(std::min(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 55;
325  }
326  // tt + no additional b jets
327  else{
328  // tt + 1 additional c jet
329  if(cJetAdditionalIds.size() == 1){
330  const int nHadronsInJet = cJetAdditionalIds.begin()->second;
331  // tt + 1 additional c jet from 1 additional c hadron
332  if(nHadronsInJet == 1) additionalJetEventId += 41;
333  // tt + 1 additional c jet from >=2 additional c hadrons
334  else additionalJetEventId += 42;
335  }
336  // tt + >=2 additional c jets
337  else if(cJetAdditionalIds.size() > 1){
338  // Check first two additional c jets (rare cases could have more)
339  const std::vector<int> orderedJetIndices = nHadronsOrderedJetIndices(cJetAdditionalIds);
340  int nHadronsInJet1 = cJetAdditionalIds.at(orderedJetIndices.at(0));
341  int nHadronsInJet2 = cJetAdditionalIds.at(orderedJetIndices.at(1));
342  // tt + >=2 additional c jets each from 1 additional c hadron
343  if(std::max(nHadronsInJet1, nHadronsInJet2) == 1) additionalJetEventId += 43;
344  // tt + >=2 additional c jets one of which from >=2 additional c hadrons
345  else if(std::min(nHadronsInJet1, nHadronsInJet2) == 1 && std::max(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 44;
346  // tt + >=2 additional c jets each from >=2 additional c hadrons
347  else if(std::min(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 45;
348  }
349  // tt + no additional c jets
350  else{
351  // tt + light jets
352  additionalJetEventId += 0;
353  }
354  }
355 
356  std::unique_ptr<int> ttbarId(new int);
357  *ttbarId = additionalJetEventId;
358  iEvent.put(std::move(ttbarId), "genTtbarId");
359 }
360 
361 // ------------ method called once each job just before starting event loop ------------
362 void
364 {}
365 
366 // ------------ method called once each job just after ending the event loop ------------
367 void
369 {}
370 
371 // ------------ method called when starting to processes a run ------------
372 /*
373 void
374 GenTtbarCategorizer::beginRun(edm::Run const&, edm::EventSetup const&)
375 {
376 }
377 */
378 
379 // ------------ method called when ending the processing of a run ------------
380 /*
381 void
382 GenTtbarCategorizer::endRun(edm::Run const&, edm::EventSetup const&)
383 {
384 }
385 */
386 
387 // ------------ method called when starting to processes a luminosity block ------------
388 /*
389 void
390 GenTtbarCategorizer::beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
391 {
392 }
393 */
394 
395 // ------------ method called when ending the processing of a luminosity block ------------
396 /*
397 void
398 GenTtbarCategorizer::endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
399 {
400 }
401 */
402 
403 // ------------ method returns a vector of jet indices from the given map, sorted by N hadrons in descending order ------------
404 std::vector<int> GenTtbarCategorizer::nHadronsOrderedJetIndices(const std::map<int, int>& m_jetIndex)
405 {
406  const int nElements = m_jetIndex.size();
407  std::vector<std::pair<int, int> > v_jetNhadIndexPair;
408  v_jetNhadIndexPair.reserve(nElements);
409  for(std::map<int, int>::const_iterator it = m_jetIndex.begin(); it != m_jetIndex.end(); ++it) {
410  const int jetIndex = it->first;
411  const int nHadrons = it->second;
412  v_jetNhadIndexPair.push_back( std::pair<int, int>(nHadrons, jetIndex) );
413  }
414  // Sorting the vector of pairs by their key value
415  std::sort(v_jetNhadIndexPair.begin(), v_jetNhadIndexPair.end(), std::greater<std::pair<int, int> >());
416  // Building the vector of indices in the proper order
417  std::vector<int> v_orderedJetIndices;
418  v_orderedJetIndices.reserve(nElements);
419  for(std::vector<std::pair<int, int> >::const_iterator it = v_jetNhadIndexPair.begin(); it != v_jetNhadIndexPair.end(); ++it) {
420  v_orderedJetIndices.push_back(it->second);
421  }
422 
423  return v_orderedJetIndices;
424 }
425 
426 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
427 void
429  //The following says we do not know what parameters are allowed so do no validation
430  // Please change this to state exactly what you do use, even if it is no parameters
432  desc.setUnknown();
433  descriptions.addDefault(desc);
434 }
435 
436 //define this as a plug-in
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
const edm::EDGetTokenT< std::vector< int > > genBHadFromTopWeakDecayToken_
const edm::EDGetTokenT< std::vector< int > > genBHadJetIndexToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
const edm::EDGetTokenT< std::vector< int > > genCHadFromTopWeakDecayToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< GenJet > GenJetCollection
collection of GenJet objects
const edm::EDGetTokenT< std::vector< int > > genCHadFlavourToken_
std::vector< int > nHadronsOrderedJetIndices(const std::map< int, int > &m_jetIndex)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< std::vector< int > > genBHadFlavourToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
const edm::EDGetTokenT< std::vector< std::vector< int > > > genBHadPlusMothersIndicesToken_
const edm::EDGetTokenT< std::vector< int > > genBHadLeptonHadronIndexToken_
const edm::EDGetTokenT< std::vector< int > > genBHadIndexToken_
const edm::EDGetTokenT< std::vector< int > > genBHadLeptonViaTauToken_
fixed size matrix
HLT enums.
const edm::EDGetTokenT< std::vector< reco::GenParticle > > genBHadPlusMothersToken_
GenTtbarCategorizer(const edm::ParameterSet &)
void produce(edm::Event &, const edm::EventSetup &) override
const edm::EDGetTokenT< std::vector< int > > genCHadJetIndexToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< std::vector< int > > genCHadBHadronIdToken_
const edm::EDGetTokenT< reco::GenJetCollection > genJetsToken_
def move(src, dest)
Definition: eostools.py:510