CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoHI/HiTracking/plugins/HITrackingRegionProducer.h

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 
00012 #include "DataFormats/Common/interface/DetSetVector.h"    
00013 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00014 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00015 
00016 #include "TMath.h"
00017 
00018 class HITrackingRegionProducer : public TrackingRegionProducer {
00019   
00020  public:
00021   
00022   HITrackingRegionProducer(const edm::ParameterSet& cfg) { 
00023     
00024     edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
00025     
00026     thePtMin            = regionPSet.getParameter<double>("ptMin");
00027     theOriginRadius     = regionPSet.getParameter<double>("originRadius");
00028     theOriginHalfLength = regionPSet.getParameter<double>("originHalfLength");
00029     double xPos         = regionPSet.getParameter<double>("originXPos");
00030     double yPos         = regionPSet.getParameter<double>("originYPos");
00031     double zPos         = regionPSet.getParameter<double>("originZPos");
00032     double xDir         = regionPSet.getParameter<double>("directionXCoord");
00033     double yDir         = regionPSet.getParameter<double>("directionYCoord");
00034     double zDir         = regionPSet.getParameter<double>("directionZCoord");
00035     thePrecise          = regionPSet.getParameter<bool>("precise"); 
00036     theSiPixelRecHits   = regionPSet.getParameter<std::string>("siPixelRecHits");
00037     theOrigin = GlobalPoint(xPos,yPos,zPos);
00038     theDirection = GlobalVector(xDir, yDir, zDir);
00039   }   
00040   
00041   virtual ~HITrackingRegionProducer(){}
00042   
00043   int estimateMultiplicity
00044     (const edm::Event& ev, const edm::EventSetup& es) const
00045   {
00046     //rechits
00047     edm::Handle<SiPixelRecHitCollection> recHitColl;
00048     ev.getByLabel(theSiPixelRecHits, recHitColl);
00049     
00050     int numRecHits = 0;
00051     //FIXME: this can be optimized quite a bit by looping only on the per-det 'items' of DetSetVector
00052     for(SiPixelRecHitCollection::const_iterator recHitIdIterator = recHitColl->begin(), recHitIdIteratorEnd = recHitColl->end();
00053         recHitIdIterator != recHitIdIteratorEnd; recHitIdIterator++) {
00054       SiPixelRecHitCollection::DetSet hits = *recHitIdIterator;
00055       DetId detId = DetId(hits.detId());   // Get the Detid object
00056       unsigned int detType=detId.det();    // det type, tracker=1
00057       unsigned int subid=detId.subdetId(); //subdetector type, barrel=1, fpix=2
00058       PXBDetId pdetId = PXBDetId(detId);
00059       unsigned int layer=0;
00060       layer=pdetId.layer();
00061       if(detType==1 && subid==1 && layer==1) {
00062         numRecHits += hits.size();
00063       }
00064     }
00065     return numRecHits;
00066   }
00067   
00068   virtual std::vector<TrackingRegion* > regions(const edm::Event& ev, const edm::EventSetup& es) const {
00069     
00070     int estMult = estimateMultiplicity(ev, es);
00071     
00072     // fit from MC information
00073     float aa = 1.90935e-04;
00074     float bb = -2.90167e-01;
00075     float cc = 3.86125e+02;
00076     
00077     float estTracks = aa*estMult*estMult+bb*estMult+cc;
00078     
00079     LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]";
00080     LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]";
00081     LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]";
00082     
00083     float regTracking = 400.;  //if we have more tracks -> regional tracking
00084     float etaB = 10.;
00085     float phiB = TMath::Pi()/2.;
00086     
00087     float decEta = estTracks/600.;
00088     etaB = 2.5/decEta;
00089     
00090     if(estTracks>regTracking) {
00091       LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]";
00092       LogTrace("heavyIonHLTVertexing")<<"  [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]";
00093       LogTrace("heavyIonHLTVertexing")<<"  [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]";
00094       LogTrace("heavyIonHLTVertexing")<<"  [Regional Tracking: factor of decrease: " << decEta*2. << "]";  // 2:from phi
00095     }
00096     
00097     // tracking region selection
00098     std::vector<TrackingRegion* > result;
00099     if(estTracks>regTracking) {  // regional tracking
00100       result.push_back( new RectangularEtaPhiTrackingRegion(theDirection, theOrigin, thePtMin, theOriginRadius, theOriginHalfLength, etaB, phiB, 0, thePrecise) );
00101     }
00102     else {                       // global tracking
00103       LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]";
00104       result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, theOriginHalfLength, thePrecise) );
00105     }
00106     return 
00107       result;
00108   }
00109   
00110  private:
00111   std::string theSiPixelRecHits;
00112   double thePtMin; 
00113   GlobalPoint theOrigin;
00114   double theOriginRadius; 
00115   double theOriginHalfLength; 
00116   bool thePrecise;
00117   GlobalVector theDirection;
00118 };
00119 
00120 #endif