CMS 3D CMS Logo

Public Member Functions | Private Types | Private Attributes

TrackWithVertexSelector Class Reference

#include <TrackWithVertexSelector.h>

List of all members.

Public Member Functions

bool operator() (const reco::Track &t, const edm::Event &iEvent) const
bool operator() (const reco::Track &t, const reco::VertexCollection &vtxs) const
bool testTrack (const reco::Track &t) const
bool testVertices (const reco::Track &t, const reco::VertexCollection &vtxs) const
 TrackWithVertexSelector (const edm::ParameterSet &)
 ~TrackWithVertexSelector ()

Private Types

typedef math::XYZPoint Point

Private Attributes

double d0Max_
double dzMax_
double etaMax_
double etaMin_
double normalizedChi2_
uint32_t numberOfLostHits_
uint32_t numberOfValidHits_
uint32_t numberOfValidPixelHits_
uint32_t nVertices_
double ptErrorCut_
double ptMax_
double ptMin_
std::string quality_
double rhoVtx_
edm::InputTag vertexTag_
bool vtxFallback_
double zetaVtx_

Detailed Description

Definition at line 16 of file TrackWithVertexSelector.h.


Member Typedef Documentation

Definition at line 39 of file TrackWithVertexSelector.h.


Constructor & Destructor Documentation

TrackWithVertexSelector::TrackWithVertexSelector ( const edm::ParameterSet iConfig) [explicit]

Definition at line 6 of file TrackWithVertexSelector.cc.

                                                                               :
  numberOfValidHits_(iConfig.getParameter<uint32_t>("numberOfValidHits")),
  numberOfValidPixelHits_(iConfig.getParameter<uint32_t>("numberOfValidPixelHits")),
  numberOfLostHits_(iConfig.getParameter<uint32_t>("numberOfLostHits")),
  normalizedChi2_(iConfig.getParameter<double>("normalizedChi2")),
  ptMin_(iConfig.getParameter<double>("ptMin")),
  ptMax_(iConfig.getParameter<double>("ptMax")),
  etaMin_(iConfig.getParameter<double>("etaMin")),
  etaMax_(iConfig.getParameter<double>("etaMax")),
  dzMax_(iConfig.getParameter<double>("dzMax")),
  d0Max_(iConfig.getParameter<double>("d0Max")),
  ptErrorCut_(iConfig.getParameter<double>("ptErrorCut")),
  quality_(iConfig.getParameter<std::string>("quality")),
  nVertices_(iConfig.getParameter<bool>("useVtx") ? iConfig.getParameter<uint32_t>("nVertices") : 0),
  vertexTag_(iConfig.getParameter<edm::InputTag>("vertexTag")),
  vtxFallback_(iConfig.getParameter<bool>("vtxFallback")),
  zetaVtx_(iConfig.getParameter<double>("zetaVtx")),
  rhoVtx_(iConfig.getParameter<double>("rhoVtx")) {
} 
TrackWithVertexSelector::~TrackWithVertexSelector ( )

Definition at line 26 of file TrackWithVertexSelector.cc.

{  }

Member Function Documentation

bool TrackWithVertexSelector::operator() ( const reco::Track t,
const edm::Event iEvent 
) const

Definition at line 65 of file TrackWithVertexSelector.cc.

References edm::Event::getByLabel(), nVertices_, testTrack(), testVertices(), and vertexTag_.

                                                                                      {
  if (!testTrack(t)) return false;
  if (nVertices_ == 0) return true;
  edm::Handle<reco::VertexCollection> hVtx;
  evt.getByLabel(vertexTag_, hVtx);
  return testVertices(t, *hVtx);
} 
bool TrackWithVertexSelector::operator() ( const reco::Track t,
const reco::VertexCollection vtxs 
) const

Definition at line 73 of file TrackWithVertexSelector.cc.

References nVertices_, testTrack(), and testVertices().

                                                                                                   {
  if (!testTrack(t)) return false;
  if (nVertices_ == 0) return true;
  return testVertices(t, vtxs);
}
bool TrackWithVertexSelector::testTrack ( const reco::Track t) const
bool TrackWithVertexSelector::testVertices ( const reco::Track t,
const reco::VertexCollection vtxs 
) const

Definition at line 47 of file TrackWithVertexSelector.cc.

References abs, reco::TrackBase::dxy(), reco::TrackBase::dz(), nVertices_, convertSQLiteXML::ok, rhoVtx_, reco::TrackBase::vertex(), vtxFallback_, and zetaVtx_.

Referenced by operator()().

                                                                                                     {
  bool ok = false;
  if (vtxs.size() > 0) {
    unsigned int tested = 1;
    for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end();
         it != ed; ++it) {
      if ((std::abs(t.dxy(it->position())) < rhoVtx_) && 
          (std::abs(t.dz(it->position())) < zetaVtx_)) {
        ok = true; break; 
      }
      if (tested++ >= nVertices_) break;
    }
  } else if (vtxFallback_) {
    return ( (std::abs(t.vertex().z()) < 15.9) && (t.vertex().Rho() < 0.2) );
  }
  return ok;
} 

Member Data Documentation

Definition at line 30 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 30 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 29 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 29 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 28 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 27 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 25 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 26 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 34 of file TrackWithVertexSelector.h.

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

Definition at line 31 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 29 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 29 of file TrackWithVertexSelector.h.

Referenced by testTrack().

std::string TrackWithVertexSelector::quality_ [private]

Definition at line 32 of file TrackWithVertexSelector.h.

Referenced by testTrack().

Definition at line 37 of file TrackWithVertexSelector.h.

Referenced by testVertices().

Definition at line 35 of file TrackWithVertexSelector.h.

Referenced by operator()().

Definition at line 36 of file TrackWithVertexSelector.h.

Referenced by testVertices().

Definition at line 37 of file TrackWithVertexSelector.h.

Referenced by testVertices().