CMS 3D CMS Logo

JetFlavourClustering.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: JetFlavourClustering
4 // Class: JetFlavourClustering
5 //
72 //
73 // Original Author: Dinko Ferencek
74 // Created: Wed Nov 6 00:49:55 CET 2013
75 //
76 //
77 
78 
79 // system include files
80 #include <memory>
81 
82 // user include files
85 
88 
90 
99 
100 #include "fastjet/JetDefinition.hh"
101 #include "fastjet/ClusterSequence.hh"
102 #include "fastjet/Selector.hh"
103 #include "fastjet/PseudoJet.hh"
104 
105 //
106 // constants, enums and typedefs
107 //
108 typedef std::shared_ptr<fastjet::ClusterSequence> ClusterSequencePtr;
109 typedef std::shared_ptr<fastjet::JetDefinition> JetDefPtr;
110 
111 //
112 // class declaration
113 //
114 class GhostInfo : public fastjet::PseudoJet::UserInfoBase{
115  public:
116  GhostInfo(const bool & isHadron,
117  const bool & isbHadron,
118  const bool & isParton,
119  const bool & isLepton,
121  m_particleRef(particleRef)
122  {
123  m_type = 0;
124  if( isHadron ) m_type |= ( 1 << 0 );
125  if( isbHadron ) m_type |= ( 1 << 1 );
126  if( isParton ) m_type |= ( 1 << 2 );
127  if( isLepton ) m_type |= ( 1 << 3 );
128  }
129 
130  const bool isHadron() const { return ( m_type & (1 << 0) ); }
131  const bool isbHadron() const { return ( m_type & (1 << 1) ); }
132  const bool isParton() const { return ( m_type & (1 << 2) ); }
133  const bool isLepton() const { return ( m_type & (1 << 3) ); }
134  const reco::GenParticleRef & particleRef() const { return m_particleRef; }
135 
136  protected:
138  int m_type;
139 };
140 
142  public:
143  explicit JetFlavourClustering(const edm::ParameterSet&);
144  ~JetFlavourClustering() override;
145 
146  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
147 
148  private:
149  void produce(edm::Event&, const edm::EventSetup&) override;
150 
151  void insertGhosts(const edm::Handle<reco::GenParticleRefVector>& particles,
152  const double ghostRescaling,
153  const bool isHadron, const bool isbHadron, const bool isParton, const bool isLepton,
154  std::vector<fastjet::PseudoJet>& constituents);
155 
156  void matchReclusteredJets(const edm::Handle<edm::View<reco::Jet> >& jets,
157  const std::vector<fastjet::PseudoJet>& matchedJets,
158  std::vector<int>& matchedIndices);
159  void matchGroomedJets(const edm::Handle<edm::View<reco::Jet> >& jets,
160  const edm::Handle<edm::View<reco::Jet> >& matchedJets,
161  std::vector<int>& matchedIndices);
162  void matchSubjets(const std::vector<int>& groomedIndices,
163  const edm::Handle<edm::View<reco::Jet> >& groomedJets,
164  const edm::Handle<edm::View<reco::Jet> >& subjets,
165  std::vector<std::vector<int> >& matchedIndices);
166 
167  void setFlavours(const reco::GenParticleRefVector& clusteredbHadrons,
168  const reco::GenParticleRefVector& clusteredcHadrons,
169  const reco::GenParticleRefVector& clusteredPartons,
170  int& hadronFlavour,
171  int& partonFlavour);
172 
173  void assignToSubjets(const reco::GenParticleRefVector& clusteredParticles,
174  const edm::Handle<edm::View<reco::Jet> >& subjets,
175  const std::vector<int>& subjetIndices,
176  std::vector<reco::GenParticleRefVector>& assignedParticles);
177 
178  // ----------member data ---------------------------
179  const edm::EDGetTokenT<edm::View<reco::Jet> > jetsToken_; // Input jet collection
180  edm::EDGetTokenT<edm::View<reco::Jet> > groomedJetsToken_; // Input groomed jet collection
186 
188  const double rParam_;
189  const double jetPtMin_;
190  const double ghostRescaling_;
191  const double relPtTolerance_;
193  const bool useSubjets_;
194  const bool useLeptons_;
195 
198 };
199 
200 //
201 // static data member definitions
202 //
203 
204 //
205 // constructors and destructor
206 //
208 
209  jetsToken_(consumes<edm::View<reco::Jet> >( iConfig.getParameter<edm::InputTag>("jets")) ),
210  bHadronsToken_(consumes<reco::GenParticleRefVector>( iConfig.getParameter<edm::InputTag>("bHadrons") )),
211  cHadronsToken_(consumes<reco::GenParticleRefVector>( iConfig.getParameter<edm::InputTag>("cHadrons") )),
212  partonsToken_(consumes<reco::GenParticleRefVector>( iConfig.getParameter<edm::InputTag>("partons") )),
213  jetAlgorithm_(iConfig.getParameter<std::string>("jetAlgorithm")),
214  rParam_(iConfig.getParameter<double>("rParam")),
215  jetPtMin_(0.), // hardcoded to 0. since we simply want to recluster all input jets which already had some PtMin applied
216  ghostRescaling_(iConfig.exists("ghostRescaling") ? iConfig.getParameter<double>("ghostRescaling") : 1e-18),
217  relPtTolerance_(iConfig.exists("relPtTolerance") ? iConfig.getParameter<double>("relPtTolerance") : 1e-03), // 0.1% relative difference in Pt should be sufficient to detect possible misconfigurations
218  hadronFlavourHasPriority_(iConfig.getParameter<bool>("hadronFlavourHasPriority")),
219  useSubjets_(iConfig.exists("groomedJets") && iConfig.exists("subjets")),
220  useLeptons_(iConfig.exists("leptons"))
221 
222 {
223  // register your products
224  produces<reco::JetFlavourInfoMatchingCollection>();
225  if( useSubjets_ )
226  produces<reco::JetFlavourInfoMatchingCollection>("SubJets");
227 
228  // set jet algorithm
229  if (jetAlgorithm_=="Kt")
230  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::kt_algorithm, rParam_) );
231  else if (jetAlgorithm_=="CambridgeAachen")
232  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::cambridge_algorithm, rParam_) );
233  else if (jetAlgorithm_=="AntiKt")
234  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm, rParam_) );
235  else
236  throw cms::Exception("InvalidJetAlgorithm") << "Jet clustering algorithm is invalid: " << jetAlgorithm_ << ", use CambridgeAachen | Kt | AntiKt" << std::endl;
237 
238  if (useSubjets_) {
239  groomedJetsToken_=consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("groomedJets"));
240  subjetsToken_=consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("subjets"));
241  }
242  if ( useLeptons_ ) {
243  leptonsToken_=consumes<reco::GenParticleRefVector>( iConfig.getParameter<edm::InputTag>("leptons") );
244  }
245 
246 }
247 
248 
250 {
251 
252  // do anything here that needs to be done at desctruction time
253  // (e.g. close files, deallocate resources etc.)
254 
255 }
256 
257 
258 //
259 // member functions
260 //
261 
262 // ------------ method called to produce the data ------------
263 void
265 {
267  iEvent.getByToken(jetsToken_, jets);
268 
269  edm::Handle<edm::View<reco::Jet> > groomedJets;
271  if( useSubjets_ )
272  {
273  iEvent.getByToken(groomedJetsToken_, groomedJets);
274  iEvent.getByToken(subjetsToken_, subjets);
275  }
276 
278  iEvent.getByToken(bHadronsToken_, bHadrons);
279 
281  iEvent.getByToken(cHadronsToken_, cHadrons);
282 
284  iEvent.getByToken(partonsToken_, partons);
285 
287  if( useLeptons_ )
288  iEvent.getByToken(leptonsToken_, leptons);
289 
290  auto jetFlavourInfos = std::make_unique<reco::JetFlavourInfoMatchingCollection>(reco::JetRefBaseProd(jets));
291  std::unique_ptr<reco::JetFlavourInfoMatchingCollection> subjetFlavourInfos;
292  if( useSubjets_ )
293  subjetFlavourInfos = std::make_unique<reco::JetFlavourInfoMatchingCollection>(reco::JetRefBaseProd(subjets));
294 
295  // vector of constituents for reclustering jets and "ghosts"
296  std::vector<fastjet::PseudoJet> fjInputs;
297  unsigned int reserve = jets->size()*128 + bHadrons->size() + cHadrons->size() + partons->size();
298  if (useLeptons_) reserve += leptons->size();
299  fjInputs.reserve(reserve);
300  // loop over all input jets and collect all their constituents
301  for(edm::View<reco::Jet>::const_iterator it = jets->begin(); it != jets->end(); ++it)
302  {
303  std::vector<edm::Ptr<reco::Candidate> > constituents = it->getJetConstituents();
304  std::vector<edm::Ptr<reco::Candidate> >::const_iterator m;
305  for( m = constituents.begin(); m != constituents.end(); ++m )
306  {
307  reco::CandidatePtr constit = *m;
308  if(!constit.isNonnull() || !constit.isAvailable()) {
309  edm::LogError("MissingJetConstituent") << "Jet constituent required for jet reclustering is missing. Reclustered jets are not guaranteed to reproduce the original jets!";
310  continue;
311  }
312  if(constit->pt() == 0)
313  {
314  edm::LogWarning("NullTransverseMomentum") << "dropping input candidate with pt=0";
315  continue;
316  }
317  fjInputs.push_back(fastjet::PseudoJet(constit->px(),constit->py(),constit->pz(),constit->energy()));
318  }
319  }
320  // insert "ghost" b hadrons in the vector of constituents
321  insertGhosts(bHadrons, ghostRescaling_, true, true, false, false, fjInputs);
322  // insert "ghost" c hadrons in the vector of constituents
323  insertGhosts(cHadrons, ghostRescaling_, true, false, false, false, fjInputs);
324  // insert "ghost" partons in the vector of constituents
325  insertGhosts(partons, ghostRescaling_, false, false, true, false, fjInputs);
326  // if used, insert "ghost" leptons in the vector of constituents
327  if( useLeptons_ )
328  insertGhosts(leptons, ghostRescaling_, false, false, false, true, fjInputs);
329 
330  // define jet clustering sequence
331  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs, *fjJetDefinition_ ) );
332  // recluster jet constituents and inserted "ghosts"
333  std::vector<fastjet::PseudoJet> inclusiveJets = fastjet::sorted_by_pt( fjClusterSeq_->inclusive_jets(jetPtMin_) );
334 
335  if( inclusiveJets.size() < jets->size() )
336  edm::LogError("TooFewReclusteredJets") << "There are fewer reclustered (" << inclusiveJets.size() << ") than original jets (" << jets->size() << "). Please check that the jet algorithm and jet size match those used for the original jet collection.";
337 
338  // match reclustered and original jets
339  std::vector<int> reclusteredIndices;
340  matchReclusteredJets(jets,inclusiveJets,reclusteredIndices);
341 
342  // match groomed and original jets
343  std::vector<int> groomedIndices;
344  if( useSubjets_ )
345  {
346  if( groomedJets->size() > jets->size() )
347  edm::LogError("TooManyGroomedJets") << "There are more groomed (" << groomedJets->size() << ") than original jets (" << jets->size() << "). Please check that the two jet collections belong to each other.";
348 
349  matchGroomedJets(jets,groomedJets,groomedIndices);
350  }
351 
352  // match subjets and original jets
353  std::vector<std::vector<int> > subjetIndices;
354  if( useSubjets_ )
355  {
356  matchSubjets(groomedIndices,groomedJets,subjets,subjetIndices);
357  }
358 
359  // determine jet flavour
360  for(size_t i=0; i<jets->size(); ++i)
361  {
362  reco::GenParticleRefVector clusteredbHadrons;
363  reco::GenParticleRefVector clusteredcHadrons;
364  reco::GenParticleRefVector clusteredPartons;
365  reco::GenParticleRefVector clusteredLeptons;
366 
367  // if matching reclustered to original jets failed
368  if( reclusteredIndices.at(i) < 0 )
369  {
370  // set an empty JetFlavourInfo for this jet
371  (*jetFlavourInfos)[jets->refAt(i)] = reco::JetFlavourInfo(clusteredbHadrons, clusteredcHadrons, clusteredPartons, clusteredLeptons, 0, 0);
372  }
373  else if( jets->at(i).pt() == 0 )
374  {
375  edm::LogWarning("NullTransverseMomentum") << "The original jet " << i << " has Pt=0. This is not expected so the jet will be skipped.";
376 
377  // set an empty JetFlavourInfo for this jet
378  (*jetFlavourInfos)[jets->refAt(i)] = reco::JetFlavourInfo(clusteredbHadrons, clusteredcHadrons, clusteredPartons, clusteredLeptons, 0, 0);
379 
380  // if subjets are used
381  if( useSubjets_ && !subjetIndices.at(i).empty() )
382  {
383  // loop over subjets
384  for(size_t sj=0; sj<subjetIndices.at(i).size(); ++sj)
385  {
386  // set an empty JetFlavourInfo for this subjet
387  (*subjetFlavourInfos)[subjets->refAt(subjetIndices.at(i).at(sj))] = reco::JetFlavourInfo(reco::GenParticleRefVector(), reco::GenParticleRefVector(), reco::GenParticleRefVector(), reco::GenParticleRefVector(), 0, 0);
388  }
389  }
390  }
391  else
392  {
393  // since the "ghosts" are extremely soft, the configuration and ordering of the reclustered and original jets should in principle stay the same
394  if( ( std::abs( inclusiveJets.at(reclusteredIndices.at(i)).pt() - jets->at(i).pt() ) / jets->at(i).pt() ) > relPtTolerance_ )
395  {
396  if( jets->at(i).pt() < 10. ) // special handling for low-Pt jets (Pt<10 GeV)
397  edm::LogWarning("JetPtMismatchAtLowPt") << "The reclustered and original jet " << i << " have different Pt's (" << inclusiveJets.at(reclusteredIndices.at(i)).pt() << " vs " << jets->at(i).pt() << " GeV, respectively).\n"
398  << "Please check that the jet algorithm and jet size match those used for the original jet collection and also make sure the original jets are uncorrected. In addition, make sure you are not using CaloJets which are presently not supported.\n"
399  << "Since the mismatch is at low Pt (Pt<10 GeV), it is ignored and only a warning is issued.\n"
400  << "\nIn extremely rare instances the mismatch could be caused by a difference in the machine precision in which case make sure the original jet collection is produced and reclustering is performed in the same job.";
401  else
402  edm::LogError("JetPtMismatch") << "The reclustered and original jet " << i << " have different Pt's (" << inclusiveJets.at(reclusteredIndices.at(i)).pt() << " vs " << jets->at(i).pt() << " GeV, respectively).\n"
403  << "Please check that the jet algorithm and jet size match those used for the original jet collection and also make sure the original jets are uncorrected. In addition, make sure you are not using CaloJets which are presently not supported.\n"
404  << "\nIn extremely rare instances the mismatch could be caused by a difference in the machine precision in which case make sure the original jet collection is produced and reclustering is performed in the same job.";
405  }
406 
407  // get jet constituents (sorted by Pt)
408  std::vector<fastjet::PseudoJet> constituents = fastjet::sorted_by_pt( inclusiveJets.at(reclusteredIndices.at(i)).constituents() );
409 
410  // loop over jet constituents and try to find "ghosts"
411  for(std::vector<fastjet::PseudoJet>::const_iterator it = constituents.begin(); it != constituents.end(); ++it)
412  {
413  if( !it->has_user_info() ) continue; // skip if not a "ghost"
414 
415  // "ghost" hadron
416  if( it->user_info<GhostInfo>().isHadron() )
417  {
418  // "ghost" b hadron
419  if( it->user_info<GhostInfo>().isbHadron() )
420  clusteredbHadrons.push_back(it->user_info<GhostInfo>().particleRef());
421  // "ghost" c hadron
422  else
423  clusteredcHadrons.push_back(it->user_info<GhostInfo>().particleRef());
424  }
425  // "ghost" parton
426  else if( it->user_info<GhostInfo>().isParton() )
427  clusteredPartons.push_back(it->user_info<GhostInfo>().particleRef());
428  // "ghost" lepton
429  else if( it->user_info<GhostInfo>().isLepton() )
430  clusteredLeptons.push_back(it->user_info<GhostInfo>().particleRef());
431  }
432 
433  int hadronFlavour = 0; // default hadron flavour set to 0 (= undefined)
434  int partonFlavour = 0; // default parton flavour set to 0 (= undefined)
435 
436  // set hadron- and parton-based flavours
437  setFlavours(clusteredbHadrons, clusteredcHadrons, clusteredPartons, hadronFlavour, partonFlavour);
438 
439  // set the JetFlavourInfo for this jet
440  (*jetFlavourInfos)[jets->refAt(i)] = reco::JetFlavourInfo(clusteredbHadrons, clusteredcHadrons, clusteredPartons, clusteredLeptons, hadronFlavour, partonFlavour);
441  }
442 
443  // if subjets are used, determine their flavour
444  if( useSubjets_ )
445  {
446  if( subjetIndices.at(i).empty() ) continue; // continue if the original jet does not have subjets assigned
447 
448  // define vectors of GenParticleRefVectors for hadrons and partons assigned to different subjets
449  std::vector<reco::GenParticleRefVector> assignedbHadrons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
450  std::vector<reco::GenParticleRefVector> assignedcHadrons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
451  std::vector<reco::GenParticleRefVector> assignedPartons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
452  std::vector<reco::GenParticleRefVector> assignedLeptons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
453 
454  // loop over clustered b hadrons and assign them to different subjets based on smallest dR
455  assignToSubjets(clusteredbHadrons, subjets, subjetIndices.at(i), assignedbHadrons);
456  // loop over clustered c hadrons and assign them to different subjets based on smallest dR
457  assignToSubjets(clusteredcHadrons, subjets, subjetIndices.at(i), assignedcHadrons);
458  // loop over clustered partons and assign them to different subjets based on smallest dR
459  assignToSubjets(clusteredPartons, subjets, subjetIndices.at(i), assignedPartons);
460  // if used, loop over clustered leptons and assign them to different subjets based on smallest dR
461  if( useLeptons_ )
462  assignToSubjets(clusteredLeptons, subjets, subjetIndices.at(i), assignedLeptons);
463 
464  // loop over subjets and determine their flavour
465  for(size_t sj=0; sj<subjetIndices.at(i).size(); ++sj)
466  {
467  int subjetHadronFlavour = 0; // default hadron flavour set to 0 (= undefined)
468  int subjetPartonFlavour = 0; // default parton flavour set to 0 (= undefined)
469 
470  // set hadron- and parton-based flavours
471  setFlavours(assignedbHadrons.at(sj), assignedcHadrons.at(sj), assignedPartons.at(sj), subjetHadronFlavour, subjetPartonFlavour);
472 
473  // set the JetFlavourInfo for this subjet
474  (*subjetFlavourInfos)[subjets->refAt(subjetIndices.at(i).at(sj))] = reco::JetFlavourInfo(assignedbHadrons.at(sj), assignedcHadrons.at(sj), assignedPartons.at(sj), assignedLeptons.at(sj), subjetHadronFlavour, subjetPartonFlavour);
475  }
476  }
477  }
478 
479  //deallocate only at the end of the event processing
480  fjClusterSeq_.reset();
481 
482  // put jet flavour infos in the event
483  iEvent.put(std::move(jetFlavourInfos));
484  // put subjet flavour infos in the event
485  if( useSubjets_ )
486  iEvent.put(std::move(subjetFlavourInfos), "SubJets" );
487 }
488 
489 // ------------ method that inserts "ghost" particles in the vector of jet constituents ------------
490 void
492  const double ghostRescaling,
493  const bool isHadron, const bool isbHadron, const bool isParton, const bool isLepton,
494  std::vector<fastjet::PseudoJet>& constituents)
495 {
496  // insert "ghost" particles in the vector of jet constituents
497  for(reco::GenParticleRefVector::const_iterator it = particles->begin(); it != particles->end(); ++it)
498  {
499  if((*it)->pt() == 0)
500  {
501  edm::LogInfo("NullTransverseMomentum") << "dropping input ghost candidate with pt=0";
502  continue;
503  }
504  fastjet::PseudoJet p((*it)->px(),(*it)->py(),(*it)->pz(),(*it)->energy());
505  p*=ghostRescaling; // rescale particle momentum
506  p.set_user_info(new GhostInfo(isHadron, isbHadron, isParton, isLepton, *it));
507  constituents.push_back(p);
508  }
509 }
510 
511 // ------------ method that matches reclustered and original jets based on minimum dR ------------
512 void
514  const std::vector<fastjet::PseudoJet>& reclusteredJets,
515  std::vector<int>& matchedIndices)
516 {
517  std::vector<bool> matchedLocks(reclusteredJets.size(),false);
518 
519  for(size_t j=0; j<jets->size(); ++j)
520  {
521  double matchedDR2 = 1e9;
522  int matchedIdx = -1;
523 
524  for(size_t rj=0; rj<reclusteredJets.size(); ++rj)
525  {
526  if( matchedLocks.at(rj) ) continue; // skip jets that have already been matched
527 
528  double tempDR2 = reco::deltaR2( jets->at(j).rapidity(), jets->at(j).phi(), reclusteredJets.at(rj).rapidity(), reclusteredJets.at(rj).phi_std() );
529  if( tempDR2 < matchedDR2 )
530  {
531  matchedDR2 = tempDR2;
532  matchedIdx = rj;
533  }
534  }
535 
536  if( matchedIdx>=0 )
537  {
538  if ( matchedDR2 > rParam_*rParam_ )
539  {
540  edm::LogError("JetMatchingFailed") << "Matched reclustered jet " << matchedIdx << " and original jet " << j <<" are separated by dR=" << sqrt(matchedDR2) << " which is greater than the jet size R=" << rParam_ << ".\n"
541  << "This is not expected so please check that the jet algorithm and jet size match those used for the original jet collection.";
542  }
543  else
544  matchedLocks.at(matchedIdx) = true;
545  }
546  else
547  edm::LogError("JetMatchingFailed") << "Matching reclustered to original jets failed. Please check that the jet algorithm and jet size match those used for the original jet collection.";
548 
549  matchedIndices.push_back(matchedIdx);
550  }
551 }
552 
553 // ------------ method that matches groomed and original jets based on minimum dR ------------
554 void
556  const edm::Handle<edm::View<reco::Jet> >& groomedJets,
557  std::vector<int>& matchedIndices)
558 {
559  std::vector<bool> jetLocks(jets->size(),false);
560  std::vector<int> jetIndices;
561 
562  for(size_t gj=0; gj<groomedJets->size(); ++gj)
563  {
564  double matchedDR2 = 1e9;
565  int matchedIdx = -1;
566 
567  if( groomedJets->at(gj).pt()>0. ) // skip pathological cases of groomed jets with Pt=0
568  {
569  for(size_t j=0; j<jets->size(); ++j)
570  {
571  if( jetLocks.at(j) ) continue; // skip jets that have already been matched
572 
573  double tempDR2 = reco::deltaR2( jets->at(j).rapidity(), jets->at(j).phi(), groomedJets->at(gj).rapidity(), groomedJets->at(gj).phi() );
574  if( tempDR2 < matchedDR2 )
575  {
576  matchedDR2 = tempDR2;
577  matchedIdx = j;
578  }
579  }
580  }
581 
582  if( matchedIdx>=0 )
583  {
584  if ( matchedDR2 > rParam_*rParam_ )
585  {
586  edm::LogWarning("MatchedJetsFarApart") << "Matched groomed jet " << gj << " and original jet " << matchedIdx <<" are separated by dR=" << sqrt(matchedDR2) << " which is greater than the jet size R=" << rParam_ << ".\n"
587  << "This is not expected so the matching of these two jets has been discarded. Please check that the two jet collections belong to each other.";
588  matchedIdx = -1;
589  }
590  else
591  jetLocks.at(matchedIdx) = true;
592  }
593  jetIndices.push_back(matchedIdx);
594  }
595 
596  for(size_t j=0; j<jets->size(); ++j)
597  {
598  std::vector<int>::iterator matchedIndex = std::find( jetIndices.begin(), jetIndices.end(), j );
599 
600  matchedIndices.push_back( matchedIndex != jetIndices.end() ? std::distance(jetIndices.begin(),matchedIndex) : -1 );
601  }
602 }
603 
604 // ------------ method that matches subjets and original jets ------------
605 void
606 JetFlavourClustering::matchSubjets(const std::vector<int>& groomedIndices,
607  const edm::Handle<edm::View<reco::Jet> >& groomedJets,
608  const edm::Handle<edm::View<reco::Jet> >& subjets,
609  std::vector<std::vector<int> >& matchedIndices)
610 {
611  for(size_t g=0; g<groomedIndices.size(); ++g)
612  {
613  std::vector<int> subjetIndices;
614 
615  if( groomedIndices.at(g)>=0 )
616  {
617  for(size_t s=0; s<groomedJets->at(groomedIndices.at(g)).numberOfDaughters(); ++s)
618  {
619  const edm::Ptr<reco::Candidate> & subjet = groomedJets->at(groomedIndices.at(g)).daughterPtr(s);
620 
621  for(size_t sj=0; sj<subjets->size(); ++sj)
622  {
623  if( subjet == edm::Ptr<reco::Candidate>(subjets->ptrAt(sj)) )
624  {
625  subjetIndices.push_back(sj);
626  break;
627  }
628  }
629  }
630 
631  if( subjetIndices.empty() )
632  edm::LogError("SubjetMatchingFailed") << "Matching subjets to original jets failed. Please check that the groomed jet and subjet collections belong to each other.";
633 
634  matchedIndices.push_back(subjetIndices);
635  }
636  else
637  matchedIndices.push_back(subjetIndices);
638  }
639 }
640 
641 // ------------ method that sets hadron- and parton-based flavours ------------
642 void
644  const reco::GenParticleRefVector& clusteredcHadrons,
645  const reco::GenParticleRefVector& clusteredPartons,
646  int& hadronFlavour,
647  int& partonFlavour)
648 {
649  reco::GenParticleRef hardestParton;
650  reco::GenParticleRef hardestLightParton;
651  reco::GenParticleRef flavourParton;
652 
653  // loop over clustered partons (already sorted by Pt)
654  for(reco::GenParticleRefVector::const_iterator it = clusteredPartons.begin(); it != clusteredPartons.end(); ++it)
655  {
656  // hardest parton
657  if( hardestParton.isNull() )
658  hardestParton = (*it);
659  // hardest light-flavour parton
660  if( hardestLightParton.isNull() )
661  {
662  if( CandMCTagUtils::isLightParton( *(*it) ) )
663  hardestLightParton = (*it);
664  }
665  // c flavour
666  if( flavourParton.isNull() && ( std::abs( (*it)->pdgId() ) == 4 ) )
667  flavourParton = (*it);
668  // b flavour gets priority
669  if( std::abs( (*it)->pdgId() ) == 5 )
670  {
671  if( flavourParton.isNull() )
672  flavourParton = (*it);
673  else if( std::abs( flavourParton->pdgId() ) != 5 )
674  flavourParton = (*it);
675  }
676  }
677 
678  // set hadron-based flavour
679  if( !clusteredbHadrons.empty() )
680  hadronFlavour = 5;
681  else if( !clusteredcHadrons.empty() && clusteredbHadrons.empty() )
682  hadronFlavour = 4;
683  // set parton-based flavour
684  if( flavourParton.isNull() )
685  {
686  if( hardestParton.isNonnull() ) partonFlavour = hardestParton->pdgId();
687  }
688  else
689  partonFlavour = flavourParton->pdgId();
690 
691  // if enabled, check for conflicts between hadron- and parton-based flavours and give priority to the hadron-based flavour
693  {
694  if( hadronFlavour==0 && (std::abs(partonFlavour)==4 || std::abs(partonFlavour)==5) )
695  partonFlavour = ( hardestLightParton.isNonnull() ? hardestLightParton->pdgId() : 0 );
696  else if( hadronFlavour!=0 && std::abs(partonFlavour)!=hadronFlavour )
697  partonFlavour = hadronFlavour;
698  }
699 }
700 
701 // ------------ method that assigns clustered particles to subjets ------------
702 void
704  const edm::Handle<edm::View<reco::Jet> >& subjets,
705  const std::vector<int>& subjetIndices,
706  std::vector<reco::GenParticleRefVector>& assignedParticles)
707 {
708  // loop over clustered particles and assign them to different subjets based on smallest dR
709  for(reco::GenParticleRefVector::const_iterator it = clusteredParticles.begin(); it != clusteredParticles.end(); ++it)
710  {
711  std::vector<double> dR2toSubjets;
712 
713  for(size_t sj=0; sj<subjetIndices.size(); ++sj)
714  dR2toSubjets.push_back( reco::deltaR2( (*it)->rapidity(), (*it)->phi(), subjets->at(subjetIndices.at(sj)).rapidity(), subjets->at(subjetIndices.at(sj)).phi() ) );
715 
716  // find the closest subjet
717  int closestSubjetIdx = std::distance( dR2toSubjets.begin(), std::min_element(dR2toSubjets.begin(), dR2toSubjets.end()) );
718 
719  assignedParticles.at(closestSubjetIdx).push_back( *it );
720  }
721 }
722 
723 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
724 void
726  //The following says we do not know what parameters are allowed so do no validation
727  // Please change this to state exactly what you do use, even if it is no parameters
729  desc.setUnknown();
730  descriptions.addDefault(desc);
731 }
732 
733 //define this as a plug-in
const reco::GenParticleRef & particleRef() const
edm::EDGetTokenT< reco::GenParticleRefVector > leptonsToken_
std::shared_ptr< fastjet::JetDefinition > JetDefPtr
T getParameter(std::string const &) const
bool isLightParton(const reco::Candidate &c)
Definition: CandMCTag.cc:63
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
void assignToSubjets(const reco::GenParticleRefVector &clusteredParticles, const edm::Handle< edm::View< reco::Jet > > &subjets, const std::vector< int > &subjetIndices, std::vector< reco::GenParticleRefVector > &assignedParticles)
Clusters hadrons, partons, and jet contituents to determine the jet flavour.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void insertGhosts(const edm::Handle< reco::GenParticleRefVector > &particles, const double ghostRescaling, const bool isHadron, const bool isbHadron, const bool isParton, const bool isLepton, std::vector< fastjet::PseudoJet > &constituents)
GhostInfo(const bool &isHadron, const bool &isbHadron, const bool &isParton, const bool &isLepton, const reco::GenParticleRef &particleRef)
edm::EDGetTokenT< edm::View< reco::Jet > > subjetsToken_
bool isAvailable() const
Definition: Ptr.h:258
const reco::GenParticleRef m_particleRef
void setFlavours(const reco::GenParticleRefVector &clusteredbHadrons, const reco::GenParticleRefVector &clusteredcHadrons, const reco::GenParticleRefVector &clusteredPartons, int &hadronFlavour, int &partonFlavour)
const edm::EDGetTokenT< reco::GenParticleRefVector > bHadronsToken_
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:19
void matchGroomedJets(const edm::Handle< edm::View< reco::Jet > > &jets, const edm::Handle< edm::View< reco::Jet > > &matchedJets, std::vector< int > &matchedIndices)
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:253
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:104
const edm::EDGetTokenT< reco::GenParticleRefVector > partonsToken_
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
const bool isParton() const
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:248
const bool isLepton() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::RefToBaseProd< reco::Jet > JetRefBaseProd
Definition: JetCollection.h:14
const edm::EDGetTokenT< reco::GenParticleRefVector > cHadronsToken_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
Definition: Jet.py:1
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
Class storing the jet flavour information.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
JetFlavourClustering(const edm::ParameterSet &)
bool isNull() const
Checks for null.
Definition: Ref.h:248
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
void matchReclusteredJets(const edm::Handle< edm::View< reco::Jet > > &jets, const std::vector< fastjet::PseudoJet > &matchedJets, std::vector< int > &matchedIndices)
const std::string jetAlgorithm_
edm::RefVector< GenParticleCollection > GenParticleRefVector
vector of reference to GenParticle in the same collection
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
const edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool isParton(const reco::Candidate &c)
Definition: CandMCTag.cc:48
hadronFlavour
Definition: jets_cff.py:476
void produce(edm::Event &, const edm::EventSetup &) override
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
partonFlavour
Definition: jets_cff.py:475
const bool isbHadron() const
edm::EDGetTokenT< edm::View< reco::Jet > > groomedJetsToken_
const bool isHadron() const
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
void matchSubjets(const std::vector< int > &groomedIndices, const edm::Handle< edm::View< reco::Jet > > &groomedJets, const edm::Handle< edm::View< reco::Jet > > &subjets, std::vector< std::vector< int > > &matchedIndices)
ClusterSequencePtr fjClusterSeq_
def move(src, dest)
Definition: eostools.py:511