CMS 3D CMS Logo

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