CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
JetFlavourClustering Class Reference

Clusters hadrons, partons, and jet contituents to determine the jet flavour. More...

#include <PhysicsTools/JetMCAlgos/plugins/JetFlavourClustering.cc>

Inheritance diagram for JetFlavourClustering:
edm::stream::EDProducer<> edm::stream::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 JetFlavourClustering (const edm::ParameterSet &)
 
 ~JetFlavourClustering ()
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
- Public Member Functions inherited from edm::stream::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducerBase ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Static Public Member Functions

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

Private Member Functions

void assignToSubjets (const reco::GenParticleRefVector &clusteredParticles, const edm::Handle< edm::View< reco::Jet > > &subjets, const std::vector< int > &subjetIndices, std::vector< reco::GenParticleRefVector > &assignedParticles)
 
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)
 
void matchGroomedJets (const edm::Handle< edm::View< reco::Jet > > &jets, const edm::Handle< edm::View< reco::Jet > > &matchedJets, std::vector< int > &matchedIndices)
 
void matchReclusteredJets (const edm::Handle< edm::View< reco::Jet > > &jets, const std::vector< fastjet::PseudoJet > &matchedJets, std::vector< int > &matchedIndices)
 
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)
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 
void setFlavours (const reco::GenParticleRefVector &clusteredbHadrons, const reco::GenParticleRefVector &clusteredcHadrons, const reco::GenParticleRefVector &clusteredPartons, int &hadronFlavour, int &partonFlavour)
 

Private Attributes

const edm::EDGetTokenT
< reco::GenParticleRefVector
bHadronsToken_
 
const edm::EDGetTokenT
< reco::GenParticleRefVector
cHadronsToken_
 
ClusterSequencePtr fjClusterSeq_
 
JetDefPtr fjJetDefinition_
 
const double ghostRescaling_
 
const edm::EDGetTokenT
< edm::View< reco::Jet > > 
groomedJetsToken_
 
const bool hadronFlavourHasPriority_
 
const std::string jetAlgorithm_
 
const double jetPtMin_
 
const edm::EDGetTokenT
< edm::View< reco::Jet > > 
jetsToken_
 
const edm::EDGetTokenT
< reco::GenParticleRefVector
leptonsToken_
 
const edm::EDGetTokenT
< reco::GenParticleRefVector
partonsToken_
 
const double relPtTolerance_
 
const double rParam_
 
const edm::EDGetTokenT
< edm::View< reco::Jet > > 
subjetsToken_
 
const bool useLeptons_
 
const bool useSubjets_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T...> CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T...> HasAbility
 
typedef
CacheTypes::LuminosityBlockCache 
LuminosityBlockCache
 
typedef
LuminosityBlockContextT
< LuminosityBlockCache,
RunCache, GlobalCache
LuminosityBlockContext
 
typedef
CacheTypes::LuminosityBlockSummaryCache 
LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache,
GlobalCache
RunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDProducerBase
typedef EDProducerAdaptorBase ModuleType
 
- Public Types inherited from edm::ProducerBase
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

Clusters hadrons, partons, and jet contituents to determine the jet flavour.

This producer clusters hadrons, partons and jet contituents to determine the jet flavour. The jet flavour information is stored in the event as an AssociationVector which associates an object of type JetFlavourInfo to each of the jets.

The producer takes as input jets and hadron and partons selected by the HadronAndPartonSelector producer. The hadron and parton four-momenta are rescaled by a very small number (default rescale factor is 10e-18) which turns them into the so-called "ghosts". The "ghost" hadrons and partons are clustered together with all of the jet constituents. It is important to use the same clustering algorithm and jet size as for the original input jet collection. Since the "ghost" hadrons and partons are extremely soft, the resulting jet collection will be practically identical to the original one but now with "ghost" hadrons and partons clustered inside jets. The jet flavour is determined based on the "ghost" hadrons clustered inside a jet:

To further assign a more specific flavour to light-flavour jets, "ghost" partons are used:

In rare instances a conflict between the hadron- and parton-based flavours can occur. In such cases it is possible to keep both flavours or to give priority to the hadron-based flavour. This is controlled by the 'hadronFlavourHasPriority' switch. The priority is given to the hadron-based flavour as follows:

The producer is also capable of assigning the flavour to subjets of fat jets, in which case it produces an additional AssociationVector providing the flavour information for subjets. In order to assign the flavour to subjets, three input jet collections are required:

