#include <HITrackingRegionForPrimaryVtxProducer.h>
Public Member Functions | |
int | estimateMultiplicity (const edm::Event &ev, const edm::EventSetup &es) const |
HITrackingRegionForPrimaryVtxProducer (const edm::ParameterSet &cfg) | |
virtual std::vector < TrackingRegion * > | regions (const edm::Event &ev, const edm::EventSetup &es) const |
virtual | ~HITrackingRegionForPrimaryVtxProducer () |
Private Attributes | |
bool | doVariablePtMin |
edm::InputTag | theBeamSpotTag |
GlobalVector | theDirection |
double | theFixedError |
double | theNSigmaZ |
double | theOriginRadius |
bool | thePrecise |
double | thePtMin |
double | theSigmaZVertex |
edm::InputTag | theSiPixelRecHits |
bool | theUseFixedError |
bool | theUseFoundVertices |
edm::InputTag | vertexCollName |
Definition at line 26 of file HITrackingRegionForPrimaryVtxProducer.h.
HITrackingRegionForPrimaryVtxProducer::HITrackingRegionForPrimaryVtxProducer | ( | const edm::ParameterSet & | cfg | ) | [inline] |
Definition at line 30 of file HITrackingRegionForPrimaryVtxProducer.h.
References doVariablePtMin, edm::ParameterSet::getParameter(), theBeamSpotTag, theDirection, theFixedError, theNSigmaZ, theOriginRadius, thePrecise, thePtMin, theSigmaZVertex, theSiPixelRecHits, theUseFixedError, theUseFoundVertices, and vertexCollName.
{ edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet"); thePtMin = regionPSet.getParameter<double>("ptMin"); theOriginRadius = regionPSet.getParameter<double>("originRadius"); theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ"); theBeamSpotTag = regionPSet.getParameter<edm::InputTag>("beamSpot"); thePrecise = regionPSet.getParameter<bool>("precise"); theSiPixelRecHits = regionPSet.getParameter<edm::InputTag>("siPixelRecHits"); doVariablePtMin = regionPSet.getParameter<bool>("doVariablePtMin"); double xDir = regionPSet.getParameter<double>("directionXCoord"); double yDir = regionPSet.getParameter<double>("directionYCoord"); double zDir = regionPSet.getParameter<double>("directionZCoord"); theDirection = GlobalVector(xDir, yDir, zDir); // for using vertex instead of beamspot theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex"); theFixedError = regionPSet.getParameter<double>("fixedError"); theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices"); theUseFixedError = regionPSet.getParameter<bool>("useFixedError"); vertexCollName = regionPSet.getParameter<edm::InputTag>("VertexCollection"); }
virtual HITrackingRegionForPrimaryVtxProducer::~HITrackingRegionForPrimaryVtxProducer | ( | ) | [inline, virtual] |
Definition at line 53 of file HITrackingRegionForPrimaryVtxProducer.h.
{}
int HITrackingRegionForPrimaryVtxProducer::estimateMultiplicity | ( | const edm::Event & | ev, |
const edm::EventSetup & | es | ||
) | const [inline] |
Definition at line 56 of file HITrackingRegionForPrimaryVtxProducer.h.
References edmNew::copyDetSetRange(), edm::Event::getByLabel(), TrackerLayerIdAccessor::pixelBarrelLayer(), and theSiPixelRecHits.
Referenced by regions().
{ //rechits edm::Handle<SiPixelRecHitCollection> recHitColl; ev.getByLabel(theSiPixelRecHits, recHitColl); std::vector<const TrackingRecHit*> theChosenHits; TrackerLayerIdAccessor acc; edmNew::copyDetSetRange(*recHitColl,theChosenHits,acc.pixelBarrelLayer(1)); return theChosenHits.size(); }
virtual std::vector<TrackingRegion* > HITrackingRegionForPrimaryVtxProducer::regions | ( | const edm::Event & | ev, |
const edm::EventSetup & | es | ||
) | const [inline, virtual] |
Implements TrackingRegionProducer.
Definition at line 69 of file HITrackingRegionForPrimaryVtxProducer.h.
References doVariablePtMin, estimateMultiplicity(), edm::Event::getByLabel(), edm::HandleBase::isValid(), LogTrace, Pi, query::result, reco::BeamSpot::sigmaZ(), theBeamSpotTag, theDirection, theFixedError, theNSigmaZ, theOriginRadius, thePrecise, thePtMin, theSigmaZVertex, theUseFixedError, theUseFoundVertices, GoodVertex_cfg::vertexCollection, vertexCollName, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().
{ int estMult = estimateMultiplicity(ev, es); // from MC relating first layer pixel hits to "findable" sim tracks with pt>1 GeV float cc = -38.6447; float bb = 0.0581765; float aa = 1.34306e-06; float estTracks = aa*estMult*estMult+bb*estMult+cc; LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]"; LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]"; LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]"; float regTracking = 60.; //if we have more tracks -> regional tracking float etaB = 10.; float phiB = TMath::Pi()/2.; float decEta = estTracks/90.; etaB = 2.5/decEta; if(estTracks>regTracking) { LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]"; LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]"; LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]"; LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: factor of decrease: " << decEta*2. << "]"; // 2:from phi } float minpt = thePtMin; float varPtCutoff = 1500; //cutoff if(doVariablePtMin && estMult < varPtCutoff) { minpt = 0.075; if(estMult > 0) minpt += estMult * (thePtMin - 0.075)/varPtCutoff; // lower ptMin linearly with pixel hit multiplicity } // tracking region selection std::vector<TrackingRegion* > result; double halflength; GlobalPoint origin; edm::Handle<reco::BeamSpot> bsHandle; ev.getByLabel( theBeamSpotTag, bsHandle); if(bsHandle.isValid()) { const reco::BeamSpot & bs = *bsHandle; origin=GlobalPoint(bs.x0(), bs.y0(), bs.z0()); halflength=theNSigmaZ*bs.sigmaZ(); if(theUseFoundVertices) { edm::Handle<reco::VertexCollection> vertexCollection; ev.getByLabel(vertexCollName,vertexCollection); for(reco::VertexCollection::const_iterator iV=vertexCollection->begin(); iV != vertexCollection->end() ; iV++) { if (iV->isFake() || !iV->isValid()) continue; origin = GlobalPoint(bs.x0(),bs.y0(),iV->z()); halflength = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex); } } if(estTracks>regTracking) { // regional tracking result.push_back( new RectangularEtaPhiTrackingRegion(theDirection, origin, thePtMin, theOriginRadius, halflength, etaB, phiB, 0, thePrecise) ); } else { // global tracking LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]"; result.push_back( new GlobalTrackingRegion(minpt, origin, theOriginRadius, halflength, thePrecise) ); } } return result; }
bool HITrackingRegionForPrimaryVtxProducer::doVariablePtMin [private] |
Definition at line 150 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
Definition at line 146 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
Definition at line 148 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
double HITrackingRegionForPrimaryVtxProducer::theFixedError [private] |
Definition at line 153 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
double HITrackingRegionForPrimaryVtxProducer::theNSigmaZ [private] |
Definition at line 145 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
double HITrackingRegionForPrimaryVtxProducer::theOriginRadius [private] |
Definition at line 144 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
bool HITrackingRegionForPrimaryVtxProducer::thePrecise [private] |
Definition at line 147 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
double HITrackingRegionForPrimaryVtxProducer::thePtMin [private] |
Definition at line 143 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
double HITrackingRegionForPrimaryVtxProducer::theSigmaZVertex [private] |
Definition at line 152 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
Definition at line 149 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by estimateMultiplicity(), and HITrackingRegionForPrimaryVtxProducer().
bool HITrackingRegionForPrimaryVtxProducer::theUseFixedError [private] |
Definition at line 155 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
bool HITrackingRegionForPrimaryVtxProducer::theUseFoundVertices [private] |
Definition at line 154 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().
Definition at line 156 of file HITrackingRegionForPrimaryVtxProducer.h.
Referenced by HITrackingRegionForPrimaryVtxProducer(), and regions().