CMS 3D CMS Logo

List of all members | Public Member Functions | Static 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)
 
std::vector< std::unique_ptr< TrackingRegion > > regions (const edm::Event &ev, const edm::EventSetup &es) const override
 
 ~HITrackingRegionForPrimaryVtxProducer () override
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

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< SiPixelRecHitCollectiontheSiPixelRecHitsToken
 
bool theUseFixedError
 
bool theUseFoundVertices
 
edm::InputTag vertexCollName
 
edm::EDGetTokenT< reco::VertexCollectionvertexCollToken
 

Detailed Description

Definition at line 29 of file HITrackingRegionForPrimaryVtxProducer.h.

Constructor & Destructor Documentation

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

Definition at line 33 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.

33  {
34 
35  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
36  thePtMin = regionPSet.getParameter<double>("ptMin");
37  theOriginRadius = regionPSet.getParameter<double>("originRadius");
38  theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
39  theBeamSpotTag = regionPSet.getParameter<edm::InputTag>("beamSpot");
41  thePrecise = regionPSet.getParameter<bool>("precise");
42  theSiPixelRecHits = regionPSet.getParameter<edm::InputTag>("siPixelRecHits");
44  doVariablePtMin = regionPSet.getParameter<bool>("doVariablePtMin");
45  double xDir = regionPSet.getParameter<double>("directionXCoord");
46  double yDir = regionPSet.getParameter<double>("directionYCoord");
47  double zDir = regionPSet.getParameter<double>("directionZCoord");
48  theDirection = GlobalVector(xDir, yDir, zDir);
49 
50  // for using vertex instead of beamspot
51  theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
52  theFixedError = regionPSet.getParameter<double>("fixedError");
53  theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
54  theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
55  vertexCollName = regionPSet.getParameter<edm::InputTag>("VertexCollection");
57  }
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
HITrackingRegionForPrimaryVtxProducer::~HITrackingRegionForPrimaryVtxProducer ( )
inlineoverride

Definition at line 59 of file HITrackingRegionForPrimaryVtxProducer.h.

59 {}

Member Function Documentation

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

Definition at line 90 of file HITrackingRegionForPrimaryVtxProducer.h.

References edmNew::copyDetSetRange(), edm::EventSetup::get(), edm::Event::getByToken(), TrackerTopology::pxbDetIdLayerComparator(), and theSiPixelRecHitsToken.

Referenced by fillDescriptions(), and regions().

91  {
92  //rechits
94  ev.getByToken(theSiPixelRecHitsToken, recHitColl);
95 
97  es.get<TrackerTopologyRcd>().get(httopo);
98 
99  std::vector<const TrackingRecHit*> theChosenHits;
100  edmNew::copyDetSetRange(*recHitColl,theChosenHits, httopo->pxbDetIdLayerComparator(1));
101  return theChosenHits.size();
102 
103  }
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:517
T get() const
Definition: EventSetup.h:71
std::pair< DetId, SameLayerComparator > pxbDetIdLayerComparator(uint32_t layer) const
static void HITrackingRegionForPrimaryVtxProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
inlinestatic

Definition at line 61 of file HITrackingRegionForPrimaryVtxProducer.h.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and estimateMultiplicity().

61  {
63 
64  desc.add<double>("ptMin", 0.7);
65  desc.add<bool>("doVariablePtMin", true);
66  desc.add<double>("originRadius", 0.2);
67  desc.add<double>("nSigmaZ", 3.0);
68  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
69  desc.add<bool>("precise", true);
70  desc.add<bool>("useMultipleScattering", false);
71  desc.add<bool>("useFakeVertices", false);
72  desc.add<edm::InputTag>("siPixelRecHits", edm::InputTag("siPixelRecHits"));
73  desc.add<double>("directionXCoord", 1.0);
74  desc.add<double>("directionYCoord", 1.0);
75  desc.add<double>("directionZCoord", 0.0);
76  desc.add<bool>("useFoundVertices", true);
77  desc.add<edm::InputTag>("VertexCollection", edm::InputTag("hiPixelClusterVertex"));
78  desc.add<bool>("useFixedError", true);
79  desc.add<double>("fixedError", 3.0);
80  desc.add<double>("sigmaZVertex", 3.0);
81 
82  // Only for backwards-compatibility
84  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
85 
86  descriptions.add("hiTrackingRegionFromClusterVtx", descRegion);
87  }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector<std::unique_ptr<TrackingRegion> > HITrackingRegionForPrimaryVtxProducer::regions ( const edm::Event ev,
const edm::EventSetup es 
) const
inlineoverridevirtual