The "ghost" hadrons and partons clustered inside a fat jet are assigned to the closest subjet in the rapidity-phi space. Once hadrons and partons have been assigned to subjets, the subjet flavour is determined in the same way as for jets. The reason for requiring three jet collections as input in order to determine the subjet flavour is to avoid possible inconsistencies between the fat jet and subjet flavours (such as a non-b fat jet having a b subjet and vice versa) as well as the fact that re-clustering the constituents of groomed fat jets will generally result in a jet collection different from the input groomed fat jets. Also note that "ghost" particles generally cannot be clustered inside subjets in the same way this is done for fat jets. This is because some of the jet grooming techniques could reject such very soft particle. So instead, the "ghost" particles are assigned to the closest subjet.

Finally, "ghost" leptons can also be clustered inside jets but they are not used in any way to determine the jet flavour. This functionality is optional and is potentially useful to identify jets from hadronic taus.

For more details, please refer to https://twiki.cern.ch/twiki/bin/view/CMSPublic/SWGuideBTagMCTools

Definition at line 141 of file JetFlavourClustering.cc.

Constructor & Destructor Documentation

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

Definition at line 207 of file JetFlavourClustering.cc.

References Exception, fjJetDefinition_, jetAlgorithm_, rParam_, and useSubjets_.

207  :
208 
210  groomedJetsToken_(mayConsume<edm::View<reco::Jet> >( iConfig.exists("groomedJets") ? iConfig.getParameter<edm::InputTag>("groomedJets") : edm::InputTag() )),
211  subjetsToken_(mayConsume<edm::View<reco::Jet> >( iConfig.exists("subjets") ? iConfig.getParameter<edm::InputTag>("subjets") : edm::InputTag() )),
212  bHadronsToken_(consumes<reco::GenParticleRefVector>( iConfig.getParameter<edm::InputTag>("bHadrons") )),
213  cHadronsToken_(consumes<reco::GenParticleRefVector>( iConfig.getParameter<edm::InputTag>("cHadrons") )),
214  partonsToken_(consumes<reco::GenParticleRefVector>( iConfig.getParameter<edm::InputTag>("partons") )),
215  leptonsToken_(mayConsume<reco::GenParticleRefVector>( iConfig.exists("leptons") ? iConfig.getParameter<edm::InputTag>("leptons") : edm::InputTag() )),
216  jetAlgorithm_(iConfig.getParameter<std::string>("jetAlgorithm")),
217  rParam_(iConfig.getParameter<double>("rParam")),
218  jetPtMin_(0.), // hardcoded to 0. since we simply want to recluster all input jets which already had some PtMin applied
219  ghostRescaling_(iConfig.exists("ghostRescaling") ? iConfig.getParameter<double>("ghostRescaling") : 1e-18),
220  relPtTolerance_(iConfig.exists("relPtTolerance") ? iConfig.getParameter<double>("relPtTolerance") : 1e-03), // 0.1% relative difference in Pt should be sufficient to detect possible misconfigurations
221  hadronFlavourHasPriority_(iConfig.getParameter<bool>("hadronFlavourHasPriority")),
222  useSubjets_(iConfig.exists("groomedJets") && iConfig.exists("subjets")),
223  useLeptons_(iConfig.exists("leptons"))
224 
225 {
226  // register your products
227  produces<reco::JetFlavourInfoMatchingCollection>();
228  if( useSubjets_ )
229  produces<reco::JetFlavourInfoMatchingCollection>("SubJets");
230 
231  // set jet algorithm
232  if (jetAlgorithm_=="Kt")
233  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::kt_algorithm, rParam_) );
234  else if (jetAlgorithm_=="CambridgeAachen")
235  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::cambridge_algorithm, rParam_) );
236  else if (jetAlgorithm_=="AntiKt")
237  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm, rParam_) );
238  else
239  throw cms::Exception("InvalidJetAlgorithm") << "Jet clustering algorithm is invalid: " << jetAlgorithm_ << ", use CambridgeAachen | Kt | AntiKt" << std::endl;
240 }
T getParameter(std::string const &) const
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
bool exists(std::string const &parameterName) const
checks if a parameter exists
const edm::EDGetTokenT< reco::GenParticleRefVector > bHadronsToken_
const edm::EDGetTokenT< reco::GenParticleRefVector > partonsToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EDGetTokenT< reco::GenParticleRefVector > cHadronsToken_
const std::string jetAlgorithm_
const edm::EDGetTokenT< edm::View< reco::Jet > > groomedJetsToken_
const edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken_
const edm::EDGetTokenT< reco::GenParticleRefVector > leptonsToken_
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
const edm::EDGetTokenT< edm::View< reco::Jet > > subjetsToken_
JetFlavourClustering::~JetFlavourClustering ( )

Definition at line 243 of file JetFlavourClustering.cc.

