CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
JetPartonMatcher Class Reference
Inheritance diagram for JetPartonMatcher:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 JetPartonMatcher (const edm::ParameterSet &)
 
 ~JetPartonMatcher ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- 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 ()
 

Private Member Functions

int fillAlgoritDefinition (const Jet &)
 
int fillPhysicsDefinition (const Jet &)
 
virtual void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

double coneSizeToAssociate
 
bool doPriority
 
edm::EDGetTokenT< edm::View
< reco::Jet > > 
m_jetsSrcToken
 
edm::EDGetTokenT
< GenParticleRefVector
m_ParticleSrcToken
 
Handle< GenParticleRefVectorparticles
 
bool physDefinition
 
vector< int > priorityList
 
int theHardest
 
int theHeaviest
 
int theNearest2
 
int theNearest3
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- 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

Definition at line 94 of file JetPartonMatcher.cc.

Constructor & Destructor Documentation

JetPartonMatcher::JetPartonMatcher ( const edm::ParameterSet iConfig)

Definition at line 123 of file JetPartonMatcher.cc.

References coneSizeToAssociate, doPriority, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), m_jetsSrcToken, m_ParticleSrcToken, and priorityList.

123  :
124  theHeaviest(0),
125  theNearest2(0),
126  theNearest3(0),
127  theHardest(0)
128 {
129  produces<JetMatchedPartonsCollection>();
130  m_jetsSrcToken = consumes<edm::View <reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
131  m_ParticleSrcToken = consumes<GenParticleRefVector>(iConfig.getParameter<edm::InputTag>("partons"));
132  coneSizeToAssociate = iConfig.getParameter<double>("coneSizeToAssociate");
133  if ( iConfig.exists("doPriority") ) {
134  doPriority = iConfig.getParameter<bool>("doPriority");
135  priorityList = iConfig.getParameter<vector<int> >("priorityList");
136  } else {
137  doPriority = false;
138  priorityList.clear();
139  }
140 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< GenParticleRefVector > m_ParticleSrcToken
vector< int > priorityList
edm::EDGetTokenT< edm::View< reco::Jet > > m_jetsSrcToken
JetPartonMatcher::~JetPartonMatcher ( )

Definition at line 144 of file JetPartonMatcher.cc.

145 {
146 }

Member Function Documentation

int JetPartonMatcher::fillAlgoritDefinition ( const Jet theJet)
private

Definition at line 216 of file JetPartonMatcher.cc.

References funct::abs(), coneSizeToAssociate, HLT_25ns14e33_v1_cff::DeltaR, doPriority, spr::find(), visualization-live-secondInstance_cfg::m, reco::Candidate::p4(), reco::LeafCandidate::p4(), particles, reco::Candidate::pdgId(), priorityList, reco::Candidate::pt(), reco::Candidate::status(), theHardest, theHeaviest, and theNearest2.

Referenced by produce().

216  {
217 
218  int tempParticle = -1;
219  int tempPartonHighestPt = -1;
220  int tempNearest = -1;
221  float maxPt = 0;
222  float minDr = 1000;
223  bool foundPriority = false;
224 
225  // Loop over the particles in question until we find a priority
226  // "hit", or if we find none in the priority list (or we don't want
227  // to consider priority), then loop through all particles and fill
228  // standard definition.
229 
230  //
231  // Matching:
232  //
233  // 1) First try to match by hand. The "priority list" is given
234  // by the user. The algorithm finds any particles in that
235  // "priority list" that are within the specified cone size.
236  // If it finds one, it counts the object as associating to that
237  // particle.
238  // NOTE! The objects in the priority list are given in order
239  // of priority. So if the user specifies:
240  // 6, 24, 21,
241  // then it will first look for top quarks, then W bosons, then gluons.
242  // 2) If no priority items are found, do the default "standard"
243  // matching.
244  for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
245  const Candidate & aParton = *(particles->at(m).get());
246 
247  // "Priority" behavoir:
248  // Associate to the first particle found in the priority list, regardless
249  // of delta R.
250  if ( doPriority ) {
251  int ipdgId = aParton.pdgId();
252  vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
253  // Check to see if this particle is in our priority list
254  if ( ipriority != priorityList.end() ) {
255  // This particle is on our priority list. If it matches,
256  // we break, since we take in order of priority, not deltaR
257  double dist = DeltaR( theJet.p4(), aParton.p4() );
258  if( dist <= coneSizeToAssociate ) {
259  tempParticle = m;
260  foundPriority = true;
261  }
262  }
263  }
264  // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
265  // is turned off.
266  else {
267  foundPriority = false;
268  }
269 
270  // Default behavior:
271  // Look for partons before the color string to associate.
272  // Do this if we don't want to do priority matching, or if
273  // we didn't find a priority object.
274 
275  if( !foundPriority && aParton.status() != 3 ) {
276  double dist = DeltaR( theJet.p4(), aParton.p4() );
277  if( dist <= coneSizeToAssociate ) {
278  if( dist < minDr ) {
279  minDr = dist;
280  tempNearest = m;
281  }
282  if( tempParticle == -1 && ( abs( aParton.pdgId() ) == 4 ) ) tempParticle = m;
283  if( abs( aParton.pdgId() ) == 5 ) tempParticle = m;
284  if( aParton.pt() > maxPt ) {
285  maxPt = aParton.pt();
286  tempPartonHighestPt = m;
287  }
288  }
289  }
290  }
291 
292  if ( foundPriority ) {
293  theHeaviest = tempParticle; // The priority-matched particle
294  theHardest = -1; // set the rest to -1
295  theNearest2 = -1; // " "
296  } else {
297  theHeaviest = tempParticle;
298  theHardest = tempPartonHighestPt;
299  theNearest2 = tempNearest;
300  if ( tempParticle == -1 ) tempParticle = tempPartonHighestPt;
301  }
302  return tempParticle;
303 }
virtual double pt() const =0
transverse momentum
virtual int status() const =0
status word
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
Handle< GenParticleRefVector > particles
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int pdgId() const =0
PDG identifier.
vector< int > priorityList
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
int JetPartonMatcher::fillPhysicsDefinition ( const Jet theJet)
private

