CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFTauQualityCutWrapper.cc
Go to the documentation of this file.
2 
3 // ****************************************************
4 // ****** Isolation filters ************************
5 // ****************************************************
6 
7 using namespace reco;
8 
9 void
10 PFTauQualityCutWrapper::isolationChargedObjects(const PFTau& pfTau, const Vertex& pv, std::vector<reco::LeafCandidate>& output)
11 {
12  if( isoQCuts.useTracksInsteadOfPF )
13  {
15  isoQCuts.minTrackPt,
16  isoQCuts.minTrackPixelHits,
17  isoQCuts.minTrackHits,
18  isoQCuts.maxTransverseImpactParameter,
19  isoQCuts.maxTrackChi2,
20  isoQCuts.maxDeltaZ,
21  pv,
22  pv.position().z() ); //????
23 
24  size_t nTracks = result.size();
25  for(size_t iTrack = 0; iTrack < nTracks; ++iTrack)
26  {
27  // this sucks
28  int charge = result[iTrack]->charge();
29  math::XYZVector p3 = result[iTrack]->momentum();
30  reco::Particle::LorentzVector p4(p3.R(), p3.x(), p3.y(), p3.z());
31  output.push_back(reco::LeafCandidate(charge, p4));
32  }
33  } else
34  {
36  isoQCuts.minTrackPt,
37  isoQCuts.minTrackPixelHits,
38  isoQCuts.minTrackHits,
39  isoQCuts.maxTransverseImpactParameter,
40  isoQCuts.maxTrackChi2,
41  10000., // isoQCuts.maxDeltaZ,
42  pv,
43  pv.position().z() ); //????
44 
45 
46  for(unsigned int i=0;i<preresult.size();++i)
47  if(preresult.at(i)->trackRef().isNonnull()) {
48  //get the vertex weight and require to be >50%
49  float w = pv.trackWeight(preresult.at(i)->trackRef());
50  if(w>0.0)
51  output.push_back(reco::LeafCandidate(preresult.at(i)->charge(), preresult[i]->p4()));
52 
53  }
54 
55 
56  }
57 
58 }
59 
60 
61 void
62 PFTauQualityCutWrapper::isolationPUObjects(const PFTau& pfTau, const Vertex& pv, std::vector<reco::LeafCandidate>& output)
63 {
64 
66  isoQCuts.minGammaEt,0,0,1000.,100000,
67  10000., // isoQCuts.maxDeltaZ,
68  pv,
69  pv.position().z() ); //????
70 
71 
72  for(unsigned int i=0;i<preresult.size();++i)
73  if(preresult.at(i)->trackRef().isNonnull()) {
74  //get the vertex weight and require to be >50%
75  float w = pv.trackWeight(preresult.at(i)->trackRef());
76  if(w==0.0)
77  output.push_back(reco::LeafCandidate(preresult.at(i)->charge(), preresult[i]->p4()));
78 
79  }
80 
81 }
82 
83 
84 
85 void
86 PFTauQualityCutWrapper::isolationGammaObjects(const PFTau& pfTau, std::vector<reco::LeafCandidate>& output)
87 {
89 
90  size_t nGammas = result.size();
91  for(size_t iGamma = 0; iGamma < nGammas; ++iGamma)
92  {
93  output.push_back(reco::LeafCandidate(result[iGamma]->charge(), result[iGamma]->p4()));
94  }
95 
96 }
97 
98 // ****************************************************
99 // ****** Signal region filters *********************
100 // ****************************************************
101 
102 void
103 PFTauQualityCutWrapper::signalChargedObjects(const PFTau& pfTau, const Vertex& pv, std::vector<reco::LeafCandidate>& output)
104 {
105  if( signalQCuts.useTracksInsteadOfPF )
106  {
108  signalQCuts.minTrackPt,
109  signalQCuts.minTrackPixelHits,
110  signalQCuts.minTrackHits,
111  signalQCuts.maxTransverseImpactParameter,
112  signalQCuts.maxTrackChi2,
113  signalQCuts.maxDeltaZ,
114  pv,
115  pv.position().z() ); //????
116 
117  size_t nTracks = result.size();
118  for(size_t iTrack = 0; iTrack < nTracks; ++iTrack)
119  {
120  // this sucks
121  int charge = result[iTrack]->charge();
122  math::XYZVector p3 = result[iTrack]->momentum();
123  reco::Particle::LorentzVector p4(p3.R(), p3.x(), p3.y(), p3.z());
124  output.push_back(reco::LeafCandidate(charge, p4));
125  }
126  } else
127  {
128 
129  //First perform basic filtering without vertex dz
131 
132  signalQCuts.minTrackPt,
133  signalQCuts.minTrackPixelHits,
134  signalQCuts.minTrackHits,
135  signalQCuts.maxTransverseImpactParameter,
136  signalQCuts.maxTrackChi2,
137  10000., // signalQCuts.maxDeltaZ,
138  pv,
139  pv.position().z() ); //????
140 
141  //Now check the vertex association
142  for(unsigned int i=0;i<preresult.size();++i)
143  if(preresult.at(i)->trackRef().isNonnull()) {
144  //get the vertex weight and require to be >50%
145  float w = pv.trackWeight(preresult.at(i)->trackRef());
146  if(w>0.0)
147  output.push_back(reco::LeafCandidate(preresult.at(i)->charge(), preresult[i]->p4()));
148 
149  }
150 
151 
152  }
153 
154 }
155 
156 void
157 PFTauQualityCutWrapper::signalGammaObjects(const PFTau& pfTau, std::vector<reco::LeafCandidate>& output)
158 {
160 
161  size_t nGammas = result.size();
162  for(size_t iGamma = 0; iGamma < nGammas; ++iGamma)
163  {
164  output.push_back(reco::LeafCandidate(result[iGamma]->charge(), result[iGamma]->p4()));
165  }
166 
167 }
int i
Definition: DBlmapReader.cc:9
unsigned int nGammas(const GenJet &jet)
void isolationGammaObjects(const reco::PFTau &, std::vector< reco::LeafCandidate > &)
retrieve filtered isolation gamma objects from the pfTau
void signalChargedObjects(const reco::PFTau &, const reco::Vertex &, std::vector< reco::LeafCandidate > &)
retrieve filtered signal charged objects from the pfTau
virtual const reco::TrackRefVector & isolationTracks() const
Definition: BaseTau.cc:30
reco::PFCandidateRefVector filteredPFGammaCands(reco::PFCandidateRefVector theInitialPFCands, double GammaCand_EcalclusMinEt)
Definition: TauTagTools.cc:171
void isolationChargedObjects(const reco::PFTau &, const reco::Vertex &, std::vector< reco::LeafCandidate > &)
retrieve filtered isolation charged objects from the pfTau
reco::TrackRefVector filteredTracks(reco::TrackRefVector theInitialTracks, double tkminPt, int tkminPixelHitsn, int tkminTrackerHitsn, double tkmaxipt, double tkmaxChi2, reco::Vertex pV)
Definition: TauTagTools.cc:77
const PFCandidateRefVector & signalPFChargedHadrCands() const
Charged hadrons in signal region.
Definition: PFTau.cc:75
const Point & position() const
position
Definition: Vertex.h:93
double charge(const std::vector< uint8_t > &Ampls)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
void signalGammaObjects(const reco::PFTau &, std::vector< reco::LeafCandidate > &)
retrieve filtered signal gamma objects from the pfTau
void isolationPUObjects(const reco::PFTau &, const reco::Vertex &, std::vector< reco::LeafCandidate > &)
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
double p4[4]
Definition: TauolaWrapper.h:92
const PFCandidateRefVector & isolationPFGammaCands() const
Gamma candidates in isolation region.
Definition: PFTau.cc:88
tuple result
Definition: query.py:137
virtual const reco::TrackRefVector & signalTracks() const
Definition: BaseTau.cc:28
reco::PFCandidateRefVector filteredPFChargedHadrCands(reco::PFCandidateRefVector theInitialPFCands, double ChargedHadrCand_tkminPt, int ChargedHadrCand_tkminPixelHitsn, int ChargedHadrCand_tkminTrackerHitsn, double ChargedHadrCand_tkmaxipt, double ChargedHadrCand_tkmaxChi2, reco::Vertex pV)
Definition: TauTagTools.cc:119
const PFCandidateRefVector & signalPFGammaCands() const
Gamma candidates in signal region.
Definition: PFTau.cc:79
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
value_type const at(size_type idx) const
Retrieve an element of the RefVector.
Definition: RefVector.h:76
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
T w() const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
const PFCandidateRefVector & isolationPFChargedHadrCands() const
Charged candidates in isolation region.
Definition: PFTau.cc:84
double p3[4]
Definition: TauolaWrapper.h:91