244 {
245 
246  // do anything here that needs to be done at desctruction time
247  // (e.g. close files, deallocate resources etc.)
248 
249 }

Member Function Documentation

void JetFlavourClustering::assignToSubjets ( const reco::GenParticleRefVector clusteredParticles,
const edm::Handle< edm::View< reco::Jet > > &  subjets,
const std::vector< int > &  subjetIndices,
std::vector< reco::GenParticleRefVector > &  assignedParticles 
)
private

Definition at line 687 of file JetFlavourClustering.cc.

References edm::RefVector< C, T, F >::begin(), reco::deltaR2(), HLT_25ns14e33_v1_cff::distance, edm::RefVector< C, T, F >::end(), and phi.

Referenced by produce().

691 {
692  // loop over clustered particles and assign them to different subjets based on smallest dR
693  for(reco::GenParticleRefVector::const_iterator it = clusteredParticles.begin(); it != clusteredParticles.end(); ++it)
694  {
695  std::vector<double> dR2toSubjets;
696 
697  for(size_t sj=0; sj<subjetIndices.size(); ++sj)
698  dR2toSubjets.push_back( reco::deltaR2( (*it)->rapidity(), (*it)->phi(), subjets->at(subjetIndices.at(sj)).rapidity(), subjets->at(subjetIndices.at(sj)).phi() ) );
699 
700  // find the closest subjet
701  int closestSubjetIdx = std::distance( dR2toSubjets.begin(), std::min_element(dR2toSubjets.begin(), dR2toSubjets.end()) );
702 
703  assignedParticles.at(closestSubjetIdx).push_back( *it );
704  }
705 }
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
const_reference at(size_type pos) const
void JetFlavourClustering::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 709 of file JetFlavourClustering.cc.

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

709  {
710  //The following says we do not know what parameters are allowed so do no validation
711  // Please change this to state exactly what you do use, even if it is no parameters
713  desc.setUnknown();
714  descriptions.addDefault(desc);
715 }
void addDefault(ParameterSetDescription const &psetDescription)
void JetFlavourClustering::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 
)
private

Definition at line 475 of file JetFlavourClustering.cc.

References AlCaHLTBitMon_ParallelJobs::p.

Referenced by produce().

479 {
480  // insert "ghost" particles in the vector of jet constituents
481  for(reco::GenParticleRefVector::const_iterator it = particles->begin(); it != particles->end(); ++it)
482  {
483  if((*it)->pt() == 0)
484  {
485  edm::LogInfo("NullTransverseMomentum") << "dropping input ghost candidate with pt=0";
486  continue;
487  }
488  fastjet::PseudoJet p((*it)->px(),(*it)->py(),(*it)->pz(),(*it)->energy());
489  p*=ghostRescaling; // rescale particle momentum
490  p.set_user_info(new GhostInfo(isHadron, isbHadron, isParton, isLepton, *it));
491  constituents.push_back(p);
492  }
493 }
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:19
bool isParton(const reco::Candidate &c)
Definition: CandMCTag.cc:48
void JetFlavourClustering::matchGroomedJets ( const edm::Handle< edm::View< reco::Jet > > &  jets,
const edm::Handle< edm::View< reco::Jet > > &  matchedJets,
std::vector< int > &  matchedIndices 
)
private

Definition at line 539 of file JetFlavourClustering.cc.

References reco::deltaR2(), HLT_25ns14e33_v1_cff::distance, spr::find(), j, fwrapper::jets, rParam_, and mathSSE::sqrt().

Referenced by produce().

542 {
543  std::vector<bool> jetLocks(jets->size(),false);
544  std::vector<int> jetIndices;
545 
546  for(size_t gj=0; gj<groomedJets->size(); ++gj)
547  {
548  double matchedDR2 = 1e9;
549  int matchedIdx = -1;
550 
551  if( groomedJets->at(gj).pt()>0. ) // skip pathological cases of groomed jets with Pt=0
552  {
553  for(size_t j=0; j<jets->size(); ++j)
554  {
555  if( jetLocks.at(j) ) continue; // skip jets that have already been matched
556 
557  double tempDR2 = reco::deltaR2( jets->at(j).rapidity(), jets->at(j).phi(), groomedJets->at(gj).rapidity(), groomedJets->at(gj).phi() );
558  if( tempDR2 < matchedDR2 )
559  {
560  matchedDR2 = tempDR2;
561  matchedIdx = j;
562  }
563  }
564  }
565 
566  if( matchedIdx>=0 )
567  {
568  if ( matchedDR2 > rParam_*rParam_ )
569  {
570  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"
571  << "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.";
572  matchedIdx = -1;
573  }
574  else
575  jetLocks.at(matchedIdx) = true;
576  }
577  jetIndices.push_back(matchedIdx);
578  }
579 
580  for(size_t j=0; j<jets->size(); ++j)
581  {
582  std::vector<int>::iterator matchedIndex = std::find( jetIndices.begin(), jetIndices.end(), j );
583 
584  matchedIndices.push_back( matchedIndex != jetIndices.end() ? std::distance(jetIndices.begin(),matchedIndex) : -1 );
585  }
586 }
size_type size() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
T sqrt(T t)
Definition: SSEVec.h:48
int j
Definition: DBlmapReader.cc:9
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
const_reference at(size_type pos) const
void JetFlavourClustering::matchReclusteredJets ( const edm::Handle< edm::View< reco::Jet > > &  jets,
const std::vector< fastjet::PseudoJet > &  matchedJets,
std::vector< int > &  matchedIndices 
)
private