Definition at line 327 of file JetPartonMatcher.cc.

References funct::abs(), coneSizeToAssociate, HLT_25ns14e33_v1_cff::DeltaR, doPriority, spr::find(), reco::flavour(), visualization-live-secondInstance_cfg::m, reco::Candidate::numberOfDaughters(), reco::Candidate::p4(), reco::LeafCandidate::p4(), particles, benchmark_cfg::pdgId, reco::Candidate::pdgId(), priorityList, reco::Candidate::status(), theHardest, theHeaviest, theNearest2, and theNearest3.

Referenced by produce().

327  {
328 
329  float TheBiggerConeSize = 0.7; // In HepMC it's 0.3 --> it's a mistake: value has to be 0.7
330  int tempParticle = -1;
331  int nInTheCone = 0;
332  int tempNearest = -1;
333  float minDr = 1000;
334  bool foundPriority = false;
335 
336  vector<const reco::Candidate *> theContaminations;
337  theContaminations.clear();
338 
339  // Loop over the particles in question until we find a priority
340  // "hit", or if we find none in the priority list (or we don't want
341  // to consider priority), then loop through all particles and fill
342  // standard definition.
343 
344  //
345  // Matching:
346  //
347  // 1) First try to match by hand. The "priority list" is given
348  // by the user. The algorithm finds any particles in that
349  // "priority list" that are within the specified cone size.
350  // If it finds one, it counts the object as associating to that
351  // particle.
352  // NOTE! The objects in the priority list are given in order
353  // of priority. So if the user specifies:
354  // 6, 24, 21,
355  // then it will first look for top quarks, then W bosons, then gluons.
356  // 2) If no priority items are found, do the default "standard"
357  // matching.
358  for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
359 
360  const Candidate & aParticle = *(particles->at(m).get());
361 
362  // "Priority" behavoir:
363  // Associate to the first particle found in the priority list, regardless
364  // of delta R.
365  if ( doPriority ) {
366  int ipdgId = aParticle.pdgId();
367  vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
368  // Check to see if this particle is in our priority list
369  if ( ipriority != priorityList.end() ) {
370  // This particle is on our priority list. If it matches,
371  // we break, since we take in order of priority, not deltaR
372  double dist = DeltaR( theJet.p4(), aParticle.p4() );
373  if( dist <= coneSizeToAssociate ) {
374  tempParticle = m;
375  foundPriority = true;
376  }
377  }
378  }
379  // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
380  // is turned off.
381  else {
382  foundPriority = false;
383  }
384 
385  // Default behavior:
386  if ( !foundPriority ) {
387 
388  // skipping all particle but udscbg (is this correct/enough?!?!)
389  bool isAParton = false;
390  int flavour = abs(aParticle.pdgId());
391  if(flavour == 1 ||
392  flavour == 2 ||
393  flavour == 3 ||
394  flavour == 4 ||
395  flavour == 5 ||
396  flavour == 21 ) isAParton = true;
397  if(!isAParton) continue;
398  double dist = DeltaR( theJet.p4(), aParticle.p4() );
399  if( aParticle.status() == 3 && dist < minDr ) {
400  minDr = dist;
401  tempNearest = m;
402  }
403  if( aParticle.status() == 3 && dist <= coneSizeToAssociate ) {
404  //cout << "particle in small cone=" << aParticle.pdgId() << endl;
405  tempParticle = m;
406  nInTheCone++;
407  }
408  // Look for heavy partons in TheBiggerConeSize now
409  if( aParticle.numberOfDaughters() > 0 && aParticle.status() != 3 ) {
410  if( flavour == 1 ||
411  flavour == 2 ||
412  flavour == 3 ||
413  flavour == 21 ) continue;
414  if( dist < TheBiggerConeSize ) theContaminations.push_back( &aParticle );
415  }
416  }
417  }
418 
419 
420  // Here's the default behavior for assignment if there is no priority.
421  if ( !foundPriority ) {
422  theNearest3 = tempNearest;
423 
424  if(nInTheCone != 1) return -1; // rejected --> only one initialParton requested
425  if(theContaminations.size() == 0 ) return tempParticle; //no contamination
426  int initialPartonFlavour = abs( (particles->at(tempParticle).get()) ->pdgId() );
427 
428  vector<const Candidate *>::const_iterator itCont = theContaminations.begin();
429  for( ; itCont != theContaminations.end(); itCont++ ) {
430  int contaminatingFlavour = abs( (*itCont)->pdgId() );
431  if( (*itCont)->numberOfMothers()>0 && (*itCont)->mother(0) == particles->at(tempParticle).get() ) continue; // mother is the initialParton --> OK
432  if( initialPartonFlavour == 4 ) {
433  if( contaminatingFlavour == 4 ) continue; // keep association --> the initialParton is a c --> the contaminated parton is a c
434  tempParticle = -1; // all the other cases reject!
435  return tempParticle;
436  }
437  }
438  }
439  // If there is priority, then just set the heaviest to priority, the rest are -1.
440  else {
441  theHeaviest = tempParticle; // Set the heaviest to tempParticle
442  theNearest2 = -1; // Set the rest to -1
443  theNearest3 = -1; // " "
444  theHardest = -1; // " "
445  }
446 
447  return tempParticle;
448 }
virtual int status() const =0
status word
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
virtual size_type numberOfDaughters() const =0
number of daughters
Handle< GenParticleRefVector > particles
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int pdgId() const =0
PDG identifier.
vector< int > priorityList
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
void JetPartonMatcher::produce ( edm::Event iEvent,
const edm::EventSetup iEs 
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 150 of file JetPartonMatcher.cc.

References reco::LeafCandidate::eta(), fillAlgoritDefinition(), fillPhysicsDefinition(), edm::Event::getByToken(), j, visualization-live-secondInstance_cfg::m, m_jetsSrcToken, m_ParticleSrcToken, particles, reco::LeafCandidate::pdgId(), reco::LeafCandidate::phi(), reco::LeafCandidate::pt(), edm::Event::put(), reco::LeafCandidate::status(), theHeaviest, theNearest2, and theNearest3.

151 {
153  iEvent.getByToken(m_jetsSrcToken, jets_h );
155 
156  edm::LogVerbatim("JetPartonMatcher") << "=== Partons size:" << particles->size();
157 
158  for( size_t m = 0; m != particles->size(); ++ m ) {
159  const GenParticle & aParton = *(particles->at(m).get());
160  edm::LogVerbatim("JetPartonMatcher") << aParton.status() << " " <<
161  aParton.pdgId() << " " <<
162  aParton.pt() << " " <<
163  aParton.eta() << " " <<
164  aParton.phi() << endl;
165  }
166 
167  auto_ptr<reco::JetMatchedPartonsCollection> jetMatchedPartons( new JetMatchedPartonsCollection(reco::JetRefBaseProd(jets_h)));
168 
169  for (size_t j = 0; j < jets_h->size(); j++) {
170 
171  const int theMappedPartonAlg = fillAlgoritDefinition( (*jets_h)[j] );
172  const int theMappedPartonPhy = fillPhysicsDefinition( (*jets_h)[j] );
173 
174  GenParticleRef pHV;
175  GenParticleRef pN2;
176  GenParticleRef pN3;
177  GenParticleRef pPH;
178  GenParticleRef pAL;
179 
180  if(theHeaviest>=0) pHV = particles->at( theHeaviest );
181  if(theNearest2>=0) pN2 = particles->at( theNearest2 );
182  if(theNearest3>=0) pN3 = particles->at( theNearest3 );
183  if(theMappedPartonPhy>=0) pPH = particles->at( theMappedPartonPhy );
184  if(theMappedPartonAlg>=0) pAL = particles->at( theMappedPartonAlg );
185 
186  (*jetMatchedPartons)[jets_h->refAt(j)]=MatchedPartons(pHV,pN2,pN3,pPH,pAL);
187  }
188 
189  iEvent.put( jetMatchedPartons );
190 
191 }
virtual int pdgId() const
PDG identifier.
int fillPhysicsDefinition(const Jet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual int status() const
status word
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
int fillAlgoritDefinition(const Jet &)
edm::EDGetTokenT< GenParticleRefVector > m_ParticleSrcToken
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
Handle< GenParticleRefVector > particles
int j
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< edm::View< reco::Jet > > m_jetsSrcToken
virtual double phi() const
momentum azimuthal angle

Member Data Documentation

double JetPartonMatcher::coneSizeToAssociate
private
bool JetPartonMatcher::doPriority
private
edm::EDGetTokenT<edm::View <reco::Jet> > JetPartonMatcher::m_jetsSrcToken
private

Definition at line 112 of file JetPartonMatcher.cc.

Referenced by JetPartonMatcher(), and produce().

edm::EDGetTokenT<GenParticleRefVector> JetPartonMatcher::m_ParticleSrcToken
private

Definition at line 113 of file JetPartonMatcher.cc.

Referenced by JetPartonMatcher(), and produce().

Handle<GenParticleRefVector> JetPartonMatcher::particles
private

Definition at line 110 of file JetPartonMatcher.cc.

Referenced by fillAlgoritDefinition(), fillPhysicsDefinition(), and produce().

bool JetPartonMatcher::physDefinition
private

Definition at line 115 of file JetPartonMatcher.cc.

vector<int> JetPartonMatcher::priorityList
private
int JetPartonMatcher::theHardest
private

Definition at line 108 of file JetPartonMatcher.cc.

Referenced by fillAlgoritDefinition(), and fillPhysicsDefinition().

int JetPartonMatcher::theHeaviest
private

Definition at line 105 of file JetPartonMatcher.cc.

Referenced by fillAlgoritDefinition(), fillPhysicsDefinition(), and produce().

int JetPartonMatcher::theNearest2
private

Definition at line 106 of file JetPartonMatcher.cc.

Referenced by fillAlgoritDefinition(), fillPhysicsDefinition(), and produce().

int JetPartonMatcher::theNearest3
private

Definition at line 107 of file JetPartonMatcher.cc.

Referenced by fillPhysicsDefinition(), and produce().