test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecoTauClusterVariables.h
Go to the documentation of this file.
1 
10 namespace
11 {
12 
13  const std::vector<reco::PFCandidatePtr>& getPFGammas(const reco::PFTau& tau, bool signal = true) {
14  if (signal)
15  return tau.signalPFGammaCands();
16  return tau.isolationPFGammaCands();
17  }
18 
19  bool isInside(float photon_pt, float deta, float dphi){
20 
21  double stripEtaAssociationDistance_0p95_p0 = 0.197077;
22  double stripEtaAssociationDistance_0p95_p1 = 0.658701;
23 
24  double stripPhiAssociationDistance_0p95_p0 = 0.352476;
25  double stripPhiAssociationDistance_0p95_p1 = 0.707716;
26 
27  if(photon_pt==0){return false;}
28 
29  if(
30  (dphi < std::min(0.3, std::max(0.05, stripPhiAssociationDistance_0p95_p0*std::pow(photon_pt, -(stripPhiAssociationDistance_0p95_p1))))) && \
31  (deta < std::min(0.15, std::max(0.05, stripEtaAssociationDistance_0p95_p0*std::pow(photon_pt, -(stripEtaAssociationDistance_0p95_p1)))))
32  ){
33  return true;
34  }
35 
36  return false;
37 
38  }
39 
40  float tau_leadTrackChi2(const reco::PFTau& tau){
41 
42  // leading charged hadron PFCand in signal cone
43 
44  Float_t LeadingTracknormalizedChi2 = 0;
45 
46  const reco::PFCandidatePtr& leadingPFCharged = tau.leadPFChargedHadrCand() ;
47  if ( leadingPFCharged.isNonnull() ) {
48  reco::TrackRef tref = leadingPFCharged -> trackRef();
49  if ( tref.isNonnull() ) {
50  LeadingTracknormalizedChi2 = (float)(tref -> normalizedChi2());
51  }
52  }
53 
54  return LeadingTracknormalizedChi2;
55  }
56 
57  float tau_Eratio(const reco::PFTau& tau){
58 
59  std::vector<reco::PFCandidatePtr> constsignal = tau.signalPFCands();
60  Float_t EcalEnInSignalPFCands = 0;
61  Float_t HcalEnInSignalPFCands = 0;
62 
63  typedef std::vector <reco::PFCandidatePtr>::iterator constituents_iterator;
64  for(constituents_iterator it=constsignal.begin(); it != constsignal.end(); ++it) {
65  reco::PFCandidatePtr & icand = *it;
66  EcalEnInSignalPFCands += icand -> ecalEnergy();
67  HcalEnInSignalPFCands += icand -> hcalEnergy();
68  }
69 
70  Float_t total = EcalEnInSignalPFCands + HcalEnInSignalPFCands;
71  if(total==0) return -1;
72  else return EcalEnInSignalPFCands/total;
73  }
74 
75  float pt_weighted_dx(const reco::PFTau& tau, int mode = 0, int var = 0, int decaymode = -1){
76 
77  float sum_pt = 0.;
78  float sum_dx_pt = 0.;
79  float signalrad = std::max(0.05, std::min(0.1, 3./std::max(1., tau.pt())));
80  int is3prong = (decaymode==10);
81 
82  auto& cands = getPFGammas(tau, mode < 2);
83 
84  for (auto& cand : cands) {
85  // only look at electrons/photons with pT > 0.5
86  if (cand->pt() < 0.5)
87  continue;
88 
89  float dr = reco::deltaR(*cand, tau);
90  float deta = std::abs(cand->eta() - tau.eta());
91  float dphi = std::abs(reco::deltaPhi(cand->phi(), tau.phi()));
92  float pt = cand->pt();
93 
94  bool flag = isInside(pt, deta, dphi);
95 
96  if(is3prong==0){
97  if (mode == 2 || (mode == 0 && dr < signalrad) || (mode == 1 && dr > signalrad)) {
98  sum_pt += pt;
99  if (var == 0)
100  sum_dx_pt += pt * dr;
101  else if (var == 1)
102  sum_dx_pt += pt * deta;
103  else if (var == 2)
104  sum_dx_pt += pt * dphi;
105  }
106  }else if(is3prong==1){
107 
108  if( (mode==2 && flag==false) || (mode==1 && flag==true) || mode==0){
109  sum_pt += pt;
110 
111  if (var == 0)
112  sum_dx_pt += pt * dr;
113  else if (var == 1)
114  sum_dx_pt += pt * deta;
115  else if (var == 2)
116  sum_dx_pt += pt * dphi;
117  }
118  }
119  }
120 
121  if (sum_pt > 0.)
122  return sum_dx_pt/sum_pt;
123  return 0.;
124  }
125 
126  float tau_pt_weighted_dr_signal(const reco::PFTau& tau, int dm) {
127  return pt_weighted_dx(tau, 0, 0, dm);
128  }
129 
130  float tau_pt_weighted_deta_strip(const reco::PFTau& tau, int dm) {
131  if(dm==10){
132  return pt_weighted_dx(tau, 2, 1, dm);
133  }else{
134  return pt_weighted_dx(tau, 1, 1, dm);
135  }
136  }
137 
138  float tau_pt_weighted_dphi_strip(const reco::PFTau& tau, int dm) {
139  if(dm==10){
140  return pt_weighted_dx(tau, 2, 2, dm);
141  }else{
142  return pt_weighted_dx(tau, 1, 2, dm);
143  }
144  }
145 
146  float tau_pt_weighted_dr_iso(const reco::PFTau& tau, int dm) {
147  return pt_weighted_dx(tau, 2, 0, dm);
148  }
149 
150  unsigned int tau_n_photons_total(const reco::PFTau& tau) {
151  unsigned int n_photons = 0;
152  for (auto& cand : tau.signalPFGammaCands()) {
153  if (cand->pt() > 0.5)
154  ++n_photons;
155  }
156  for (auto& cand : tau.isolationPFGammaCands()) {
157  if (cand->pt() > 0.5)
158  ++n_photons;
159  }
160  return n_photons;
161  }
162 
163 }
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
const PFCandidatePtr & leadPFChargedHadrCand() const
Definition: PFTau.cc:67
const std::vector< reco::PFCandidatePtr > & signalPFGammaCands() const
Gamma candidates in signal region.
Definition: PFTau.cc:84
virtual double phi() const final
momentum azimuthal angle
const std::vector< reco::PFCandidatePtr > & signalPFCands() const
PFCandidates in signal region.
Definition: PFTau.cc:78
const std::vector< reco::PFCandidatePtr > & isolationPFGammaCands() const
Gamma candidates in isolation region.
Definition: PFTau.cc:93
tuple dm
Definition: symbols.py:65
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
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:169
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
virtual double eta() const final
momentum pseudorapidity
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
virtual double pt() const final
transverse momentum