Definition at line 497 of file JetFlavourClustering.cc.

References reco::deltaR2(), j, fwrapper::jets, rParam_, and mathSSE::sqrt().

Referenced by produce().

500 {
501  std::vector<bool> matchedLocks(reclusteredJets.size(),false);
502 
503  for(size_t j=0; j<jets->size(); ++j)
504  {
505  double matchedDR2 = 1e9;
506  int matchedIdx = -1;
507 
508  for(size_t rj=0; rj<reclusteredJets.size(); ++rj)
509  {
510  if( matchedLocks.at(rj) ) continue; // skip jets that have already been matched
511 
512  double tempDR2 = reco::deltaR2( jets->at(j).rapidity(), jets->at(j).phi(), reclusteredJets.at(rj).rapidity(), reclusteredJets.at(rj).phi_std() );
513  if( tempDR2 < matchedDR2 )
514  {
515  matchedDR2 = tempDR2;
516  matchedIdx = rj;
517  }
518  }
519 
520  if( matchedIdx>=0 )
521  {
522  if ( matchedDR2 > rParam_*rParam_ )
523  {
524  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"
525  << "This is not expected so please check that the jet algorithm and jet size match those used for the original jet collection.";
526  }
527  else
528  matchedLocks.at(matchedIdx) = true;
529  }
530  else
531  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.";
532 
533  matchedIndices.push_back(matchedIdx);
534  }
535 }
size_type size() const
T sqrt(T t)
Definition: SSEVec.h:48
int j
Definition: DBlmapReader.cc:9
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
const_reference at(size_type pos) const
void JetFlavourClustering::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 
)
private

Definition at line 590 of file JetFlavourClustering.cc.

References g, and alignCSCRings::s.

Referenced by produce().

594 {
595  for(size_t g=0; g<groomedIndices.size(); ++g)
596  {
597  std::vector<int> subjetIndices;
598 
599  if( groomedIndices.at(g)>=0 )
600  {
601  for(size_t s=0; s<groomedJets->at(groomedIndices.at(g)).numberOfDaughters(); ++s)
602  {
603  const edm::Ptr<reco::Candidate> & subjet = groomedJets->at(groomedIndices.at(g)).daughterPtr(s);
604 
605  for(size_t sj=0; sj<subjets->size(); ++sj)
606  {
607  if( subjet == edm::Ptr<reco::Candidate>(subjets->ptrAt(sj)) )
608  {
609  subjetIndices.push_back(sj);
610  break;
611  }
612  }
613  }
614 
615  if( subjetIndices.size() == 0 )
616  edm::LogError("SubjetMatchingFailed") << "Matching subjets to original jets failed. Please check that the groomed jet and subjet collections belong to each other.";
617 
618  matchedIndices.push_back(subjetIndices);
619  }
620  else
621  matchedIndices.push_back(subjetIndices);
622  }
623 }
Ptr< value_type > ptrAt(size_type i) const
size_type size() const
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_reference at(size_type pos) const
void JetFlavourClustering::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
privatevirtual

Implements edm::stream::EDProducerBase.

Definition at line 258 of file JetFlavourClustering.cc.

References funct::abs(), assignToSubjets(), bHadronsToken_, cHadronsToken_, fjClusterSeq_, fjJetDefinition_, edm::Event::getByToken(), ghostRescaling_, groomedJetsToken_, i, insertGhosts(), GhostInfo::isbHadron(), GhostInfo::isHadron(), GhostInfo::isLepton(), GhostInfo::isParton(), jetPtMin_, fwrapper::jets, jetsToken_, HLT_25ns14e33_v1_cff::leptons, leptonsToken_, visualization-live-secondInstance_cfg::m, matchGroomedJets(), matchReclusteredJets(), matchSubjets(), GhostInfo::particleRef(), partonsToken_, EnergyCorrector::pt, edm::RefVector< C, T, F >::push_back(), edm::Event::put(), relPtTolerance_, setFlavours(), subjetsToken_, useLeptons_, and useSubjets_.

