CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions
reco::tau::qcuts Namespace Reference

Functions

bool AND (const TrackBaseRef &track, const RecoTauQualityCuts::TrackQCutFuncCollection &cuts)
 
bool AND_cand (const PFCandidate &cand, const RecoTauQualityCuts::CandQCutFuncCollection &cuts)
 
bool etMin_cand (const PFCandidate &cand, double cut)
 
bool mapAndCutByType (const PFCandidate &cand, const RecoTauQualityCuts::CandQCutFuncMap &funcMap)
 
bool minTrackVertexWeight (const TrackBaseRef &track, const reco::VertexRef *pv, double cut)
 
bool minTrackVertexWeight_cand (const PFCandidate &cand, const reco::VertexRef *pv, double cut)
 
bool ptMin (const TrackBaseRef &track, double cut)
 
bool ptMin_cand (const PFCandidate &cand, double cut)
 
bool trkChi2 (const TrackBaseRef &track, double cut)
 
bool trkChi2_cand (const PFCandidate &cand, double cut)
 
bool trkLongitudinalImpactParameter (const TrackBaseRef &track, const reco::VertexRef *pv, double cut)
 
bool trkLongitudinalImpactParameter_cand (const PFCandidate &cand, const reco::VertexRef *pv, double cut)
 
bool trkLongitudinalImpactParameterWrtTrack (const TrackBaseRef &track, const reco::TrackBaseRef *leadTrack, const reco::VertexRef *pv, double cut)
 DZ cut, with respect to the current lead rack. More...
 
bool trkLongitudinalImpactParameterWrtTrack_cand (const PFCandidate &cand, const reco::TrackBaseRef *leadTrack, const reco::VertexRef *pv, double cut)
 
bool trkPixelHits (const TrackBaseRef &track, int cut)
 
bool trkPixelHits_cand (const PFCandidate &cand, int cut)
 
bool trkTrackerHits (const TrackBaseRef &track, int cut)
 
bool trkTrackerHits_cand (const PFCandidate &cand, int cut)
 
bool trkTransverseImpactParameter (const TrackBaseRef &track, const reco::VertexRef *pv, double cut)
 
bool trkTransverseImpactParameter_cand (const PFCandidate &cand, const reco::VertexRef *pv, double cut)
 

Function Documentation

bool reco::tau::qcuts::AND ( const TrackBaseRef &  track,
const RecoTauQualityCuts::TrackQCutFuncCollection &  cuts 
)

Definition at line 192 of file RecoTauQualityCuts.cc.

References statics::func.

Referenced by egHLT::EgObjTrigCut< T >::pass(), egHLT::EgEvtTrigCut< T >::pass(), and reco::tau::RecoTauQualityCuts::RecoTauQualityCuts().

193 {
194  BOOST_FOREACH( const RecoTauQualityCuts::TrackQCutFunc& func, cuts ) {
195  if ( !func(track) ) return false;
196  }
197  return true;
198 }
string func
Definition: statics.py:48
bool reco::tau::qcuts::AND_cand ( const PFCandidate &  cand,
const RecoTauQualityCuts::CandQCutFuncCollection &  cuts 
)

Definition at line 200 of file RecoTauQualityCuts.cc.

References statics::func.

Referenced by mapAndCutByType().

201 {
202  BOOST_FOREACH( const RecoTauQualityCuts::CandQCutFunc& func, cuts ) {
203  if ( !func(cand) ) return false;
204  }
205  return true;
206 }
string func
Definition: statics.py:48
bool reco::tau::qcuts::etMin_cand ( const PFCandidate &  cand,
double  cut 
)

Definition at line 36 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, and reco::LeafCandidate::et().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts().

37 {
38  //std::cout << "<etMin_cand>: Et = " << cand.et() << ", cut = " << cut << std::endl;
39  return (cand.et() > cut);
40 }
bool reco::tau::qcuts::mapAndCutByType ( const PFCandidate &  cand,
const RecoTauQualityCuts::CandQCutFuncMap &  funcMap 
)

Definition at line 209 of file RecoTauQualityCuts.cc.

References AND_cand(), hpstanc_transforms::cuts, and reco::PFCandidate::particleId().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts().

