CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
HITrackingRegionForPrimaryVtxProducer Class Reference

#include <HITrackingRegionForPrimaryVtxProducer.h>

Inheritance diagram for HITrackingRegionForPrimaryVtxProducer:
TrackingRegionProducer

Public Member Functions

int estimateMultiplicity (const edm::Event &ev, const edm::EventSetup &es) const
 
 HITrackingRegionForPrimaryVtxProducer (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
virtual std::vector
< TrackingRegion * > 
regions (const edm::Event &ev, const edm::EventSetup &es) const
 
virtual ~HITrackingRegionForPrimaryVtxProducer ()
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Private Attributes

bool doVariablePtMin
 
edm::InputTag theBeamSpotTag
 
edm::EDGetTokenT< reco::BeamSpottheBeamSpotToken
 
GlobalVector theDirection
 
double theFixedError
 
double theNSigmaZ
 
double theOriginRadius
 
bool thePrecise
 
double thePtMin
 
double theSigmaZVertex
 
edm::InputTag theSiPixelRecHits
 
edm::EDGetTokenT
< SiPixelRecHitCollection
theSiPixelRecHitsToken
 
bool theUseFixedError
 
bool theUseFoundVertices
 
edm::InputTag vertexCollName
 
edm::EDGetTokenT
< reco::VertexCollection
vertexCollToken
 

Detailed Description

Definition at line 27 of file HITrackingRegionForPrimaryVtxProducer.h.

Constructor & Destructor Documentation

HITrackingRegionForPrimaryVtxProducer::HITrackingRegionForPrimaryVtxProducer ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 31 of file HITrackingRegionForPrimaryVtxProducer.h.

References doVariablePtMin, edm::ParameterSet::getParameter(), theBeamSpotTag, theBeamSpotToken, theDirection, theFixedError, theNSigmaZ, theOriginRadius, thePrecise, thePtMin, theSigmaZVertex, theSiPixelRecHits, theSiPixelRecHitsToken, theUseFixedError, theUseFoundVertices, vertexCollName, and vertexCollToken.

31  {
32 
33  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
34  thePtMin = regionPSet.getParameter<double>("ptMin");
35  theOriginRadius = regionPSet.getParameter<double>("originRadius");
36  theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
37  theBeamSpotTag = regionPSet.getParameter<edm::InputTag>("beamSpot");
39  thePrecise = regionPSet.getParameter<bool>("precise");
40  theSiPixelRecHits = regionPSet.getParameter<edm::InputTag>("siPixelRecHits");
42  doVariablePtMin = regionPSet.getParameter<bool>("doVariablePtMin");
43  double xDir = regionPSet.getParameter<double>("directionXCoord");
44  double yDir = regionPSet.getParameter<double>("directionYCoord");
45  double zDir = regionPSet.getParameter<double>("directionZCoord");
46  theDirection = GlobalVector(xDir, yDir, zDir);
47 
48  // for using vertex instead of beamspot
49  theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
50  theFixedError = regionPSet.getParameter<double>("fixedError");
51  theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
52  theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
53  vertexCollName = regionPSet.getParameter<edm::InputTag>("VertexCollection");
55  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< SiPixelRecHitCollection > theSiPixelRecHitsToken
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > vertexCollToken
Global3DVector GlobalVector
Definition: GlobalVector.h:10
virtual HITrackingRegionForPrimaryVtxProducer::~HITrackingRegionForPrimaryVtxProducer ( )
inlinevirtual

Definition at line 57 of file HITrackingRegionForPrimaryVtxProducer.h.

57 {}

Member Function Documentation

int HITrackingRegionForPrimaryVtxProducer::estimateMultiplicity ( const edm::Event ev,
const edm::EventSetup es 
) const
inline

Definition at line 60 of file HITrackingRegionForPrimaryVtxProducer.h.

References edmNew::copyDetSetRange(), edm::Event::getByToken(), TrackerLayerIdAccessor::pixelBarrelLayer(), and theSiPixelRecHitsToken.

Referenced by regions().

61  {
62  //rechits
64  ev.getByToken(theSiPixelRecHitsToken, recHitColl);
65 
66  std::vector<const TrackingRecHit*> theChosenHits;
68  edmNew::copyDetSetRange(*recHitColl,theChosenHits,acc.pixelBarrelLayer(1));
69  return theChosenHits.size();
70 
71  }
void copyDetSetRange(DSTV const &dstv, std::vector< T const * > &v, std::pair< A, B > const &sel)
edm::EDGetTokenT< SiPixelRecHitCollection > theSiPixelRecHitsToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::pair< DetId, DetIdPXBSameLayerComparator > pixelBarrelLayer(int layer)
virtual std::vector<TrackingRegion* > HITrackingRegionForPrimaryVtxProducer::regions ( const edm::Event ev,
const edm::EventSetup es 
) const
inlinevirtual

Implements TrackingRegionProducer.

Definition at line 73 of file HITrackingRegionForPrimaryVtxProducer.h.

References doVariablePtMin, estimateMultiplicity(), edm::Event::getByToken(), edm::HandleBase::isValid(), LogTrace, Pi, query::result, reco::BeamSpot::sigmaZ(), theBeamSpotToken, theDirection, theFixedError, theNSigmaZ, theOriginRadius, thePrecise, thePtMin, theSigmaZVertex, theUseFixedError, theUseFoundVertices, GoodVertex_cfg::vertexCollection, vertexCollToken, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().

73  {
74 
75  int estMult = estimateMultiplicity(ev, es);
76 
77  // from MC relating first layer pixel hits to "findable" sim tracks with pt>1 GeV
78  float cc = -38.6447;
79  float bb = 0.0581765;
80  float aa = 1.34306e-06;
81 
82  float estTracks = aa*estMult*estMult+bb*estMult+cc;
83 
84  LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]";
85  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]";
86  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]";
87 
88  float regTracking = 60.; //if we have more tracks -> regional tracking
89  float etaB = 10.;
90  float phiB = TMath::Pi()/2.;
91 
92  float decEta = estTracks/90.;
93  etaB = 2.5/decEta;
94 
95  if(estTracks>regTracking) {
96  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]";
97  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]";
98  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]";
99  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: factor of decrease: " << decEta*2. << "]"; // 2:from phi
100  }
101 
102  float minpt = thePtMin;
103  float varPtCutoff = 1500; //cutoff
104  if(doVariablePtMin && estMult < varPtCutoff) {
105  minpt = 0.075;
106  if(estMult > 0) minpt += estMult * (thePtMin - 0.075)/varPtCutoff; // lower ptMin linearly with pixel hit multiplicity
107  }
108 
109  // tracking region selection
110  std::vector<TrackingRegion* > result;
111  double halflength;
112  GlobalPoint origin;
114  ev.getByToken(theBeamSpotToken, bsHandle);
115  if(bsHandle.isValid()) {
116  const reco::BeamSpot & bs = *bsHandle;
117  origin=GlobalPoint(bs.x0(), bs.y0(), bs.z0());
118  halflength=theNSigmaZ*bs.sigmaZ();
119 
121  {
123  ev.getByToken(vertexCollToken,vertexCollection);
124 
125  for(reco::VertexCollection::const_iterator iV=vertexCollection->begin();
126  iV != vertexCollection->end() ; iV++) {
127  if (iV->isFake() || !iV->isValid()) continue;
128  origin = GlobalPoint(bs.x0(),bs.y0(),iV->z());
129  halflength = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex);
130  }
131  }
132 
133  if(estTracks>regTracking) { // regional tracking
134  result.push_back(
135  new RectangularEtaPhiTrackingRegion(theDirection, origin, thePtMin, theOriginRadius, halflength, etaB, phiB, 0, thePrecise) );
136  }
137  else { // global tracking
138  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]";
139  result.push_back(
140  new GlobalTrackingRegion(minpt, origin, theOriginRadius, halflength, thePrecise) );
141  }
142  }
143  return result;
144  }
const double Pi
double z0() const
z coordinate
Definition: BeamSpot.h:68
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
tuple vertexCollection
tuple result
Definition: query.py:137
bool isValid() const
Definition: HandleBase.h:76
#define LogTrace(id)
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
edm::EDGetTokenT< reco::VertexCollection > vertexCollToken
double y0() const
y coordinate
Definition: BeamSpot.h:66
int estimateMultiplicity(const edm::Event &ev, const edm::EventSetup &es) const
double x0() const
x coordinate
Definition: BeamSpot.h:64

