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)
 
virtual std::vector< std::unique_ptr< TrackingRegion > > regions (const edm::Event &ev, const edm::EventSetup &es) const
 
virtual ~HITrackingRegionForPrimaryVtxProducer ()
 
- 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 28 of file HITrackingRegionForPrimaryVtxProducer.h.

Constructor & Destructor Documentation

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

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

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

58 {}

Member Function Documentation

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

Definition at line 89 of file HITrackingRegionForPrimaryVtxProducer.h.

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

Referenced by fillDescriptions(), and regions().

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

Definition at line 60 of file HITrackingRegionForPrimaryVtxProducer.h.

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

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

Implements TrackingRegionProducer.

Definition at line 104 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().

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