CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HIProtoTrackFilter.cc
Go to the documentation of this file.
2 
8 
11 
14 
18 
19 using namespace std;
20 using namespace edm;
21 
22 /*****************************************************************************/
24 theTIPMax( ps.getParameter<double>("tipMax") ),
25 theChi2Max( ps.getParameter<double>("chi2") ),
26 thePtMin( ps.getParameter<double>("ptMin") ),
27 doVariablePtMin( ps.getParameter<bool>("doVariablePtMin") ),
28 theBeamSpotTag( ps.getParameter<InputTag>("beamSpot")),
29 theSiPixelRecHits( ps.getParameter<InputTag>("siPixelRecHits")),
30 theBeamSpot(0),
31 theVariablePtMin(0)
32 {
33 }
34 
35 /*****************************************************************************/
37 theTIPMax( ps.getParameter<double>("tipMax") ),
38 theChi2Max( ps.getParameter<double>("chi2") ),
39 thePtMin( ps.getParameter<double>("ptMin") ),
40 doVariablePtMin( ps.getParameter<bool>("doVariablePtMin") ),
41 theBeamSpotTag( ps.getParameter<InputTag>("beamSpot")),
42 theSiPixelRecHits( ps.getParameter<InputTag>("siPixelRecHits")),
43 theBeamSpot(0),
44 theVariablePtMin(0)
45 {
46 }
47 
48 /*****************************************************************************/
50 { }
51 
52 /*****************************************************************************/
53 bool HIProtoTrackFilter::operator() (const reco::Track* track,const PixelTrackFilter::Hits & recHits) const
54 {
55 
56  if (!track) return false;
57 
58  float minpt = thePtMin;
60 
61  if (track->chi2() > theChi2Max || track->pt() < minpt) return false;
62 
63  math::XYZPoint vtxPoint(0.0,0.0,0.0);
64 
65  if(theBeamSpot)
66  vtxPoint = theBeamSpot->position();
67 
68  double d0=0.0;
69  d0 = -1.*track->dxy(vtxPoint);
70 
71  if (theTIPMax>0 && fabs(d0)>theTIPMax) return false;
72 
73  return true;
74 }
75 
76 /*****************************************************************************/
78 {
79 
80  // Get the beam spot
82  ev.getByLabel( theBeamSpotTag, bsHandle);
83  theBeamSpot = bsHandle.product();
84 
85  if(theBeamSpot) {
86  LogInfo("HeavyIonVertexing")
87  << "[HIProtoTrackFilter] Proto track selection based on beamspot"
88  << "\n (x,y,z) = (" << theBeamSpot->x0() << "," << theBeamSpot->y0() << "," << theBeamSpot->z0() << ")";
89  } else {
90  LogError("HeavyIonVertexing") // this can be made a warning when operator() is fixed
91  << "No beamspot found with tag '" << theBeamSpotTag << "'";
92  }
93 
94  // Estimate multiplicity
96  ev.getByLabel(theSiPixelRecHits, recHitColl);
97 
98  vector<const TrackingRecHit*> theChosenHits;
100  edmNew::copyDetSetRange(*recHitColl,theChosenHits,acc.pixelBarrelLayer(1));
101  float estMult = theChosenHits.size();
102 
104 
105  // parameterize ptMin such that a roughly constant number of selected prototracks passed are to vertexing
106  float varPtCutoff = 1500; //cutoff for variable ptMin
107  if(estMult < varPtCutoff) {
108  theVariablePtMin = 0.075;
109  if(estMult > 0) theVariablePtMin = (13. - (varPtCutoff/estMult) )/12.;
110  if(theVariablePtMin<0.075) theVariablePtMin = 0.075; // don't lower the cut past 75 MeV
111  }
112 
113  LogTrace("heavyIonHLTVertexing")<<" [HIProtoTrackFilter: theVariablePtMin: " << theVariablePtMin << "]";
114 
115 
116  return;
117 
118 }
double z0() const
z coordinate
Definition: BeamSpot.h:69
void copyDetSetRange(DSTV const &dstv, std::vector< T const * > &v, std::pair< A, B > const &sel)
virtual void update(edm::Event &ev)
std::vector< const TrackingRecHit * > Hits
edm::InputTag theBeamSpotTag
edm::InputTag theSiPixelRecHits
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:107
double pt() const
track transverse momentum
Definition: TrackBase.h:131
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#define LogTrace(id)
virtual bool operator()(const reco::Track *, const PixelTrackFilter::Hits &hits) const
const reco::BeamSpot * theBeamSpot
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
T const * product() const
Definition: Handle.h:74
HIProtoTrackFilter(const edm::ParameterSet &ps, const edm::EventSetup &es)
std::pair< DetId, DetIdPXBSameLayerComparator > pixelBarrelLayer(int layer)
double y0() const
y coordinate
Definition: BeamSpot.h:67
const Point & position() const
position
Definition: BeamSpot.h:63
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:121
double x0() const
x coordinate
Definition: BeamSpot.h:65