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  theAllowEmpty = regionPSet.getParameter<bool>("allowEmpty");
45  theMinOriginR = regionPSet.getParameter<double>("minOriginR");
46  theMaxPtMin = regionPSet.getParameter<double>("maxPtMin");
47  theMinHalfLength = regionPSet.getParameter<double>("minHalfLength");
48  theScalingStart = regionPSet.getParameter<double>("scalingStartNPix");
49  theScalingEnd = regionPSet.getParameter<double>("scalingEndNPix");
50  edm::InputTag pixelClustersForScaling = regionPSet.getParameter<edm::InputTag>("pixelClustersForScaling");
53 
54  if (theUseMS) {
55  token_msmaker = iC.esConsumes();
56  }
57  }
58 
60 
61  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
63 
64  desc.add<bool>("precise", true);
65  desc.add<bool>("useMultipleScattering", false);
66  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
67  desc.add<bool>("useFixedError", true);
68  desc.add<double>("originRadius", 0.2);
69  desc.add<double>("sigmaZVertex", 3.0);
70  desc.add<double>("fixedError", 0.2);
71  desc.add<edm::InputTag>("VertexCollection", edm::InputTag("firstStepPrimaryVertices"));
72  desc.add<double>("ptMin", 0.9);
73  desc.add<bool>("useFoundVertices", true);
74  desc.add<bool>("useFakeVertices", false);
75  desc.add<int>("maxNVertices", -1)->setComment("-1 for all vertices");
76  desc.add<double>("nSigmaZ", 4.0);
77  desc.add<edm::InputTag>("pixelClustersForScaling", edm::InputTag("siPixelClusters"));
78  desc.add<bool>("originRScaling4BigEvts", false);
79  desc.add<bool>("ptMinScaling4BigEvts", false);
80  desc.add<bool>("halfLengthScaling4BigEvts", false);
81  desc.add<bool>("allowEmpty", false);
82  desc.add<double>("minOriginR", 0);
83  desc.add<double>("maxPtMin", 1000);
84  desc.add<double>("minHalfLength", 0);
85  desc.add<double>("scalingStartNPix", 0.0);
86  desc.add<double>("scalingEndNPix", 1.0);
87 
88  // Only for backwards-compatibility
90  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
91 
92  descriptions.add("globalTrackingRegionWithVertices", descRegion);
93  }
94 
95  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev,
96  const edm::EventSetup& es) const override {
97  std::vector<std::unique_ptr<TrackingRegion> > result;
98 
99  GlobalPoint theOrigin;
101  ev.getByToken(token_beamSpot, bsHandle);
102  double bsSigmaZ;
103  if (bsHandle.isValid()) {
104  const reco::BeamSpot& bs = *bsHandle;
105  bsSigmaZ = theNSigmaZ * bs.sigmaZ();
106  theOrigin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
107  } else {
108  throw cms::Exception("Seeding") << "ERROR: input beamSpot is not valid in GlobalTrackingRegionWithVertices";
109  }
110 
111  const MultipleScatteringParametrisationMaker* msmaker = nullptr;
112  if (theUseMS) {
113  msmaker = &es.getData(token_msmaker);
114  }
115 
116  if (theUseFoundVertices) {
118  ev.getByToken(token_vertex, vertexCollection);
119 
122  if (doScaling)
123  ev.getByToken(token_pc, pixelClusterDSV);
124 
125  for (reco::VertexCollection::const_iterator iV = vertexCollection->begin(); iV != vertexCollection->end(); iV++) {
126  if (!iV->isValid())
127  continue;
128  if (iV->isFake() && !(theUseFakeVertices && theUseFixedError))
129  continue;
130  GlobalPoint theOrigin_ = GlobalPoint(iV->x(), iV->y(), iV->z());
131 
132  //scaling origin radius, half length, min pt for high-occupancy HI events to keep timing reasonable
133  if (doScaling) {
134  //Use the unscaled radius unless one of the two conditions below is met
135  double scaledOriginRadius = theOriginRadius;
136  double scaledHalfLength = theFixedError;
137  double scaledPtMin = thePtMin;
138 
139  //calculate nPixels (adapted from TkSeedGenerator/src/ClusterChecker.cc)
140  double nPix = 0;
141 
142  const edmNew::DetSetVector<SiPixelCluster>& input = *pixelClusterDSV;
143  nPix = input.dataSize();
144 
145  //first condition is for high occupancy, second makes sure we won't divide by zero or a negative number
146  if ((nPix > theScalingEnd) || ((theScalingEnd - theScalingStart) <= 0)) {
147  if (theOriginRScaling)
148  scaledOriginRadius = theMinOriginR; // sets parameters to minimum value from PSet
150  scaledHalfLength = theMinHalfLength;
151  if (thePtMinScaling)
152  scaledPtMin = theMaxPtMin;
153  }
154  //second condition - scale radius linearly by Npix in the region from ScalingStart to ScalingEnd
155  else if ((nPix <= theScalingEnd) && (nPix > theScalingStart)) {
156  float slopeFactor = (nPix - theScalingStart) / (theScalingEnd - theScalingStart);
157  if (theOriginRScaling)
158  scaledOriginRadius = theOriginRadius - (theOriginRadius - theMinOriginR) * slopeFactor;
160  scaledHalfLength = theFixedError - (theFixedError - theMinHalfLength) * slopeFactor;
161  if (thePtMinScaling)
162  scaledPtMin = thePtMin - (thePtMin - theMaxPtMin) * slopeFactor;
163  }
164  //if region has 0 size, return 'result' empty, otherwise make a tracking region
165  if (scaledOriginRadius != 0 && scaledHalfLength != 0) {
166  result.push_back(std::make_unique<GlobalTrackingRegion>(
167  scaledPtMin, theOrigin_, scaledOriginRadius, scaledHalfLength, thePrecise, theUseMS, msmaker));
168  }
169  } //end of region scaling code, pp behavior below
170 
171  else {
172  double theOriginHalfLength_ = (theUseFixedError ? theFixedError : (iV->zError()) * theSigmaZVertex);
173  result.push_back(std::make_unique<GlobalTrackingRegion>(
174  thePtMin, theOrigin_, theOriginRadius, theOriginHalfLength_, thePrecise, theUseMS, msmaker));
175  if (theMaxNVertices >= 0 && result.size() >= static_cast<unsigned>(theMaxNVertices))
176  break;
177  }
178  }
179 
181  result.push_back(std::make_unique<GlobalTrackingRegion>(
182  thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS, msmaker));
183  }
184  } else {
185  result.push_back(std::make_unique<GlobalTrackingRegion>(
186  thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS, msmaker));
187  }
188 
189  return result;
190  }
191 
192 private:
193  double thePtMin;
195  double theNSigmaZ;
197 
202  bool theUseMS;
203 
210 
211  //HI-related variables
217  double theMaxPtMin;
222 };
223 
224 #endif
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
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:50
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)