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