CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
GenTtbarCategorizer Class Reference

#include <TopQuarkAnalysis/TopTools/plugins/GenTtbarCategorizer.cc>

Inheritance diagram for GenTtbarCategorizer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 GenTtbarCategorizer (const edm::ParameterSet &)
 
 ~GenTtbarCategorizer () override
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducer () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
virtual ~ProducerBase () noexcept(false)
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

void beginJob () override
 
void endJob () override
 
std::vector< int > nHadronsOrderedJetIndices (const std::map< int, int > &m_jetIndex)
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

const edm::EDGetTokenT< std::vector< int > > genBHadFlavourToken_
 
const edm::EDGetTokenT< std::vector< int > > genBHadFromTopWeakDecayToken_
 
const edm::EDGetTokenT< std::vector< int > > genBHadIndexToken_
 
const edm::EDGetTokenT< std::vector< int > > genBHadJetIndexToken_
 
const edm::EDGetTokenT< std::vector< int > > genBHadLeptonHadronIndexToken_
 
const edm::EDGetTokenT< std::vector< int > > genBHadLeptonViaTauToken_
 
const edm::EDGetTokenT< std::vector< std::vector< int > > > genBHadPlusMothersIndicesToken_
 
const edm::EDGetTokenT< std::vector< reco::GenParticle > > genBHadPlusMothersToken_
 
const edm::EDGetTokenT< std::vector< int > > genCHadBHadronIdToken_
 
const edm::EDGetTokenT< std::vector< int > > genCHadFlavourToken_
 
const edm::EDGetTokenT< std::vector< int > > genCHadFromTopWeakDecayToken_
 
const edm::EDGetTokenT< std::vector< int > > genCHadJetIndexToken_
 
const double genJetAbsEtaMax_
 
const double genJetPtMin_
 
const edm::EDGetTokenT< reco::GenJetCollectiongenJetsToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: Categorization of different tt+xx processes, returning unique ID for each process as e.g. tt+bb, tt+b, tt+2b, tt+cc, ...

Implementation:

The classification scheme returns an ID per event, and works as follows:

All jets in the following need to be in the acceptance as given by the config parameters |eta|, pt.
A c jet must contain at least one c hadron and should contain no b hadrons

First, b jets from top are identified, i.e. jets containing a b hadron from t->b decay
They are encoded in the ID as numberOfBjetsFromTop*100, i.e.
0xx: no b jets from top in acceptance
1xx: 1 b jet from top in acceptance
2xx: both b jets from top in acceptance

Then, b jets from W are identified, i.e. jets containing a b hadron from W->b decay
They are encoded in the ID as numberOfBjetsFromW*1000, i.e.
0xxx: no b jets from W in acceptance
1xxx: 1 b jet from W in acceptance
2xxx: 2 b jets from W in acceptance

Then, c jets from W are identified, i.e. jets containing a c hadron from W->c decay, but no b hadrons
They are encoded in the ID as numberOfCjetsFromW*10000, i.e.
0xxxx: no c jets from W in acceptance
1xxxx: 1 c jet from W in acceptance
2xxxx: 2 c jets from W in acceptance

From the remaining jets, the ID is formed based on the additional b jets (IDs x5x) and c jets (IDs x4x) in the following order:
x55: at least 2 additional b jets with at least two of them having >= 2 b hadrons in each
x54: at least 2 additional b jets with one of them having >= 2 b hadrons, the others having =1 b hadron
x53: at least 2 additional b jets with all having =1 b hadron
x52: exactly 1 additional b jet having >=2 b hadrons
x51: exactly 1 additional b jet having =1 b hadron
x45: at least 2 additional c jets with at least two of them having >= 2 c hadrons in each
x44: at least 2 additional c jets with one of them having >= 2 c hadrons, the others having =1 c hadron
x43: at least 2 additional c jets with all having =1 c hadron
x42: exactly 1 additional c jet having >=2 c hadrons
x41: exactly 1 additional c jet having =1 c hadron
x00: No additional b or c jet, i.e. only light flavour jets or no additional jets

Definition at line 75 of file GenTtbarCategorizer.cc.

Constructor & Destructor Documentation

GenTtbarCategorizer::GenTtbarCategorizer ( const edm::ParameterSet iConfig)
explicit

Definition at line 132 of file GenTtbarCategorizer.cc.

132  :
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 }
T getParameter(std::string const &) const
const edm::EDGetTokenT< std::vector< int > > genBHadFromTopWeakDecayToken_
const edm::EDGetTokenT< std::vector< int > > genBHadJetIndexToken_
const edm::EDGetTokenT< std::vector< int > > genCHadFromTopWeakDecayToken_
const edm::EDGetTokenT< std::vector< int > > genCHadFlavourToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EDGetTokenT< std::vector< int > > genBHadFlavourToken_
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_
const edm::EDGetTokenT< std::vector< reco::GenParticle > > genBHadPlusMothersToken_
const edm::EDGetTokenT< std::vector< int > > genCHadJetIndexToken_
const edm::EDGetTokenT< std::vector< int > > genCHadBHadronIdToken_
const edm::EDGetTokenT< reco::GenJetCollection > genJetsToken_
GenTtbarCategorizer::~GenTtbarCategorizer ( )
override