259 {
261  iEvent.getByToken(jetsToken_, jets);
262 
263  edm::Handle<edm::View<reco::Jet> > groomedJets;
265  if( useSubjets_ )
266  {
267  iEvent.getByToken(groomedJetsToken_, groomedJets);
268  iEvent.getByToken(subjetsToken_, subjets);
269  }
270 
272  iEvent.getByToken(bHadronsToken_, bHadrons);
273 
275  iEvent.getByToken(cHadronsToken_, cHadrons);
276 
278  iEvent.getByToken(partonsToken_, partons);
279 
281  if( useLeptons_ )
282  iEvent.getByToken(leptonsToken_, leptons);
283 
284  std::auto_ptr<reco::JetFlavourInfoMatchingCollection> jetFlavourInfos( new reco::JetFlavourInfoMatchingCollection(reco::JetRefBaseProd(jets)) );
285  std::auto_ptr<reco::JetFlavourInfoMatchingCollection> subjetFlavourInfos;
286  if( useSubjets_ )
287  subjetFlavourInfos = std::auto_ptr<reco::JetFlavourInfoMatchingCollection>( new reco::JetFlavourInfoMatchingCollection(reco::JetRefBaseProd(subjets)) );
288 
289  // vector of constituents for reclustering jets and "ghosts"
290  std::vector<fastjet::PseudoJet> fjInputs;
291  // loop over all input jets and collect all their constituents
292  for(edm::View<reco::Jet>::const_iterator it = jets->begin(); it != jets->end(); ++it)
293  {
294  std::vector<edm::Ptr<reco::Candidate> > constituents = it->getJetConstituents();
295  std::vector<edm::Ptr<reco::Candidate> >::const_iterator m;
296  for( m = constituents.begin(); m != constituents.end(); ++m )
297  {
298  reco::CandidatePtr constit = *m;
299  if(constit->pt() == 0)
300  {
301  edm::LogWarning("NullTransverseMomentum") << "dropping input candidate with pt=0";
302  continue;
303  }
304  fjInputs.push_back(fastjet::PseudoJet(constit->px(),constit->py(),constit->pz(),constit->energy()));
305  }
306  }
307  // insert "ghost" b hadrons in the vector of constituents
308  insertGhosts(bHadrons, ghostRescaling_, true, true, false, false, fjInputs);
309  // insert "ghost" c hadrons in the vector of constituents
310  insertGhosts(cHadrons, ghostRescaling_, true, false, false, false, fjInputs);
311  // insert "ghost" partons in the vector of constituents
312  insertGhosts(partons, ghostRescaling_, false, false, true, false, fjInputs);
313  // if used, insert "ghost" leptons in the vector of constituents
314  if( useLeptons_ )
315  insertGhosts(leptons, ghostRescaling_, false, false, false, true, fjInputs);
316 
317  // define jet clustering sequence
318  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs, *fjJetDefinition_ ) );
319  // recluster jet constituents and inserted "ghosts"
320  std::vector<fastjet::PseudoJet> inclusiveJets = fastjet::sorted_by_pt( fjClusterSeq_->inclusive_jets(jetPtMin_) );
321 
322  if( inclusiveJets.size() < jets->size() )
323  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.";
324 
325  // match reclustered and original jets
326  std::vector<int> reclusteredIndices;
327  matchReclusteredJets(jets,inclusiveJets,reclusteredIndices);
328 
329  // match groomed and original jets
330  std::vector<int> groomedIndices;
331  if( useSubjets_ )
332  {
333  if( groomedJets->size() > jets->size() )
334  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.";
335 
336  matchGroomedJets(jets,groomedJets,groomedIndices);
337  }
338 
339  // match subjets and original jets
340  std::vector<std::vector<int> > subjetIndices;
341  if( useSubjets_ )
342  {
343  matchSubjets(groomedIndices,groomedJets,subjets,subjetIndices);
344  }
345 
346  // determine jet flavour
347  for(size_t i=0; i<jets->size(); ++i)
348  {
349  reco::GenParticleRefVector clusteredbHadrons;
350  reco::GenParticleRefVector clusteredcHadrons;
351  reco::GenParticleRefVector clusteredPartons;
352  reco::GenParticleRefVector clusteredLeptons;
353 
354  // if matching reclustered to original jets failed
355  if( reclusteredIndices.at(i) < 0 )
356  {
357  // set an empty JetFlavourInfo for this jet
358  (*jetFlavourInfos)[jets->refAt(i)] = reco::JetFlavourInfo(clusteredbHadrons, clusteredcHadrons, clusteredPartons, clusteredLeptons, 0, 0);
359  }
360  else if( jets->at(i).pt() == 0 )
361  {
362  edm::LogWarning("NullTransverseMomentum") << "The original jet " << i << " has Pt=0. This is not expected so the jet will be skipped.";
363 
364  // set an empty JetFlavourInfo for this jet
365  (*jetFlavourInfos)[jets->refAt(i)] = reco::JetFlavourInfo(clusteredbHadrons, clusteredcHadrons, clusteredPartons, clusteredLeptons, 0, 0);
366 
367  // if subjets are used
368  if( useSubjets_ && subjetIndices.at(i).size()>0 )
369  {
370  // loop over subjets
371  for(size_t sj=0; sj<subjetIndices.at(i).size(); ++sj)
372  {
373  // set an empty JetFlavourInfo for this subjet
374  (*subjetFlavourInfos)[subjets->refAt(subjetIndices.at(i).at(sj))] = reco::JetFlavourInfo(reco::GenParticleRefVector(), reco::GenParticleRefVector(), reco::GenParticleRefVector(), reco::GenParticleRefVector(), 0, 0);
375  }
376  }
377  }
378  else
379  {
380  // since the "ghosts" are extremely soft, the configuration and ordering of the reclustered and original jets should in principle stay the same
381  if( ( std::abs( inclusiveJets.at(reclusteredIndices.at(i)).pt() - jets->at(i).pt() ) / jets->at(i).pt() ) > relPtTolerance_ )
382  {
383  if( jets->at(i).pt() < 10. ) // special handling for low-Pt jets (Pt<10 GeV)
384  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"
385  << "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"
386  << "Since the mismatch is at low Pt (Pt<10 GeV), it is ignored and only a warning is issued.\n"
387  << "\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.";
388  else
389  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"
390  << "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"
391  << "\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.";
392  }
393 
394  // get jet constituents (sorted by Pt)
395  std::vector<fastjet::PseudoJet> constituents = fastjet::sorted_by_pt( inclusiveJets.at(reclusteredIndices.at(i)).constituents() );
396 
397  // loop over jet constituents and try to find "ghosts"
398  for(std::vector<fastjet::PseudoJet>::const_iterator it = constituents.begin(); it != constituents.end(); ++it)
399  {
400  if( !it->has_user_info() ) continue; // skip if not a "ghost"
401 
402  // "ghost" hadron
403  if( it->user_info<GhostInfo>().isHadron() )
404  {
405  // "ghost" b hadron
406  if( it->user_info<GhostInfo>().isbHadron() )
407  clusteredbHadrons.push_back(it->user_info<GhostInfo>().particleRef());
408  // "ghost" c hadron
409  else
410  clusteredcHadrons.push_back(it->user_info<GhostInfo>().particleRef());
411  }
412  // "ghost" parton
413  else if( it->user_info<GhostInfo>().isParton() )
414  clusteredPartons.push_back(it->user_info<GhostInfo>().particleRef());
415  // "ghost" lepton
416  else if( it->user_info<GhostInfo>().isLepton() )
417  clusteredLeptons.push_back(it->user_info<GhostInfo>().particleRef());
418  }
419 
420  int hadronFlavour = 0; // default hadron flavour set to 0 (= undefined)
421  int partonFlavour = 0; // default parton flavour set to 0 (= undefined)
422 
423  // set hadron- and parton-based flavours
424  setFlavours(clusteredbHadrons, clusteredcHadrons, clusteredPartons, hadronFlavour, partonFlavour);
425 
426  // set the JetFlavourInfo for this jet
427  (*jetFlavourInfos)[jets->refAt(i)] = reco::JetFlavourInfo(clusteredbHadrons, clusteredcHadrons, clusteredPartons, clusteredLeptons, hadronFlavour, partonFlavour);
428  }
429 
430  // if subjets are used, determine their flavour
431  if( useSubjets_ )
432  {
433  if( subjetIndices.at(i).size()==0 ) continue; // continue if the original jet does not have subjets assigned
434 
435  // define vectors of GenParticleRefVectors for hadrons and partons assigned to different subjets
436  std::vector<reco::GenParticleRefVector> assignedbHadrons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
437  std::vector<reco::GenParticleRefVector> assignedcHadrons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
438  std::vector<reco::GenParticleRefVector> assignedPartons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
439  std::vector<reco::GenParticleRefVector> assignedLeptons(subjetIndices.at(i).size(),reco::GenParticleRefVector());
440 
441  // loop over clustered b hadrons and assign them to different subjets based on smallest dR
442  assignToSubjets(clusteredbHadrons, subjets, subjetIndices.at(i), assignedbHadrons);
443  // loop over clustered c hadrons and assign them to different subjets based on smallest dR
444  assignToSubjets(clusteredcHadrons, subjets, subjetIndices.at(i), assignedcHadrons);
445  // loop over clustered partons and assign them to different subjets based on smallest dR
446  assignToSubjets(clusteredPartons, subjets, subjetIndices.at(i), assignedPartons);
447  // if used, loop over clustered leptons and assign them to different subjets based on smallest dR
448  if( useLeptons_ )
449  assignToSubjets(clusteredLeptons, subjets, subjetIndices.at(i), assignedLeptons);
450 
451  // loop over subjets and determine their flavour
452  for(size_t sj=0; sj<subjetIndices.at(i).size(); ++sj)
453  {
454  int subjetHadronFlavour = 0; // default hadron flavour set to 0 (= undefined)
455  int subjetPartonFlavour = 0; // default parton flavour set to 0 (= undefined)
456 
457  // set hadron- and parton-based flavours
458  setFlavours(assignedbHadrons.at(sj), assignedcHadrons.at(sj), assignedPartons.at(sj), subjetHadronFlavour, subjetPartonFlavour);
459 
460  // set the JetFlavourInfo for this subjet
461  (*subjetFlavourInfos)[subjets->refAt(subjetIndices.at(i).at(sj))] = reco::JetFlavourInfo(assignedbHadrons.at(sj), assignedcHadrons.at(sj), assignedPartons.at(sj), assignedLeptons.at(sj), subjetHadronFlavour, subjetPartonFlavour);
462  }
463  }
464  }
465 
466  // put jet flavour infos in the event
467  iEvent.put( jetFlavourInfos );
468  // put subjet flavour infos in the event
469  if( useSubjets_ )
470  iEvent.put( subjetFlavourInfos, "SubJets" );
471 }
const reco::GenParticleRef & particleRef() const
int i
Definition: DBlmapReader.cc:9
void assignToSubjets(const reco::GenParticleRefVector &clusteredParticles, const edm::Handle< edm::View< reco::Jet > > &subjets, const std::vector< int > &subjetIndices, std::vector< reco::GenParticleRefVector > &assignedParticles)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
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)
void setFlavours(const reco::GenParticleRefVector &clusteredbHadrons, const reco::GenParticleRefVector &clusteredcHadrons, const reco::GenParticleRefVector &clusteredPartons, int &hadronFlavour, int &partonFlavour)
const edm::EDGetTokenT< reco::GenParticleRefVector > bHadronsToken_
void matchGroomedJets(const edm::Handle< edm::View< reco::Jet > > &jets, const edm::Handle< edm::View< reco::Jet > > &matchedJets, std::vector< int > &matchedIndices)
const edm::EDGetTokenT< reco::GenParticleRefVector > partonsToken_
const bool isParton() const
const bool isLepton() const
const edm::EDGetTokenT< reco::GenParticleRefVector > cHadronsToken_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
vector< PseudoJet > jets
Class storing the jet flavour information.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void matchReclusteredJets(const edm::Handle< edm::View< reco::Jet > > &jets, const std::vector< fastjet::PseudoJet > &matchedJets, std::vector< int > &matchedIndices)
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_
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
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:62
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_
void JetFlavourClustering::setFlavours ( const reco::GenParticleRefVector clusteredbHadrons,
const reco::GenParticleRefVector clusteredcHadrons,
const reco::GenParticleRefVector clusteredPartons,
int &  hadronFlavour,
int &  partonFlavour 
)
private

