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 (std::string const &iProcessName, std::string const &iModuleLabel, bool iPrint, std::vector< char const * > &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_
 
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_
 
edm::EDGetTokenT
< reco::GenParticleRefVector
leptonsToken_
 
const edm::EDGetTokenT
< reco::GenParticleRefVector
partonsToken_
 
const double relPtTolerance_
 
const double rParam_
 
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_, edm::ParameterSet::getParameter(), groomedJetsToken_, jetAlgorithm_, leptonsToken_, rParam_, subjetsToken_, useLeptons_, and useSubjets_.

207  :
208 
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 }
edm::EDGetTokenT< reco::GenParticleRefVector > leptonsToken_
T getParameter(std::string const &) const
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
edm::EDGetTokenT< edm::View< reco::Jet > > subjetsToken_
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 > > jetsToken_
edm::EDGetTokenT< edm::View< reco::Jet > > groomedJetsToken_
JetFlavourClustering::~JetFlavourClustering ( )

Definition at line 249 of file JetFlavourClustering.cc.

250 {
251 
252  // do anything here that needs to be done at desctruction time
253  // (e.g. close files, deallocate resources etc.)
254 
255 }

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 697 of file JetFlavourClustering.cc.

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

Referenced by produce().

701 {
702  // loop over clustered particles and assign them to different subjets based on smallest dR
703  for(reco::GenParticleRefVector::const_iterator it = clusteredParticles.begin(); it != clusteredParticles.end(); ++it)
704  {
705  std::vector<double> dR2toSubjets;
706 
707  for(size_t sj=0; sj<subjetIndices.size(); ++sj)
708  dR2toSubjets.push_back( reco::deltaR2( (*it)->rapidity(), (*it)->phi(), subjets->at(subjetIndices.at(sj)).rapidity(), subjets->at(subjetIndices.at(sj)).phi() ) );
709 
710  // find the closest subjet
711  int closestSubjetIdx = std::distance( dR2toSubjets.begin(), std::min_element(dR2toSubjets.begin(), dR2toSubjets.end()) );
712 
713  assignedParticles.at(closestSubjetIdx).push_back( *it );
714  }
715 }
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
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
const_reference at(size_type pos) const
void JetFlavourClustering::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 719 of file JetFlavourClustering.cc.

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

719  {
720  //The following says we do not know what parameters are allowed so do no validation
721  // Please change this to state exactly what you do use, even if it is no parameters
723  desc.setUnknown();
724  descriptions.addDefault(desc);
725 }
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 485 of file JetFlavourClustering.cc.

References AlCaHLTBitMon_ParallelJobs::p.

Referenced by produce().

489 {
490  // insert "ghost" particles in the vector of jet constituents
491  for(reco::GenParticleRefVector::const_iterator it = particles->begin(); it != particles->end(); ++it)
492  {
493  if((*it)->pt() == 0)
494  {
495  edm::LogInfo("NullTransverseMomentum") << "dropping input ghost candidate with pt=0";
496  continue;
497  }
498  fastjet::PseudoJet p((*it)->px(),(*it)->py(),(*it)->pz(),(*it)->energy());
499  p*=ghostRescaling; // rescale particle momentum
500  p.set_user_info(new GhostInfo(isHadron, isbHadron, isParton, isLepton, *it));
501  constituents.push_back(p);
502  }
503 }
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 549 of file JetFlavourClustering.cc.

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

Referenced by produce().

552 {
553  std::vector<bool> jetLocks(jets->size(),false);
554  std::vector<int> jetIndices;
555 
556  for(size_t gj=0; gj<groomedJets->size(); ++gj)
557  {
558  double matchedDR2 = 1e9;
559  int matchedIdx = -1;
560 
561  if( groomedJets->at(gj).pt()>0. ) // skip pathological cases of groomed jets with Pt=0
562  {
563  for(size_t j=0; j<jets->size(); ++j)
564  {
565  if( jetLocks.at(j) ) continue; // skip jets that have already been matched
566 
567  double tempDR2 = reco::deltaR2( jets->at(j).rapidity(), jets->at(j).phi(), groomedJets->at(gj).rapidity(), groomedJets->at(gj).phi() );
568  if( tempDR2 < matchedDR2 )
569  {
570  matchedDR2 = tempDR2;
571  matchedIdx = j;
572  }
573  }
574  }
575 
576  if( matchedIdx>=0 )
577  {
578  if ( matchedDR2 > rParam_*rParam_ )
579  {
580  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"
581  << "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.";
582  matchedIdx = -1;
583  }
584  else
585  jetLocks.at(matchedIdx) = true;
586  }
587  jetIndices.push_back(matchedIdx);
588  }
589 
590  for(size_t j=0; j<jets->size(); ++j)
591  {
592  std::vector<int>::iterator matchedIndex = std::find( jetIndices.begin(), jetIndices.end(), j );
593 
594  matchedIndices.push_back( matchedIndex != jetIndices.end() ? std::distance(jetIndices.begin(),matchedIndex) : -1 );
595  }
596 }
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:18
int j
Definition: DBlmapReader.cc:9
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
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 507 of file JetFlavourClustering.cc.

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

Referenced by produce().

