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 
9 
12 
16 
19 
20 using namespace std;
21 using namespace edm;
22 
23 /*****************************************************************************/
25 theTIPMax( ps.getParameter<double>("tipMax") ),
26 theChi2Max( ps.getParameter<double>("chi2") ),
27 thePtMin( ps.getParameter<double>("ptMin") ),
28 doVariablePtMin( ps.getParameter<bool>("doVariablePtMin") ),
29 theBeamSpotTag( ps.getParameter<InputTag>("beamSpot")),
30 theBeamSpotToken( iC.consumes<reco::BeamSpot>(theBeamSpotTag)),
31 theSiPixelRecHitsToken( iC.consumes<SiPixelRecHitCollection>(ps.getParameter<InputTag>("siPixelRecHits"))),
32 theBeamSpot(0),
33 theVariablePtMin(0)
34 {
35 }
36 
37 /*****************************************************************************/
39 { }
40 
41 /*****************************************************************************/
43 {
44 
45  if (!track) return false;
46 
47  float minpt = thePtMin;
49 
50  if (track->chi2() > theChi2Max || track->pt() < minpt) return false;
51 
52  math::XYZPoint vtxPoint(0.0,0.0,0.0);
53 
54  if(theBeamSpot)
55  vtxPoint = theBeamSpot->position();
56 
57  double d0=0.0;
58  d0 = -1.*track->dxy(vtxPoint);
59 
60  if (theTIPMax>0 && fabs(d0)>theTIPMax) return false;
61 
62  return true;
63 }
64 
65 /*****************************************************************************/
67 {
68 
69  // Get the beam spot
71  ev.getByToken( theBeamSpotToken, bsHandle);
72  theBeamSpot = bsHandle.product();
73 
74  if(theBeamSpot) {
75  LogInfo("HeavyIonVertexing")
76  << "[HIProtoTrackFilter] Proto track selection based on beamspot"
77  << "\n (x,y,z) = (" << theBeamSpot->x0() << "," << theBeamSpot->y0() << "," << theBeamSpot->z0() << ")";
78  } else {
79  LogError("HeavyIonVertexing") // this can be made a warning when operator() is fixed
80  << "No beamspot found with tag '" << theBeamSpotTag << "'";
81  }
82 
83  // Estimate multiplicity
85  ev.getByToken(theSiPixelRecHitsToken, recHitColl);
86 
88  es.get<TrackerTopologyRcd>().get(httopo);
89 
90  vector<const TrackingRecHit*> theChosenHits;
91  edmNew::copyDetSetRange(*recHitColl,theChosenHits, httopo->pxbDetIdLayerComparator(1));
92  float estMult = theChosenHits.size();
93 
95 
96  // parameterize ptMin such that a roughly constant number of selected prototracks passed are to vertexing
97  float varPtCutoff = 1500; //cutoff for variable ptMin
98  if(estMult < varPtCutoff) {
99  theVariablePtMin = 0.075;
100  if(estMult > 0) theVariablePtMin = (13. - (varPtCutoff/estMult) )/12.;
101  if(theVariablePtMin<0.075) theVariablePtMin = 0.075; // don't lower the cut past 75 MeV
102  }
103 
104  LogTrace("heavyIonHLTVertexing")<<" [HIProtoTrackFilter: theVariablePtMin: " << theVariablePtMin << "]";
105 
106 
107  return;
108 
109 }
double z0() const
z coordinate
Definition: BeamSpot.h:68
void copyDetSetRange(DSTV const &dstv, std::vector< T const * > &v, std::pair< A, B > const &sel)
std::vector< const TrackingRecHit * > Hits
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::InputTag theBeamSpotTag
edm::EDGetTokenT< SiPixelRecHitCollection > theSiPixelRecHitsToken
bool ev
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:536
double pt() const
track transverse momentum
Definition: TrackBase.h:608
#define LogTrace(id)
virtual void update(const edm::Event &ev, const edm::EventSetup &es) override
virtual bool operator()(const reco::Track *, const PixelTrackFilter::Hits &hits) const
T const * product() const
Definition: Handle.h:81
const reco::BeamSpot * theBeamSpot
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const T & get() const
Definition: EventSetup.h:56
HIProtoTrackFilter(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
double y0() const
y coordinate
Definition: BeamSpot.h:66
const Point & position() const
position
Definition: BeamSpot.h:62
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
edm::EDGetTokenT< reco::BeamSpot > theBeamSpotToken
double x0() const
x coordinate
Definition: BeamSpot.h:64