13 const std::vector<reco::PFCandidatePtr>& getPFGammas(
const reco::PFTau&
tau,
bool signal =
true) {
19 bool isInside(
float photon_pt,
float deta,
float dphi){
21 double stripEtaAssociationDistance_0p95_p0 = 0.197077;
22 double stripEtaAssociationDistance_0p95_p1 = 0.658701;
24 double stripPhiAssociationDistance_0p95_p0 = 0.352476;
25 double stripPhiAssociationDistance_0p95_p1 = 0.707716;
27 if(photon_pt==0){
return false;}
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)))))
44 Float_t LeadingTracknormalizedChi2 = 0;
50 LeadingTracknormalizedChi2 = (float)(tref -> normalizedChi2());
54 return LeadingTracknormalizedChi2;
59 std::vector<reco::PFCandidatePtr> constsignal = tau.
signalPFCands();
60 Float_t EcalEnInSignalPFCands = 0;
61 Float_t HcalEnInSignalPFCands = 0;
63 typedef std::vector <reco::PFCandidatePtr>::iterator constituents_iterator;
64 for(constituents_iterator it=constsignal.begin(); it != constsignal.end(); ++it) {
66 EcalEnInSignalPFCands += icand -> ecalEnergy();
67 HcalEnInSignalPFCands += icand -> hcalEnergy();
70 Float_t
total = EcalEnInSignalPFCands + HcalEnInSignalPFCands;
71 if(total==0)
return -1;
72 else return EcalEnInSignalPFCands/
total;
75 float pt_weighted_dx(
const reco::PFTau& tau,
int mode = 0,
int var = 0,
int decaymode = -1){
80 int is3prong = (decaymode==10);
82 auto& cands = getPFGammas(tau,
mode < 2);
84 for (
auto& cand : cands) {
92 float pt = cand->pt();
94 bool flag = isInside(pt, deta, dphi);
97 if (
mode == 2 || (
mode == 0 && dr < signalrad) || (
mode == 1 && dr > signalrad)) {
100 sum_dx_pt += pt * dr;
102 sum_dx_pt += pt * deta;
104 sum_dx_pt += pt * dphi;
106 }
else if(is3prong==1){
108 if( (
mode==2 && flag==
false) || (
mode==1 && flag==
true) ||
mode==0){
112 sum_dx_pt += pt * dr;
114 sum_dx_pt += pt * deta;
116 sum_dx_pt += pt * dphi;
122 return sum_dx_pt/sum_pt;
126 float tau_pt_weighted_dr_signal(
const reco::PFTau& tau,
int dm) {
127 return pt_weighted_dx(tau, 0, 0, dm);
130 float tau_pt_weighted_deta_strip(
const reco::PFTau& tau,
int dm) {
132 return pt_weighted_dx(tau, 2, 1, dm);
134 return pt_weighted_dx(tau, 1, 1, dm);
138 float tau_pt_weighted_dphi_strip(
const reco::PFTau& tau,
int dm) {
140 return pt_weighted_dx(tau, 2, 2, dm);
142 return pt_weighted_dx(tau, 1, 2, dm);
146 float tau_pt_weighted_dr_iso(
const reco::PFTau& tau,
int dm) {
147 return pt_weighted_dx(tau, 2, 0, dm);
150 unsigned int tau_n_photons_total(
const reco::PFTau& tau) {
151 unsigned int n_photons = 0;
153 if (cand->pt() > 0.5)
157 if (cand->pt() > 0.5)
bool isNonnull() const
Checks for non-null.
const PFCandidatePtr & leadPFChargedHadrCand() const
const std::vector< reco::PFCandidatePtr > & signalPFGammaCands() const
Gamma candidates in signal region.
virtual double phi() const final
momentum azimuthal angle
const std::vector< reco::PFCandidatePtr > & signalPFCands() const
PFCandidates in signal region.
const std::vector< reco::PFCandidatePtr > & isolationPFGammaCands() const
Gamma candidates in isolation region.
Abs< T >::type abs(const T &t)
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
bool isNonnull() const
Checks for non-null.
double deltaPhi(double phi1, double phi2)
virtual double eta() const final
momentum pseudorapidity
Power< A, B >::type pow(const A &a, const B &b)
virtual double pt() const final
transverse momentum