Definition at line 627 of file JetFlavourClustering.cc.

References funct::abs(), edm::RefVector< C, T, F >::begin(), edm::RefVector< C, T, F >::end(), hadronFlavourHasPriority_, CandMCTagUtils::isLightParton(), edm::Ref< C, T, F >::isNonnull(), edm::Ref< C, T, F >::isNull(), and edm::RefVector< C, T, F >::size().

Referenced by produce().

632 {
633  reco::GenParticleRef hardestParton;
634  reco::GenParticleRef hardestLightParton;
635  reco::GenParticleRef flavourParton;
636 
637  // loop over clustered partons (already sorted by Pt)
638  for(reco::GenParticleRefVector::const_iterator it = clusteredPartons.begin(); it != clusteredPartons.end(); ++it)
639  {
640  // hardest parton
641  if( hardestParton.isNull() )
642  hardestParton = (*it);
643  // hardest light-flavour parton
644  if( hardestLightParton.isNull() )
645  {
646  if( CandMCTagUtils::isLightParton( *(*it) ) )
647  hardestLightParton = (*it);
648  }
649  // c flavour
650  if( flavourParton.isNull() && ( std::abs( (*it)->pdgId() ) == 4 ) )
651  flavourParton = (*it);
652  // b flavour gets priority
653  if( std::abs( (*it)->pdgId() ) == 5 )
654  {
655  if( flavourParton.isNull() )
656  flavourParton = (*it);
657  else if( std::abs( flavourParton->pdgId() ) != 5 )
658  flavourParton = (*it);
659  }
660  }
661 
662  // set hadron-based flavour
663  if( clusteredbHadrons.size()>0 )
664  hadronFlavour = 5;
665  else if( clusteredcHadrons.size()>0 && clusteredbHadrons.size()==0 )
666  hadronFlavour = 4;
667  // set parton-based flavour
668  if( flavourParton.isNull() )
669  {
670  if( hardestParton.isNonnull() ) partonFlavour = hardestParton->pdgId();
671  }
672  else
673  partonFlavour = flavourParton->pdgId();
674 
675  // if enabled, check for conflicts between hadron- and parton-based flavours and give priority to the hadron-based flavour
677  {
678  if( hadronFlavour==0 && (std::abs(partonFlavour)==4 || std::abs(partonFlavour)==5) )
679  partonFlavour = ( hardestLightParton.isNonnull() ? hardestLightParton->pdgId() : 0 );
680  else if( hadronFlavour!=0 && std::abs(partonFlavour)!=hadronFlavour )
681  partonFlavour = hadronFlavour;
682  }
683 }
bool isLightParton(const reco::Candidate &c)
Definition: CandMCTag.cc:63
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isNull() const
Checks for null.
Definition: Ref.h:249
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99

