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

56 {}

Member Function Documentation

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

Definition at line 86 of file HITrackingRegionForPrimaryVtxProducer.h.

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

Referenced by regions().

86  {
87  //rechits
89  ev.getByToken(theSiPixelRecHitsToken, recHitColl);
90 
92  es.get<TrackerTopologyRcd>().get(httopo);
93 
94  std::vector<const TrackingRecHit*> theChosenHits;
95  edmNew::copyDetSetRange(*recHitColl, theChosenHits, httopo->pxbDetIdLayerComparator(1));
96  return theChosenHits.size();
97  }
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:525
T get() const
Definition: EventSetup.h:73
std::pair< DetId, SameLayerComparator > pxbDetIdLayerComparator(uint32_t layer) const
static void HITrackingRegionForPrimaryVtxProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
inlinestatic

Definition at line 58 of file HITrackingRegionForPrimaryVtxProducer.h.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and HLT_2018_cff::InputTag.

58  {
60 
61  desc.add<double>("ptMin", 0.7);
62  desc.add<bool>("doVariablePtMin", true);
63  desc.add<double>("originRadius", 0.2);
64  desc.add<double>("nSigmaZ", 3.0);
65  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
66  desc.add<bool>("precise", true);
67  desc.add<bool>("useMultipleScattering", false);
68  desc.add<bool>("useFakeVertices", false);
69  desc.add<edm::InputTag>("siPixelRecHits", edm::InputTag("siPixelRecHits"));
70  desc.add<double>("directionXCoord", 1.0);
71  desc.add<double>("directionYCoord", 1.0);
72  desc.add<double>("directionZCoord", 0.0);
73  desc.add<bool>("useFoundVertices", true);
74  desc.add<edm::InputTag>("VertexCollection", edm::InputTag("hiPixelClusterVertex"));
75  desc.add<bool>("useFixedError", true);
76  desc.add<double>("fixedError", 3.0);
77  desc.add<double>("sigmaZVertex", 3.0);
78 
79  // Only for backwards-compatibility
81  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
82 
83  descriptions.add("hiTrackingRegionFromClusterVtx", descRegion);
84  }
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 99 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, spclusmultinvestigator_cfi::vertexCollection, vertexCollToken, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().

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

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