CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HIPixelTrackFilter.cc
Go to the documentation of this file.
2 
9 
13 
14 using namespace std;
15 using namespace edm;
16 
17 /*****************************************************************************/
20 theTIPMax( ps.getParameter<double>("tipMax") ),
21 theNSigmaTipMaxTolerance( ps.getParameter<double>("nSigmaTipMaxTolerance")),
22 theLIPMax( ps.getParameter<double>("lipMax") ),
23 theNSigmaLipMaxTolerance( ps.getParameter<double>("nSigmaLipMaxTolerance")),
24 theChi2Max( ps.getParameter<double>("chi2") ),
25 thePtMin( ps.getParameter<double>("ptMin") ),
26 useClusterShape( ps.getParameter<bool>("useClusterShape") ),
27 theVertexCollection( ps.getParameter<edm::InputTag>("VertexCollection")),
28 theVertexCollectionToken( iC.consumes<reco::VertexCollection>(theVertexCollection)),
29 theVertices(0)
30 {
31 }
32 
33 /*****************************************************************************/
35 { }
36 
37 /*****************************************************************************/
39  const TrackerTopology *tTopo) const
40 {
41 
42  if (!track) return false;
43  if (track->chi2() > theChi2Max || track->pt() < thePtMin) return false;
44 
45 
46  math::XYZPoint vtxPoint(0.0,0.0,0.0);
47  double vzErr =0.0, vxErr=0.0, vyErr=0.0;
48 
49  if(theVertices->size()>0) {
50  vtxPoint=theVertices->begin()->position();
51  vzErr=theVertices->begin()->zError();
52  vxErr=theVertices->begin()->xError();
53  vyErr=theVertices->begin()->yError();
54  } else {
55  // THINK OF SOMETHING TO DO IF THERE IS NO VERTEX
56  }
57 
58  double d0=0.0, dz=0.0, d0sigma=0.0, dzsigma=0.0;
59  d0 = -1.*track->dxy(vtxPoint);
60  dz = track->dz(vtxPoint);
61  d0sigma = sqrt(track->d0Error()*track->d0Error()+vxErr*vyErr);
62  dzsigma = sqrt(track->dzError()*track->dzError()+vzErr*vzErr);
63 
64  if (theTIPMax>0 && fabs(d0)>theTIPMax) return false;
65  if (theNSigmaTipMaxTolerance>0 && (fabs(d0)/d0sigma)>theNSigmaTipMaxTolerance) return false;
66  if (theLIPMax>0 && fabs(dz)>theLIPMax) return false;
67  if (theNSigmaLipMaxTolerance>0 && (fabs(dz)/dzsigma)>theNSigmaLipMaxTolerance) return false;
68 
69  bool ok = true;
70  if(useClusterShape) ok = ClusterShapeTrackFilter::operator() (track,recHits,tTopo);
71 
72  return ok;
73 }
74 
75 /*****************************************************************************/
77 {
79 
80  // Get reco vertices
83  theVertices = vc.product();
84 
85  if(theVertices->size()>0) {
86  LogInfo("HeavyIonVertexing")
87  << "[HIPixelTrackFilter] Pixel track selection based on best vertex"
88  << "\n vz = " << theVertices->begin()->z()
89  << "\n vz sigma = " << theVertices->begin()->zError();
90  } else {
91  LogError("HeavyIonVertexing") // this can be made a warning when operator() is fixed
92  << "No vertex found in collection '" << theVertexCollection << "'";
93  }
94 
95  return;
96 
97 }
std::vector< const TrackingRecHit * > Hits
double d0Error() const
error on d0
Definition: TrackBase.h:789
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual bool operator()(const reco::Track *, const std::vector< const TrackingRecHit * > &hits, const TrackerTopology *tTopo) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken
const reco::VertexCollection * theVertices
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:536
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:608
virtual bool operator()(const reco::Track *, const PixelTrackFilter::Hits &hits, const TrackerTopology *tTopo) const
edm::InputTag theVertexCollection
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:596
double dzError() const
error on dz
Definition: TrackBase.h:801
T const * product() const
Definition: Handle.h:81
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
HIPixelTrackFilter(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
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:578
virtual void update(const edm::Event &ev, const edm::EventSetup &es) override
void update(const edm::Event &ev, const edm::EventSetup &es) override