CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/PhysicsTools/JetMCAlgos/plugins/JetPartonMatcher.cc

Go to the documentation of this file.
00001 // 
00002 // Translation of BTag MCJetFlavour tool to identify real flavour of a jet 
00003 // work with CaloJet objects
00004 // Store Infos by Values in JetFlavour.h
00005 // Author: Attilio  
00006 // Date: 05.10.2007
00007 //
00008 //
00009 // \class JetPartonMatcher
00010 //
00011 // \brief Interface to create a specified "matched" parton to jet match based 
00012 //        on the user interpretation. 
00013 //
00014 // Algorithm for definitions:
00015 //
00016 // If the particle is matched to any item in the priority list (given by user),
00017 // then count that prioritized particle as the match. 
00018 //
00019 // If not resort to default behavior:
00020 //
00021 //
00022 // 1) Algorithmic Definition:
00023 // 
00024 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93) 
00025 // If (parton is within the cone defined by coneSizeToAssociate) then:
00026 //           if (parton is a b)                                   then associatedParton is the b
00027 //      else if (associatedParton =! b and parton is a c)         then associatedParton is the c
00028 //      else if (associatedParton =! b and associatedParton =! c) then associatedParton is the one with the highest pT
00029 // associatedParton can be -1 --> no partons in the cone
00030 // True Flavour of the jet --> flavour of the associatedParton
00031 //
00032 // ToDo: if more than one b(c) in the cone --> the selected parton is not always the one with highest pT
00033 //
00034 //
00035 // 2) Physics Definition:
00036 //
00037 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93)
00038 // if( only one initialParticle within the cone defined by coneSizeToAssociate) associatedInitialParticle is the initialParticle
00039 // TheBiggerConeSize = 0.7 --> it's hard coded!
00040 // if( a parton not coming from associatedInitialParticle in TheBiggerConeSize is a b or a c) reject the association
00041 //
00042 // associatedInitialParticle can be -1 --> no initialParticle in the cone or rejected association
00043 // True Flavour of the jet --> flavour of the associatedInitialParticle
00044 //
00045 //
00046 // Modifications:
00047 //
00048 //     09.03.2008: Sal Rappoccio.
00049 //                 Added capability for "priority" list. 
00050 // 
00051 
00052 //=======================================================================
00053 
00054 // user include files
00055 #include "FWCore/Framework/interface/EDProducer.h"
00056 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00057 #include "FWCore/ParameterSet/interface/ParameterSetfwd.h"
00058 #include "FWCore/Utilities/interface/InputTag.h"
00059  
00060 #include "FWCore/Framework/interface/Event.h"
00061 #include "FWCore/Framework/interface/EventSetup.h"
00062 #include "FWCore/Framework/interface/MakerMacros.h"
00063 #include "FWCore/Framework/interface/ESHandle.h"
00064 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00065 
00066 #include "DataFormats/JetReco/interface/Jet.h"
00067 #include "DataFormats/JetReco/interface/JetCollection.h"
00068 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00069 #include "SimDataFormats/JetMatching/interface/JetFlavour.h"
00070 #include "SimDataFormats/JetMatching/interface/JetFlavourMatching.h"
00071 #include "SimDataFormats/JetMatching/interface/MatchedPartons.h"
00072 #include "SimDataFormats/JetMatching/interface/JetMatchedPartons.h"
00073 
00074 #include "DataFormats/Common/interface/View.h"
00075 #include "DataFormats/Common/interface/Ref.h"
00076 #include "DataFormats/Common/interface/getRef.h" 
00077 //#include "DataFormats/Candidate/interface/Candidate.h"
00078 //#include "DataFormats/Candidate/interface/CandidateFwd.h"
00079 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00080 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00081 
00082 #include <memory>
00083 #include <string>
00084 #include <iostream>
00085 #include <vector>
00086 #include <Math/VectorUtil.h>
00087 #include <TMath.h>
00088 
00089 using namespace std;
00090 using namespace reco;
00091 using namespace edm;
00092 using namespace ROOT::Math::VectorUtil;
00093 
00094 class JetPartonMatcher : public edm::EDProducer 
00095 {
00096   public:
00097     JetPartonMatcher( const edm::ParameterSet & );
00098     ~JetPartonMatcher();
00099 
00100   private:
00101     virtual void produce(edm::Event&, const edm::EventSetup& );
00102 
00103     int fillAlgoritDefinition( const Jet& );
00104     int fillPhysicsDefinition( const Jet& );
00105     int theHeaviest;
00106     int theNearest2;
00107     int theNearest3;
00108     int theHardest;
00109 
00110     Handle <GenParticleRefVector> particles;
00111 
00112     edm::InputTag m_jetsSrc, m_ParticleSrc;
00113     double coneSizeToAssociate;
00114     bool physDefinition;
00115     bool doPriority;
00116     vector<int> priorityList;
00117 
00118 };
00119 
00120 //=========================================================================
00121 
00122 JetPartonMatcher::JetPartonMatcher( const edm::ParameterSet& iConfig ) :
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 }
00140 
00141 //=========================================================================
00142 
00143 JetPartonMatcher::~JetPartonMatcher() 
00144 {
00145 }
00146 
00147 // ------------ method called to produce the data  ------------
00148 
00149 void JetPartonMatcher::produce( Event& iEvent, const EventSetup& iEs ) 
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 }
00191 
00192 //
00193 // Algorithmic Definition: 
00194 // Output: define one associatedParton
00195 // Loop on all particle.
00196 // 
00197 // (Note: This part added by Salvatore Rappoccio)
00198 // [begin]
00199 // If the particle is matched to any item in the priority list (given by user),
00200 // then count that prioritized particle as the match. 
00201 //
00202 // If not resort to default behavior:
00203 // [end]
00204 //
00205 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93) 
00206 // If (parton is within the cone defined by coneSizeToAssociate) then:
00207 //           if (parton is a b)                                   then associatedParton is the b
00208 //      else if (associatedParton =! b and parton is a c)         then associatedParton is the c
00209 //      else if (associatedParton =! b and associatedParton =! c) then associatedParton is the one with the highest pT
00210 // associatedParton can be -1 --> no partons in the cone
00211 // True Flavour of the jet --> flavour of the associatedParton
00212 //
00213 // ToDo: if more than one b(c) in the cone --> the selected parton is not always the one with highest pT
00214 //
00215 int JetPartonMatcher::fillAlgoritDefinition( const Jet& theJet ) {
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     
00274     if( !foundPriority && aParton.status() != 3 ) {
00275       double dist = DeltaR( theJet.p4(), aParton.p4() );
00276       if( dist <= coneSizeToAssociate ) {
00277         if( dist < minDr ) {
00278           minDr = dist;
00279           tempNearest = m;
00280         }
00281         if( tempParticle == -1 && ( abs( aParton.pdgId() ) == 4 )  ) tempParticle = m;
00282         if(                         abs( aParton.pdgId() ) == 5    ) tempParticle = m;
00283         if( aParton.pt() > maxPt ) {
00284           maxPt = aParton.pt();
00285           tempPartonHighestPt = m;
00286         }
00287       }
00288     }
00289   }
00290 
00291   if ( foundPriority ) {
00292     theHeaviest = tempParticle; // The priority-matched particle
00293     theHardest  = -1;  //  set the rest to -1
00294     theNearest2 = -1;  // "                  "
00295   } else {
00296     theHeaviest = tempParticle;
00297     theHardest  = tempPartonHighestPt;
00298     theNearest2 = tempNearest;
00299     if ( tempParticle == -1 ) tempParticle = tempPartonHighestPt;
00300   }
00301   return tempParticle;
00302 }
00303 
00304 //
00305 // Physics Definition: 
00306 // A initialParticle is a particle with status=3
00307 // Output: define one associatedInitialParticle
00308 // Loop on all particles
00309 // 
00310 // (Note: This part added by Salvatore Rappoccio)
00311 // [begin]
00312 // If the particle is matched to any item in the priority list (given by user),
00313 // then count that prioritized particle as the match. 
00314 //
00315 // If not resort to default behavior:
00316 // [end]
00317 //
00318 // A particle is a parton if its daughter is a string(code=92) or a cluster(code=93)
00319 // if( only one initialParticle within the cone defined by coneSizeToAssociate) associatedInitialParticle is the initialParticle
00320 // TheBiggerConeSize = 0.7 --> it's hard coded!
00321 // if( a parton not coming from associatedInitialParticle in TheBiggerConeSize is a b or a c) reject the association
00322 //
00323 // associatedInitialParticle can be -1 --> no initialParticle in the cone or rejected association
00324 // True Flavour of the jet --> flavour of the associatedInitialParticle
00325 //
00326 int JetPartonMatcher::fillPhysicsDefinition( const Jet& theJet ) {
00327 
00328   float TheBiggerConeSize = 0.7; // In HepMC it's 0.3 --> it's a mistake: value has to be 0.7
00329   int tempParticle = -1;
00330   int nInTheCone = 0;
00331   int tempNearest = -1;
00332   float minDr = 1000;
00333   bool foundPriority = false;
00334 
00335   vector<const reco::Candidate *> theContaminations;
00336   theContaminations.clear();
00337 
00338   // Loop over the particles in question until we find a priority
00339   // "hit", or if we find none in the priority list (or we don't want
00340   // to consider priority), then loop through all particles and fill
00341   // standard definition. 
00342 
00343   // 
00344   // Matching:
00345   //
00346   // 1) First try to match by hand. The "priority list" is given
00347   //    by the user. The algorithm finds any particles in that 
00348   //    "priority list" that are within the specified cone size.
00349   //    If it finds one, it counts the object as associating to that
00350   //    particle.
00351   //    NOTE! The objects in the priority list are given in order
00352   //    of priority. So if the user specifies:
00353   //    6, 24, 21,
00354   //    then it will first look for top quarks, then W bosons, then gluons.
00355   // 2) If no priority items are found, do the default "standard"
00356   //    matching. 
00357   for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
00358 
00359     const Candidate & aParticle = *(particles->at(m).get());
00360 
00361     // "Priority" behavoir:
00362     // Associate to the first particle found in the priority list, regardless
00363     // of delta R. 
00364     if ( doPriority ) {
00365       int ipdgId = aParticle.pdgId();
00366       vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
00367       // Check to see if this particle is in our priority list
00368       if ( ipriority != priorityList.end() ) {
00369         // This particle is on our priority list. If it matches,
00370         // we break, since we take in order of priority, not deltaR
00371         double dist = DeltaR( theJet.p4(), aParticle.p4() );
00372         if( dist <= coneSizeToAssociate ) {
00373           tempParticle = m;
00374           foundPriority = true;
00375         }
00376       }
00377     }
00378     // Here we do not want to do priority matching. Ensure the "foundPriority" swtich
00379     // is turned off. 
00380     else {
00381       foundPriority = false;
00382     }
00383 
00384     // Default behavior: 
00385     if ( !foundPriority ) {
00386 
00387       // skipping all particle but udscbg (is this correct/enough?!?!)
00388       bool isAParton = false;
00389       int flavour = abs(aParticle.pdgId());
00390       if(flavour == 1 || 
00391          flavour == 2 ||
00392          flavour == 3 ||
00393          flavour == 4 ||
00394          flavour == 5 ||
00395          flavour == 21 ) isAParton = true;
00396       if(!isAParton) continue;
00397       double dist = DeltaR( theJet.p4(), aParticle.p4() );
00398       if( aParticle.status() == 3 && dist < minDr ) {
00399         minDr = dist;
00400         tempNearest = m;
00401       }
00402       if( aParticle.status() == 3 && dist <= coneSizeToAssociate ) {
00403         //cout << "particle in small cone=" << aParticle.pdgId() << endl;
00404         tempParticle = m;
00405         nInTheCone++;
00406       }
00407       // Look for heavy partons in TheBiggerConeSize now
00408       if( aParticle.numberOfDaughters() > 0 && aParticle.status() != 3 ) {
00409         if( flavour ==  1 ||
00410             flavour ==  2 ||
00411             flavour ==  3 ||
00412             flavour == 21 ) continue;
00413         if( dist < TheBiggerConeSize ) theContaminations.push_back( &aParticle );
00414       }
00415     }
00416   }
00417 
00418 
00419   // Here's the default behavior for assignment if there is no priority. 
00420   if ( !foundPriority ) {
00421     theNearest3 = tempNearest;
00422 
00423     if(nInTheCone != 1) return -1; // rejected --> only one initialParton requested
00424     if(theContaminations.size() == 0 ) return tempParticle; //no contamination
00425     int initialPartonFlavour = abs( (particles->at(tempParticle).get()) ->pdgId() );
00426 
00427     vector<const Candidate *>::const_iterator itCont = theContaminations.begin();
00428     for( ; itCont != theContaminations.end(); itCont++ ) {
00429       int contaminatingFlavour = abs( (*itCont)->pdgId() );
00430       if( (*itCont)->numberOfMothers()>0 && (*itCont)->mother(0) == particles->at(tempParticle).get() ) continue; // mother is the initialParton --> OK
00431       if( initialPartonFlavour == 4 ) {  
00432         if( contaminatingFlavour == 4 ) continue; // keep association --> the initialParton is a c --> the contaminated parton is a c
00433         tempParticle = -1; // all the other cases reject!
00434         return tempParticle;
00435       }
00436     } 
00437   }
00438   // If there is priority, then just set the heaviest to priority, the rest are -1. 
00439   else {
00440     theHeaviest = tempParticle; // Set the heaviest to tempParticle
00441     theNearest2 = -1; //  Set the rest to -1
00442     theNearest3 = -1; // "                  "
00443     theHardest = -1;  // "                  "
00444   }
00445 
00446   return tempParticle;   
00447 }
00448 
00449 //define this as a plug-in
00450 DEFINE_FWK_MODULE(JetPartonMatcher);
00451