CMS 3D CMS Logo

JetPartonMatcher Class Reference

Inheritance diagram for JetPartonMatcher:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 JetPartonMatcher (const edm::ParameterSet &)
 ~JetPartonMatcher ()

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< intpriorityList
int theHardest
int theHeaviest
int theNearest2
int theNearest3


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.

00122                                                                    :
00123   theHeaviest(0),
00124   theNearest2(0),
00125   theNearest3(0),
00126   theHardest(0)
00127 {
00128     produces<JetMatchedPartonsCollection>();
00129     m_jetsSrc           = iConfig.getParameter<edm::InputTag>("jets");
00130     m_ParticleSrc       = iConfig.getParameter<edm::InputTag>("partons");
00131     coneSizeToAssociate = iConfig.getParameter<double>("coneSizeToAssociate");
00132     if ( iConfig.exists("doPriority") ) {
00133       doPriority = iConfig.getParameter<bool>("doPriority");
00134       priorityList = iConfig.getParameter<vector<int> >("priorityList");
00135     } else {
00136       doPriority = false;
00137       priorityList.clear();
00138     }
00139 }

JetPartonMatcher::~JetPartonMatcher (  ) 

Definition at line 143 of file JetPartonMatcher.cc.

00144 {
00145 }


Member Function Documentation

int JetPartonMatcher::fillAlgoritDefinition ( const Jet theJet  )  [private]

Definition at line 215 of file JetPartonMatcher.cc.

References funct::abs(), coneSizeToAssociate, reco::Candidate::daughter(), dist(), doPriority, find(), m, reco::Candidate::numberOfDaughters(), reco::Particle::p4(), particles, reco::Particle::pdgId(), priorityList, reco::Particle::pt(), theHardest, theHeaviest, and theNearest2.

Referenced by produce().

00215                                                                {
00216 
00217   int tempParticle = -1;
00218   int tempPartonHighestPt = -1;
00219   int tempNearest = -1;
00220   float maxPt = 0;
00221   float minDr = 1000;
00222   bool foundPriority = false;
00223 
00224   // Loop over the particles in question until we find a priority
00225   // "hit", or if we find none in the priority list (or we don't want
00226   // to consider priority), then loop through all particles and fill
00227   // standard definition. 
00228 
00229   // 
00230   // Matching:
00231   //
00232   // 1) First try to match by hand. The "priority list" is given
00233   //    by the user. The algorithm finds any particles in that 
00234   //    "priority list" that are within the specified cone size.
00235   //    If it finds one, it counts the object as associating to that
00236   //    particle.
00237   //    NOTE! The objects in the priority list are given in order
00238   //    of priority. So if the user specifies:
00239   //    6, 24, 21,
00240   //    then it will first look for top quarks, then W bosons, then gluons.
00241   // 2) If no priority items are found, do the default "standard"
00242   //    matching. 
00243   for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
00244     const Candidate & aParton = *(particles->at(m).get());
00245 
00246     // "Priority" behavoir:
00247     // Associate to the first particle found in the priority list, regardless
00248     // of delta R. 
00249     if ( doPriority ) {
00250       int ipdgId = aParton.pdgId();
00251       vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
00252       // Check to see if this particle is in our priority list
00253       if ( ipriority != priorityList.end() ) {
00254         // This particle is on our priority list. If it matches,
00255         // we break, since we take in order of priority, not deltaR
00256         double dist = DeltaR( theJet.p4(), aParton.p4() );
00257         if( dist <= coneSizeToAssociate ) {
00258           tempParticle = m;
00259           foundPriority = true;
00260         }
00261       }
00262     }
00263     // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
00264     // is turned off. 
00265     else {
00266       foundPriority = false;
00267     }
00268 
00269     // Default behavior: 
00270     // Look for partons before the color string to associate. 
00271     // Do this if we don't want to do priority matching, or if
00272     // we didn't find a priority object. 
00273     if( !foundPriority && aParton.numberOfDaughters() > 0  && ( aParton.daughter(0)->pdgId() == 91 || aParton.daughter(0)->pdgId() == 92 ) ) {
00274       double dist = DeltaR( theJet.p4(), aParton.p4() );
00275       if( dist <= coneSizeToAssociate ) {
00276         if( dist < minDr ) {
00277           minDr = dist;
00278           tempNearest = m;
00279         }
00280         if( tempParticle == -1 && ( abs( aParton.pdgId() ) == 4 )  ) tempParticle = m;
00281         if(                         abs( aParton.pdgId() ) == 5    ) tempParticle = m;
00282         if( aParton.pt() > maxPt ) {
00283           maxPt = aParton.pt();
00284           tempPartonHighestPt = m;
00285         }
00286       }
00287     }
00288   }
00289 
00290   if ( foundPriority ) {
00291     theHeaviest = tempParticle; // The priority-matched particle
00292     theHardest  = -1;  //  set the rest to -1
00293     theNearest2 = -1;  // "                  "
00294   } else {
00295     theHeaviest = tempParticle;
00296     theHardest  = tempPartonHighestPt;
00297     theNearest2 = tempNearest;
00298     if ( tempParticle == -1 ) tempParticle = tempPartonHighestPt;
00299   }
00300   return tempParticle;
00301 }