Implements TrackingRegionProducer.

Definition at line 105 of file HITrackingRegionForPrimaryVtxProducer.h.

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

105  {
106 
107  int estMult = estimateMultiplicity(ev, es);
108 
109  // from MC relating first layer pixel hits to "findable" sim tracks with pt>1 GeV
110  float cc = -38.6447;
111  float bb = 0.0581765;
112  float aa = 1.34306e-06;
113 
114  float estTracks = aa*estMult*estMult+bb*estMult+cc;
115 
116  LogTrace("heavyIonHLTVertexing")<<"[HIVertexing]";
117  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: hits in the 1. layer:" << estMult << "]";
118  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: estimated number of tracks:" << estTracks << "]";
119 
120  float regTracking = 60.; //if we have more tracks -> regional tracking
121  float etaB = 10.;
122  float phiB = TMath::Pi()/2.;
123 
124  float decEta = estTracks/90.;
125  etaB = 2.5/decEta;
126 
127  if(estTracks>regTracking) {
128  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Regional Tracking]";
129  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: eta range: -" << etaB << ", "<< etaB <<"]";
130  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: phi range: -" << phiB << ", "<< phiB <<"]";
131  LogTrace("heavyIonHLTVertexing")<<" [Regional Tracking: factor of decrease: " << decEta*2. << "]"; // 2:from phi
132  }
133 
134  float minpt = thePtMin;
135  float varPtCutoff = 1500; //cutoff
136  if(doVariablePtMin && estMult < varPtCutoff) {
137  minpt = 0.075;
138  if(estMult > 0) minpt += estMult * (thePtMin - 0.075)/varPtCutoff; // lower ptMin linearly with pixel hit multiplicity
139  }
140 
141  // tracking region selection
142  std::vector<std::unique_ptr<TrackingRegion> > result;
143  double halflength;
144  GlobalPoint origin;
146  ev.getByToken(theBeamSpotToken, bsHandle);
147  if(bsHandle.isValid()) {
148  const reco::BeamSpot & bs = *bsHandle;
149  origin=GlobalPoint(bs.x0(), bs.y0(), bs.z0());
150  halflength=theNSigmaZ*bs.sigmaZ();
151 
153  {
155  ev.getByToken(vertexCollToken,vertexCollection);
156 
157  for(reco::VertexCollection::const_iterator iV=vertexCollection->begin();
158  iV != vertexCollection->end() ; iV++) {
159  if (iV->isFake() || !iV->isValid()) continue;
160  origin = GlobalPoint(bs.x0(),bs.y0(),iV->z());
161  halflength = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex);
162  }
163  }
164 
165  if(estTracks>regTracking) { // regional tracking
166  result.push_back(
167  std::make_unique<RectangularEtaPhiTrackingRegion>(theDirection, origin, thePtMin, theOriginRadius, halflength, etaB, phiB, RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever, thePrecise) );
168  }
169  else { // global tracking
170  LogTrace("heavyIonHLTVertexing")<<" [HIVertexing: Global Tracking]";
171  result.push_back(
172  std::make_unique<GlobalTrackingRegion>(minpt, origin, theOriginRadius, halflength, thePrecise) );
173  }
174  }
175  return result;
176  }
const double Pi
double z0() const
z coordinate
Definition: BeamSpot.h:68
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
bool isValid() const
Definition: HandleBase.h:74
#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