CMS 3D CMS Logo

PFRecoTauClusterVariables.h
Go to the documentation of this file.
1 
12  public:
18  float tau_leadTrackChi2(const reco::PFTau& tau) const {
19  float LeadingTracknormalizedChi2 = 0;
20  const reco::PFCandidatePtr& leadingPFCharged = tau.leadPFChargedHadrCand() ;
21  if (leadingPFCharged.isNonnull()) {
22  reco::TrackRef tref = leadingPFCharged -> trackRef();
23  if (tref.isNonnull()) {
24  LeadingTracknormalizedChi2 = (float)(tref -> normalizedChi2());
25  }
26  }
27  return LeadingTracknormalizedChi2;
28  }
30  float tau_Eratio(const reco::PFTau& tau) const {
31  std::vector<reco::PFCandidatePtr> constsignal = tau.signalPFCands();
32  float EcalEnInSignalPFCands = 0;
33  float HcalEnInSignalPFCands = 0;
34  typedef std::vector <reco::PFCandidatePtr>::iterator constituents_iterator;
35  for(constituents_iterator it=constsignal.begin(); it != constsignal.end(); ++it) {
36  reco::PFCandidatePtr & icand = *it;
37  EcalEnInSignalPFCands += icand -> ecalEnergy();
38  HcalEnInSignalPFCands += icand -> hcalEnergy();
39  }
40  float total = EcalEnInSignalPFCands + HcalEnInSignalPFCands;
41  if(total==0){
42  return -1;
43  }
44  return EcalEnInSignalPFCands/total;
45  }
46  float tau_Eratio(const pat::Tau& tau) const {
47  float EcalEnInSignalCands = tau.ecalEnergy();
48  float HcalEnInSignalCands = tau.hcalEnergy();
49  float total = EcalEnInSignalCands + HcalEnInSignalCands;
50  if(total == 0){
51  return -1;
52  }
53  return EcalEnInSignalCands/total;
54  }
57  float pt_weighted_dx(const reco::PFTau& tau, int mode = 0, int var = 0, int decaymode = -1) const {
58  float sum_pt = 0.;
59  float sum_dx_pt = 0.;
60  float signalrad = std::max(0.05, std::min(0.1, 3./std::max(1., tau.pt())));
61  int is3prong = (decaymode==10);
62  const auto& cands = getPFGammas(tau, mode < 2);
63  for (const auto& cand : cands) {
64  // only look at electrons/photons with pT > 0.5
65  if ((float)cand->pt() < 0.5){
66  continue;
67  }
68  float dr = reco::deltaR((float)cand->eta(),(float)cand->phi(),(float)tau.eta(),(float)tau.phi());
69  float deta = std::abs((float)cand->eta() - (float)tau.eta());
70  float dphi = std::abs(reco::deltaPhi((float)cand->phi(), (float)tau.phi()));
71  float pt = cand->pt();
72  bool flag = isInside(pt, deta, dphi);
73  if(is3prong==0){
74  if (mode == 2 || (mode == 0 && dr < signalrad) || (mode == 1 && dr > signalrad)) {
75  sum_pt += pt;
76  if (var == 0)
77  sum_dx_pt += pt * dr;
78  else if (var == 1)
79  sum_dx_pt += pt * deta;
80  else if (var == 2)
81  sum_dx_pt += pt * dphi;
82  }
83  }
84  else if(is3prong==1){
85  if( (mode==2 && flag==false) || (mode==1 && flag==true) || mode==0){
86  sum_pt += pt;
87  if (var == 0)
88  sum_dx_pt += pt * dr;
89  else if (var == 1)
90  sum_dx_pt += pt * deta;
91  else if (var == 2)
92  sum_dx_pt += pt * dphi;
93  }
94  }
95  }
96  if (sum_pt > 0.){
97  return sum_dx_pt/sum_pt;
98  }
99  return 0.;
100  }
101  float pt_weighted_dx(const pat::Tau& tau, int mode = 0, int var = 0, int decaymode = -1) const {
102  float sum_pt = 0.;
103  float sum_dx_pt = 0.;
104  float signalrad = std::max(0.05, std::min(0.1, 3./std::max(1., tau.pt())));
105  int is3prong = (decaymode==10);
106  const auto cands = getGammas(tau, mode < 2);
107  for (const auto& cand : cands) {
108  // only look at electrons/photons with pT > 0.5
109  if (cand->pt() < 0.5){
110  continue;
111  }
112  float dr = reco::deltaR(*cand, tau);
113  float deta = std::abs(cand->eta() - tau.eta());
114  float dphi = std::abs(reco::deltaPhi(cand->phi(), tau.phi()));
115  float pt = cand->pt();
116  bool flag = isInside(pt, deta, dphi);
117  if(is3prong==0){
118  if (mode == 2 || (mode == 0 && dr < signalrad) || (mode == 1 && dr > signalrad)) {
119  sum_pt += pt;
120  if (var == 0)
121  sum_dx_pt += pt * dr;
122  else if (var == 1)
123  sum_dx_pt += pt * deta;
124  else if (var == 2)
125  sum_dx_pt += pt * dphi;
126  }
127  }
128  else if(is3prong==1){
129  if( (mode==2 && flag==false) || (mode==1 && flag==true) || mode==0){
130  sum_pt += pt;
131  if (var == 0)
132  sum_dx_pt += pt * dr;
133  else if (var == 1)
134  sum_dx_pt += pt * deta;
135  else if (var == 2)
136  sum_dx_pt += pt * dphi;
137  }
138  }
139  }
140  if (sum_pt > 0.){
141  return sum_dx_pt/sum_pt;
142  }
143  return 0.;
144  }
147  float tau_pt_weighted_dr_signal(const reco::PFTau& tau, int dm) const {
148  return pt_weighted_dx(tau, 0, 0, dm);
149  }
150  float tau_pt_weighted_dr_signal(const pat::Tau& tau, int dm) const {
151  return pt_weighted_dx(tau, 0, 0, dm);
152  }
155  float tau_pt_weighted_deta_strip(const reco::PFTau& tau, int dm) const {
156  return pt_weighted_dx(tau, dm==10 ? 2 : 1, 1, dm);
157  }
158  float tau_pt_weighted_deta_strip(const pat::Tau& tau, int dm) const {
159  return pt_weighted_dx(tau, dm==10 ? 2 : 1, 1, dm);
160  }
163  float tau_pt_weighted_dphi_strip(const reco::PFTau& tau, int dm) const {
164  return pt_weighted_dx(tau, dm==10 ? 2 : 1, 2, dm);
165  }
166  float tau_pt_weighted_dphi_strip(const pat::Tau& tau, int dm) const {
167  return pt_weighted_dx(tau, dm==10 ? 2 : 1, 2, dm);
168  }
171  float tau_pt_weighted_dr_iso(const reco::PFTau& tau, int dm) const {
172  return pt_weighted_dx(tau, 2, 0, dm);
173  }
174  float tau_pt_weighted_dr_iso(const pat::Tau& tau, int dm) const {
175  return pt_weighted_dx(tau, 2, 0, dm);
176  }
178  unsigned int tau_n_photons_total(const reco::PFTau& tau) const {
179  unsigned int n_photons = 0;
180  for (auto& cand : tau.signalPFGammaCands()) {
181  if ((float)cand->pt() > 0.5)
182  ++n_photons;
183  }
184  for (auto& cand : tau.isolationPFGammaCands()) {
185  if ((float)cand->pt() > 0.5)
186  ++n_photons;
187  }
188  return n_photons;
189  }
190  unsigned int tau_n_photons_total(const pat::Tau& tau) const {
191  unsigned int n_photons = 0;
192  for (auto& cand : tau.signalGammaCands()) {
193  if (cand->pt() > 0.5)
194  ++n_photons;
195  }
196  for (auto& cand : tau.isolationGammaCands()) {
197  if (cand->pt() > 0.5)
198  ++n_photons;
199  }
200  return n_photons;
201  }
202 
203  private:
205  const std::vector<reco::PFCandidatePtr>& getPFGammas(const reco::PFTau& tau, bool signal = true) const {
206  if (signal){
207  return tau.signalPFGammaCands();
208  }
209  return tau.isolationPFGammaCands();
210  }
211  reco::CandidatePtrVector getGammas(const pat::Tau& tau, bool signal = true) const {
212  if(signal){
213  return tau.signalGammaCands();
214  }
215  return tau.isolationGammaCands();
216  }
218  bool isInside(float photon_pt, float deta, float dphi) const {
219  const double stripEtaAssociationDistance_0p95_p0 = 0.197077;
220  const double stripEtaAssociationDistance_0p95_p1 = 0.658701;
221  const double stripPhiAssociationDistance_0p95_p0 = 0.352476;
222  const double stripPhiAssociationDistance_0p95_p1 = 0.707716;
223  if(photon_pt==0){
224  return false;
225  }
226  if((dphi<0.3 && dphi<std::max(0.05, stripPhiAssociationDistance_0p95_p0*std::pow(photon_pt, -stripPhiAssociationDistance_0p95_p1))) && \
227  (deta<0.15 && deta<std::max(0.05, stripEtaAssociationDistance_0p95_p0*std::pow(photon_pt, -stripEtaAssociationDistance_0p95_p1)))){
228  return true;
229  }
230  return false;
231  }
232 };
float tau_pt_weighted_dr_signal(const pat::Tau &tau, int dm) const
virtual double pt() const final
transverse momentum
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
float tau_Eratio(const reco::PFTau &tau) const
return ratio of energy in ECAL over sum of energy in ECAL and HCAL
const PFCandidatePtr & leadPFChargedHadrCand() const
Definition: PFTau.cc:67
float tau_pt_weighted_deta_strip(const pat::Tau &tau, int dm) const
virtual double eta() const final
momentum pseudorapidity
const std::vector< reco::PFCandidatePtr > & getPFGammas(const reco::PFTau &tau, bool signal=true) const
return pf photon candidates that are associated to signal
reco::CandidatePtrVector getGammas(const pat::Tau &tau, bool signal=true) const
float tau_pt_weighted_dphi_strip(const pat::Tau &tau, int dm) const
const std::vector< reco::PFCandidatePtr > & signalPFGammaCands() const
Gamma candidates in signal region.
Definition: PFTau.cc:84
~TauIdMVAAuxiliaries()
default destructor
TauIdMVAAuxiliaries()
default constructor
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
const std::vector< reco::PFCandidatePtr > & isolationPFGammaCands() const
Gamma candidates in isolation region.
Definition: PFTau.cc:93
bool isInside(float photon_pt, float deta, float dphi) const
decide if photon candidate is inside the cone to be associated to the tau signal
virtual double phi() const final
momentum azimuthal angle
float tau_pt_weighted_dphi_strip(const reco::PFTau &tau, int dm) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
T min(T a, T b)
Definition: MathUtil.h:58
float tau_pt_weighted_dr_signal(const reco::PFTau &tau, int dm) const
float tau_leadTrackChi2(const reco::PFTau &tau) const
return chi2 of the leading track ==> deprecated? <==
Analysis-level tau class.
Definition: Tau.h:55
float tau_pt_weighted_dr_iso(const reco::PFTau &tau, int dm) const
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:169
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
float tau_pt_weighted_deta_strip(const reco::PFTau &tau, int dm) const
float pt_weighted_dx(const reco::PFTau &tau, int mode=0, int var=0, int decaymode=-1) const
float pt_weighted_dx(const pat::Tau &tau, int mode=0, int var=0, int decaymode=-1) const
reco::CandidatePtrVector isolationGammaCands() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
float tau_Eratio(const pat::Tau &tau) const
float tau_pt_weighted_dr_iso(const pat::Tau &tau, int dm) const
unsigned int tau_n_photons_total(const pat::Tau &tau) const
unsigned int tau_n_photons_total(const reco::PFTau &tau) const
return total number of pf photon candidates with pT>500 MeV, which are associated to signal ...