CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
RecoTrackSelectorBase Class Reference

#include <RecoTrackSelectorBase.h>

Inheritance diagram for RecoTrackSelectorBase:
RecoTrackRefSelector RecoTrackSelector RecoTrackViewRefSelector

Public Member Functions

void init (const edm::Event &event, const edm::EventSetup &es)
 
bool operator() (const reco::TrackRef &tref) const
 
bool operator() (const reco::Track &t) const
 
bool operator() (const reco::Track &t, const reco::Track::Point &vertex) const
 
 RecoTrackSelectorBase ()
 
 RecoTrackSelectorBase (const edm::ParameterSet &cfg)
 
 RecoTrackSelectorBase (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 

Private Attributes

std::vector< reco::TrackBase::TrackAlgorithmalgorithm_
 
std::vector< reco::TrackBase::TrackAlgorithmalgorithmMask_
 
edm::EDGetTokenT< reco::BeamSpotbsSrcToken_
 
double lip_
 
double maxChi2_
 
double maxRapidity_
 
double meanPhi_
 
int min3DLayer_
 
int minHit_
 
int minLayer_
 
int minPixelHit_
 
double minRapidity_
 
std::vector< reco::TrackBase::TrackAlgorithmoriginalAlgorithm_
 
double ptMin_
 
std::vector< reco::TrackBase::TrackQualityquality_
 
double rangePhi_
 
double tip_
 
bool usePV_
 
reco::Track::Point vertex_
 
edm::EDGetTokenT< reco::VertexCollectionvertexToken_
 

Detailed Description

Definition at line 16 of file RecoTrackSelectorBase.h.

Constructor & Destructor Documentation

RecoTrackSelectorBase::RecoTrackSelectorBase ( )
inline

Definition at line 18 of file RecoTrackSelectorBase.h.

18 {}
RecoTrackSelectorBase::RecoTrackSelectorBase ( const edm::ParameterSet cfg)
inline

Definition at line 19 of file RecoTrackSelectorBase.h.

References reco::TrackBase::algoByName(), electronCleaner_cfi::algorithm, algorithm_, algorithmMask_, Exception, edm::ParameterSet::getParameter(), M_PI, trackingParticleSelector_cfi::maxPhi, trackingParticleSelector_cfi::minPhi, originalAlgorithm_, jets_cff::quality, quality_, reco::TrackBase::qualityByName(), and AlCaHLTBitMon_QueryRunRegistry::string.

19  :
20  ptMin_(cfg.getParameter<double>("ptMin")),
21  minRapidity_(cfg.getParameter<double>("minRapidity")),
22  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
23  meanPhi_((cfg.getParameter<double>("minPhi")+cfg.getParameter<double>("maxPhi"))/2.),
24  rangePhi_((cfg.getParameter<double>("maxPhi")-cfg.getParameter<double>("minPhi"))/2.),
25  tip_(cfg.getParameter<double>("tip")),
26  lip_(cfg.getParameter<double>("lip")),
27  maxChi2_(cfg.getParameter<double>("maxChi2")),
28  minHit_(cfg.getParameter<int>("minHit")),
29  minPixelHit_(cfg.getParameter<int>("minPixelHit")),
30  minLayer_(cfg.getParameter<int>("minLayer")),
31  min3DLayer_(cfg.getParameter<int>("min3DLayer")),
32  usePV_(false) {
33  const auto minPhi = cfg.getParameter<double>("minPhi");
34  const auto maxPhi = cfg.getParameter<double>("maxPhi");
35  if(minPhi >= maxPhi) {
36  throw cms::Exception("Configuration") << "RecoTrackSelectorPhase: minPhi (" << minPhi << ") must be smaller than maxPhi (" << maxPhi << "). The range is constructed from minPhi to maxPhi around their average.";
37  }
38  if(minPhi >= M_PI) {
39  throw cms::Exception("Configuration") << "RecoTrackSelectorPhase: minPhi (" << minPhi << ") must be smaller than PI. The range is constructed from minPhi to maxPhi around their average.";
40  }
41  if(maxPhi <= -M_PI) {
42  throw cms::Exception("Configuration") << "RecoTrackSelectorPhase: maxPhi (" << maxPhi << ") must be larger than -PI. The range is constructed from minPhi to maxPhi around their average.";
43  }
44 
45  for(const std::string& quality: cfg.getParameter<std::vector<std::string> >("quality"))
47  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithm"))
49  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("originalAlgorithm"))
51  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithmMaskContains"))
53  }
T getParameter(std::string const &) const
std::vector< reco::TrackBase::TrackAlgorithm > algorithmMask_
std::vector< reco::TrackBase::TrackQuality > quality_
#define M_PI
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
std::vector< reco::TrackBase::TrackAlgorithm > algorithm_
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:137
std::vector< reco::TrackBase::TrackAlgorithm > originalAlgorithm_
RecoTrackSelectorBase::RecoTrackSelectorBase ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)
inline