Member Data Documentation

const edm::EDGetTokenT<reco::GenParticleRefVector> JetFlavourClustering::bHadronsToken_
private

Definition at line 182 of file JetFlavourClustering.cc.

Referenced by produce().

const edm::EDGetTokenT<reco::GenParticleRefVector> JetFlavourClustering::cHadronsToken_
private

Definition at line 183 of file JetFlavourClustering.cc.

Referenced by produce().

ClusterSequencePtr JetFlavourClustering::fjClusterSeq_
private

Definition at line 196 of file JetFlavourClustering.cc.

Referenced by produce().

JetDefPtr JetFlavourClustering::fjJetDefinition_
private

Definition at line 197 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering(), and produce().

const double JetFlavourClustering::ghostRescaling_
private

Definition at line 190 of file JetFlavourClustering.cc.

Referenced by produce().

const edm::EDGetTokenT<edm::View<reco::Jet> > JetFlavourClustering::groomedJetsToken_
private

Definition at line 180 of file JetFlavourClustering.cc.

Referenced by produce().

const bool JetFlavourClustering::hadronFlavourHasPriority_
private

Definition at line 192 of file JetFlavourClustering.cc.

Referenced by setFlavours().

const std::string JetFlavourClustering::jetAlgorithm_
private

Definition at line 187 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering().