int JetPartonMatcher::fillPhysicsDefinition ( const Jet theJet  )  [private]

Definition at line 325 of file JetPartonMatcher.cc.

References funct::abs(), coneSizeToAssociate, reco::Candidate::daughter(), dist(), doPriority, find(), flavour, m, reco::Candidate::numberOfDaughters(), reco::Particle::p4(), particles, reco::Particle::pdgId(), priorityList, reco::Particle::status(), theHardest, theHeaviest, theNearest2, and theNearest3.

Referenced by produce().

00325                                                                {
00326 
00327   float TheBiggerConeSize = 0.7; // In HepMC it's 0.3 --> it's a mistake: value has to be 0.7
00328   int tempParticle = -1;
00329   int nInTheCone = 0;
00330   int tempNearest = -1;
00331   float minDr = 1000;
00332   bool foundPriority = false;
00333 
00334   vector<const reco::Candidate *> theContaminations;
00335   theContaminations.clear();
00336 
00337   // Loop over the particles in question until we find a priority
00338   // "hit", or if we find none in the priority list (or we don't want
00339   // to consider priority), then loop through all particles and fill
00340   // standard definition. 
00341 
00342   // 
00343   // Matching:
00344   //
00345   // 1) First try to match by hand. The "priority list" is given
00346   //    by the user. The algorithm finds any particles in that 
00347   //    "priority list" that are within the specified cone size.
00348   //    If it finds one, it counts the object as associating to that
00349   //    particle.
00350   //    NOTE! The objects in the priority list are given in order
00351   //    of priority. So if the user specifies:
00352   //    6, 24, 21,
00353   //    then it will first look for top quarks, then W bosons, then gluons.
00354   // 2) If no priority items are found, do the default "standard"
00355   //    matching. 
00356   for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
00357 
00358     const Candidate & aParticle = *(particles->at(m).get());
00359 
00360     // "Priority" behavoir:
00361     // Associate to the first particle found in the priority list, regardless
00362     // of delta R. 
00363     if ( doPriority ) {
00364       int ipdgId = aParticle.pdgId();
00365       vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
00366       // Check to see if this particle is in our priority list
00367       if ( ipriority != priorityList.end() ) {
00368         // This particle is on our priority list. If it matches,
00369         // we break, since we take in order of priority, not deltaR
00370         double dist = DeltaR( theJet.p4(), aParticle.p4() );
00371         if( dist <= coneSizeToAssociate ) {
00372           tempParticle = m;
00373           foundPriority = true;
00374         }
00375       }
00376     }
00377     // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
00378     // is turned off. 
00379     else {
00380       foundPriority = false;
00381     }
00382 
00383     // Default behavior: 
00384     if ( !foundPriority ) {
00385 
00386       // skipping all particle but udscbg (is this correct/enough?!?!)
00387       bool isAParton = false;
00388       int flavour = abs(aParticle.pdgId());
00389       if(flavour == 1 || 
00390          flavour == 2 ||
00391          flavour == 3 ||
00392          flavour == 4 ||
00393          flavour == 5 ||
00394          flavour == 21 ) isAParton = true;
00395       if(!isAParton) continue;
00396       double dist = DeltaR( theJet.p4(), aParticle.p4() );
00397       if( aParticle.status() == 3 && dist < minDr ) {
00398         minDr = dist;
00399         tempNearest = m;
00400       }
00401       if( aParticle.status() == 3 && dist <= coneSizeToAssociate ) {
00402         //cout << "particle in small cone=" << aParticle.pdgId() << endl;
00403         tempParticle = m;
00404         nInTheCone++;
00405       }
00406       // Look for heavy partons in TheBiggerConeSize now
00407       if( aParticle.numberOfDaughters() > 0  && ( aParticle.daughter(0)->pdgId() == 91 || aParticle.daughter(0)->pdgId() == 92 ) ) {
00408         if( flavour ==  1 ||
00409             flavour ==  2 ||
00410             flavour ==  3 ||
00411             flavour == 21 ) continue;
00412         if( dist < TheBiggerConeSize ) theContaminations.push_back( &aParticle );
00413       }
00414     }
00415   }
00416 
00417 
00418   // Here's the default behavior for assignment if there is no priority. 
00419   if ( !foundPriority ) {
00420     theNearest3 = tempNearest;
00421 
00422     if(nInTheCone != 1) return -1; // rejected --> only one initialParton requested
00423     if(theContaminations.size() == 0 ) return tempParticle; //no contamination
00424     int initialPartonFlavour = abs( (particles->at(tempParticle).get()) ->pdgId() );
00425 
00426     vector<const Candidate *>::const_iterator itCont = theContaminations.begin();
00427     for( ; itCont != theContaminations.end(); itCont++ ) {
00428       int contaminatingFlavour = abs( (*itCont)->pdgId() );
00429       if( (*itCont)->numberOfMothers()>0 && (*itCont)->mother(0) == particles->at(tempParticle).get() ) continue; // mother is the initialParton --> OK
00430       if( initialPartonFlavour == 4 ) {  
00431         if( contaminatingFlavour == 4 ) continue; // keep association --> the initialParton is a c --> the contaminated parton is a c
00432         tempParticle = -1; // all the other cases reject!
00433         return tempParticle;
00434       }
00435     } 
00436   }
00437   // If there is priority, then just set the heaviest to priority, the rest are -1. 
00438   else {
00439     theHeaviest = tempParticle; // Set the heaviest to tempParticle
00440     theNearest2 = -1; //  Set the rest to -1
00441     theNearest3 = -1; // "                  "
00442     theHardest = -1;  // "                  "
00443   }
00444 
00445   return tempParticle;   
00446 }

