85 virtual void endJob()
override;
133 genJetPtMin_(iConfig.getParameter<double>(
"genJetPtMin")),
134 genJetAbsEtaMax_(iConfig.getParameter<double>(
"genJetAbsEtaMax")),
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")))
149 produces<int>(
"genTtbarId");
212 std::map<int, int> bJetFromTopIds;
214 std::map<int, int> bJetFromWIds;
216 std::map<int, int> cJetFromWIds;
218 std::map<int, int> bJetAdditionalIds;
220 std::map<int, int> cJetAdditionalIds;
224 for(
size_t hadronId = 0; hadronId < genBHadIndex->size(); ++hadronId) {
226 const int jetIndex = genBHadJetIndex->at(hadronId);
228 if(jetIndex < 0)
continue;
233 const int flavour = genBHadFlavour->at(hadronId);
236 if(bJetFromTopIds.count(jetIndex) < 1) bJetFromTopIds[jetIndex] = 1;
237 else bJetFromTopIds[jetIndex]++;
242 if(bJetFromWIds.count(jetIndex) < 1) bJetFromWIds[jetIndex] = 1;
243 else bJetFromWIds[jetIndex]++;
247 if(bJetAdditionalIds.count(jetIndex) < 1) bJetAdditionalIds[jetIndex] = 1;
248 else bJetAdditionalIds[jetIndex]++;
252 for(std::map<int, int>::iterator it = bJetFromWIds.begin(); it != bJetFromWIds.end(); ) {
254 if(bJetFromTopIds.count(it->first) > 0) bJetFromWIds.erase(it++);
259 for(std::map<int, int>::iterator it = bJetAdditionalIds.begin(); it != bJetAdditionalIds.end(); ) {
261 if(bJetFromTopIds.count(it->first) > 0) bJetAdditionalIds.erase(it++);
263 else if(bJetFromWIds.count(it->first) > 0) bJetAdditionalIds.erase(it++);
268 for(
size_t hadronId = 0; hadronId < genCHadJetIndex->size(); ++hadronId) {
270 const int jetIndex = genCHadJetIndex->at(hadronId);
272 if(jetIndex < 0)
continue;
274 if(genCHadBHadronId->at(hadronId) >= 0)
continue;
279 if(bJetFromTopIds.count(jetIndex) > 0)
continue;
280 if(bJetFromWIds.count(jetIndex) > 0)
continue;
281 if(bJetAdditionalIds.count(jetIndex) > 0)
continue;
283 const int flavour = genCHadFlavour->at(hadronId);
286 if(cJetFromWIds.count(jetIndex) < 1) cJetFromWIds[jetIndex] = 1;
287 else cJetFromWIds[jetIndex]++;
291 if(cJetAdditionalIds.count(jetIndex) < 1) cJetAdditionalIds[jetIndex] = 1;
292 else cJetAdditionalIds[jetIndex]++;
296 for(std::map<int, int>::iterator it = cJetAdditionalIds.begin(); it != cJetAdditionalIds.end(); ) {
298 if(cJetFromWIds.count(it->first) > 0) cJetAdditionalIds.erase(it++);
304 int additionalJetEventId = bJetFromTopIds.size()*100 + bJetFromWIds.size()*1000 + cJetFromWIds.size()*10000;
306 if(bJetAdditionalIds.size() == 1){
307 const int nHadronsInJet = bJetAdditionalIds.begin()->second;
309 if(nHadronsInJet == 1) additionalJetEventId += 51;
311 else additionalJetEventId += 52;
314 else if(bJetAdditionalIds.size() > 1){
317 int nHadronsInJet1 = bJetAdditionalIds.at(orderedJetIndices.at(0));
318 int nHadronsInJet2 = bJetAdditionalIds.at(orderedJetIndices.at(1));
320 if(
std::max(nHadronsInJet1, nHadronsInJet2) == 1) additionalJetEventId += 53;
322 else if(
std::min(nHadronsInJet1, nHadronsInJet2) == 1 &&
std::max(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 54;
324 else if(
std::min(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 55;
329 if(cJetAdditionalIds.size() == 1){
330 const int nHadronsInJet = cJetAdditionalIds.begin()->second;
332 if(nHadronsInJet == 1) additionalJetEventId += 41;
334 else additionalJetEventId += 42;
337 else if(cJetAdditionalIds.size() > 1){
340 int nHadronsInJet1 = cJetAdditionalIds.at(orderedJetIndices.at(0));
341 int nHadronsInJet2 = cJetAdditionalIds.at(orderedJetIndices.at(1));
343 if(
std::max(nHadronsInJet1, nHadronsInJet2) == 1) additionalJetEventId += 43;
345 else if(
std::min(nHadronsInJet1, nHadronsInJet2) == 1 &&
std::max(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 44;
347 else if(
std::min(nHadronsInJet1, nHadronsInJet2) > 1) additionalJetEventId += 45;
352 additionalJetEventId += 0;
356 std::unique_ptr<int> ttbarId(
new int);
357 *ttbarId = additionalJetEventId;
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) );
415 std::sort(v_jetNhadIndexPair.begin(), v_jetNhadIndexPair.end(), std::greater<std::pair<int, int> >());
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);
423 return v_orderedJetIndices;
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
const edm::EDGetTokenT< std::vector< int > > genBHadFromTopWeakDecayToken_
const edm::EDGetTokenT< std::vector< int > > genBHadJetIndexToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
const edm::EDGetTokenT< std::vector< int > > genCHadFromTopWeakDecayToken_
#define DEFINE_FWK_MODULE(type)
std::vector< GenJet > GenJetCollection
collection of GenJet objects
const double genJetPtMin_
const edm::EDGetTokenT< std::vector< int > > genCHadFlavourToken_
const double genJetAbsEtaMax_
std::vector< int > nHadronsOrderedJetIndices(const std::map< int, int > &m_jetIndex)
void addDefault(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< std::vector< int > > genBHadFlavourToken_
Abs< T >::type abs(const T &t)
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_
virtual void beginJob() override
GenTtbarCategorizer(const edm::ParameterSet &)
virtual 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_
virtual void endJob() override
const edm::EDGetTokenT< reco::GenJetCollection > genJetsToken_
int flavour(const Candidate &part)