CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetPartonMatching.cc
Go to the documentation of this file.
3 #include <Math/VectorUtil.h>
4 
5 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<reco::GenJet>& j,
6  const int algorithm = totalMinDist, const bool useMaxDist = true,
7  const bool useDeltaR = true, const double maxDist = 0.3)
8  : partons(p), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
9 {
10  std::vector<const reco::Candidate*> js;
11  for(unsigned int i=0; i<j.size(); ++i) js.push_back( &(j[i]) );
12  jets = js;
13  calculate();
14 }
15 
16 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<reco::CaloJet>& j,
17  const int algorithm = totalMinDist, const bool useMaxDist = true,
18  const bool useDeltaR = true, const double maxDist = 0.3)
19  : partons(p), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
20 {
21  std::vector<const reco::Candidate*> js;
22  for(unsigned int i = 0; i < j.size(); ++i) js.push_back( &(j[i]) );
23  jets = js;
24  calculate();
25 }
26 
27 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<pat::Jet>& j,
28  const int algorithm = totalMinDist, const bool useMaxDist = true,
29  const bool useDeltaR = true, const double maxDist = 0.3)
30  : partons(p), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
31 {
32  std::vector<const reco::Candidate*> js;
33  for(unsigned int i = 0; i < j.size(); ++i) js.push_back( &(j[i]) );
34  jets = js;
35  calculate();
36 }
37 
38 JetPartonMatching::JetPartonMatching(const std::vector<const reco::Candidate*>& p, const std::vector<const reco::Candidate*>& j,
39  const int algorithm = totalMinDist, const bool useMaxDist = true,
40  const bool useDeltaR = true, const double maxDist = 0.3)
41  : partons(p), jets(j), algorithm_(algorithm), useMaxDist_(useMaxDist), useDeltaR_(useDeltaR), maxDist_(maxDist)
42 {
43  calculate();
44 }
45 
46 void
48 {
49  // use maximal distance between objects
50  // in case of unambiguousOnly algorithmm
52 
53  // check if there are empty partons in
54  // the vector, which happpens if the
55  // event is not ttbar or the decay is
56  // not as expected
57  bool emptyParton=false;
58  for(unsigned int ip=0; ip<partons.size(); ++ip){
59  if( partons[ip]->pdgId() ==0 ){
60  emptyParton=true;
61  break;
62  }
63  }
64 
65  // switch algorithm, default is to match
66  // on the minimal sum of the distance
67  // (if jets or a parton is empty fill match with blanks)
68  if( jets.empty() || emptyParton ) {
69  MatchingCollection dummyMatch;
70  for(unsigned int ip=0; ip<partons.size(); ++ip)
71  dummyMatch.push_back( std::make_pair(ip, -1) );
72  matching.push_back( dummyMatch );
73  }
74  else {
75  switch(algorithm_) {
76 
77  case totalMinDist:
79  break;
80 
81  case minSumDist:
83  break;
84 
85  case ptOrderedMinDist:
87  break;
88 
89  case unambiguousOnly:
91  break;
92 
93  default:
95  }
96  }
97 
99  sumDeltaE .clear();
100  sumDeltaPt.clear();
101  sumDeltaR .clear();
102  for(unsigned int comb = 0; comb < matching.size(); ++comb) {
104  std::sort(match.begin(), match.end());
105  matching[comb] = match;
106 
107  int nUnmatchedPartons = partons.size();
108  for(unsigned int part=0; part<partons.size(); ++part)
109  if(getMatchForParton(part,comb)>=0) --nUnmatchedPartons;
110 
111  double sumDE = -999.;
112  double sumDPt = -999.;
113  double sumDR = -999.;
114  if(nUnmatchedPartons==0){
115  sumDE = 0;
116  sumDPt = 0;
117  sumDR = 0;
118  for(unsigned int i=0; i<match.size(); ++i){
119  sumDE += fabs(partons[match[i].first]->energy() - jets[match[i].second]->energy());
120  sumDPt += fabs(partons[match[i].first]->pt() - jets[match[i].second]->pt());
121  sumDR += distance(partons[match[i].first]->p4(), jets[match[i].second]->p4());
122  }
123  }
124 
125  numberOfUnmatchedPartons.push_back( nUnmatchedPartons );
126  sumDeltaE .push_back( sumDE );
127  sumDeltaPt.push_back( sumDPt );
128  sumDeltaR .push_back( sumDR );
129  }
130 }
131 
132 double
134 {
135  // calculate the distance between two lorentz vectors
136  // using DeltaR(eta, phi) or normal space angle(theta, phi)
137  if(useDeltaR_) return ROOT::Math::VectorUtil::DeltaR(v1, v2);
138  return ROOT::Math::VectorUtil::Angle(v1, v2);
139 }
140 
141 void
143 {
144  // match parton to jet with shortest distance
145  // starting with the shortest distance available
146  // apply some outlier rejection if desired
147 
148  // prepare vector of pairs with distances between
149  // all partons to all jets in the input vectors
150  std::vector< std::pair<double, unsigned int> > distances;
151  for(unsigned int ip=0; ip<partons.size(); ++ip){
152  for(unsigned int ij=0; ij<jets.size(); ++ij){
153  double dist = distance(jets[ij]->p4(), partons[ip]->p4());
154  distances.push_back(std::pair<double, unsigned int>(dist, ip*jets.size()+ij));
155  }
156  }
157  std::sort(distances.begin(), distances.end());
158 
160 
161  while(match.size() < partons.size()){
162  unsigned int partonIndex = distances[0].second/jets.size();
163  int jetIndex = distances[0].second-jets.size()*partonIndex;
164 
165  // use primitive outlier rejection if desired
166  if(useMaxDist_&& distances[0].first>maxDist_) jetIndex = -1;
167 
168  // prevent underflow in case of too few jets
169  if( distances.empty() )
170  match.push_back(std::make_pair(partonIndex, -1));
171  else
172  match.push_back(std::make_pair(partonIndex, jetIndex));
173 
174  // remove all values for the matched parton
175  // and the matched jet
176  for(unsigned int a=0; a<distances.size(); ++a){
177  unsigned int pIndex = distances[a].second/jets.size();
178  int jIndex = distances[a].second-jets.size()*pIndex;
179  if((pIndex == partonIndex) || (jIndex == jetIndex)){
180  distances.erase(distances.begin()+a, distances.begin()+a+1);
181  --a;
182  }
183  }
184  }
185 
186  matching.clear();
187  matching.push_back( match );
188  return;
189 }
190 
191 void
192 JetPartonMatching::minSumDist_recursion(const unsigned int ip, std::vector<unsigned int> & jetIndices,
193  std::vector<bool> & usedJets, std::vector<std::pair<double, MatchingCollection> > & distMatchVec)
194 {
195  // build up jet combinations recursively
196  if(ip<partons.size()){
197  for(unsigned int ij=0; ij<jets.size(); ++ij){
198  if(usedJets[ij]) continue;
199  usedJets[ij] = true;
200  jetIndices[ip] = ij;
201  minSumDist_recursion(ip+1, jetIndices, usedJets, distMatchVec);
202  usedJets[ij] = false;
203  }
204  return;
205  }
206 
207  // calculate sumDist for each completed combination
208  double sumDist = 0;
210  for(unsigned int ip=0; ip<partons.size(); ++ip){
211  double dist = distance(partons[ip]->p4(), jets[jetIndices[ip]]->p4());
212  if(useMaxDist_ && dist > maxDist_) return; // outlier rejection
213  sumDist += distance(partons[ip]->p4(), jets[jetIndices[ip]]->p4());
214  match.push_back(std::make_pair(ip, jetIndices[ip]));
215  }
216 
217  distMatchVec.push_back( std::make_pair(sumDist, match) );
218  return;
219 }
220 
222 {
223  // match partons to jets with minimal sum of
224  // the distances between all partons and jets
225 
226  std::vector<std::pair<double, MatchingCollection> > distMatchVec;
227 
228  std::vector<bool> usedJets;
229  for(unsigned int i=0; i<jets.size(); ++i){
230  usedJets.push_back(false);
231  }
232 
233  std::vector<unsigned int> jetIndices;
234  jetIndices.reserve(partons.size());
235 
236  minSumDist_recursion(0, jetIndices, usedJets, distMatchVec);
237 
238  std::sort(distMatchVec.begin(), distMatchVec.end());
239 
240  matching.clear();
241 
242  if(distMatchVec.empty()) {
243  MatchingCollection dummyMatch;
244  for(unsigned int ip=0; ip<partons.size(); ++ip)
245  dummyMatch.push_back(std::make_pair(ip, -1));
246  matching.push_back( dummyMatch );
247  }
248  else
249  for(unsigned int i=0; i<distMatchVec.size(); ++i)
250  matching.push_back( distMatchVec[i].second );
251 
252  return;
253 }
254 
255 void
257 {
258  // match partons to jets with minimal sum of
259  // the distances between all partons and jets
260  // order partons in pt first
261  std::vector<std::pair <double, unsigned int> > ptOrderedPartons;
262 
263  for(unsigned int ip=0; ip<partons.size(); ++ip)
264  ptOrderedPartons.push_back(std::make_pair(partons[ip]->pt(), ip));
265 
266  std::sort(ptOrderedPartons.begin(), ptOrderedPartons.end());
267  std::reverse(ptOrderedPartons.begin(), ptOrderedPartons.end());
268 
269  std::vector<unsigned int> jetIndices;
270  for(unsigned int ij=0; ij<jets.size(); ++ij) jetIndices.push_back(ij);
271 
273 
274  for(unsigned int ip=0; ip<ptOrderedPartons.size(); ++ip){
275  double minDist = 999.;
276  int ijMin = -1;
277 
278  for(unsigned int ij=0; ij<jetIndices.size(); ++ij){
279  double dist = distance(partons[ptOrderedPartons[ip].second]->p4(), jets[jetIndices[ij]]->p4());
280  if(dist < minDist){
281  if(!useMaxDist_ || dist <= maxDist_){
282  minDist = dist;
283  ijMin = ij;
284  }
285  }
286  }
287 
288  if(ijMin >= 0){
289  match.push_back( std::make_pair(ptOrderedPartons[ip].second, jetIndices[ijMin]) );
290  jetIndices.erase(jetIndices.begin() + ijMin, jetIndices.begin() + ijMin + 1);
291  }
292  else
293  match.push_back( std::make_pair(ptOrderedPartons[ip].second, -1) );
294  }
295 
296  matching.clear();
297  matching.push_back( match );
298  return;
299 }
300 
301 void
303 {
304  // match partons to jets, only accept event
305  // if there are no ambiguities
306  std::vector<bool> jetMatched;
307  for(unsigned int ij=0; ij<jets.size(); ++ij) jetMatched.push_back(false);
308 
310 
311  for(unsigned int ip=0; ip<partons.size(); ++ip){
312  int iMatch = -1;
313  for(unsigned int ij=0; ij<jets.size(); ++ij){
314  double dist = distance(partons[ip]->p4(), jets[ij]->p4());
315  if(dist <= maxDist_){
316  if(!jetMatched[ij]){ // new match for jet
317  jetMatched[ij] = true;
318  if(iMatch == -1) // new match for parton and jet
319  iMatch = ij;
320  else // parton already matched: ambiguity!
321  iMatch = -2;
322  }
323  else // jet already matched: ambiguity!
324  iMatch = -2;
325  }
326  }
327  match.push_back(std::make_pair(ip, iMatch));
328  }
329 
330  matching.clear();
331  matching.push_back( match );
332  return;
333 }
334 
335 int
336 JetPartonMatching::getMatchForParton(const unsigned int part, const unsigned int comb)
337 {
338  // return index of the matched jet for a given parton
339  // (if arguments for parton index and combinatoric index
340  // are in the valid range)
341  if(comb >= matching.size()) return -9;
342  if(part >= matching[comb].size()) return -9;
343  return (matching[comb])[part].second;
344 }
345 
346 std::vector<int>
348 {
349  // return a vector with the indices of the matched jets
350  // (ordered according to the vector of partons)
351  std::vector<int> jetIndices;
352  for(unsigned int part=0; part<partons.size(); ++part)
353  jetIndices.push_back( getMatchForParton(part, comb) );
354  return jetIndices;
355 }
356 
357 double
358 JetPartonMatching::getDistanceForParton(const unsigned int part, const unsigned int comb)
359 {
360  // get the distance between parton and its best matched jet
361  if(getMatchForParton(part, comb) < 0) return -999.;
362  return distance( jets[getMatchForParton(part,comb)]->p4(), partons[part]->p4() );
363 }
364 
365 double
367 {
368  // get sum of distances between partons and matched jets
369  double sumDists = 0.;
370  for(unsigned int part=0; part<partons.size(); ++part){
371  double dist = getDistanceForParton(part, comb);
372  if(dist < 0.) return -999.;
373  sumDists += dist;
374  }
375  return sumDists;
376 }
377 
378 void
380 {
381  //report using MessageLogger
382  edm::LogInfo log("JetPartonMatching");
383  log << "++++++++++++++++++++++++++++++++++++++++++++++ \n";
384  log << " algorithm : ";
385  switch(algorithm_) {
386  case totalMinDist : log << "totalMinDist "; break;
387  case minSumDist : log << "minSumDist "; break;
388  case ptOrderedMinDist : log << "ptOrderedMinDist"; break;
389  case unambiguousOnly : log << "unambiguousOnly "; break;
390  default : log << "UNKNOWN ";
391  }
392  log << "\n";
393  log << " useDeltaR : ";
394  switch(useDeltaR_) {
395  case false : log << "false"; break;
396  case true : log << "true ";
397  }
398  log << "\n";
399  log << " useMaxDist: ";
400  switch(useMaxDist_) {
401  case false : log << "false"; break;
402  case true : log << "true ";
403  }
404  log << " maxDist: " << maxDist_ << "\n";
405  log << " number of partons / jets: " << partons.size() << " / " << jets.size() << "\n";
406  log << " number of available combinations: " << getNumberOfAvailableCombinations() << "\n";
407  for(unsigned int comb = 0; comb < matching.size(); ++comb) {
408  log << " -------------------------------------------- \n";
409  log << " ind. of matched jets:";
410  for(unsigned int part = 0; part < partons.size(); ++part)
411  log << std::setw(4) << getMatchForParton(part, comb);
412  log << "\n";
413  log << " sumDeltaR : " << getSumDeltaR(comb) << "\n";
414  log << " sumDeltaPt / sumDeltaE: " << getSumDeltaPt(comb) << " / " << getSumDeltaE(comb);
415  log << "\n";
416  }
417  log << "++++++++++++++++++++++++++++++++++++++++++++++";
418 }
std::vector< const reco::Candidate * > partons
int i
Definition: DBlmapReader.cc:9
std::vector< std::pair< unsigned int, int > > MatchingCollection
std::vector< double > sumDeltaE
unsigned int getNumberOfAvailableCombinations()
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double p4[4]
Definition: TauolaWrapper.h:92
dictionary comb
vector< PseudoJet > jets
std::vector< int > getMatchesForPartons(const unsigned int comb=0)
int j
Definition: DBlmapReader.cc:9
double getSumDeltaE(const unsigned int comb=0)
std::vector< unsigned int > numberOfUnmatchedPartons
part
Definition: HCALResponse.h:20
void minSumDist_recursion(const unsigned int, std::vector< unsigned int > &, std::vector< bool > &, std::vector< std::pair< double, MatchingCollection > > &)
int getMatchForParton(const unsigned int part, const unsigned int comb=0)
double distance(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &)
std::vector< const reco::Candidate * > jets
std::vector< MatchingCollection > matching
std::vector< double > sumDeltaR
double a
Definition: hdecay.h:121
double getSumDeltaR(const unsigned int comb=0)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< double > sumDeltaPt
double getSumDistances(const unsigned int comb=0)
tuple size
Write out results.
double getSumDeltaPt(const unsigned int comb=0)
double getDistanceForParton(const unsigned int part, const unsigned int comb=0)
tuple log
Definition: cmsBatch.py:347