CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/PhysicsTools/JetMCAlgos/plugins/JetFlavourIdentifier.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 JetFlavourIdentifier
00010 //
00011 // \brief Interface to pull out the proper flavour identifier from a 
00012 //        jet->parton matching collection.
00013 //
00014 // In detail, the definitions are as follows:
00015 //
00016 // Definitions:
00017 // The default behavior is that the "definition" is NULL_DEF,
00018 // so the software will fall back to either "physics" or "algorithmic" definition
00019 // as per the "physDefinition" switch.
00020 //
00021 // If the user specifies "definition", then that definition is taken. However,
00022 // if the requested definition is not defined, the software reverts back to
00023 // either "physics" or "algorithmic" based on the "physDefinition" switch.
00024 // For example, if the user specifies "heaviest" as a flavor ID, and there
00025 // are no bottom, charm, or top quarks in the event that match to the jet,
00026 // then the software will fall back to the "physics" or "algorithmic" definition.
00027 //
00028 // Modifications:
00029 //
00030 //     09.03.2008: Sal Rappoccio.
00031 //                 Added capability for all methods of identification, not just
00032 //                 "physics" or "algorithmic". If the requested method does not exist
00033 //                 (i.e. is unphysical), the "physics" or "algorithmic" definition
00034 //                 is defaulted to. 
00035 
00036 //=======================================================================
00037 
00038 // user include files
00039 #include "FWCore/Framework/interface/EDProducer.h"
00040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00041 #include "FWCore/ParameterSet/interface/ParameterSetfwd.h"
00042 #include "FWCore/Utilities/interface/InputTag.h"
00043  
00044 #include "FWCore/Framework/interface/Event.h"
00045 #include "FWCore/Framework/interface/EventSetup.h"
00046 #include "FWCore/Framework/interface/MakerMacros.h"
00047 #include "FWCore/Framework/interface/ESHandle.h"
00048 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00049 
00050 //#include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00051 
00052 
00053 #include "DataFormats/JetReco/interface/Jet.h"
00054 #include "DataFormats/JetReco/interface/JetCollection.h"
00055 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00056 #include "SimDataFormats/JetMatching/interface/JetFlavour.h"
00057 #include "SimDataFormats/JetMatching/interface/JetFlavourMatching.h"
00058 #include "SimDataFormats/JetMatching/interface/MatchedPartons.h"
00059 #include "SimDataFormats/JetMatching/interface/JetMatchedPartons.h"
00060 
00061 #include "DataFormats/Common/interface/Ref.h"
00062 #include "DataFormats/Candidate/interface/Candidate.h"
00063 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00064 
00065 #include "DataFormats/Math/interface/Point3D.h"
00066 #include "DataFormats/Math/interface/LorentzVector.h"
00067 
00068 #include <memory>
00069 #include <string>
00070 #include <iostream>
00071 #include <vector>
00072 #include <Math/VectorUtil.h>
00073 #include <TMath.h>
00074 
00075 using namespace std;
00076 using namespace reco;
00077 using namespace edm;
00078 using namespace ROOT::Math::VectorUtil;
00079 
00080 namespace reco { namespace modules {
00081 
00082 //--------------------------------------------------------------------------
00083 // 
00084 //--------------------------------------------------------------------------
00085 class JetFlavourIdentifier : public edm::EDProducer 
00086 {
00087   public:
00088   enum DEFINITION_T { PHYSICS=0, ALGO, NEAREST_STATUS2, NEAREST_STATUS3, HEAVIEST, 
00089                       N_DEFINITIONS,
00090                       NULL_DEF};
00091   
00092     JetFlavourIdentifier( const edm::ParameterSet & );
00093     ~JetFlavourIdentifier();
00094 
00095   private:
00096     virtual void produce(edm::Event&, const edm::EventSetup& );
00097 
00098     JetFlavour::Leptons findLeptons(const GenParticleRef &);
00099     std::vector<const reco::Candidate*> findCandidates(const reco::Candidate*, int);
00100     void fillLeptons(const std::vector<const reco::Candidate*> &, JetFlavour::Leptons &, int, int);
00101     static int heaviestFlavour(int);
00102 
00103     Handle<JetMatchedPartonsCollection> theTagByRef;
00104     InputTag sourceByRefer_;
00105     bool physDefinition;
00106     bool leptonInfo_;
00107     DEFINITION_T definition;
00108     math::XYZTLorentzVector thePartonLV;
00109 
00110 };
00111 } }
00112 using reco::modules::JetFlavourIdentifier;
00113 
00114 //=========================================================================
00115 
00116 JetFlavourIdentifier::JetFlavourIdentifier( const edm::ParameterSet& iConfig )
00117 {
00118     produces<JetFlavourMatchingCollection>();
00119     sourceByRefer_ = iConfig.getParameter<InputTag>("srcByReference");
00120     physDefinition = iConfig.getParameter<bool>("physicsDefinition");
00121     leptonInfo_ = iConfig.exists("leptonInfo") ? iConfig.getParameter<bool>("leptonInfo") : true;
00122     // If we have a definition of which parton to identify, use it,
00123     // otherwise we default to the "old" behavior of either "physics" or "algorithmic".
00124     // Furthermore, if the specified definition is not sensible for the given jet,
00125     // then the "physDefinition" switch is used to identify the flavour of the jet. 
00126     if ( iConfig.exists("definition") ) {
00127       definition = static_cast<DEFINITION_T>( iConfig.getParameter<int>("definition") );
00128     } else {
00129       definition = NULL_DEF;
00130     }
00131 }
00132 
00133 //=========================================================================
00134 
00135 JetFlavourIdentifier::~JetFlavourIdentifier() 
00136 {
00137 }
00138 
00139 // ------------ method called to produce the data  ------------
00140 
00141 void JetFlavourIdentifier::produce( Event& iEvent, const EventSetup& iEs ) 
00142 {
00143   // Get the JetMatchedPartons
00144   iEvent.getByLabel (sourceByRefer_, theTagByRef);
00145 
00146   // Create a JetFlavourMatchingCollection
00147   JetFlavourMatchingCollection *jfmc;
00148   if (!theTagByRef->empty()) {
00149     RefToBase<Jet> jj = theTagByRef->begin()->first;
00150     jfmc = new JetFlavourMatchingCollection(RefToBaseProd<Jet>(jj));
00151   } else {
00152     jfmc = new JetFlavourMatchingCollection();
00153   }
00154   auto_ptr<reco::JetFlavourMatchingCollection> jetFlavMatching(jfmc);
00155 
00156   // Loop over the matched partons and see which match. 
00157   for ( JetMatchedPartonsCollection::const_iterator j  = theTagByRef->begin();
00158                                                     j != theTagByRef->end();
00159                                                     j ++ ) {
00160 
00161 
00162     // Consider this match. 
00163     const MatchedPartons aMatch = (*j).second;
00164 
00165     // This will hold the 4-vector, vertex, flavour and the leptonian decays (0: no lepton, xyz: x leptons in layer 2, y in layer 1 and z in layer 0) of the requested object.
00166     math::XYZTLorentzVector thePartonLorentzVector(0,0,0,0);
00167     math::XYZPoint          thePartonVertex(0,0,0);
00168     int                     thePartonFlavour = 0;
00169     JetFlavour::Leptons     theLeptons;
00170 
00171     // get the partons based on which definition to use. 
00172     switch (definition) {
00173       case PHYSICS: {
00174         const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
00175         if (aPartPhy.isNonnull()) {
00176                 thePartonLorentzVector = aPartPhy.get()->p4();
00177                 thePartonVertex        = aPartPhy.get()->vertex();
00178                 thePartonFlavour       = aPartPhy.get()->pdgId();
00179           if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
00180         }
00181         break;
00182       }
00183       case ALGO: {
00184         const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
00185         if (aPartAlg.isNonnull()) {
00186                 thePartonLorentzVector = aPartAlg.get()->p4();
00187                 thePartonVertex        = aPartAlg.get()->vertex();
00188                 thePartonFlavour       = aPartAlg.get()->pdgId();
00189           if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
00190         }
00191         break;
00192       }
00193       case NEAREST_STATUS2 : {
00194         const GenParticleRef aPartN2 = aMatch.nearest_status2();
00195         if (aPartN2.isNonnull()) {
00196                 thePartonLorentzVector = aPartN2.get()->p4();
00197                 thePartonVertex        = aPartN2.get()->vertex();
00198                 thePartonFlavour       = aPartN2.get()->pdgId();
00199           if (leptonInfo_) theLeptons = findLeptons(aPartN2);
00200         }
00201         break;
00202       }
00203       case NEAREST_STATUS3: {
00204         const GenParticleRef aPartN3 = aMatch.nearest_status3();
00205         if (aPartN3.isNonnull()) {
00206                 thePartonLorentzVector = aPartN3.get()->p4();
00207                 thePartonVertex        = aPartN3.get()->vertex();
00208                 thePartonFlavour       = aPartN3.get()->pdgId();
00209           if (leptonInfo_) theLeptons = findLeptons(aPartN3);
00210         }
00211         break;
00212       }
00213       case HEAVIEST: {
00214         const GenParticleRef aPartHeaviest = aMatch.heaviest();
00215         if (aPartHeaviest.isNonnull()) {
00216                 thePartonLorentzVector = aPartHeaviest.get()->p4();
00217                 thePartonVertex        = aPartHeaviest.get()->vertex();
00218                 thePartonFlavour       = aPartHeaviest.get()->pdgId();
00219           if (leptonInfo_) theLeptons = findLeptons(aPartHeaviest);
00220         }
00221         break;
00222       }
00223       // Default case is backwards-compatible
00224       default:{
00225         if (physDefinition) {
00226           const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
00227           if (aPartPhy.isNonnull()) {
00228             thePartonLorentzVector = aPartPhy.get()->p4();
00229             thePartonVertex        = aPartPhy.get()->vertex();
00230             thePartonFlavour       = aPartPhy.get()->pdgId();
00231             if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
00232           }
00233         } else {
00234           const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
00235           if (aPartAlg.isNonnull()) {
00236             thePartonLorentzVector = aPartAlg.get()->p4();
00237             thePartonVertex        = aPartAlg.get()->vertex();
00238             thePartonFlavour       = aPartAlg.get()->pdgId();
00239             if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
00240           }
00241         }
00242       } break;
00243     } // end switch on definition
00244 
00245     // Now make sure we have a match. If the user specified "heaviest", for instance,
00246     // and there is no b- or c-quarks, then fall back to the "physDefinition" switch.
00247     
00248     if (thePartonFlavour == 0) {
00249       if (physDefinition) {
00250         const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
00251         if (aPartPhy.isNonnull()) {
00252           thePartonLorentzVector = aPartPhy.get()->p4();
00253           thePartonVertex        = aPartPhy.get()->vertex();
00254           thePartonFlavour       = aPartPhy.get()->pdgId();
00255           if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
00256         }
00257       } else {
00258         const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
00259         if (aPartAlg.isNonnull()) {
00260           thePartonLorentzVector = aPartAlg.get()->p4();
00261           thePartonVertex        = aPartAlg.get()->vertex();
00262           thePartonFlavour       = aPartAlg.get()->pdgId();
00263           if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
00264         }
00265       }
00266     }
00267 
00268 /*
00269      std::cout << "Leptons of " <<thePartonFlavour << " Jet: " << std::endl;
00270      std::cout << "  electrons: " <<theLeptons.electron << std::endl;
00271      std::cout << "  muons    : " <<theLeptons.muon << std::endl;
00272      std::cout << "  tau      : " <<theLeptons.tau << std::endl;
00273 */
00274     // Add the jet->flavour match to the map. 
00275     (*jetFlavMatching)[(*j).first] = JetFlavour(thePartonLorentzVector, thePartonVertex, thePartonFlavour, theLeptons);
00276   }// end loop over jets
00277 
00278 
00279   // Put the object into the event. 
00280   iEvent.put(  jetFlavMatching );
00281 
00282 }
00283 
00284 JetFlavour::Leptons JetFlavourIdentifier::findLeptons(const GenParticleRef &parton)
00285 {
00286   JetFlavour::Leptons theLeptons;
00287 
00288   thePartonLV = parton->p4();
00289 
00291   const reco::Candidate *mcstring = parton->daughter(0);
00292   int partonFlavour = std::abs(parton->pdgId());
00293 //  std::cout << "parton DeltaR: " << DeltaR(thePartonLV, parton->p4()) << std::endl;
00294 
00296   std::vector<const reco::Candidate*> candidates = findCandidates(mcstring, partonFlavour);
00297 //  std::cout << "Candidates are:" << std::endl;
00298 //  for(unsigned int j = 0; j < candidates.size(); j++) std::cout << "   --> " << candidates[j]->pdgId() << std::endl;
00299 
00301   fillLeptons(candidates, theLeptons, 1, partonFlavour);
00302 
00303   return theLeptons;
00304 }
00305 
00306 std::vector<const reco::Candidate*> JetFlavourIdentifier::findCandidates(const reco::Candidate *cand, int partonFlavour)
00307 {
00308   std::vector<const reco::Candidate*> cands;
00309   
00310   for(unsigned int i = 0; i < cand->numberOfDaughters(); i++) {
00311 /*
00312     std::cout << "DeltaR - " << partonFlavour << " ";
00313     if (DeltaR(thePartonLV, cand->daughter(i)->p4()) > 0.7) std::cout << "(";
00314     std::cout << cand->daughter(i)->pdgId() << ": " << DeltaR(thePartonLV, cand->daughter(i)->p4());
00315     if (DeltaR(thePartonLV, cand->daughter(i)->p4()) > 0.7) std::cout << ")";
00316     std::cout << std::endl;
00317 */
00318 
00319     if (DeltaR(thePartonLV, cand->daughter(i)->p4()) < 0.7) {
00320       int pdgId = std::abs(cand->daughter(i)->pdgId());
00321       int flavour = heaviestFlavour(pdgId);
00322       if (flavour == partonFlavour ||
00323           (flavour >= 10 && partonFlavour >= 10)) {
00324 //        std::cout << "<------- " << std::endl;
00325         std::vector<const reco::Candidate*> newcands = findCandidates(cand->daughter(i), partonFlavour);
00326 //        std::cout << " ------->" << std::endl;
00327         std::copy(newcands.begin(), newcands.end(), std::back_inserter(cands));
00328       }
00329       if (partonFlavour >= 10)
00330         cands.push_back(cand->daughter(i));
00331     }
00332   }
00333 
00334   if (cands.empty() && std::abs(cand->pdgId()) > 110 &&
00335       !(partonFlavour >= 4 && partonFlavour < 10 &&
00336         heaviestFlavour(cand->pdgId()) < 4))
00337     cands.push_back(cand);
00338 
00339   return cands;
00340 }
00341 
00342 void JetFlavourIdentifier::fillLeptons(const std::vector<const reco::Candidate*> &cands, JetFlavour::Leptons &leptons, int rank, int flavour)
00343 {
00344   for(unsigned int j = 0; j < cands.size(); j++) {
00345     for(unsigned int i = 0; i < cands[j]->numberOfDaughters(); i++) {
00346       int pdgId = std::abs(cands[j]->daughter(i)->pdgId());
00347 
00348 //      for(int z = 1; z <= rank; z *= 10) std::cout << " ------ ";
00349 //      std::cout << pdgId << std::endl;
00350 
00352       if (pdgId == 12)
00353         leptons.electron += rank;
00354       else if (pdgId == 14)
00355         leptons.muon += rank;
00356       else if (pdgId == 16)
00357         leptons.tau += rank;
00358       else {
00359         int heaviest = heaviestFlavour(pdgId);
00360         int heaviest_ = heaviest < 10 ? heaviest : 0;
00361         if (!heaviest || (flavour < 4 ? (heaviest_ < 4) : (heaviest >= 4))) {
00362           std::vector<const reco::Candidate*> newcands = findCandidates(cands[j]->daughter(i), heaviest);
00363           if (pdgId <= 110) newcands.push_back(cands[j]->daughter(i));
00364           fillLeptons(newcands, leptons, rank * 10, std::max(heaviest_, flavour));
00365         }
00366       }
00367     }
00368   }
00369 }
00370 
00371 int JetFlavourIdentifier::heaviestFlavour(int pdgId)
00372 {
00373   int flavour = 0;
00374 
00375   pdgId = std::abs(pdgId) % 100000;
00376   if (pdgId > 110) {
00377     while(pdgId % 10 > 0 && pdgId % 10 < 6) {
00378       pdgId /= 10;
00379       if (pdgId % 10 > flavour)
00380         flavour = pdgId % 10;
00381     }
00382   } else
00383     flavour = pdgId;
00384 
00385   return flavour;
00386 }
00387 
00388 
00389 //define this as a plug-in
00390 DEFINE_FWK_MODULE(JetFlavourIdentifier);