Definition at line 55 of file RecoTrackSelectorBase.h.

References bsSrcToken_, edm::ConsumesCollector::consumes(), edm::ParameterSet::getParameter(), usePV_, and vertexToken_.

55  :
57  usePV_ = cfg.getParameter<bool>("usePV");
59  if (usePV_)
61  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
edm::EDGetTokenT< reco::BeamSpot > bsSrcToken_

Member Function Documentation

void RecoTrackSelectorBase::init ( const edm::Event event,
const edm::EventSetup es 
)
inline

Definition at line 63 of file RecoTrackSelectorBase.h.

References ecalDrivenElectronSeedsParameters_cff::beamSpot, bsSrcToken_, reco::BeamSpot::position(), usePV_, vertex_, and vertexToken_.

Referenced by RecoTrackViewRefSelector::select(), RecoTrackRefSelector::select(), and RecoTrackSelector::select().

63  {
65  event.getByToken(bsSrcToken_,beamSpot);
66  vertex_ = beamSpot->position();
67  if (!usePV_) return;
69  event.getByToken(vertexToken_, hVtx);
70  if (hVtx->empty()) return;
71  vertex_ = (*hVtx)[0].position();
72  }
reco::Track::Point vertex_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
edm::EDGetTokenT< reco::BeamSpot > bsSrcToken_
const Point & position() const
position
Definition: BeamSpot.h:62
bool RecoTrackSelectorBase::operator() ( const reco::TrackRef tref) const
inline

Definition at line 74 of file RecoTrackSelectorBase.h.

74  {
75  return (*this)(*tref);
76  }
bool RecoTrackSelectorBase::operator() ( const reco::Track t) const
inline

Definition at line 78 of file RecoTrackSelectorBase.h.

References lumiQTWidget::t, and vertex_.

78  {
79  return (*this)(t, vertex_);
80  }
reco::Track::Point vertex_
bool RecoTrackSelectorBase::operator() ( const reco::Track t,
const reco::Track::Point vertex 
) const
inline

Definition at line 82 of file RecoTrackSelectorBase.h.

References patPFMETCorrections_cff::algo, reco::TrackBase::algo(), reco::TrackBase::algoMask(), algorithm_, algorithmMask_, hiPixelPairStep_cff::deltaPhi, reco::TrackBase::dsz(), reco::TrackBase::dxy(), reco::TrackBase::eta(), spr::find(), reco::TrackBase::hitPattern(), mps_fire::i, lip_, maxChi2_, maxRapidity_, meanPhi_, min3DLayer_, minHit_, minLayer_, minPixelHit_, minRapidity_, reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), reco::HitPattern::numberOfValidStripLayersWithMonoAndStereo(), reco::TrackBase::originalAlgo(), originalAlgorithm_, reco::TrackBase::phi(), reco::HitPattern::pixelLayersWithMeasurement(), reco::TrackBase::pt(), ptMin_, reco::TrackBase::quality(), quality_, rangePhi_, tip_, and reco::HitPattern::trackerLayersWithMeasurement().

