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