210 {
211  // Find the cuts that for this particle type
212  RecoTauQualityCuts::CandQCutFuncMap::const_iterator cuts = funcMap.find(cand.particleId());
213  // Return false if we dont' know how to deal with this particle type
214  if ( cuts == funcMap.end() ) return false;
215  return AND_cand(cand, cuts->second); // Otherwise AND all the cuts
216 }
bool AND_cand(const PFCandidate &cand, const RecoTauQualityCuts::CandQCutFuncCollection &cuts)
bool reco::tau::qcuts::minTrackVertexWeight ( const TrackBaseRef &  track,
const reco::VertexRef pv,
double  cut 
)

Definition at line 148 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, and edm::Ref< C, T, F >::isNull().

Referenced by minTrackVertexWeight_cand(), and reco::tau::RecoTauQualityCuts::RecoTauQualityCuts().

149 {
150  if ( pv->isNull() ) {
151  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
152  "RecoTauQualityCuts is invalid. - minTrackVertexWeight";
153  return false;
154  }
155  //std::cout << "<minTrackVertexWeight>:" << std::endl;
156  //std::cout << " track: Pt = " << track->pt() << ", eta = " << track->eta() << ", phi = " << track->phi() << std::endl;
157  //std::cout << " vertex: x = " << (*pv)->position().x() << ", y = " << (*pv)->position().y() << ", z = " << (*pv)->position().z() << std::endl;
158  //std::cout << "--> trackWeight = " << (*pv)->trackWeight(track) << " (cut = " << cut << ")" << std::endl;
159  return ((*pv)->trackWeight(track) >= cut);
160 }
bool isNull() const
Checks for null.
Definition: Ref.h:247
bool reco::tau::qcuts::minTrackVertexWeight_cand ( const PFCandidate &  cand,
const reco::VertexRef pv,
double  cut 
)

Definition at line 162 of file RecoTauQualityCuts.cc.

References minTrackVertexWeight().

163 {
164  auto track = getTrackRef(cand);
165  if ( track.isNonnull() ) {
166  return minTrackVertexWeight(track, pv, cut);
167  } else {
168  //std::cout << "<minTrackVertexWeight_cand>: weight = N/A, cut = " << cut << std::endl;
169  return false;
170  }
171 }
bool minTrackVertexWeight(const TrackBaseRef &track, const reco::VertexRef *pv, double cut)
bool reco::tau::qcuts::ptMin ( const TrackBaseRef &  track,
double  cut 
)

Definition at line 24 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, and reco::TrackBase::pt().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts().

25 {
26  //std::cout << "<ptMin>: Pt = " << track->pt() << ", cut = " << cut << std::endl;
27  return (track->pt() > cut);
28 }
bool reco::tau::qcuts::ptMin_cand ( const PFCandidate &  cand,
double  cut 
)

Definition at line 30 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, and reco::LeafCandidate::pt().

31 {
32  //std::cout << "<ptMin_cand>: Pt = " << cand.pt() << ", cut = " << cut << std::endl;
33  return (cand.pt() > cut);
34 }
bool reco::tau::qcuts::trkChi2 ( const TrackBaseRef &  track,
double  cut 
)

Definition at line 173 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, and reco::TrackBase::normalizedChi2().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts(), and trkChi2_cand().

174 {
175  //std::cout << "<trkChi2>: chi^2 = " << track->normalizedChi2() << ", cut = " << cut << std::endl;
176  return (track->normalizedChi2() <= cut);
177 }
bool reco::tau::qcuts::trkChi2_cand ( const PFCandidate &  cand,
double  cut 
)

Definition at line 179 of file RecoTauQualityCuts.cc.

References trkChi2().

180 {
181  auto track = getTrackRef(cand);
182  if ( track.isNonnull() ) {
183  //std::cout << "<trkChi2_cand>: chi^2 = " << track->normalizedChi2() << ", cut = " << cut << std::endl;
184  return trkChi2(track, cut);
185  } else {
186  //std::cout << "<trkChi2_cand>: chi^2 = N/A, cut = " << cut << std::endl;
187  return false;
188  }
189 }
bool trkChi2(const TrackBaseRef &track, double cut)
bool reco::tau::qcuts::trkLongitudinalImpactParameter ( const TrackBaseRef &  track,
const reco::VertexRef pv,
double  cut 
)