Definition at line 153 of file GenTtbarCategorizer.cc.

154 {}

Member Function Documentation

void GenTtbarCategorizer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 363 of file GenTtbarCategorizer.cc.

364 {}
void GenTtbarCategorizer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 368 of file GenTtbarCategorizer.cc.

369 {}
void GenTtbarCategorizer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 428 of file GenTtbarCategorizer.cc.

References edm::ConfigurationDescriptions::addDefault(), DEFINE_FWK_MODULE, and edm::ParameterSetDescription::setUnknown().

428  {
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 }
void addDefault(ParameterSetDescription const &psetDescription)
std::vector< int > GenTtbarCategorizer::nHadronsOrderedJetIndices ( const std::map< int, int > &  m_jetIndex)
private

Definition at line 404 of file GenTtbarCategorizer.cc.

Referenced by produce().

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 }
void GenTtbarCategorizer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 163 of file GenTtbarCategorizer.cc.

References funct::abs(), GenHFHadronMatcher_cff::flavour, ttHFGenFilter_cfi::genBHadFlavour, genBHadFlavourToken_, ttHFGenFilter_cfi::genBHadFromTopWeakDecay, genBHadFromTopWeakDecayToken_, ttHFGenFilter_cfi::genBHadIndex, genBHadIndexToken_, GenTtbarCategorizer_cfi::genBHadJetIndex, genBHadJetIndexToken_, GenTtbarCategorizer_cfi::genBHadLeptonHadronIndex, genBHadLeptonHadronIndexToken_, GenTtbarCategorizer_cfi::genBHadLeptonViaTau, genBHadLeptonViaTauToken_, ttHFGenFilter_cfi::genBHadPlusMothers, ttHFGenFilter_cfi::genBHadPlusMothersIndices, genBHadPlusMothersIndicesToken_, genBHadPlusMothersToken_, GenTtbarCategorizer_cfi::genCHadBHadronId, genCHadBHadronIdToken_, GenTtbarCategorizer_cfi::genCHadFlavour, genCHadFlavourToken_, GenTtbarCategorizer_cfi::genCHadFromTopWeakDecay, genCHadFromTopWeakDecayToken_, GenTtbarCategorizer_cfi::genCHadJetIndex, genCHadJetIndexToken_, genJetAbsEtaMax_, genJetPtMin_, slimmedGenJetsFlavourInfos_cfi::genJets, genJetsToken_, edm::Event::getByToken(), SiStripPI::max, min(), eostools::move(), nHadronsOrderedJetIndices(), and edm::Event::put().

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 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
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:508
const edm::EDGetTokenT< std::vector< int > > genCHadFromTopWeakDecayToken_
const edm::EDGetTokenT< std::vector< int > > genCHadFlavourToken_
std::vector< int > nHadronsOrderedJetIndices(const std::map< int, int > &m_jetIndex)
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_
const edm::EDGetTokenT< std::vector< reco::GenParticle > > genBHadPlusMothersToken_
const edm::EDGetTokenT< std::vector< int > > genCHadJetIndexToken_
const edm::EDGetTokenT< std::vector< int > > genCHadBHadronIdToken_
const edm::EDGetTokenT< reco::GenJetCollection > genJetsToken_
def move(src, dest)
Definition: eostools.py:510

Member Data Documentation

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genBHadFlavourToken_
private

Definition at line 104 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genBHadFromTopWeakDecayToken_
private

Definition at line 105 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genBHadIndexToken_
private

Definition at line 108 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genBHadJetIndexToken_
private

Definition at line 103 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genBHadLeptonHadronIndexToken_
private

Definition at line 109 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genBHadLeptonViaTauToken_
private

Definition at line 110 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<std::vector<int> > > GenTtbarCategorizer::genBHadPlusMothersIndicesToken_
private

Definition at line 107 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<reco::GenParticle> > GenTtbarCategorizer::genBHadPlusMothersToken_
private

Definition at line 106 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genCHadBHadronIdToken_
private

Definition at line 115 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genCHadFlavourToken_
private

Definition at line 113 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genCHadFromTopWeakDecayToken_
private

Definition at line 114 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<std::vector<int> > GenTtbarCategorizer::genCHadJetIndexToken_
private

Definition at line 112 of file GenTtbarCategorizer.cc.

Referenced by produce().

const double GenTtbarCategorizer::genJetAbsEtaMax_
private

Definition at line 98 of file GenTtbarCategorizer.cc.

Referenced by produce().

const double GenTtbarCategorizer::genJetPtMin_
private

Definition at line 97 of file GenTtbarCategorizer.cc.

Referenced by produce().

const edm::EDGetTokenT<reco::GenJetCollection> GenTtbarCategorizer::genJetsToken_
private

Definition at line 101 of file GenTtbarCategorizer.cc.

Referenced by produce().