510 {
511  std::vector<bool> matchedLocks(reclusteredJets.size(),false);
512 
513  for(size_t j=0; j<jets->size(); ++j)
514  {
515  double matchedDR2 = 1e9;
516  int matchedIdx = -1;
517 
518  for(size_t rj=0; rj<reclusteredJets.size(); ++rj)
519  {
520  if( matchedLocks.at(rj) ) continue; // skip jets that have already been matched
521 
522  double tempDR2 = reco::deltaR2( jets->at(j).rapidity(), jets->at(j).phi(), reclusteredJets.at(rj).rapidity(), reclusteredJets.at(rj).phi_std() );
523  if( tempDR2 < matchedDR2 )
524  {
525  matchedDR2 = tempDR2;
526  matchedIdx = rj;
527  }
528  }
529 
530  if( matchedIdx>=0 )
531  {
532  if ( matchedDR2 > rParam_*rParam_ )
533  {
534  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"
535  << "This is not expected so please check that the jet algorithm and jet size match those used for the original jet collection.";
536  }
537  else
538  matchedLocks.at(matchedIdx) = true;
539  }
540  else
541  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.";
542 
543  matchedIndices.push_back(matchedIdx);
544  }
545 }
size_type size() const
T sqrt(T t)
Definition: SSEVec.h:18
int j
Definition: DBlmapReader.cc:9
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
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 600 of file JetFlavourClustering.cc.

References g, and alignCSCRings::s.

Referenced by produce().

604 {
605  for(size_t g=0; g<groomedIndices.size(); ++g)
606  {
607  std::vector<int> subjetIndices;
608 
609  if( groomedIndices.at(g)>=0 )
610  {
611  for(size_t s=0; s<groomedJets->at(groomedIndices.at(g)).numberOfDaughters(); ++s)
612  {
613  const edm::Ptr<reco::Candidate> & subjet = groomedJets->at(groomedIndices.at(g)).daughterPtr(s);
614 
615  for(size_t sj=0; sj<subjets->size(); ++sj)
616  {
617  if( subjet == edm::Ptr<reco::Candidate>(subjets->ptrAt(sj)) )
618  {
619  subjetIndices.push_back(sj);
620  break;
621  }
622  }
623  }
624 
625  if( subjetIndices.size() == 0 )
626  edm::LogError("SubjetMatchingFailed") << "Matching subjets to original jets failed. Please check that the groomed jet and subjet collections belong to each other.";
627 
628  matchedIndices.push_back(subjetIndices);
629  }
630  else
631  matchedIndices.push_back(subjetIndices);
632  }
633 }
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 264 of file JetFlavourClustering.cc.

References funct::abs(), assignToSubjets(), bHadronsToken_, cHadronsToken_, fjClusterSeq_, fjJetDefinition_, edm::Event::getByToken(), ghostRescaling_, groomedJetsToken_, i, insertGhosts(), edm::Ptr< T >::isAvailable(), GhostInfo::isbHadron(), GhostInfo::isHadron(), GhostInfo::isLepton(), edm::Ptr< T >::isNonnull(), GhostInfo::isParton(), jetPtMin_, fwrapper::jets, jetsToken_, HLT_25ns10e33_v2_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_.

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

Definition at line 637 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().

642 {
643  reco::GenParticleRef hardestParton;
644  reco::GenParticleRef hardestLightParton;
645  reco::GenParticleRef flavourParton;
646 
647  // loop over clustered partons (already sorted by Pt)
648  for(reco::GenParticleRefVector::const_iterator it = clusteredPartons.begin(); it != clusteredPartons.end(); ++it)
649  {
650  // hardest parton
651  if( hardestParton.isNull() )
652  hardestParton = (*it);
653  // hardest light-flavour parton
654  if( hardestLightParton.isNull() )
655  {
656  if( CandMCTagUtils::isLightParton( *(*it) ) )
657  hardestLightParton = (*it);
658  }
659  // c flavour
660  if( flavourParton.isNull() && ( std::abs( (*it)->pdgId() ) == 4 ) )
661  flavourParton = (*it);
662  // b flavour gets priority
663  if( std::abs( (*it)->pdgId() ) == 5 )
664  {
665  if( flavourParton.isNull() )
666  flavourParton = (*it);
667  else if( std::abs( flavourParton->pdgId() ) != 5 )
668  flavourParton = (*it);
669  }
670  }
671 
672  // set hadron-based flavour
673  if( clusteredbHadrons.size()>0 )
674  hadronFlavour = 5;
675  else if( clusteredcHadrons.size()>0 && clusteredbHadrons.size()==0 )
676  hadronFlavour = 4;
677  // set parton-based flavour
678  if( flavourParton.isNull() )
679  {
680  if( hardestParton.isNonnull() ) partonFlavour = hardestParton->pdgId();
681  }
682  else
683  partonFlavour = flavourParton->pdgId();
684 
685  // if enabled, check for conflicts between hadron- and parton-based flavours and give priority to the hadron-based flavour
687  {
688  if( hadronFlavour==0 && (std::abs(partonFlavour)==4 || std::abs(partonFlavour)==5) )
689  partonFlavour = ( hardestLightParton.isNonnull() ? hardestLightParton->pdgId() : 0 );
690  else if( hadronFlavour!=0 && std::abs(partonFlavour)!=hadronFlavour )
691  partonFlavour = hadronFlavour;
692  }
693 }
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().

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

Definition at line 180 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering(), and 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().

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

Definition at line 185 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering(), and 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
edm::EDGetTokenT<edm::View<reco::Jet> > JetFlavourClustering::subjetsToken_
private

Definition at line 181 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering(), and produce().

const bool JetFlavourClustering::useLeptons_
private

Definition at line 194 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering(), and produce().

const bool JetFlavourClustering::useSubjets_
private

Definition at line 193 of file JetFlavourClustering.cc.

Referenced by JetFlavourClustering(), and produce().