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::ProductRegistryHelper

Public Member Functions

 JetPartonMatcher (const edm::ParameterSet &)
 
 ~JetPartonMatcher ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

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

Private Attributes

double coneSizeToAssociate
 
bool doPriority
 
edm::InputTag m_jetsSrc
 
edm::InputTag m_ParticleSrc
 
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
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- 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::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 94 of file JetPartonMatcher.cc.

Constructor & Destructor Documentation

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

Definition at line 122 of file JetPartonMatcher.cc.

References coneSizeToAssociate, doPriority, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), m_jetsSrc, m_ParticleSrc, and priorityList.

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

Definition at line 143 of file JetPartonMatcher.cc.

144 {
145 }

Member Function Documentation

int JetPartonMatcher::fillAlgoritDefinition ( const Jet theJet)
private

Definition at line 215 of file JetPartonMatcher.cc.

References abs, coneSizeToAssociate, doPriority, spr::find(), 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().

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

Definition at line 326 of file JetPartonMatcher.cc.

References abs, coneSizeToAssociate, doPriority, spr::find(), reco::flavour(), 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().

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

Implements edm::EDProducer.

Definition at line 149 of file JetPartonMatcher.cc.

References reco::LeafCandidate::eta(), fillAlgoritDefinition(), fillPhysicsDefinition(), edm::Event::getByLabel(), j, m, m_jetsSrc, m_ParticleSrc, particles, reco::LeafCandidate::pdgId(), reco::LeafCandidate::phi(), reco::LeafCandidate::pt(), edm::Event::put(), reco::LeafCandidate::status(), theHeaviest, theNearest2, and theNearest3.

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

Member Data Documentation

double JetPartonMatcher::coneSizeToAssociate
private
bool JetPartonMatcher::doPriority
private
edm::InputTag JetPartonMatcher::m_jetsSrc
private

Definition at line 112 of file JetPartonMatcher.cc.

Referenced by JetPartonMatcher(), and produce().

edm::InputTag JetPartonMatcher::m_ParticleSrc
private

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