Definition at line 105 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, reco::TrackBase::dz(), and edm::Ref< C, T, F >::isNull().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts(), and trkLongitudinalImpactParameter_cand().

106 {
107  if ( pv->isNull() ) {
108  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
109  "RecoTauQualityCuts is invalid. - trkLongitudinalImpactParameter";
110  return false;
111  }
112  //std::cout << "<trkLongitudinalImpactParameter>:" << std::endl;
113  //std::cout << " track: Pt = " << track->pt() << ", eta = " << track->eta() << ", phi = " << track->phi() << std::endl;
114  //std::cout << " vertex: x = " << (*pv)->position().x() << ", y = " << (*pv)->position().y() << ", z = " << (*pv)->position().z() << std::endl;
115  //std::cout << "--> dz = " << std::fabs(track->dz((*pv)->position())) << " (cut = " << cut << ")" << std::endl;
116  return (std::fabs(track->dz((*pv)->position())) <= cut);
117 }
bool isNull() const
Checks for null.
Definition: Ref.h:247
bool reco::tau::qcuts::trkLongitudinalImpactParameter_cand ( const PFCandidate &  cand,
const reco::VertexRef pv,
double  cut 
)

Definition at line 119 of file RecoTauQualityCuts.cc.

References trkLongitudinalImpactParameter().

120 {
121  auto track = getTrackRef(cand);
122  if ( track.isNonnull() ) {
123  return trkLongitudinalImpactParameter(track, pv, cut);
124  } else {
125  //std::cout << "<trkLongitudinalImpactParameter_cand>: dZ = N/A, cut = " << cut << std::endl;
126  return false;
127  }
128 }
bool trkLongitudinalImpactParameter(const TrackBaseRef &track, const reco::VertexRef *pv, double cut)
bool reco::tau::qcuts::trkLongitudinalImpactParameterWrtTrack ( const TrackBaseRef &  track,
const reco::TrackBaseRef leadTrack,
const reco::VertexRef pv,
double  cut 
)

DZ cut, with respect to the current lead rack.

Definition at line 131 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, and reco::TrackBase::dz().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts(), and trkLongitudinalImpactParameterWrtTrack_cand().

132 {
133  if ( leadTrack->isNull()) {
134  edm::LogError("QCutsNoValidLeadTrack") << "Lead track Ref in " <<
135  "RecoTauQualityCuts is invalid. - trkLongitudinalImpactParameterWrtTrack";
136  return false;
137  }
138  return (std::fabs(track->dz((*pv)->position()) - (*leadTrack)->dz((*pv)->position())) <= cut);
139 }
bool reco::tau::qcuts::trkLongitudinalImpactParameterWrtTrack_cand ( const PFCandidate &  cand,
const reco::TrackBaseRef leadTrack,
const reco::VertexRef pv,
double  cut 
)

Definition at line 141 of file RecoTauQualityCuts.cc.

References trkLongitudinalImpactParameterWrtTrack().

142 {
143  auto track = getTrackRef(cand);
144  if ( track.isNonnull() ) return trkLongitudinalImpactParameterWrtTrack(track, leadTrack, pv, cut);
145  else return false;
146 }
bool trkLongitudinalImpactParameterWrtTrack(const TrackBaseRef &track, const reco::TrackBaseRef *leadTrack, const reco::VertexRef *pv, double cut)
DZ cut, with respect to the current lead rack.
bool reco::tau::qcuts::trkPixelHits ( const TrackBaseRef &  track,
int  cut 
)

Definition at line 42 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, reco::TrackBase::hitPattern(), and reco::HitPattern::numberOfValidPixelHits().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts(), and trkPixelHits_cand().

43 {
44  // For some reason, the number of hits is signed
45  //std::cout << "<trkPixelHits>: #Pxl hits = " << track->hitPattern().numberOfValidPixelHits() << ", cut = " << cut << std::endl;
46  return (track->hitPattern().numberOfValidPixelHits() >= cut);
47 }
bool reco::tau::qcuts::trkPixelHits_cand ( const PFCandidate &  cand,
int  cut 
)

