CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/TopQuarkAnalysis/TopTools/src/JetPartonMatching.cc

Go to the documentation of this file.
00001 #include "TopQuarkAnalysis/TopTools/interface/JetPartonMatching.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include <Math/VectorUtil.h>
00004 
00005 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<reco::GenJet>& j,
00006                                      const int algorithm = totalMinDist, const bool useMaxDist = true, 
00007                                      const bool useDeltaR = true, const double maxDist = 0.3)
00008   : partons(p), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
00009 {
00010   std::vector<const reco::Candidate*> js;
00011   for(unsigned int i=0; i<j.size(); ++i) js.push_back( &(j[i]) );
00012   jets = js;
00013   calculate();
00014 }
00015 
00016 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<reco::CaloJet>& j,
00017                                      const int algorithm = totalMinDist, const bool useMaxDist = true,
00018                                      const bool useDeltaR = true, const double maxDist = 0.3)
00019   : partons(p), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
00020 {
00021   std::vector<const reco::Candidate*> js;
00022   for(unsigned int i = 0; i < j.size(); ++i) js.push_back( &(j[i]) );
00023   jets = js; 
00024   calculate(); 
00025 }
00026 
00027 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<pat::Jet>& j,
00028                                      const int algorithm = totalMinDist, const bool useMaxDist = true,
00029                                      const bool useDeltaR = true, const double maxDist = 0.3)
00030   : partons(p), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
00031 {
00032   std::vector<const reco::Candidate*> js;
00033   for(unsigned int i = 0; i < j.size(); ++i) js.push_back( &(j[i]) );
00034   jets = js; 
00035   calculate(); 
00036 }
00037 
00038 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<const reco::Candidate*>& j,
00039                                      const int algorithm = totalMinDist, const bool useMaxDist = true,
00040                                      const bool useDeltaR = true, const double maxDist = 0.3)
00041   : partons(p), jets(j), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
00042 {
00043   calculate();
00044 }
00045 
00046 void 
00047 JetPartonMatching::calculate()
00048 {
00049   // use maximal distance between objects 
00050   // in case of unambiguousOnly algorithmm
00051   if(algorithm_==unambiguousOnly) useMaxDist_=true;
00052 
00053   // check if there are empty partons in
00054   // the vector, which happpens if the 
00055   // event is not ttbar or the decay is 
00056   // not as expected
00057   bool emptyParton=false;
00058   for(unsigned int ip=0; ip<partons.size(); ++ip){
00059     if( partons[ip]->pdgId() ==0 ){
00060       emptyParton=true;
00061       break;
00062     }
00063   }
00064 
00065   // switch algorithm, default is to match
00066   // on the minimal sum of the distance 
00067   // (if jets or a parton is empty fill match with blanks)
00068   if( jets.empty() || emptyParton ) {
00069     MatchingCollection dummyMatch;
00070     for(unsigned int ip=0; ip<partons.size(); ++ip)
00071       dummyMatch.push_back( std::make_pair(ip, -1) );
00072     matching.push_back( dummyMatch );
00073   }
00074   else {
00075     switch(algorithm_) {
00076       
00077     case totalMinDist: 
00078       matchingTotalMinDist();    
00079       break;
00080       
00081     case minSumDist: 
00082       matchingMinSumDist();
00083       break;
00084       
00085     case ptOrderedMinDist: 
00086       matchingPtOrderedMinDist();
00087       break;
00088       
00089     case unambiguousOnly:
00090       matchingUnambiguousOnly();
00091       break;
00092       
00093     default:
00094       matchingMinSumDist();
00095     }
00096   }
00097 
00098   numberOfUnmatchedPartons.clear();
00099   sumDeltaE .clear();
00100   sumDeltaPt.clear();
00101   sumDeltaR .clear();
00102   for(unsigned int comb = 0; comb < matching.size(); ++comb) {
00103     MatchingCollection match = matching[comb];
00104     std::sort(match.begin(), match.end());
00105     matching[comb] = match;
00106     
00107     int nUnmatchedPartons = partons.size();
00108     for(unsigned int part=0; part<partons.size(); ++part)
00109       if(getMatchForParton(part,comb)>=0) --nUnmatchedPartons;
00110 
00111     double sumDE  = -999.;
00112     double sumDPt = -999.;
00113     double sumDR  = -999.;
00114     if(nUnmatchedPartons==0){
00115       sumDE  = 0;
00116       sumDPt = 0;
00117       sumDR  = 0;
00118       for(unsigned int i=0; i<match.size(); ++i){
00119         sumDE  += fabs(partons[match[i].first]->energy() - jets[match[i].second]->energy());
00120         sumDPt += fabs(partons[match[i].first]->pt()     - jets[match[i].second]->pt());
00121         sumDR  += distance(partons[match[i].first]->p4(), jets[match[i].second]->p4());
00122       }
00123     }
00124 
00125     numberOfUnmatchedPartons.push_back( nUnmatchedPartons );
00126     sumDeltaE .push_back( sumDE  );
00127     sumDeltaPt.push_back( sumDPt );
00128     sumDeltaR .push_back( sumDR  );
00129   }
00130 }
00131 
00132 double 
00133 JetPartonMatching::distance(const math::XYZTLorentzVector& v1, const math::XYZTLorentzVector& v2)
00134 {
00135   // calculate the distance between two lorentz vectors 
00136   // using DeltaR(eta, phi) or normal space angle(theta, phi)
00137   if(useDeltaR_) return ROOT::Math::VectorUtil::DeltaR(v1, v2);
00138   return ROOT::Math::VectorUtil::Angle(v1, v2);
00139 }
00140 
00141 void 
00142 JetPartonMatching::matchingTotalMinDist()
00143 {
00144   // match parton to jet with shortest distance
00145   // starting with the shortest distance available
00146   // apply some outlier rejection if desired
00147 
00148   // prepare vector of pairs with distances between
00149   // all partons to all jets in the input vectors
00150   std::vector< std::pair<double, unsigned int> > distances;
00151   for(unsigned int ip=0; ip<partons.size(); ++ip){
00152     for(unsigned int ij=0; ij<jets.size(); ++ij){ 
00153       double dist = distance(jets[ij]->p4(), partons[ip]->p4());
00154       distances.push_back(std::pair<double, unsigned int>(dist, ip*jets.size()+ij));
00155     }
00156   }
00157   std::sort(distances.begin(), distances.end());
00158 
00159   MatchingCollection match;
00160 
00161   while(match.size() < partons.size()){
00162     unsigned int partonIndex = distances[0].second/jets.size();
00163     int jetIndex = distances[0].second-jets.size()*partonIndex;
00164     
00165     // use primitive outlier rejection if desired
00166     if(useMaxDist_&& distances[0].first>maxDist_) jetIndex = -1;
00167 
00168     // prevent underflow in case of too few jets
00169     if( distances.empty() )
00170       match.push_back(std::make_pair(partonIndex, -1));
00171     else
00172       match.push_back(std::make_pair(partonIndex, jetIndex));
00173     
00174     // remove all values for the matched parton 
00175     // and the matched jet
00176     for(unsigned int a=0; a<distances.size(); ++a){
00177       unsigned int pIndex = distances[a].second/jets.size();
00178       int jIndex = distances[a].second-jets.size()*pIndex;
00179       if((pIndex == partonIndex) || (jIndex == jetIndex)){
00180         distances.erase(distances.begin()+a, distances.begin()+a+1); 
00181         --a;
00182       }
00183     }
00184   }
00185 
00186   matching.clear();
00187   matching.push_back( match );
00188   return;
00189 }
00190 
00191 void 
00192 JetPartonMatching::minSumDist_recursion(const unsigned int ip, std::vector<unsigned int> & jetIndices,
00193                                         std::vector<bool> & usedJets, std::vector<std::pair<double, MatchingCollection> > & distMatchVec)
00194 {
00195   // build up jet combinations recursively
00196   if(ip<partons.size()){
00197     for(unsigned int ij=0; ij<jets.size(); ++ij){
00198       if(usedJets[ij]) continue;
00199       usedJets[ij] = true;
00200       jetIndices[ip] = ij;
00201       minSumDist_recursion(ip+1, jetIndices, usedJets, distMatchVec);
00202       usedJets[ij] = false;
00203     }
00204     return;
00205   }
00206 
00207   // calculate sumDist for each completed combination
00208   double sumDist = 0;
00209   MatchingCollection match;
00210   for(unsigned int ip=0; ip<partons.size(); ++ip){
00211     double dist  = distance(partons[ip]->p4(), jets[jetIndices[ip]]->p4());
00212     if(useMaxDist_ && dist > maxDist_) return; // outlier rejection
00213     sumDist += distance(partons[ip]->p4(), jets[jetIndices[ip]]->p4());
00214     match.push_back(std::make_pair(ip, jetIndices[ip]));
00215   }
00216 
00217   distMatchVec.push_back( std::make_pair(sumDist, match)  );
00218   return;
00219 }
00220 
00221 void JetPartonMatching::matchingMinSumDist()
00222 {
00223   // match partons to jets with minimal sum of
00224   // the distances between all partons and jets
00225 
00226   std::vector<std::pair<double, MatchingCollection> > distMatchVec;
00227 
00228   std::vector<bool> usedJets;
00229   for(unsigned int i=0; i<jets.size(); ++i){
00230     usedJets.push_back(false);
00231   }
00232 
00233   std::vector<unsigned int> jetIndices;
00234   jetIndices.reserve(partons.size());
00235 
00236   minSumDist_recursion(0, jetIndices, usedJets, distMatchVec);
00237 
00238   std::sort(distMatchVec.begin(), distMatchVec.end());
00239 
00240   matching.clear();
00241 
00242   if(distMatchVec.empty()) {
00243     MatchingCollection dummyMatch;
00244     for(unsigned int ip=0; ip<partons.size(); ++ip)
00245       dummyMatch.push_back(std::make_pair(ip, -1));
00246     matching.push_back( dummyMatch );
00247   }
00248   else
00249     for(unsigned int i=0; i<distMatchVec.size(); ++i)
00250       matching.push_back( distMatchVec[i].second );
00251 
00252   return;
00253 }
00254 
00255 void 
00256 JetPartonMatching::matchingPtOrderedMinDist()
00257 {
00258   // match partons to jets with minimal sum of
00259   // the distances between all partons and jets
00260   // order partons in pt first
00261   std::vector<std::pair <double, unsigned int> > ptOrderedPartons;
00262 
00263   for(unsigned int ip=0; ip<partons.size(); ++ip)
00264     ptOrderedPartons.push_back(std::make_pair(partons[ip]->pt(), ip));
00265 
00266   std::sort(ptOrderedPartons.begin(), ptOrderedPartons.end());
00267   std::reverse(ptOrderedPartons.begin(), ptOrderedPartons.end());
00268 
00269   std::vector<unsigned int> jetIndices;
00270   for(unsigned int ij=0; ij<jets.size(); ++ij) jetIndices.push_back(ij);
00271 
00272   MatchingCollection match;
00273 
00274   for(unsigned int ip=0; ip<ptOrderedPartons.size(); ++ip){
00275     double minDist = 999.;
00276     int ijMin = -1;
00277 
00278     for(unsigned int ij=0; ij<jetIndices.size(); ++ij){
00279       double dist = distance(partons[ptOrderedPartons[ip].second]->p4(), jets[jetIndices[ij]]->p4());
00280       if(dist < minDist){
00281         if(!useMaxDist_ || dist <= maxDist_){
00282           minDist = dist;
00283           ijMin = ij;
00284         }
00285       }
00286     }
00287     
00288     if(ijMin >= 0){
00289       match.push_back( std::make_pair(ptOrderedPartons[ip].second, jetIndices[ijMin]) );
00290       jetIndices.erase(jetIndices.begin() + ijMin, jetIndices.begin() + ijMin + 1);
00291     }
00292     else
00293       match.push_back( std::make_pair(ptOrderedPartons[ip].second, -1) );
00294   }
00295 
00296   matching.clear();
00297   matching.push_back( match );
00298   return;
00299 }
00300 
00301 void 
00302 JetPartonMatching::matchingUnambiguousOnly()
00303 {
00304   // match partons to jets, only accept event 
00305   // if there are no ambiguities
00306   std::vector<bool> jetMatched;
00307   for(unsigned int ij=0; ij<jets.size(); ++ij) jetMatched.push_back(false);
00308 
00309   MatchingCollection match;
00310   
00311   for(unsigned int ip=0; ip<partons.size(); ++ip){
00312     int iMatch = -1;
00313     for(unsigned int ij=0; ij<jets.size(); ++ij){
00314       double dist = distance(partons[ip]->p4(), jets[ij]->p4());
00315       if(dist <= maxDist_){
00316         if(!jetMatched[ij]){ // new match for jet
00317           jetMatched[ij] = true;
00318           if(iMatch == -1) // new match for parton and jet
00319             iMatch = ij;
00320           else // parton already matched: ambiguity!
00321             iMatch = -2;
00322         }
00323         else // jet already matched: ambiguity!
00324           iMatch = -2;
00325       }
00326     }
00327     match.push_back(std::make_pair(ip, iMatch));
00328   }
00329 
00330   matching.clear();
00331   matching.push_back( match );
00332   return;
00333 }
00334 
00335 int
00336 JetPartonMatching::getMatchForParton(const unsigned int part, const unsigned int comb)
00337 {
00338   // return index of the matched jet for a given parton
00339   // (if arguments for parton index and combinatoric index
00340   // are in the valid range)
00341   if(comb >= matching.size()) return -9;
00342   if(part >= matching[comb].size()) return -9;
00343   return (matching[comb])[part].second;
00344 }
00345 
00346 std::vector<int>
00347 JetPartonMatching::getMatchesForPartons(const unsigned int comb)
00348 {
00349   // return a vector with the indices of the matched jets
00350   // (ordered according to the vector of partons)
00351   std::vector<int> jetIndices;
00352   for(unsigned int part=0; part<partons.size(); ++part)
00353     jetIndices.push_back( getMatchForParton(part, comb) );
00354   return jetIndices;
00355 }
00356 
00357 double 
00358 JetPartonMatching::getDistanceForParton(const unsigned int part, const unsigned int comb)
00359 {
00360   // get the distance between parton and its best matched jet
00361   if(getMatchForParton(part, comb) < 0) return -999.;
00362   return distance( jets[getMatchForParton(part,comb)]->p4(), partons[part]->p4() );
00363 }
00364 
00365 double  
00366 JetPartonMatching::getSumDistances(const unsigned int comb)
00367 {
00368   // get sum of distances between partons and matched jets
00369   double sumDists = 0.;
00370   for(unsigned int part=0; part<partons.size(); ++part){
00371     double dist = getDistanceForParton(part, comb);
00372     if(dist < 0.) return -999.;
00373     sumDists += dist;
00374   }
00375   return sumDists;
00376 }
00377 
00378 void
00379 JetPartonMatching::print()
00380 {
00381   //report using MessageLogger
00382   edm::LogInfo log("JetPartonMatching");
00383   log << "++++++++++++++++++++++++++++++++++++++++++++++ \n";
00384   log << " algorithm : ";
00385   switch(algorithm_) {
00386   case totalMinDist     : log << "totalMinDist    "; break;
00387   case minSumDist       : log << "minSumDist      "; break;
00388   case ptOrderedMinDist : log << "ptOrderedMinDist"; break;
00389   case unambiguousOnly  : log << "unambiguousOnly "; break;
00390   default               : log << "UNKNOWN         ";
00391   }
00392   log << "\n";
00393   log << " useDeltaR : ";
00394   switch(useDeltaR_) {
00395   case false : log << "false"; break;
00396   case true  : log << "true ";
00397   }
00398   log << "\n";
00399   log << " useMaxDist: ";
00400   switch(useMaxDist_) {
00401   case false : log << "false"; break;
00402   case true  : log << "true ";
00403   }
00404   log << "      maxDist: " << maxDist_ << "\n";
00405   log << " number of partons / jets: " << partons.size() << " / " << jets.size() << "\n";
00406   log << " number of available combinations: " << getNumberOfAvailableCombinations() << "\n";
00407   for(unsigned int comb = 0; comb < matching.size(); ++comb) {
00408     log << " -------------------------------------------- \n";
00409     log << " ind. of matched jets:";
00410     for(unsigned int part = 0; part < partons.size(); ++part)
00411       log << std::setw(4) << getMatchForParton(part, comb);
00412     log << "\n";
00413     log << " sumDeltaR             : " << getSumDeltaR(comb) << "\n";
00414     log << " sumDeltaPt / sumDeltaE: " << getSumDeltaPt(comb) << " / "  << getSumDeltaE(comb);
00415     log << "\n";
00416   }
00417   log << "++++++++++++++++++++++++++++++++++++++++++++++";
00418 }