CMS 3D CMS Logo

GlobalTrackingRegionWithVerticesProducer.h
Go to the documentation of this file.
1 #ifndef RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H
2 #define RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H
3 
17 
19 public:
21  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
22 
23  thePtMin = regionPSet.getParameter<double>("ptMin");
24  theOriginRadius = regionPSet.getParameter<double>("originRadius");
25  theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
26  token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
27  thePrecise = regionPSet.getParameter<bool>("precise");
28  theUseMS = regionPSet.getParameter<bool>("useMultipleScattering");
29 
30  theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
31  theFixedError = regionPSet.getParameter<double>("fixedError");
32 
33  theMaxNVertices = regionPSet.getParameter<int>("maxNVertices");
34 
35  theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
36  theUseFakeVertices = regionPSet.getParameter<bool>("useFakeVertices");
37  theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
38  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("VertexCollection"));
39 
40  //information for Heavy ion region scaling
41  theOriginRScaling = regionPSet.getParameter<bool>("originRScaling4BigEvts");
42  thePtMinScaling = regionPSet.getParameter<bool>("ptMinScaling4BigEvts");
43  theHalfLengthScaling = regionPSet.getParameter<bool>("halfLengthScaling4BigEvts");
44  theMinOriginR = regionPSet.getParameter<double>("minOriginR");
45  theMaxPtMin = regionPSet.getParameter<double>("maxPtMin");
46  theMinHalfLength = regionPSet.getParameter<double>("minHalfLength");
47  theScalingStart = regionPSet.getParameter<double>("scalingStartNPix");
48  theScalingEnd = regionPSet.getParameter<double>("scalingEndNPix");
49  edm::InputTag pixelClustersForScaling = regionPSet.getParameter<edm::InputTag>("pixelClustersForScaling");
52 
53  if (theUseMS) {
54  token_msmaker = iC.esConsumes();
55  }
56  }
57 
59 
60  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
62 
63  desc.add<bool>("precise", true);
64  desc.add<bool>("useMultipleScattering", false);
65  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
66  desc.add<bool>("useFixedError", true);
67  desc.add<double>("originRadius", 0.2);
68  desc.add<double>("sigmaZVertex", 3.0);
69  desc.add<double>("fixedError", 0.2);
70  desc.add<edm::InputTag>("VertexCollection", edm::InputTag("firstStepPrimaryVertices"));
71  desc.add<double>("ptMin", 0.9);
72  desc.add<bool>("useFoundVertices", true);
73  desc.add<bool>("useFakeVertices", false);
74  desc.add<int>("maxNVertices", -1)->setComment("-1 for all vertices");
75  desc.add<double>("nSigmaZ", 4.0);
76  desc.add<edm::InputTag>("pixelClustersForScaling", edm::InputTag("siPixelClusters"));
77  desc.add<bool>("originRScaling4BigEvts", false);
78  desc.add<bool>("ptMinScaling4BigEvts", false);
79  desc.add<bool>("halfLengthScaling4BigEvts", false);
80  desc.add<double>("minOriginR", 0);
81  desc.add<double>("maxPtMin", 1000);
82  desc.add<double>("minHalfLength", 0);
83  desc.add<double>("scalingStartNPix", 0.0);
84  desc.add<double>("scalingEndNPix", 1.0);
85 
86  // Only for backwards-compatibility
88  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
89 
90  descriptions.add("globalTrackingRegionWithVertices", descRegion);
91  }
92 
93  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev,
94  const edm::EventSetup& es) const override {
95  std::vector<std::unique_ptr<TrackingRegion> > result;
96 
97  GlobalPoint theOrigin;
99  ev.getByToken(token_beamSpot, bsHandle);
100  double bsSigmaZ;
101  if (bsHandle.isValid()) {
102  const reco::BeamSpot& bs = *bsHandle;
103  bsSigmaZ = theNSigmaZ * bs.sigmaZ();
104  theOrigin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
105  } else {
106  throw cms::Exception("Seeding") << "ERROR: input beamSpot is not valid in GlobalTrackingRegionWithVertices";
107  }
108 
109  const MultipleScatteringParametrisationMaker* msmaker = nullptr;
110  if (theUseMS) {
111  msmaker = &es.getData(token_msmaker);
112  }
113 
114  if (theUseFoundVertices) {
116  ev.getByToken(token_vertex, vertexCollection);
117 
120  if (doScaling)
121  ev.getByToken(token_pc, pixelClusterDSV);
122 
123  for (reco::VertexCollection::const_iterator iV = vertexCollection->begin(); iV != vertexCollection->end(); iV++) {
124  if (!iV->isValid())
125  continue;
126  if (iV->isFake() && !(theUseFakeVertices && theUseFixedError))
127  continue;
128  GlobalPoint theOrigin_ = GlobalPoint(iV->x(), iV->y(), iV->z());
129 
130  //scaling origin radius, half length, min pt for high-occupancy HI events to keep timing reasonable
131  if (doScaling) {
132  //Use the unscaled radius unless one of the two conditions below is met
133  double scaledOriginRadius = theOriginRadius;
134  double scaledHalfLength = theFixedError;
135  double scaledPtMin = thePtMin;
136 
137  //calculate nPixels (adapted from TkSeedGenerator/src/ClusterChecker.cc)
138  double nPix = 0;
139 
140  const edmNew::DetSetVector<SiPixelCluster>& input = *pixelClusterDSV;
141  nPix = input.dataSize();
142 
143  //first condition is for high occupancy, second makes sure we won't divide by zero or a negative number
144  if ((nPix > theScalingEnd) || ((theScalingEnd - theScalingStart) <= 0)) {
145  if (theOriginRScaling)
146  scaledOriginRadius = theMinOriginR; // sets parameters to minimum value from PSet
148  scaledHalfLength = theMinHalfLength;
149  if (thePtMinScaling)
150  scaledPtMin = theMaxPtMin;
151  }
152  //second condition - scale radius linearly by Npix in the region from ScalingStart to ScalingEnd
153  else if ((nPix <= theScalingEnd) && (nPix > theScalingStart)) {
154  float slopeFactor = (nPix - theScalingStart) / (theScalingEnd - theScalingStart);
155  if (theOriginRScaling)
156  scaledOriginRadius = theOriginRadius - (theOriginRadius - theMinOriginR) * slopeFactor;
158  scaledHalfLength = theFixedError - (theFixedError - theMinHalfLength) * slopeFactor;
159  if (thePtMinScaling)
160  scaledPtMin = thePtMin - (thePtMin - theMaxPtMin) * slopeFactor;
161  }
162  //if region has 0 size, return 'result' empty, otherwise make a tracking region
163  if (scaledOriginRadius != 0 && scaledHalfLength != 0) {
164  result.push_back(std::make_unique<GlobalTrackingRegion>(
165  scaledPtMin, theOrigin_, scaledOriginRadius, scaledHalfLength, thePrecise, theUseMS, msmaker));
166  }
167  } //end of region scaling code, pp behavior below
168 
169  else {
170  double theOriginHalfLength_ = (theUseFixedError ? theFixedError : (iV->zError()) * theSigmaZVertex);
171  result.push_back(std::make_unique<GlobalTrackingRegion>(
172  thePtMin, theOrigin_, theOriginRadius, theOriginHalfLength_, thePrecise, theUseMS, msmaker));
173  if (theMaxNVertices >= 0 && result.size() >= static_cast<unsigned>(theMaxNVertices))
174  break;
175  }
176  }
177 
179  result.push_back(std::make_unique<GlobalTrackingRegion>(
180  thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS, msmaker));
181  }
182  } else {
183  result.push_back(std::make_unique<GlobalTrackingRegion>(
184  thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS, msmaker));
185  }
186 
187  return result;
188  }
189 
190 private:
191  double thePtMin;
193  double theNSigmaZ;
195 
200  bool theUseMS;
201 
208 
209  //HI-related variables
214  double theMaxPtMin;
219 };
220 
221 #endif
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::EDGetTokenT< reco::VertexCollection > token_vertex
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > token_pc
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &ev, const edm::EventSetup &es) const override
static std::string const input
Definition: EdmProvDump.cc:47
bool getData(T &iHolder) const
Definition: EventSetup.h:122
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
GlobalTrackingRegionWithVerticesProducer(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)