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 198 of file RecoTauQualityCuts.cc.

References statics::func.

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

199 {
200  BOOST_FOREACH( const RecoTauQualityCuts::TrackQCutFunc& func, cuts ) {
201  if ( !func(track) ) return false;
202  }
203  return true;
204 }
string func
Definition: statics.py:48
bool reco::tau::qcuts::AND_cand ( const PFCandidate &  cand,
const RecoTauQualityCuts::CandQCutFuncCollection &  cuts 
)

Definition at line 206 of file RecoTauQualityCuts.cc.

References statics::func.

Referenced by mapAndCutByType().

207 {
208  BOOST_FOREACH( const RecoTauQualityCuts::CandQCutFunc& func, cuts ) {
209  if ( !func(cand) ) return false;
210  }
211  return true;
212 }
string func
Definition: statics.py:48
bool reco::tau::qcuts::etMin_cand ( const PFCandidate &  cand,
double  cut 
)

Definition at line 42 of file RecoTauQualityCuts.cc.

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

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

Definition at line 215 of file RecoTauQualityCuts.cc.

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

216 {
217  // Find the cuts that for this particle type
218  RecoTauQualityCuts::CandQCutFuncMap::const_iterator cuts = funcMap.find(cand.particleId());
219  // Return false if we dont' know how to deal with this particle type
220  if ( cuts == funcMap.end() ) return false;
221  return AND_cand(cand, cuts->second); // Otherwise AND all the cuts
222 }
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 154 of file RecoTauQualityCuts.cc.

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

Referenced by minTrackVertexWeight_cand().

155 {
156  if ( pv->isNull() ) {
157  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
158  "RecoTauQualityCuts is invalid. - minTrackVertexWeight";
159  return false;
160  }
161  //std::cout << "<minTrackVertexWeight>:" << std::endl;
162  //std::cout << " track: Pt = " << track->pt() << ", eta = " << track->eta() << ", phi = " << track->phi() << std::endl;
163  //std::cout << " vertex: x = " << (*pv)->position().x() << ", y = " << (*pv)->position().y() << ", z = " << (*pv)->position().z() << std::endl;
164  //std::cout << "--> trackWeight = " << (*pv)->trackWeight(track) << " (cut = " << cut << ")" << std::endl;
165  return ((*pv)->trackWeight(track) >= cut);
166 }
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 168 of file RecoTauQualityCuts.cc.

References minTrackVertexWeight().

169 {
170  auto track = getTrackRef(cand);
171  if ( track.isNonnull() ) {
172  return minTrackVertexWeight(track, pv, cut);
173  } else {
174  //std::cout << "<minTrackVertexWeight_cand>: weight = N/A, cut = " << cut << std::endl;
175  return false;
176  }
177 }
bool minTrackVertexWeight(const TrackBaseRef &track, const reco::VertexRef *pv, double cut)
bool reco::tau::qcuts::ptMin ( const TrackBaseRef &  track,
double  cut 
)

Definition at line 30 of file RecoTauQualityCuts.cc.

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

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

Definition at line 36 of file RecoTauQualityCuts.cc.

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

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

Definition at line 179 of file RecoTauQualityCuts.cc.

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

Referenced by trkChi2_cand().

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

Definition at line 185 of file RecoTauQualityCuts.cc.

References trkChi2().

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

Definition at line 111 of file RecoTauQualityCuts.cc.

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

Referenced by trkLongitudinalImpactParameter_cand().

112 {
113  if ( pv->isNull() ) {
114  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
115  "RecoTauQualityCuts is invalid. - trkLongitudinalImpactParameter";
116  return false;
117  }
118  //std::cout << "<trkLongitudinalImpactParameter>:" << std::endl;
119  //std::cout << " track: Pt = " << track->pt() << ", eta = " << track->eta() << ", phi = " << track->phi() << std::endl;
120  //std::cout << " vertex: x = " << (*pv)->position().x() << ", y = " << (*pv)->position().y() << ", z = " << (*pv)->position().z() << std::endl;
121  //std::cout << "--> dz = " << std::fabs(track->dz((*pv)->position())) << " (cut = " << cut << ")" << std::endl;
122  return (std::fabs(track->dz((*pv)->position())) <= cut);
123 }
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 125 of file RecoTauQualityCuts.cc.

References trkLongitudinalImpactParameter().

126 {
127  auto track = getTrackRef(cand);
128  if ( track.isNonnull() ) {
129  return trkLongitudinalImpactParameter(track, pv, cut);
130  } else {
131  //std::cout << "<trkLongitudinalImpactParameter_cand>: dZ = N/A, cut = " << cut << std::endl;
132  return false;
133  }
134 }
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 137 of file RecoTauQualityCuts.cc.

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

Referenced by trkLongitudinalImpactParameterWrtTrack_cand().

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

Definition at line 147 of file RecoTauQualityCuts.cc.

References trkLongitudinalImpactParameterWrtTrack().

148 {
149  auto track = getTrackRef(cand);
150  if ( track.isNonnull() ) return trkLongitudinalImpactParameterWrtTrack(track, leadTrack, pv, cut);
151  else return false;
152 }
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 48 of file RecoTauQualityCuts.cc.

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

Referenced by trkPixelHits_cand().

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

Definition at line 55 of file RecoTauQualityCuts.cc.

References trkPixelHits().

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

Definition at line 68 of file RecoTauQualityCuts.cc.

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

Referenced by trkTrackerHits_cand().

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

Definition at line 74 of file RecoTauQualityCuts.cc.

References trkTrackerHits().

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

Definition at line 86 of file RecoTauQualityCuts.cc.

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

Referenced by trkTransverseImpactParameter_cand().

87 {
88  if ( pv->isNull() ) {
89  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
90  "RecoTauQualityCuts is invalid. - trkTransverseImpactParameter";
91  return false;
92  }
93  //std::cout << "<trkTransverseImpactParameter>:" << std::endl;
94  //std::cout << " track: Pt = " << track->pt() << ", eta = " << track->eta() << ", phi = " << track->phi() << std::endl;
95  //std::cout << " vertex: x = " << (*pv)->position().x() << ", y = " << (*pv)->position().y() << ", z = " << (*pv)->position().z() << std::endl;
96  //std::cout << "--> dxy = " << std::fabs(track->dxy((*pv)->position())) << " (cut = " << cut << ")" << std::endl;
97  return (std::fabs(track->dxy((*pv)->position())) <= cut);
98 }
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 100 of file RecoTauQualityCuts.cc.

References trkTransverseImpactParameter().

101 {
102  auto track = getTrackRef(cand);
103  if ( track.isNonnull() ) {
104  return trkTransverseImpactParameter(track, pv, cut);
105  } else {
106  //std::cout << "<trkTransverseImpactParameter_cand>: dXY = N/A, cut = " << cut << std::endl;
107  return false;
108  }
109 }
bool trkTransverseImpactParameter(const TrackBaseRef &track, const reco::VertexRef *pv, double cut)