void JetPartonMatcher::produce ( edm::Event iEvent,
const edm::EventSetup iEs 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 149 of file JetPartonMatcher.cc.

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

00150 {
00151   edm::Handle <edm::View <reco::Jet> > jets_h;
00152   iEvent.getByLabel(m_jetsSrc,     jets_h    );
00153   iEvent.getByLabel(m_ParticleSrc, particles );
00154 
00155   edm::LogVerbatim("JetPartonMatcher") << "=== Partons size:" << particles->size();
00156 
00157   for( size_t m = 0; m != particles->size(); ++ m ) {
00158     const GenParticle & aParton = *(particles->at(m).get());
00159     edm::LogVerbatim("JetPartonMatcher") <<  aParton.status() << " " <<
00160                                              aParton.pdgId()  << " " <<
00161                                              aParton.pt()     << " " << 
00162                                              aParton.eta()    << " " <<
00163                                              aParton.phi()    << endl;
00164   }
00165 
00166   auto_ptr<reco::JetMatchedPartonsCollection> jetMatchedPartons( new JetMatchedPartonsCollection(reco::JetRefBaseProd(jets_h)));
00167   
00168   for (size_t j = 0; j < jets_h->size(); j++) {
00169 
00170     const int theMappedPartonAlg = fillAlgoritDefinition( (*jets_h)[j] );
00171     const int theMappedPartonPhy = fillPhysicsDefinition( (*jets_h)[j] );
00172 
00173     GenParticleRef pHV;
00174     GenParticleRef pN2;
00175     GenParticleRef pN3;
00176     GenParticleRef pPH;
00177     GenParticleRef pAL;
00178 
00179     if(theHeaviest>=0)        pHV = particles->at( theHeaviest        );
00180     if(theNearest2>=0)        pN2 = particles->at( theNearest2        );
00181     if(theNearest3>=0)        pN3 = particles->at( theNearest3        );
00182     if(theMappedPartonPhy>=0) pPH = particles->at( theMappedPartonPhy );
00183     if(theMappedPartonAlg>=0) pAL = particles->at( theMappedPartonAlg );
00184 
00185     (*jetMatchedPartons)[jets_h->refAt(j)]=MatchedPartons(pHV,pN2,pN3,pPH,pAL);
00186   }
00187   
00188   iEvent.put(  jetMatchedPartons );
00189 
00190 }


Member Data Documentation

double JetPartonMatcher::coneSizeToAssociate [private]

Definition at line 113 of file JetPartonMatcher.cc.

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

bool JetPartonMatcher::doPriority [private]

Definition at line 115 of file JetPartonMatcher.cc.

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

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]

Definition at line 116 of file JetPartonMatcher.cc.

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

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().


The documentation for this class was generated from the following file:
Generated on Tue Jun 9 18:26:04 2009 for CMSSW by  doxygen 1.5.4