CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 #ifndef RecoHI_HiTracking_HITrackingRegionForPrimaryVtxProducer_H 
00002 #define RecoHI_HiTracking_HITrackingRegionForPrimaryVtxProducer _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/BeamSpot/interface/BeamSpot.h"
00013 #include "FWCore/Utilities/interface/InputTag.h"
00014 
00015 #include "DataFormats/VertexReco/interface/Vertex.h"
00016 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00017 
00018 #include "Geometry/TrackerGeometryBuilder/interface/TrackerLayerIdAccessor.h"    
00019 #include "DataFormats/Common/interface/DetSetAlgorithm.h"
00020 
00021 #include "DataFormats/Common/interface/DetSetVector.h"    
00022 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00023 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00024 
00025 #include "TMath.h"
00026 
00027 class HITrackingRegionForPrimaryVtxProducer : public TrackingRegionProducer {
00028   
00029  public:
00030   
00031   HITrackingRegionForPrimaryVtxProducer(const edm::ParameterSet& cfg) { 
00032     
00033     edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
00034     thePtMin            = regionPSet.getParameter<double>("ptMin");
00035     theOriginRadius     = regionPSet.getParameter<double>("originRadius");
00036     theNSigmaZ          = regionPSet.getParameter<double>("nSigmaZ");
00037     theBeamSpotTag      = regionPSet.getParameter<edm::InputTag>("beamSpot");
00038     thePrecise          = regionPSet.getParameter<bool>("precise"); 
00039     theSiPixelRecHits   = regionPSet.getParameter<edm::InputTag>("siPixelRecHits");  
00040     doVariablePtMin     = regionPSet.getParameter<bool>("doVariablePtMin"); 
00041     double xDir         = regionPSet.getParameter<double>("directionXCoord");
00042     double yDir         = regionPSet.getParameter<double>("directionYCoord");
00043     double zDir         = regionPSet.getParameter<double>("directionZCoord");
00044     theDirection = GlobalVector(xDir, yDir, zDir);
00045 
00046     // for using vertex instead of beamspot
00047     theSigmaZVertex     = regionPSet.getParameter<double>("sigmaZVertex");
00048     theFixedError       = regionPSet.getParameter<double>("fixedError");
00049     theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
00050     theUseFixedError    = regionPSet.getParameter<bool>("useFixedError");
00051     vertexCollName      = regionPSet.getParameter<edm::InputTag>("VertexCollection");
00052   }   
00053   
00054   virtual ~HITrackingRegionForPrimaryVtxProducer(){}
00055   
00056   int estimateMultiplicity
00057     (const edm::Event& ev, const edm::EventSetup& es) const
00058     {
00059       //rechits
00060       edm::Handle<SiPixelRecHitCollection> recHitColl;
00061       ev.getByLabel(theSiPixelRecHits, recHitColl);
00062       
00063       std::vector<const TrackingRecHit*> theChosenHits;          
00064       TrackerLayerIdAccessor acc;        
00065       edmNew::copyDetSetRange(*recHitColl,theChosenHits,acc.pixelBarrelLayer(1));        
00066       return theChosenHits.size();       
00067       
00068     }
00069   
00070   virtual std::vector<TrackingRegion* > regions(const edm::Event& ev, const edm::EventSetup& es) const {
00071     
00072     int estMult = estimateMultiplicity(ev, es);
00073     
00074     // from MC relating first layer pixel hits to "findable" sim tracks with pt>1 GeV
00075     float cc = -38.6447;
00076     float bb = 0.0581765;
00077     float aa = 1.34306e-06;
00078 
00079     float estTracks = aa*estMult*estMult+bb*estMult+cc;
00080     
00081     LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]";
00082     LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]";
00083     LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]";
00084     
00085     float regTracking = 60.;  //if we have more tracks -> regional tracking
00086     float etaB = 10.;
00087     float phiB = TMath::Pi()/2.;
00088     
00089     float decEta = estTracks/90.;
00090     etaB = 2.5/decEta;
00091     
00092     if(estTracks>regTracking) {
00093       LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]";
00094       LogTrace("heavyIonHLTVertexing")<<"  [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]";
00095       LogTrace("heavyIonHLTVertexing")<<"  [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]";
00096       LogTrace("heavyIonHLTVertexing")<<"  [Regional Tracking: factor of decrease: " << decEta*2. << "]";  // 2:from phi
00097     }
00098     
00099     float minpt = thePtMin;
00100     float varPtCutoff = 1500; //cutoff
00101     if(doVariablePtMin && estMult < varPtCutoff) {
00102       minpt = 0.075;
00103       if(estMult > 0) minpt += estMult * (thePtMin - 0.075)/varPtCutoff; // lower ptMin linearly with pixel hit multiplicity
00104     }
00105     
00106     // tracking region selection
00107     std::vector<TrackingRegion* > result;
00108     double halflength;
00109     GlobalPoint origin;
00110     edm::Handle<reco::BeamSpot> bsHandle;
00111     ev.getByLabel( theBeamSpotTag, bsHandle);
00112     if(bsHandle.isValid()) {
00113       const reco::BeamSpot & bs = *bsHandle; 
00114       origin=GlobalPoint(bs.x0(), bs.y0(), bs.z0()); 
00115       halflength=theNSigmaZ*bs.sigmaZ();
00116       
00117     if(theUseFoundVertices)
00118     {
00119       edm::Handle<reco::VertexCollection> vertexCollection;
00120       ev.getByLabel(vertexCollName,vertexCollection);
00121 
00122       for(reco::VertexCollection::const_iterator iV=vertexCollection->begin(); 
00123           iV != vertexCollection->end() ; iV++) {
00124           if (iV->isFake() || !iV->isValid()) continue;
00125           origin     = GlobalPoint(bs.x0(),bs.y0(),iV->z());
00126           halflength = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex); 
00127       }
00128     }
00129 
00130       if(estTracks>regTracking) {  // regional tracking
00131         result.push_back( 
00132           new RectangularEtaPhiTrackingRegion(theDirection, origin, thePtMin, theOriginRadius, halflength, etaB, phiB, 0, thePrecise) );
00133       }
00134       else {                       // global tracking
00135         LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]";
00136         result.push_back( 
00137           new GlobalTrackingRegion(minpt, origin, theOriginRadius, halflength, thePrecise) );
00138       }
00139     } 
00140     return result;
00141   }
00142   
00143  private:
00144   double thePtMin; 
00145   double theOriginRadius; 
00146   double theNSigmaZ;
00147   edm::InputTag theBeamSpotTag; 
00148   bool thePrecise;
00149   GlobalVector theDirection;
00150   edm::InputTag theSiPixelRecHits;
00151   bool doVariablePtMin;
00152 
00153   double theSigmaZVertex;
00154   double theFixedError;  
00155   bool theUseFoundVertices;
00156   bool theUseFixedError;
00157   edm::InputTag vertexCollName;
00158 
00159 
00160 };
00161 
00162 #endif