Go to the documentation of this file.00001 #ifndef RecoHI_HiTracking_HITrackingRegionProducer_H
00002 #define RecoHI_HiTracking_HITrackingRegionProducer_H
00003
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005
00006 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00007 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
00008 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "FWCore/Framework/interface/Event.h"
00011 #include "FWCore/Framework/interface/ESHandle.h"
00012 #include "DataFormats/Common/interface/DetSetVector.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00014 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00015 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00016
00017 #include "TMath.h"
00018
00019 class HITrackingRegionProducer : public TrackingRegionProducer {
00020
00021 public:
00022
00023 HITrackingRegionProducer(const edm::ParameterSet& cfg) {
00024
00025 edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
00026
00027 thePtMin = regionPSet.getParameter<double>("ptMin");
00028 theOriginRadius = regionPSet.getParameter<double>("originRadius");
00029 theOriginHalfLength = regionPSet.getParameter<double>("originHalfLength");
00030 double xPos = regionPSet.getParameter<double>("originXPos");
00031 double yPos = regionPSet.getParameter<double>("originYPos");
00032 double zPos = regionPSet.getParameter<double>("originZPos");
00033 double xDir = regionPSet.getParameter<double>("directionXCoord");
00034 double yDir = regionPSet.getParameter<double>("directionYCoord");
00035 double zDir = regionPSet.getParameter<double>("directionZCoord");
00036 thePrecise = regionPSet.getParameter<bool>("precise");
00037 theSiPixelRecHits = regionPSet.getParameter<std::string>("siPixelRecHits");
00038 theOrigin = GlobalPoint(xPos,yPos,zPos);
00039 theDirection = GlobalVector(xDir, yDir, zDir);
00040 }
00041
00042 virtual ~HITrackingRegionProducer(){}
00043
00044 int estimateMultiplicity
00045 (const edm::Event& ev, const edm::EventSetup& es) const
00046 {
00047
00048 edm::Handle<SiPixelRecHitCollection> recHitColl;
00049 ev.getByLabel(theSiPixelRecHits, recHitColl);
00050
00051
00052 edm::ESHandle<TrackerTopology> tTopo;
00053 es.get<IdealGeometryRecord>().get(tTopo);
00054
00055 int numRecHits = 0;
00056
00057 for(SiPixelRecHitCollection::const_iterator recHitIdIterator = recHitColl->begin(), recHitIdIteratorEnd = recHitColl->end();
00058 recHitIdIterator != recHitIdIteratorEnd; recHitIdIterator++) {
00059 SiPixelRecHitCollection::DetSet hits = *recHitIdIterator;
00060 DetId detId = DetId(hits.detId());
00061 unsigned int detType=detId.det();
00062 unsigned int subid=detId.subdetId();
00063
00064 unsigned int layer=0;
00065 layer=tTopo->pxbLayer(detId);
00066 if(detType==1 && subid==1 && layer==1) {
00067 numRecHits += hits.size();
00068 }
00069 }
00070 return numRecHits;
00071 }
00072
00073 virtual std::vector<TrackingRegion* > regions(const edm::Event& ev, const edm::EventSetup& es) const {
00074
00075 int estMult = estimateMultiplicity(ev, es);
00076
00077
00078 float aa = 1.90935e-04;
00079 float bb = -2.90167e-01;
00080 float cc = 3.86125e+02;
00081
00082 float estTracks = aa*estMult*estMult+bb*estMult+cc;
00083
00084 LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]";
00085 LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]";
00086 LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]";
00087
00088 float regTracking = 400.;
00089 float etaB = 10.;
00090 float phiB = TMath::Pi()/2.;
00091
00092 float decEta = estTracks/600.;
00093 etaB = 2.5/decEta;
00094
00095 if(estTracks>regTracking) {
00096 LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]";
00097 LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]";
00098 LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]";
00099 LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: factor of decrease: " << decEta*2. << "]";
00100 }
00101
00102
00103 std::vector<TrackingRegion* > result;
00104 if(estTracks>regTracking) {
00105 result.push_back( new RectangularEtaPhiTrackingRegion(theDirection, theOrigin, thePtMin, theOriginRadius, theOriginHalfLength, etaB, phiB, 0, thePrecise) );
00106 }
00107 else {
00108 LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]";
00109 result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, theOriginHalfLength, thePrecise) );
00110 }
00111 return
00112 result;
00113 }
00114
00115 private:
00116 std::string theSiPixelRecHits;
00117 double thePtMin;
00118 GlobalPoint theOrigin;
00119 double theOriginRadius;
00120 double theOriginHalfLength;
00121 bool thePrecise;
00122 GlobalVector theDirection;
00123 };
00124
00125 #endif