Definition at line 49 of file RecoTauQualityCuts.cc.

References trkPixelHits().

50 {
51  // For some reason, the number of hits is signed
52  auto track = getTrackRef(cand);
53  if ( track.isNonnull() ) {
54  //std::cout << "<trkPixelHits_cand>: #Pxl hits = " << trkPixelHits(track, cut) << ", cut = " << cut << std::endl;
55  return trkPixelHits(track, cut);
56  } else {
57  //std::cout << "<trkPixelHits_cand>: #Pxl hits = N/A, cut = " << cut << std::endl;
58  return false;
59  }
60 }
bool trkPixelHits(const TrackBaseRef &track, int cut)
bool reco::tau::qcuts::trkTrackerHits ( const TrackBaseRef &  track,
int  cut 
)

Definition at line 62 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, reco::TrackBase::hitPattern(), and reco::HitPattern::numberOfValidHits().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts(), and trkTrackerHits_cand().

63 {
64  //std::cout << "<trkTrackerHits>: #Trk hits = " << track->hitPattern().numberOfValidHits() << ", cut = " << cut << std::endl;
65  return (track->hitPattern().numberOfValidHits() >= cut);
66 }
bool reco::tau::qcuts::trkTrackerHits_cand ( const PFCandidate &  cand,
int  cut 
)

Definition at line 68 of file RecoTauQualityCuts.cc.

References trkTrackerHits().

69 {
70  auto track = getTrackRef(cand);
71  if ( track.isNonnull() ) {
72  //std::cout << "<trkTrackerHits>: #Trk hits = " << track->hitPattern().numberOfValidHits() << ", cut = " << cut << std::endl;
73  return trkTrackerHits(track, cut);
74  } else {
75  //std::cout << "<trkTrackerHits>: #Trk hits = N/A, cut = " << cut << std::endl;
76  return false;
77  }
78 }
bool trkTrackerHits(const TrackBaseRef &track, int cut)
bool reco::tau::qcuts::trkTransverseImpactParameter ( const TrackBaseRef &  track,
const reco::VertexRef pv,
double  cut 
)

Definition at line 80 of file RecoTauQualityCuts.cc.

References GOODCOLL_filter_cfg::cut, reco::TrackBase::dxy(), and edm::Ref< C, T, F >::isNull().

Referenced by reco::tau::RecoTauQualityCuts::RecoTauQualityCuts(), and trkTransverseImpactParameter_cand().

81 {
82  if ( pv->isNull() ) {
83  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
84  "RecoTauQualityCuts is invalid. - trkTransverseImpactParameter";
85  return false;
86  }
87  //std::cout << "<trkTransverseImpactParameter>:" << std::endl;
88  //std::cout << " track: Pt = " << track->pt() << ", eta = " << track->eta() << ", phi = " << track->phi() << std::endl;
89  //std::cout << " vertex: x = " << (*pv)->position().x() << ", y = " << (*pv)->position().y() << ", z = " << (*pv)->position().z() << std::endl;
90  //std::cout << "--> dxy = " << std::fabs(track->dxy((*pv)->position())) << " (cut = " << cut << ")" << std::endl;
91  return (std::fabs(track->dxy((*pv)->position())) <= cut);
92 }
bool isNull() const
Checks for null.
Definition: Ref.h:247
bool reco::tau::qcuts::trkTransverseImpactParameter_cand ( const PFCandidate &  cand,
const reco::VertexRef pv,
double  cut 
)

Definition at line 94 of file RecoTauQualityCuts.cc.

References trkTransverseImpactParameter().

95 {
96  auto track = getTrackRef(cand);
97  if ( track.isNonnull() ) {
98  return trkTransverseImpactParameter(track, pv, cut);
99  } else {
100  //std::cout << "<trkTransverseImpactParameter_cand>: dXY = N/A, cut = " << cut << std::endl;
101  return false;
102  }
103 }
bool trkTransverseImpactParameter(const TrackBaseRef &track, const reco::VertexRef *pv, double cut)