Member Data Documentation

bool HITrackingRegionForPrimaryVtxProducer::doVariablePtMin
private
edm::InputTag HITrackingRegionForPrimaryVtxProducer::theBeamSpotTag
private
edm::EDGetTokenT<reco::BeamSpot> HITrackingRegionForPrimaryVtxProducer::theBeamSpotToken
private
GlobalVector HITrackingRegionForPrimaryVtxProducer::theDirection
private
double HITrackingRegionForPrimaryVtxProducer::theFixedError
private
double HITrackingRegionForPrimaryVtxProducer::theNSigmaZ
private
double HITrackingRegionForPrimaryVtxProducer::theOriginRadius
private
bool HITrackingRegionForPrimaryVtxProducer::thePrecise
private
double HITrackingRegionForPrimaryVtxProducer::thePtMin
private
double HITrackingRegionForPrimaryVtxProducer::theSigmaZVertex
private
edm::InputTag HITrackingRegionForPrimaryVtxProducer::theSiPixelRecHits
private
edm::EDGetTokenT<SiPixelRecHitCollection> HITrackingRegionForPrimaryVtxProducer::theSiPixelRecHitsToken
private
bool HITrackingRegionForPrimaryVtxProducer::theUseFixedError
private
bool HITrackingRegionForPrimaryVtxProducer::theUseFoundVertices
private
edm::InputTag HITrackingRegionForPrimaryVtxProducer::vertexCollName
private
edm::EDGetTokenT<reco::VertexCollection> HITrackingRegionForPrimaryVtxProducer::vertexCollToken
private