CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecoTauEnergyAlgorithmPlugin.cc
Go to the documentation of this file.
1 /*
2  * =============================================================================
3  * Filename: RecoTauEnergyAlgorithmPlugin.cc
4  *
5  * Description: Determine best estimate for tau energy
6  * for tau candidates reconstructed in different decay modes
7  *
8  * Created: 04/09/2013 11:40:00
9  *
10  * Authors: Christian Veelken (LLR)
11  *
12  * =============================================================================
13  */
14 
16 
27 
30 
31 #include <vector>
32 #include <cmath>
33 
34 namespace reco { namespace tau {
35 
37 {
38  public:
39 
42  void operator()(PFTau&) const;
43  virtual void beginJob(edm::EDProducer*);
44  virtual void beginEvent();
45  virtual void endEvent();
46 
47  private:
48 
51  double dRaddPhoton_;
52  double minGammaEt_;
53 
55 };
56 
58  : RecoTauModifierPlugin(cfg, std::move(iC)),
59  dRaddNeutralHadron_(cfg.getParameter<double>("dRaddNeutralHadron")),
60  minNeutralHadronEt_(cfg.getParameter<double>("minNeutralHadronEt")),
61  dRaddPhoton_(cfg.getParameter<double>("dRaddPhoton")),
62  minGammaEt_(cfg.getParameter<double>("minGammaEt"))
63 {
64  verbosity_ = ( cfg.exists("verbosity") ) ?
65  cfg.getParameter<int>("verbosity") : 0;
66 }
67 
69 {}
70 
72 {}
73 
75 {}
76 
77 namespace
78 {
79  double getTrackPerr2(const reco::Track& track)
80  {
81  double trackPerr = track.p()*(track.ptError()/track.pt());
82  return trackPerr*trackPerr;
83  }
84 
85  void updateTauP4(PFTau& tau, double sf, const reco::Candidate::LorentzVector& addP4)
86  {
87  // preserve tau candidate mass when adding extra neutral energy
88  double tauPx_modified = tau.px() + sf*addP4.px();
89  double tauPy_modified = tau.py() + sf*addP4.py();
90  double tauPz_modified = tau.pz() + sf*addP4.pz();
91  double tauMass = tau.mass();
92  double tauEn_modified = sqrt(tauPx_modified*tauPx_modified + tauPy_modified*tauPy_modified + tauPz_modified*tauPz_modified + tauMass*tauMass);
93  reco::Candidate::LorentzVector tauP4_modified(tauPx_modified, tauPy_modified, tauPz_modified, tauEn_modified);
94  tau.setP4(tauP4_modified);
95  }
96 
97  void killTau(PFTau& tau)
98  {
99  reco::Candidate::LorentzVector tauP4_modified(0.,0.,0.,0.);
100  tau.setP4(tauP4_modified);
101  tau.setStatus(-1);
102  }
103 }
104 
106 {
107  if ( verbosity_ ) {
108  std::cout << "<PFRecoTauEnergyAlgorithmPlugin::operator()>:" << std::endl;
109  std::cout << "tau: Pt = " << tau.pt() << ", eta = " << tau.eta() << ", phi = " << tau.phi() << " (En = " << tau.energy() << ", decayMode = " << tau.decayMode() << ")" << std::endl;
110  }
111 
112  // Add high Pt PFNeutralHadrons and PFGammas that are not "used" by tau decay mode object
113  std::vector<reco::PFCandidatePtr> addNeutrals;
114  reco::Candidate::LorentzVector addNeutralsSumP4;
115  std::vector<reco::PFCandidatePtr> jetConstituents = tau.jetRef()->getPFConstituents();
116  for ( std::vector<reco::PFCandidatePtr>::const_iterator jetConstituent = jetConstituents.begin();
117  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
118  reco::PFCandidate::ParticleType jetConstituentType = (*jetConstituent)->particleId();
119  if ( !((jetConstituentType == reco::PFCandidate::h0 && (*jetConstituent)->et() > minNeutralHadronEt_) ||
120  (jetConstituentType == reco::PFCandidate::gamma && (*jetConstituent)->et() > minGammaEt_ )) ) continue;
121 
122  bool isSignalPFCand = false;
123  const std::vector<reco::PFCandidatePtr>& signalPFCands = tau.signalPFCands();
124  for ( std::vector<reco::PFCandidatePtr>::const_iterator signalPFCand = signalPFCands.begin();
125  signalPFCand != signalPFCands.end(); ++signalPFCand ) {
126  if ( (*jetConstituent) == (*signalPFCand) ) isSignalPFCand = true;
127  }
128  if ( isSignalPFCand ) continue;
129 
130  double dR = deltaR((*jetConstituent)->p4(), tau.p4());
131  double dRadd = -1.;
132  if ( jetConstituentType == reco::PFCandidate::h0 ) dRadd = dRaddNeutralHadron_;
133  else if ( jetConstituentType == reco::PFCandidate::gamma ) dRadd = dRaddPhoton_;
134  if ( dR < dRadd ) {
135  addNeutrals.push_back(*jetConstituent);
136  addNeutralsSumP4 += (*jetConstituent)->p4();
137  }
138  }
139  if ( verbosity_ ) {
140  std::cout << "addNeutralsSumP4: En = " << addNeutralsSumP4.energy() << std::endl;
141  }
142 
143  unsigned numNonPFCandTracks = 0;
144  double nonPFCandTracksSumP = 0.;
145  double nonPFCandTracksSumPerr2 = 0.;
146  const std::vector<PFRecoTauChargedHadron>& chargedHadrons = tau.signalTauChargedHadronCandidates();
147  for ( std::vector<PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
148  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
149  if ( chargedHadron->algoIs(PFRecoTauChargedHadron::kTrack) ) {
150  ++numNonPFCandTracks;
151  const edm::Ptr<Track>& chargedHadronTrack = chargedHadron->getTrack();
152  nonPFCandTracksSumP += chargedHadronTrack->p();
153  nonPFCandTracksSumPerr2 += getTrackPerr2(*chargedHadronTrack);
154  }
155  }
156  if ( verbosity_ ) {
157  std::cout << "nonPFCandTracksSumP = " << nonPFCandTracksSumP << " +/- " << sqrt(nonPFCandTracksSumPerr2)
158  << " (numNonPFCandTracks = " << numNonPFCandTracks << ")" << std::endl;
159  }
160 
161  if ( numNonPFCandTracks == 0 ) {
162  // This is the easy case:
163  // All tau energy is taken from PFCandidates reconstructed by PFlow algorithm
164  // and there is no issue with double-counting of energy.
165  if ( verbosity_ ) {
166  std::cout << "easy case: all tracks are associated to PFCandidates --> leaving tau momentum untouched." << std::endl;
167  }
168  updateTauP4(tau, 1., addNeutralsSumP4);
169  return;
170  } else {
171  // This is the difficult case:
172  // The tau energy needs to be computed for an arbitrary mix of charged and neutral PFCandidates plus reco::Tracks.
173  // We need to make sure not to double-count energy deposited by reco::Track in ECAL and/or HCAL as neutral PFCandidates.
174 
175  // Check if we have enough energy in collection of PFNeutralHadrons and PFGammas that are not "used" by tau decay mode object
176  // to balance track momenta:
177  if ( nonPFCandTracksSumP < addNeutralsSumP4.energy() ) {
178  double scaleFactor = 1. - nonPFCandTracksSumP/addNeutralsSumP4.energy();
179  if ( !(scaleFactor >= 0. && scaleFactor <= 1.) ) {
180  edm::LogWarning("PFRecoTauEnergyAlgorithmPlugin::operator()")
181  << "Failed to compute tau energy --> killing tau candidate !!" << std::endl;
182  killTau(tau);
183  return;
184  }
185  if ( verbosity_ ) {
186  std::cout << "case (2): addNeutralsSumEn > nonPFCandTracksSumP --> adjusting tau momentum." << std::endl;
187  }
188  updateTauP4(tau, scaleFactor, addNeutralsSumP4);
189  return;
190  }
191 
192  // Determine which neutral PFCandidates are close to PFChargedHadrons
193  // and have been merged into ChargedHadrons
194  std::vector<reco::PFCandidatePtr> mergedNeutrals;
195  reco::Candidate::LorentzVector mergedNeutralsSumP4;
196  for ( std::vector<PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
197  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
198  if ( chargedHadron->algoIs(PFRecoTauChargedHadron::kTrack) ) {
199  const std::vector<reco::PFCandidatePtr>& neutralPFCands = chargedHadron->getNeutralPFCandidates();
200  for ( std::vector<reco::PFCandidatePtr>::const_iterator neutralPFCand = neutralPFCands.begin();
201  neutralPFCand != neutralPFCands.end(); ++neutralPFCand ) {
202  mergedNeutrals.push_back(*neutralPFCand);
203  mergedNeutralsSumP4 += (*neutralPFCand)->p4();
204  }
205  }
206  }
207  if ( verbosity_ ) {
208  std::cout << "mergedNeutralsSumP4: En = " << mergedNeutralsSumP4.energy() << std::endl;
209  }
210 
211  // Check if track momenta are balanced by sum of PFNeutralHadrons and PFGammas that are not "used" by tau decay mode object
212  // plus neutral PFCandidates close to PFChargedHadrons:
213  if ( nonPFCandTracksSumP < (addNeutralsSumP4.energy() + mergedNeutralsSumP4.energy()) ) {
214  double scaleFactor = ((addNeutralsSumP4.energy() + mergedNeutralsSumP4.energy()) - nonPFCandTracksSumP)/mergedNeutralsSumP4.energy();
215  if ( !(scaleFactor >= 0. && scaleFactor <= 1.) ) {
216  edm::LogWarning("PFRecoTauEnergyAlgorithmPlugin::operator()")
217  << "Failed to compute tau energy --> killing tau candidate !!" << std::endl;
218  killTau(tau);
219  return;
220  }
222  size_t numChargedHadrons = chargedHadrons.size();
223  for ( size_t iChargedHadron = 0; iChargedHadron < numChargedHadrons; ++iChargedHadron ) {
224  const PFRecoTauChargedHadron& chargedHadron = chargedHadrons[iChargedHadron];
225  if ( chargedHadron.getNeutralPFCandidates().size() >= 1 ) {
226  PFRecoTauChargedHadron chargedHadron_modified = chargedHadron;
227  setChargedHadronP4(chargedHadron_modified, scaleFactor);
228  tau.signalTauChargedHadronCandidates_[iChargedHadron] = chargedHadron_modified;
229  diffP4 += (chargedHadron.p4() - chargedHadron_modified.p4());
230  }
231  }
232  if ( verbosity_ ) {
233  std::cout << "case (3): (addNeutralsSumEn + mergedNeutralsSumEn) > nonPFCandTracksSumP --> adjusting tau momentum." << std::endl;
234  }
235  updateTauP4(tau, -1., diffP4);
236  return;
237  }
238 
239  // Determine energy sum of all PFNeutralHadrons interpreted as ChargedHadrons with missing track
240  unsigned numChargedHadronNeutrals = 0;
241  std::vector<reco::PFCandidatePtr> chargedHadronNeutrals;
242  reco::Candidate::LorentzVector chargedHadronNeutralsSumP4;
243  for ( std::vector<PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
244  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
245  if ( chargedHadron->algoIs(PFRecoTauChargedHadron::kPFNeutralHadron) ) {
246  ++numChargedHadronNeutrals;
247  chargedHadronNeutrals.push_back(chargedHadron->getChargedPFCandidate());
248  chargedHadronNeutralsSumP4 += chargedHadron->getChargedPFCandidate()->p4();
249  }
250  }
251  if ( verbosity_ ) {
252  std::cout << "chargedHadronNeutralsSumP4: En = " << chargedHadronNeutralsSumP4.energy()
253  << " (numChargedHadronNeutrals = " << numChargedHadronNeutrals << ")" << std::endl;
254  }
255 
256  // Check if sum of PFNeutralHadrons and PFGammas that are not "used" by tau decay mode object
257  // plus neutral PFCandidates close to PFChargedHadrons plus PFNeutralHadrons interpreted as ChargedHadrons with missing track balances track momenta
258  if ( nonPFCandTracksSumP < (addNeutralsSumP4.energy() + mergedNeutralsSumP4.energy() + chargedHadronNeutralsSumP4.energy()) ) {
259  double scaleFactor = ((addNeutralsSumP4.energy() + mergedNeutralsSumP4.energy() + chargedHadronNeutralsSumP4.energy()) - nonPFCandTracksSumP)/chargedHadronNeutralsSumP4.energy();
260  if ( !(scaleFactor >= 0. && scaleFactor <= 1.) ) {
261  edm::LogWarning("PFRecoTauEnergyAlgorithmPlugin::operator()")
262  << "Failed to compute tau energy --> killing tau candidate !!" << std::endl;
263  killTau(tau);
264  return;
265  }
267  size_t numChargedHadrons = chargedHadrons.size();
268  for ( size_t iChargedHadron = 0; iChargedHadron < numChargedHadrons; ++iChargedHadron ) {
269  const PFRecoTauChargedHadron& chargedHadron = chargedHadrons[iChargedHadron];
270  if ( chargedHadron.algoIs(PFRecoTauChargedHadron::kPFNeutralHadron) ) {
271  PFRecoTauChargedHadron chargedHadron_modified = chargedHadron;
272  chargedHadron_modified.neutralPFCandidates_.clear();
273  const PFCandidatePtr& chargedPFCand = chargedHadron.getChargedPFCandidate();
274  double chargedHadronPx_modified = scaleFactor*chargedPFCand->px();
275  double chargedHadronPy_modified = scaleFactor*chargedPFCand->py();
276  double chargedHadronPz_modified = scaleFactor*chargedPFCand->pz();
277  reco::Candidate::LorentzVector chargedHadronP4_modified = compChargedHadronP4fromPxPyPz(chargedHadronPx_modified, chargedHadronPy_modified, chargedHadronPz_modified);
278  chargedHadron_modified.setP4(chargedHadronP4_modified);
279  tau.signalTauChargedHadronCandidates_[iChargedHadron] = chargedHadron_modified;
280  diffP4 += (chargedHadron.p4() - chargedHadron_modified.p4());
281  }
282  }
283  if ( verbosity_ ) {
284  std::cout << "case (4): (addNeutralsSumEn + mergedNeutralsSumEn + chargedHadronNeutralsSumEn) > nonPFCandTracksSumP --> adjusting momenta of tau and chargedHadrons." << std::endl;
285  }
286  updateTauP4(tau, -1., diffP4);
287  return;
288  } else {
289  double allTracksSumP = 0.;
290  double allTracksSumPerr2 = 0.;
291  const std::vector<PFRecoTauChargedHadron> chargedHadrons = tau.signalTauChargedHadronCandidates();
292  for ( std::vector<PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
293  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
294  if ( chargedHadron->algoIs(PFRecoTauChargedHadron::kChargedPFCandidate) || chargedHadron->algoIs(PFRecoTauChargedHadron::kTrack) ) {
295  const edm::Ptr<Track>& chargedHadronTrack = chargedHadron->getTrack();
296  if ( chargedHadronTrack.isNonnull() ) {
297  allTracksSumP += chargedHadronTrack->p();
298  allTracksSumPerr2 += getTrackPerr2(*chargedHadronTrack);
299  } else {
300  edm::LogWarning("PFRecoTauEnergyAlgorithmPlugin::operator()")
301  << "PFRecoTauChargedHadron has no associated reco::Track !!" << std::endl;
302  if ( verbosity_ ) {
303  chargedHadron->print();
304  }
305  }
306  }
307  }
308  if ( verbosity_ ) {
309  std::cout << "allTracksSumP = " << allTracksSumP << " +/- " << sqrt(allTracksSumPerr2) << std::endl;
310  }
311  double allNeutralsSumEn = 0.;
312  const std::vector<reco::PFCandidatePtr>& signalPFCands = tau.signalPFCands();
313  for ( std::vector<reco::PFCandidatePtr>::const_iterator signalPFCand = signalPFCands.begin();
314  signalPFCand != signalPFCands.end(); ++signalPFCand ) {
315  if ( verbosity_ ) {
316  std::cout << "PFCandidate #" << signalPFCand->id() << ":" << signalPFCand->key() << ":"
317  << " Pt = " << (*signalPFCand)->pt() << ", eta = " << (*signalPFCand)->eta() << ", phi = " << (*signalPFCand)->phi() << std::endl;
318  std::cout << "calorimeter energy:"
319  << " ECAL = " << (*signalPFCand)->ecalEnergy() << ","
320  << " HCAL = " << (*signalPFCand)->hcalEnergy() << ","
321  << " HO = " << (*signalPFCand)->hoEnergy() << std::endl;
322  }
323  if ( edm::isFinite((*signalPFCand)->ecalEnergy()) ) allNeutralsSumEn += (*signalPFCand)->ecalEnergy();
324  if ( edm::isFinite((*signalPFCand)->hcalEnergy()) ) allNeutralsSumEn += (*signalPFCand)->hcalEnergy();
325  if ( edm::isFinite((*signalPFCand)->hoEnergy()) ) allNeutralsSumEn += (*signalPFCand)->hoEnergy();
326  }
327  allNeutralsSumEn += addNeutralsSumP4.energy();
328  if ( allNeutralsSumEn < 0. ) allNeutralsSumEn = 0.;
329  if ( verbosity_ ) {
330  std::cout << "allNeutralsSumEn = " << allNeutralsSumEn << std::endl;
331  }
332  if ( allNeutralsSumEn > allTracksSumP ) {
333  // Adjust momenta of neutral PFCandidates merged into ChargedHadrons
334  size_t numChargedHadrons = chargedHadrons.size();
335  for ( size_t iChargedHadron = 0; iChargedHadron < numChargedHadrons; ++iChargedHadron ) {
336  const PFRecoTauChargedHadron& chargedHadron = chargedHadrons[iChargedHadron];
338  PFRecoTauChargedHadron chargedHadron_modified = chargedHadron;
339  chargedHadron_modified.neutralPFCandidates_.clear();
340  chargedHadron_modified.setP4(chargedHadron.getChargedPFCandidate()->p4());
341  if ( verbosity_ ) {
342  std::cout << "chargedHadron #" << iChargedHadron << ": changing En = " << chargedHadron.energy() << " to " << chargedHadron_modified.energy() << std::endl;
343  }
344  tau.signalTauChargedHadronCandidates_[iChargedHadron] = chargedHadron_modified;
345  } else if ( chargedHadron.algoIs(PFRecoTauChargedHadron::kTrack) ) {
346  PFRecoTauChargedHadron chargedHadron_modified = chargedHadron;
347  chargedHadron_modified.neutralPFCandidates_.clear();
348  reco::Candidate::LorentzVector chargedHadronP4_modified(0.,0.,0.,0.);
349  if ( (chargedHadron.getTrack()).isNonnull() ) {
350  const Track& chargedHadronTrack = *(chargedHadron.getTrack());
351  double chargedHadronPx_modified = chargedHadronTrack.px();
352  double chargedHadronPy_modified = chargedHadronTrack.py();
353  double chargedHadronPz_modified = chargedHadronTrack.pz();
354  chargedHadronP4_modified = compChargedHadronP4fromPxPyPz(chargedHadronPx_modified, chargedHadronPy_modified, chargedHadronPz_modified);
355  } else {
356  edm::LogWarning("PFRecoTauEnergyAlgorithmPlugin::operator()")
357  << "PFRecoTauChargedHadron has no associated reco::Track !!" << std::endl;
358  if ( verbosity_ ) {
359  chargedHadron.print();
360  }
361  }
362  chargedHadron_modified.setP4(chargedHadronP4_modified);
363  if ( verbosity_ ) {
364  std::cout << "chargedHadron #" << iChargedHadron << ": changing En = " << chargedHadron.energy() << " to " << chargedHadron_modified.energy() << std::endl;
365  }
366  tau.signalTauChargedHadronCandidates_[iChargedHadron] = chargedHadron_modified;
367  }
368  }
369  double scaleFactor = allNeutralsSumEn/tau.energy();
370  if ( verbosity_ ) {
371  std::cout << "case (5): allNeutralsSumEn > allTracksSumP --> adjusting momenta of tau and chargedHadrons." << std::endl;
372  }
373  updateTauP4(tau, scaleFactor - 1., tau.p4());
374  return;
375  } else {
376  if ( numChargedHadronNeutrals == 0 && tau.signalPiZeroCandidates().size() == 0 ) {
377  // Adjust momenta of ChargedHadrons build from reco::Tracks to match sum of energy deposits in ECAL + HCAL + HO
378  size_t numChargedHadrons = chargedHadrons.size();
379  for ( size_t iChargedHadron = 0; iChargedHadron < numChargedHadrons; ++iChargedHadron ) {
380  const PFRecoTauChargedHadron& chargedHadron = chargedHadrons[iChargedHadron];
382  PFRecoTauChargedHadron chargedHadron_modified = chargedHadron;
383  chargedHadron_modified.neutralPFCandidates_.clear();
384  reco::Candidate::LorentzVector chargedHadronP4_modified(0.,0.,0.,0.);
385  const edm::Ptr<Track>& chargedHadronTrack = chargedHadron.getTrack();
386  if ( chargedHadronTrack.isNonnull() ) {
387  double trackP = chargedHadronTrack->p();
388  double trackPerr2 = getTrackPerr2(*chargedHadronTrack);
389  if ( verbosity_ ) {
390  std::cout << "trackP = " << trackP << " +/- " << sqrt(trackPerr2) << std::endl;
391  }
392  // CV: adjust track momenta such that difference beeen (measuredTrackP - adjustedTrackP)/sigmaMeasuredTrackP is minimal
393  // (expression derived using Mathematica)
394  double trackP_modified =
395  (trackP*(allTracksSumPerr2 - trackPerr2)
396  + trackPerr2*(allNeutralsSumEn - (allTracksSumP - trackP)))/
397  allTracksSumPerr2;
398  // CV: trackP_modified may actually become negative in case sum of energy deposits in ECAL + HCAL + HO is small
399  // and one of the tracks has a significantly larger momentum uncertainty than the other tracks.
400  // In this case set track momentum to small positive value.
401  if ( trackP_modified < 1.e-1 ) trackP_modified = 1.e-1;
402  if ( verbosity_ ) {
403  std::cout << "trackP (modified) = " << trackP_modified << std::endl;
404  }
405  double scaleFactor = trackP_modified/trackP;
406  if ( !(scaleFactor >= 0. && scaleFactor <= 1.) ) {
407  edm::LogWarning("PFRecoTauEnergyAlgorithmPlugin::operator()")
408  << "Failed to compute tau energy --> killing tau candidate !!" << std::endl;
409  killTau(tau);
410  return;
411  }
412  double chargedHadronPx_modified = scaleFactor*chargedHadronTrack->px();
413  double chargedHadronPy_modified = scaleFactor*chargedHadronTrack->py();
414  double chargedHadronPz_modified = scaleFactor*chargedHadronTrack->pz();
415  chargedHadronP4_modified = compChargedHadronP4fromPxPyPz(chargedHadronPx_modified, chargedHadronPy_modified, chargedHadronPz_modified);
416  } else {
417  edm::LogWarning("PFRecoTauEnergyAlgorithmPlugin::operator()")
418  << "PFRecoTauChargedHadron has no associated reco::Track !!" << std::endl;
419  if ( verbosity_ ) {
420  chargedHadron.print();
421  }
422  }
423  chargedHadron_modified.setP4(chargedHadronP4_modified);
424  if ( verbosity_ ) {
425  std::cout << "chargedHadron #" << iChargedHadron << ": changing En = " << chargedHadron.energy() << " to " << chargedHadron_modified.energy() << std::endl;
426  }
427  tau.signalTauChargedHadronCandidates_[iChargedHadron] = chargedHadron_modified;
428  }
429  }
430  double scaleFactor = allNeutralsSumEn/tau.energy();
431  if ( verbosity_ ) {
432  std::cout << "case (6): allNeutralsSumEn < allTracksSumP --> adjusting momenta of tau and chargedHadrons." << std::endl;
433  }
434  updateTauP4(tau, scaleFactor - 1., tau.p4());
435  return;
436  } else {
437  // Interpretation of PFNeutralHadrons as ChargedHadrons with missing track and/or reconstruction of extra PiZeros
438  // is not compatible with the fact that sum of reco::Track momenta exceeds sum of energy deposits in ECAL + HCAL + HO:
439  // kill tau candidate (by setting its four-vector to zero)
440  if ( verbosity_ ) {
441  std::cout << "case (7): allNeutralsSumEn < allTracksSumP not compatible with tau decay mode hypothesis --> killing tau candidate." << std::endl;
442  }
443  killTau(tau);
444  return;
445  }
446  }
447  }
448  }
449 
450  // CV: You should never come here.
451  if ( verbosity_ ) {
452  std::cout << "undefined case: you should never come here !!" << std::endl;
453  }
454  assert(0);
455 }
456 
458 {}
459 
460 }} // end namespace reco::tau
461 
463 
double p() const
momentum vector magnitude
Definition: TrackBase.h:127
virtual double energy() const GCC11_FINAL
energy
T getParameter(std::string const &) const
ParticleType
particle types
Definition: PFCandidate.h:43
const PFJetRef & jetRef() const
Definition: PFTau.cc:52
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
std::vector< reco::PFRecoTauChargedHadron > signalTauChargedHadronCandidates_
Definition: PFTau.h:255
bool exists(std::string const &parameterName) const
checks if a parameter exists
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:131
virtual double pz() const GCC11_FINAL
z coordinate of momentum vector
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
const std::vector< reco::PFCandidatePtr > & signalPFCands() const
PFCandidates in signal region.
Definition: PFTau.cc:72
bool isFinite(T x)
hadronicDecayMode decayMode() const
Definition: PFTau.cc:176
void print(std::ostream &stream=std::cout) const
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:152
const double tauMass
const std::vector< RecoTauPiZero > & signalPiZeroCandidates() const
Retrieve the association of signal region gamma candidates into candidate PiZeros.
Definition: PFTau.cc:91
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:129
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:192
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
const TrackPtr & getTrack() const
reference to reco::Track
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:135
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
virtual float mass() const GCC11_FINAL
mass
const std::vector< PFCandidatePtr > & getNeutralPFCandidates() const
references to additional neutral PFCandidates
std::vector< PFCandidatePtr > neutralPFCandidates_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
virtual void setP4(const LorentzVector &p4) GCC11_FINAL
set 4-momentum
PFRecoTauEnergyAlgorithmPlugin(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
bool algoIs(PFRecoTauChargedHadronAlgorithm algo) const
Check whether a given algo produced this charged hadron.
const std::vector< PFRecoTauChargedHadron > & signalTauChargedHadronCandidates() const
Retrieve the association of signal region PF candidates into candidate PFRecoTauChargedHadrons.
Definition: PFTau.cc:138
tuple cout
Definition: gather_cfg.py:121
#define DEFINE_EDM_PLUGIN(factory, type, name)
virtual float pt() const GCC11_FINAL
transverse momentum
void setChargedHadronP4(reco::PFRecoTauChargedHadron &chargedHadron, double scaleFactor_neutralPFCands=1.0)
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:133
const PFCandidatePtr & getChargedPFCandidate() const
reference to &quot;charged&quot; PFCandidate (either charged PFCandidate or PFNeutralHadron) ...
reco::Candidate::LorentzVector compChargedHadronP4fromPxPyPz(double, double, double)