82  {
83 
84  bool quality_ok = true;
85  if (!quality_.empty()) {
86  quality_ok = false;
87  for (unsigned int i = 0; i<quality_.size();++i) {
88  if (t.quality(quality_[i])){
89  quality_ok = true;
90  break;
91  }
92  }
93  }
94 
95  bool algo_ok = true;
96  if (!algorithm_.empty()) {
97  if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
98  }
99  if (!originalAlgorithm_.empty() && algo_ok) {
100  if (std::find(originalAlgorithm_.begin(), originalAlgorithm_.end(), t.originalAlgo()) == originalAlgorithm_.end()) algo_ok = false;
101  }
102  if(!algorithmMask_.empty() && algo_ok) {
103  if(std::find_if(algorithmMask_.begin(), algorithmMask_.end(), [&](reco::TrackBase::TrackAlgorithm algo) -> bool { // for some reason I have to either explicitly give the return type, or use static_cast<bool>()
104  return t.algoMask()[algo];
105  }) == algorithmMask_.end()) algo_ok = false;
106  }
107 
108  const auto dphi = deltaPhi(t.phi(), meanPhi_);
109 
110  return
111  (
112  (algo_ok & quality_ok) &&
118  fabs(t.pt()) >= ptMin_ &&
119  t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
120  dphi >= -rangePhi_ && dphi <= rangePhi_ &&
121  fabs(t.dxy(vertex)) <= tip_ &&
122  fabs(t.dsz(vertex)) <= lip_ &&
124  );
125  }
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:561
int numberOfValidHits() const
Definition: HitPattern.h:824
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:501
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
double dsz() const
dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0...
Definition: TrackBase.h:603
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:520
TrackAlgorithm
track algorithm
Definition: TrackBase.h:99
TrackAlgorithm algo() const
Definition: TrackBase.h:497
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:343
double pt() const
track transverse momentum
Definition: TrackBase.h:621
std::vector< reco::TrackBase::TrackAlgorithm > algorithmMask_
std::vector< reco::TrackBase::TrackQuality > quality_
TrackAlgorithm originalAlgo() const
Definition: TrackBase.h:501
AlgoMask algoMask() const
Definition: TrackBase.h:364
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:510
std::vector< reco::TrackBase::TrackAlgorithm > algorithm_
int numberOfValidPixelHits() const
Definition: HitPattern.h:839
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:591
std::vector< reco::TrackBase::TrackAlgorithm > originalAlgorithm_

Member Data Documentation

std::vector<reco::TrackBase::TrackAlgorithm> RecoTrackSelectorBase::algorithm_
private

Definition at line 147 of file RecoTrackSelectorBase.h.

Referenced by operator()(), and RecoTrackSelectorBase().

std::vector<reco::TrackBase::TrackAlgorithm> RecoTrackSelectorBase::algorithmMask_
private

Definition at line 149 of file RecoTrackSelectorBase.h.

Referenced by operator()(), and RecoTrackSelectorBase().

edm::EDGetTokenT<reco::BeamSpot> RecoTrackSelectorBase::bsSrcToken_
private

Definition at line 143 of file RecoTrackSelectorBase.h.

Referenced by init(), and RecoTrackSelectorBase().

double RecoTrackSelectorBase::lip_
private

Definition at line 135 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::maxChi2_
private

Definition at line 136 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::maxRapidity_
private

Definition at line 131 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::meanPhi_
private

Definition at line 132 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::min3DLayer_
private

Definition at line 140 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::minHit_
private

Definition at line 137 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::minLayer_
private

Definition at line 139 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::minPixelHit_
private

Definition at line 138 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::minRapidity_
private

Definition at line 130 of file RecoTrackSelectorBase.h.

Referenced by operator()().

std::vector<reco::TrackBase::TrackAlgorithm> RecoTrackSelectorBase::originalAlgorithm_
private

Definition at line 148 of file RecoTrackSelectorBase.h.

Referenced by operator()(), and RecoTrackSelectorBase().

double RecoTrackSelectorBase::ptMin_
private

Definition at line 129 of file RecoTrackSelectorBase.h.

Referenced by operator()().

std::vector<reco::TrackBase::TrackQuality> RecoTrackSelectorBase::quality_
private

Definition at line 146 of file RecoTrackSelectorBase.h.

Referenced by operator()(), and RecoTrackSelectorBase().

double RecoTrackSelectorBase::rangePhi_
private

Definition at line 133 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::tip_
private

Definition at line 134 of file RecoTrackSelectorBase.h.

Referenced by operator()().

bool RecoTrackSelectorBase::usePV_
private

Definition at line 141 of file RecoTrackSelectorBase.h.

Referenced by init(), and RecoTrackSelectorBase().

reco::Track::Point RecoTrackSelectorBase::vertex_
private

Definition at line 151 of file RecoTrackSelectorBase.h.

Referenced by init(), and operator()().

edm::EDGetTokenT<reco::VertexCollection> RecoTrackSelectorBase::vertexToken_
private

Definition at line 144 of file RecoTrackSelectorBase.h.

Referenced by init(), and RecoTrackSelectorBase().