CMS 3D CMS Logo

HITrackingRegionProducer.h
Go to the documentation of this file.
1 #ifndef RecoHI_HiTracking_HITrackingRegionProducer_H
2 #define RecoHI_HiTracking_HITrackingRegionProducer_H
3 
6 
17 
18 #include "TMath.h"
19 
21 
22  public:
23 
25 
26  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
27 
28  thePtMin = regionPSet.getParameter<double>("ptMin");
29  theOriginRadius = regionPSet.getParameter<double>("originRadius");
30  theOriginHalfLength = regionPSet.getParameter<double>("originHalfLength");
31  double xPos = regionPSet.getParameter<double>("originXPos");
32  double yPos = regionPSet.getParameter<double>("originYPos");
33  double zPos = regionPSet.getParameter<double>("originZPos");
34  double xDir = regionPSet.getParameter<double>("directionXCoord");
35  double yDir = regionPSet.getParameter<double>("directionYCoord");
36  double zDir = regionPSet.getParameter<double>("directionZCoord");
37  thePrecise = regionPSet.getParameter<bool>("precise");
38  theSiPixelRecHits = regionPSet.getParameter<edm::InputTag>("siPixelRecHits");
40  theOrigin = GlobalPoint(xPos,yPos,zPos);
41  theDirection = GlobalVector(xDir, yDir, zDir);
42  }
43 
45 
47  (const edm::Event& ev, const edm::EventSetup& es) const
48  {
49  //rechits
51  ev.getByToken(theSiPixelRecHitsToken, recHitColl);
52 
53  //Retrieve tracker topology from geometry
55  es.get<IdealGeometryRecord>().get(tTopo);
56 
57  int numRecHits = 0;
58  //FIXME: this can be optimized quite a bit by looping only on the per-det 'items' of DetSetVector
59  for(SiPixelRecHitCollection::const_iterator recHitIdIterator = recHitColl->begin(), recHitIdIteratorEnd = recHitColl->end();
60  recHitIdIterator != recHitIdIteratorEnd; recHitIdIterator++) {
61  SiPixelRecHitCollection::DetSet hits = *recHitIdIterator;
62  DetId detId = DetId(hits.detId()); // Get the Detid object
63  unsigned int detType=detId.det(); // det type, tracker=1
64  unsigned int subid=detId.subdetId(); //subdetector type, barrel=1, fpix=2
65 
66  unsigned int layer=0;
67  layer=tTopo->pxbLayer(detId);
68  if(detType==1 && subid==1 && layer==1) {
69  numRecHits += hits.size();
70  }
71  }
72  return numRecHits;
73  }
74 
75  virtual std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev, const edm::EventSetup& es) const override {
76 
77  int estMult = estimateMultiplicity(ev, es);
78 
79  // fit from MC information
80  float aa = 1.90935e-04;
81  float bb = -2.90167e-01;
82  float cc = 3.86125e+02;
83 
84  float estTracks = aa*estMult*estMult+bb*estMult+cc;
85 
86  LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]";
87  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]";
88  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]";
89 
90  float regTracking = 400.; //if we have more tracks -> regional tracking
91  float etaB = 10.;
92  float phiB = TMath::Pi()/2.;
93 
94  float decEta = estTracks/600.;
95  etaB = 2.5/decEta;
96 
97  if(estTracks>regTracking) {
98  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]";
99  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]";
100  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]";
101  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: factor of decrease: " << decEta*2. << "]"; // 2:from phi
102  }
103 
104  // tracking region selection
105  std::vector<std::unique_ptr<TrackingRegion> > result;
106  if(estTracks>regTracking) { // regional tracking
107  result.push_back( std::make_unique<RectangularEtaPhiTrackingRegion>(theDirection, theOrigin, thePtMin, theOriginRadius, theOriginHalfLength, etaB, phiB, RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever, thePrecise) );
108  }
109  else { // global tracking
110  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]";
111  result.push_back( std::make_unique<GlobalTrackingRegion>(thePtMin, theOrigin, theOriginRadius, theOriginHalfLength, thePrecise) );
112  }
113  return
114  result;
115  }
116 
117  private:
120  double thePtMin;
126 };
127 
128 #endif
const double Pi
T getParameter(std::string const &) const
int estimateMultiplicity(const edm::Event &ev, const edm::EventSetup &es) const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
HITrackingRegionProducer(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
virtual std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &ev, const edm::EventSetup &es) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
edm::EDGetTokenT< SiPixelRecHitCollection > theSiPixelRecHitsToken
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
bool ev
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
#define LogTrace(id)
unsigned int pxbLayer(const DetId &id) const
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:56
id_type detId() const
Definition: DetSetNew.h:84
size_type size() const
Definition: DetSetNew.h:87
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
const_iterator begin(bool update=false) const
Global3DVector GlobalVector
Definition: GlobalVector.h:10