CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HITrackingRegionForPrimaryVtxProducer.h
Go to the documentation of this file.
1 #ifndef RecoHI_HiTracking_HITrackingRegionForPrimaryVtxProducer_H
2 #define RecoHI_HiTracking_HITrackingRegionForPrimaryVtxProducer _H
3 
5 
11 
14 
17 
20 
24 
25 #include "TMath.h"
26 
28 
29  public:
30 
32 
33  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
34  thePtMin = regionPSet.getParameter<double>("ptMin");
35  theOriginRadius = regionPSet.getParameter<double>("originRadius");
36  theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
37  theBeamSpotTag = regionPSet.getParameter<edm::InputTag>("beamSpot");
38  thePrecise = regionPSet.getParameter<bool>("precise");
39  theSiPixelRecHits = regionPSet.getParameter<edm::InputTag>("siPixelRecHits");
40  doVariablePtMin = regionPSet.getParameter<bool>("doVariablePtMin");
41  double xDir = regionPSet.getParameter<double>("directionXCoord");
42  double yDir = regionPSet.getParameter<double>("directionYCoord");
43  double zDir = regionPSet.getParameter<double>("directionZCoord");
44  theDirection = GlobalVector(xDir, yDir, zDir);
45 
46  // for using vertex instead of beamspot
47  theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
48  theFixedError = regionPSet.getParameter<double>("fixedError");
49  theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
50  theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
51  vertexCollName = regionPSet.getParameter<edm::InputTag>("VertexCollection");
52  }
53 
55 
57  (const edm::Event& ev, const edm::EventSetup& es) const
58  {
59  //rechits
61  ev.getByLabel(theSiPixelRecHits, recHitColl);
62 
63  std::vector<const TrackingRecHit*> theChosenHits;
65  edmNew::copyDetSetRange(*recHitColl,theChosenHits,acc.pixelBarrelLayer(1));
66  return theChosenHits.size();
67 
68  }
69 
70  virtual std::vector<TrackingRegion* > regions(const edm::Event& ev, const edm::EventSetup& es) const {
71 
72  int estMult = estimateMultiplicity(ev, es);
73 
74  // from MC relating first layer pixel hits to "findable" sim tracks with pt>1 GeV
75  float cc = -38.6447;
76  float bb = 0.0581765;
77  float aa = 1.34306e-06;
78 
79  float estTracks = aa*estMult*estMult+bb*estMult+cc;
80 
81  LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]";
82  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]";
83  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]";
84 
85  float regTracking = 60.; //if we have more tracks -> regional tracking
86  float etaB = 10.;
87  float phiB = TMath::Pi()/2.;
88 
89  float decEta = estTracks/90.;
90  etaB = 2.5/decEta;
91 
92  if(estTracks>regTracking) {
93  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]";
94  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]";
95  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]";
96  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: factor of decrease: " << decEta*2. << "]"; // 2:from phi
97  }
98 
99  float minpt = thePtMin;
100  float varPtCutoff = 1500; //cutoff
101  if(doVariablePtMin && estMult < varPtCutoff) {
102  minpt = 0.075;
103  if(estMult > 0) minpt += estMult * (thePtMin - 0.075)/varPtCutoff; // lower ptMin linearly with pixel hit multiplicity
104  }
105 
106  // tracking region selection
107  std::vector<TrackingRegion* > result;
108  double halflength;
109  GlobalPoint origin;
111  ev.getByLabel( theBeamSpotTag, bsHandle);
112  if(bsHandle.isValid()) {
113  const reco::BeamSpot & bs = *bsHandle;
114  origin=GlobalPoint(bs.x0(), bs.y0(), bs.z0());
115  halflength=theNSigmaZ*bs.sigmaZ();
116 
118  {
120  ev.getByLabel(vertexCollName,vertexCollection);
121 
122  for(reco::VertexCollection::const_iterator iV=vertexCollection->begin();
123  iV != vertexCollection->end() ; iV++) {
124  if (iV->isFake() || !iV->isValid()) continue;
125  origin = GlobalPoint(bs.x0(),bs.y0(),iV->z());
126  halflength = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex);
127  }
128  }
129 
130  if(estTracks>regTracking) { // regional tracking
131  result.push_back(
132  new RectangularEtaPhiTrackingRegion(theDirection, origin, thePtMin, theOriginRadius, halflength, etaB, phiB, 0, thePrecise) );
133  }
134  else { // global tracking
135  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]";
136  result.push_back(
137  new GlobalTrackingRegion(minpt, origin, theOriginRadius, halflength, thePrecise) );
138  }
139  }
140  return result;
141  }
142 
143  private:
144  double thePtMin;
146  double theNSigmaZ;
152 
154  double theFixedError;
158 
159 
160 };
161 
162 #endif
const double Pi
T getParameter(std::string const &) const
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)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
tuple vertexCollection
virtual std::vector< TrackingRegion * > regions(const edm::Event &ev, const edm::EventSetup &es) const
tuple result
Definition: query.py:137
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#define LogTrace(id)
double sigmaZ() const
sigma z
Definition: BeamSpot.h:81
std::pair< DetId, DetIdPXBSameLayerComparator > pixelBarrelLayer(int layer)
double y0() const
y coordinate
Definition: BeamSpot.h:67
int estimateMultiplicity(const edm::Event &ev, const edm::EventSetup &es) const
HITrackingRegionForPrimaryVtxProducer(const edm::ParameterSet &cfg)
Global3DVector GlobalVector
Definition: GlobalVector.h:10
double x0() const
x coordinate
Definition: BeamSpot.h:65