CMS 3D CMS Logo

PFRecoTauClusterVariables.cc
Go to the documentation of this file.
2 
3 namespace {
4  struct PFTau_traits{ typedef reco::PFTau Tau_t; typedef const std::vector<reco::PFCandidatePtr>& Ret_t; };
5  struct PATTau_traits{ typedef pat::Tau Tau_t; typedef reco::CandidatePtrVector Ret_t; };
6 
7  template<typename T>
8  typename T::Ret_t getGammas_T(const typename T::Tau_t& tau, bool signal) {
9  return typename T::Ret_t();
10  }
12  template<>
13  const std::vector<reco::PFCandidatePtr>& getGammas_T<PFTau_traits>(const reco::PFTau& tau, bool signal) {
14  if (signal){
15  return tau.signalPFGammaCands();
16  }
17  return tau.isolationPFGammaCands();
18  }
19 
20  template<>
21  reco::CandidatePtrVector getGammas_T<PATTau_traits>(const pat::Tau& tau, bool signal) {
22  if(signal){
23  return tau.signalGammaCands();
24  }
25  return tau.isolationGammaCands();
26  }
27 
29  bool isInside(float photon_pt, float deta, float dphi) {
30  constexpr double stripEtaAssociationDistance_0p95_p0 = 0.197077;
31  constexpr double stripEtaAssociationDistance_0p95_p1 = 0.658701;
32  constexpr double stripPhiAssociationDistance_0p95_p0 = 0.352476;
33  constexpr double stripPhiAssociationDistance_0p95_p1 = 0.707716;
34  if(photon_pt==0){
35  return false;
36  } if((dphi<0.3 && dphi<std::max(0.05, stripPhiAssociationDistance_0p95_p0*std::pow(photon_pt, -stripPhiAssociationDistance_0p95_p1))) && (deta<0.15 && deta<std::max(0.05, stripEtaAssociationDistance_0p95_p0*std::pow(photon_pt, -stripEtaAssociationDistance_0p95_p1)))){
37  return true;
38  }
39  return false;
40  }
41 }
42 
43 namespace reco { namespace tau {
45  float lead_track_chi2(const reco::PFTau& tau) {
46  float LeadingTracknormalizedChi2 = 0;
47  const reco::PFCandidatePtr& leadingPFCharged = tau.leadPFChargedHadrCand() ;
48  if (leadingPFCharged.isNonnull()) {
49  reco::TrackRef tref = leadingPFCharged -> trackRef();
50  if (tref.isNonnull()) {
51  LeadingTracknormalizedChi2 = (float)(tref -> normalizedChi2());
52  }
53  }
54  return LeadingTracknormalizedChi2;
55  }
57  float eratio(const reco::PFTau& tau) {
58  float ecal_en_in_signal_pf_cands = 0;
59  float hcal_en_in_signal_pf_cands = 0;
60  for (const auto& signal_cand : tau.signalPFCands()) {
61  ecal_en_in_signal_pf_cands += signal_cand->ecalEnergy();
62  hcal_en_in_signal_pf_cands += signal_cand->hcalEnergy();
63  }
64  float total = ecal_en_in_signal_pf_cands + hcal_en_in_signal_pf_cands;
65  if(total==0){
66  return -1;
67  }
68  return ecal_en_in_signal_pf_cands/total;
69  }
70  float eratio(const pat::Tau& tau) {
71  float ecal_en_in_signal_cands = tau.ecalEnergy();
72  float hcal_en_in_signal_cands = tau.hcalEnergy();
73  float total = ecal_en_in_signal_cands + hcal_en_in_signal_cands;
74  if(total == 0){
75  return -1;
76  }
77  return ecal_en_in_signal_cands/total;
78  }
81  template<typename T>
82  float pt_weighted_dx_T(const typename T::Tau_t& tau, int mode, int var, int decaymode) {
83  float sum_pt = 0.;
84  float sum_dx_pt = 0.;
85  float signalrad = std::max(0.05, std::min(0.1, 3./std::max(1., tau.pt())));
86  int is3prong = (decaymode==10);
87  const auto& cands = getGammas_T<T>(tau, mode < 2);
88  for (const auto& cand : cands) {
89  // only look at electrons/photons with pT > 0.5
90  if (cand->pt() < 0.5){
91  continue;
92  }
93  float dr = reco::deltaR(cand->eta(), cand->phi(), tau.eta(), tau.phi());
94  float deta = std::abs(cand->eta() - tau.eta());
95  float dphi = std::abs(reco::deltaPhi(cand->phi(), tau.phi()));
96  float pt = cand->pt();
97  bool flag = isInside(pt, deta, dphi);
98  if(is3prong==0){
99  if (mode == 2 || (mode == 0 && dr < signalrad) || (mode == 1 && dr > signalrad)) {
100  sum_pt += pt;
101  if (var == 0)
102  sum_dx_pt += pt * dr;
103  else if (var == 1)
104  sum_dx_pt += pt * deta;
105  else if (var == 2)
106  sum_dx_pt += pt * dphi;
107  }
108  }
109  else if(is3prong==1){
110  if( (mode==2 && flag==false) || (mode==1 && flag==true) || mode==0){
111  sum_pt += pt;
112  if (var == 0)
113  sum_dx_pt += pt * dr;
114  else if (var == 1)
115  sum_dx_pt += pt * deta;
116  else if (var == 2)
117  sum_dx_pt += pt * dphi;
118  }
119  }
120  }
121  if (sum_pt > 0.){
122  return sum_dx_pt/sum_pt;
123  }
124  return 0.;
125  }
126  float pt_weighted_dx(const reco::PFTau& tau, int mode, int var, int decaymode) {
127  return pt_weighted_dx_T<PFTau_traits>(tau, mode, var, decaymode);
128  }
129  float pt_weighted_dx(const pat::Tau& tau, int mode, int var, int decaymode) {
130  return pt_weighted_dx_T<PATTau_traits>(tau, mode, var, decaymode);
131  }
133  unsigned int n_photons_total(const reco::PFTau& tau) {
134  unsigned int n_photons = 0;
135  for (auto& cand : tau.signalPFGammaCands()) {
136  if (cand->pt() > 0.5)
137  ++n_photons;
138  }
139  for (auto& cand : tau.isolationPFGammaCands()) {
140  if (cand->pt() > 0.5)
141  ++n_photons;
142  }
143  return n_photons;
144  }
145  unsigned int n_photons_total(const pat::Tau& tau) {
146  unsigned int n_photons = 0;
147  for (auto& cand : tau.signalGammaCands()) {
148  if (cand->pt() > 0.5)
149  ++n_photons;
150  }
151  for (auto& cand : tau.isolationGammaCands()) {
152  if (cand->pt() > 0.5)
153  ++n_photons;
154  }
155  return n_photons;
156  }
157 
158  bool fillIsoMVARun2Inputs(float* mvaInput, const pat::Tau& tau, int mvaOpt, const std::string& nameCharged,
159  const std::string& nameNeutral, const std::string& namePu,
160  const std::string& nameOutside, const std::string& nameFootprint)
161  {
162  int tauDecayMode = tau.decayMode();
163 
164  if ( ((mvaOpt == kOldDMwoLT || mvaOpt == kOldDMwLT || mvaOpt == kDBoldDMwLT || mvaOpt == kPWoldDMwLT || mvaOpt == kDBoldDMwLTwGJ)
165  && (tauDecayMode == 0 || tauDecayMode == 1 || tauDecayMode == 2 || tauDecayMode == 10))
166  ||
167  ((mvaOpt == kNewDMwoLT || mvaOpt == kNewDMwLT || mvaOpt == kDBnewDMwLT || mvaOpt == kPWnewDMwLT || mvaOpt == kDBnewDMwLTwGJ)
168  && (tauDecayMode == 0 || tauDecayMode == 1 || tauDecayMode == 2 || tauDecayMode == 5 || tauDecayMode == 6 || tauDecayMode == 10 || tauDecayMode == 11))
169  ) {
170 
171  float chargedIsoPtSum = tau.tauID(nameCharged);
172  float neutralIsoPtSum = tau.tauID(nameNeutral);
173  float puCorrPtSum = tau.tauID(namePu);
174  float photonPtSumOutsideSignalCone = tau.tauID(nameOutside);
175  float footprintCorrection = tau.tauID(nameFootprint);
176 
177  float decayDistX = tau.flightLength().x();
178  float decayDistY = tau.flightLength().y();
179  float decayDistZ = tau.flightLength().z();
180  float decayDistMag = std::sqrt(decayDistX*decayDistX + decayDistY*decayDistY + decayDistZ*decayDistZ);
181 
182  // --- The following 5 variables differ slightly between AOD & MiniAOD
183  // because they are recomputed using packedCandidates saved in the tau
184  float nPhoton = reco::tau::n_photons_total(tau);
185  float ptWeightedDetaStrip = reco::tau::pt_weighted_deta_strip(tau, tauDecayMode);
186  float ptWeightedDphiStrip = reco::tau::pt_weighted_dphi_strip(tau, tauDecayMode);
187  float ptWeightedDrSignal = reco::tau::pt_weighted_dr_signal(tau, tauDecayMode);
188  float ptWeightedDrIsolation = reco::tau::pt_weighted_dr_iso(tau, tauDecayMode);
189  // ---
190  float leadingTrackChi2 = tau.leadingTrackNormChi2();
191  float eRatio = reco::tau::eratio(tau);
192 
193  // Difference between measured and maximally allowed Gottfried-Jackson angle
194  float gjAngleDiff = -999;
195  if ( tauDecayMode == 10 ) {
196  double mTau = 1.77682;
197  double mAOne = tau.p4().M();
198  double pAOneMag = tau.p();
199  double argumentThetaGJmax = (std::pow(mTau,2) - std::pow(mAOne,2) ) / ( 2 * mTau * pAOneMag );
200  double argumentThetaGJmeasured = ( tau.p4().px() * decayDistX + tau.p4().py() * decayDistY + tau.p4().pz() * decayDistZ ) / ( pAOneMag * decayDistMag );
201  if ( std::abs(argumentThetaGJmax) <= 1. && std::abs(argumentThetaGJmeasured) <= 1. ) {
202  double thetaGJmax = std::asin( argumentThetaGJmax );
203  double thetaGJmeasured = std::acos( argumentThetaGJmeasured );
204  gjAngleDiff = thetaGJmeasured - thetaGJmax;
205  }
206  }
207 
208  if ( mvaOpt == kOldDMwoLT || mvaOpt == kNewDMwoLT ) {
209  mvaInput[0] = std::log(std::max(1.f, (float)tau.pt()));
210  mvaInput[1] = std::abs((float)tau.eta());
211  mvaInput[2] = std::log(std::max(1.e-2f, chargedIsoPtSum));
212  mvaInput[3] = std::log(std::max(1.e-2f, neutralIsoPtSum - 0.125f*puCorrPtSum));
213  mvaInput[4] = std::log(std::max(1.e-2f, puCorrPtSum));
214  mvaInput[5] = tauDecayMode;
215  } else if ( mvaOpt == kOldDMwLT || mvaOpt == kNewDMwLT ) {
216  mvaInput[0] = std::log(std::max(1.f, (float)tau.pt()));
217  mvaInput[1] = std::abs((float)tau.eta());
218  mvaInput[2] = std::log(std::max(1.e-2f, chargedIsoPtSum));
219  mvaInput[3] = std::log(std::max(1.e-2f, neutralIsoPtSum - 0.125f*puCorrPtSum));
220  mvaInput[4] = std::log(std::max(1.e-2f, puCorrPtSum));
221  mvaInput[5] = tauDecayMode;
222  mvaInput[6] = std::copysign(+1.f, tau.dxy());
223  mvaInput[7] = std::sqrt(std::min(1.f, std::abs(tau.dxy())));
224  mvaInput[8] = std::min(10.f, std::abs(tau.dxy_Sig()));
225  mvaInput[9] = ( tau.hasSecondaryVertex() ) ? 1. : 0.;
226  mvaInput[10] = std::sqrt(decayDistMag);
227  mvaInput[11] = std::min(10.f, tau.flightLengthSig());
228  } else if ( mvaOpt == kDBoldDMwLT || mvaOpt == kDBnewDMwLT ) {
229  mvaInput[0] = std::log(std::max(1.f, (float)tau.pt()));
230  mvaInput[1] = std::abs((float)tau.eta());
231  mvaInput[2] = std::log(std::max(1.e-2f, chargedIsoPtSum));
232  mvaInput[3] = std::log(std::max(1.e-2f, neutralIsoPtSum));
233  mvaInput[4] = std::log(std::max(1.e-2f, puCorrPtSum));
234  mvaInput[5] = std::log(std::max(1.e-2f, photonPtSumOutsideSignalCone));
235  mvaInput[6] = tauDecayMode;
236  mvaInput[7] = std::min(30.f, nPhoton);
237  mvaInput[8] = std::min(0.5f, ptWeightedDetaStrip);
238  mvaInput[9] = std::min(0.5f, ptWeightedDphiStrip);
239  mvaInput[10] = std::min(0.5f, ptWeightedDrSignal);
240  mvaInput[11] = std::min(0.5f, ptWeightedDrIsolation);
241  mvaInput[12] = std::min(100.f, leadingTrackChi2);
242  mvaInput[13] = std::min(1.f, eRatio);
243  mvaInput[14] = std::copysign(+1.f, tau.dxy());
244  mvaInput[15] = std::sqrt(std::min(1.f, std::abs(tau.dxy())));
245  mvaInput[16] = std::min(10.f, std::abs(tau.dxy_Sig()));
246  mvaInput[17] = std::copysign(+1.f, tau.ip3d());
247  mvaInput[18] = std::sqrt(std::min(1.f, std::abs(tau.ip3d())));
248  mvaInput[19] = std::min(10.f, std::abs(tau.ip3d_Sig()));
249  mvaInput[20] = ( tau.hasSecondaryVertex() ) ? 1. : 0.;
250  mvaInput[21] = std::sqrt(decayDistMag);
251  mvaInput[22] = std::min(10.f, tau.flightLengthSig());
252  } else if ( mvaOpt == kPWoldDMwLT || mvaOpt == kPWnewDMwLT ) {
253  mvaInput[0] = std::log(std::max(1.f, (float)tau.pt()));
254  mvaInput[1] = std::abs((float)tau.eta());
255  mvaInput[2] = std::log(std::max(1.e-2f, chargedIsoPtSum));
256  mvaInput[3] = std::log(std::max(1.e-2f, neutralIsoPtSum));
257  mvaInput[4] = std::log(std::max(1.e-2f, footprintCorrection));
258  mvaInput[5] = std::log(std::max(1.e-2f, photonPtSumOutsideSignalCone));
259  mvaInput[6] = tauDecayMode;
260  mvaInput[7] = std::min(30.f, nPhoton);
261  mvaInput[8] = std::min(0.5f, ptWeightedDetaStrip);
262  mvaInput[9] = std::min(0.5f, ptWeightedDphiStrip);
263  mvaInput[10] = std::min(0.5f, ptWeightedDrSignal);
264  mvaInput[11] = std::min(0.5f, ptWeightedDrIsolation);
265  mvaInput[12] = std::min(100.f, leadingTrackChi2);
266  mvaInput[13] = std::min(1.f, eRatio);
267  mvaInput[14] = std::copysign(+1.f, tau.dxy());
268  mvaInput[15] = std::sqrt(std::min(1.f, std::abs(tau.dxy())));
269  mvaInput[16] = std::min(10.f, std::abs(tau.dxy_Sig()));
270  mvaInput[17] = std::copysign(+1.f, tau.ip3d());
271  mvaInput[18] = std::sqrt(std::min(1.f, std::abs(tau.ip3d())));
272  mvaInput[19] = std::min(10.f, std::abs(tau.ip3d_Sig()));
273  mvaInput[20] = ( tau.hasSecondaryVertex() ) ? 1. : 0.;
274  mvaInput[21] = std::sqrt(decayDistMag);
275  mvaInput[22] = std::min(10.f, tau.flightLengthSig());
276  } else if ( mvaOpt == kDBoldDMwLTwGJ || mvaOpt == kDBnewDMwLTwGJ ) {
277  mvaInput[0] = std::log(std::max(1.f, (float)tau.pt()));
278  mvaInput[1] = std::abs((float)tau.eta());
279  mvaInput[2] = std::log(std::max(1.e-2f, chargedIsoPtSum));
280  mvaInput[3] = std::log(std::max(1.e-2f, neutralIsoPtSum));
281  mvaInput[4] = std::log(std::max(1.e-2f, puCorrPtSum));
282  mvaInput[5] = std::log(std::max(1.e-2f, photonPtSumOutsideSignalCone));
283  mvaInput[6] = tauDecayMode;
284  mvaInput[7] = std::min(30.f, nPhoton);
285  mvaInput[8] = std::min(0.5f, ptWeightedDetaStrip);
286  mvaInput[9] = std::min(0.5f, ptWeightedDphiStrip);
287  mvaInput[10] = std::min(0.5f, ptWeightedDrSignal);
288  mvaInput[11] = std::min(0.5f, ptWeightedDrIsolation);
289  mvaInput[12] = std::min(1.f, eRatio);
290  mvaInput[13] = std::copysign(+1.f, tau.dxy());
291  mvaInput[14] = std::sqrt(std::min(1.f, std::abs(tau.dxy())));
292  mvaInput[15] = std::min(10.f, std::abs(tau.dxy_Sig()));
293  mvaInput[16] = std::copysign(+1.f, tau.ip3d());
294  mvaInput[17] = std::sqrt(std::min(1.f, std::abs(tau.ip3d())));
295  mvaInput[18] = std::min(10.f, std::abs(tau.ip3d_Sig()));
296  mvaInput[19] = ( tau.hasSecondaryVertex() ) ? 1. : 0.;
297  mvaInput[20] = std::sqrt(decayDistMag);
298  mvaInput[21] = std::min(10.f, tau.flightLengthSig());
299  mvaInput[22] = std::max(-1.f, gjAngleDiff);
300  }
301 
302  return true;
303  }
304  return false;
305  }
306 
307 }} // namespaces
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
unsigned int n_photons_total(const reco::PFTau &tau)
return total number of pf photon candidates with pT>500 MeV, which are associated to signal ...
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
double eta() const final
momentum pseudorapidity
const PFCandidatePtr & leadPFChargedHadrCand() const
Definition: PFTau.cc:67
float pt_weighted_dx(const reco::PFTau &tau, int mode=0, int var=0, int decaymode=-1)
float pt_weighted_dr_signal(const reco::PFTau &tau, int dm)
float dxy() const
Definition: Tau.h:314
const std::vector< reco::PFCandidatePtr > & signalPFGammaCands() const
Gamma candidates in signal region.
Definition: PFTau.cc:84
bool fillIsoMVARun2Inputs(float *mvaInput, const pat::Tau &tau, int mvaOpt, const std::string &nameCharged, const std::string &nameNeutral, const std::string &namePu, const std::string &nameOutside, const std::string &nameFootprint)
double pt() const final
transverse momentum
float ip3d_Sig() const
float tauID(const std::string &name) const
float hcalEnergy() const
return sum of hcal energies from signal candidates
Definition: Tau.h:336
float ecalEnergy() const
Definition: Tau.h:334
const std::vector< reco::PFCandidatePtr > & signalPFCands() const
PFCandidates in signal region.
Definition: PFTau.cc:78
reco::CandidatePtrVector signalGammaCands() const
#define constexpr
const std::vector< reco::PFCandidatePtr > & isolationPFGammaCands() const
Gamma candidates in isolation region.
Definition: PFTau.cc:93
float ip3d() const
Definition: Tau.h:327
float pt_weighted_deta_strip(const reco::PFTau &tau, int dm)
T sqrt(T t)
Definition: SSEVec.h:18
const pat::tau::TauPFEssential::Vector & flightLength() const
Definition: Tau.h:321
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float pt_weighted_dr_iso(const reco::PFTau &tau, int dm)
double f[11][100]
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
bool hasSecondaryVertex() const
Definition: Tau.h:320
T min(T a, T b)
Definition: MathUtil.h:58
Analysis-level tau class.
Definition: Tau.h:55
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
double p() const final
magnitude of momentum vector
float pt_weighted_dphi_strip(const reco::PFTau &tau, int dm)
int decayMode() const
reconstructed tau decay mode (specific to PFTau)
Definition: Tau.h:365
float eratio(const reco::PFTau &tau)
return ratio of energy in ECAL over sum of energy in ECAL and HCAL
float flightLengthSig() const
Definition: Tau.h:322
fixed size matrix
float dxy_Sig() const
float leadingTrackNormChi2() const
return normalized chi2 of leading track
Definition: Tau.h:338
float pt_weighted_dx_T(const typename T::Tau_t &tau, int mode, int var, int decaymode)
reco::CandidatePtrVector isolationGammaCands() const
float lead_track_chi2(const reco::PFTau &tau)
return chi2 of the leading track ==> deprecated? <==
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40