const double JetFlavourClustering::jetPtMin_
private

Definition at line 189 of file JetFlavourClustering.cc.

Referenced by produce().

const edm::EDGetTokenT<edm::View<reco::Jet> > JetFlavourClustering::jetsToken_
private

Definition at line 179 of file JetFlavourClustering.cc.

Referenced by produce().

const edm::EDGetTokenT<reco::GenParticleRefVector> JetFlavourClustering::leptonsToken_
private

Definition at line 185 of file JetFlavourClustering.cc.

Referenced by produce().

const edm::EDGetTokenT<reco::GenParticleRefVector> JetFlavourClustering::partonsToken_
private

Definition at line 184 of file JetFlavourClustering.cc.

Referenced by produce().

const double JetFlavourClustering::relPtTolerance_
private

Definition at line 191 of file JetFlavourClustering.cc.

Referenced by produce().

const double JetFlavourClustering::rParam_
private
const edm::EDGetTokenT<edm::View<reco::Jet> > JetFlavourClustering::subjetsToken_
private

Definition at line 181 of file JetFlavourClustering.cc.

Referenced by produce().

const bool JetFlavourClustering::useLeptons_
private

Definition at line 194 of file JetFlavourClustering.cc.

Referenced by produce().

const bool JetFlavourClustering::useSubjets_
private

Definition at line 193 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering(), and produce().