CMS 3D CMS Logo

PointSeededTrackingRegionsProducer.h
Go to the documentation of this file.
1 #ifndef PointSeededTrackingRegionsProducer_h
2 #define PointSeededTrackingRegionsProducer_h
3 
7 
20 
21 #include <TMath.h>
22 #include <TLorentzVector.h>
23 
49 public:
51 
53  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
54 
55  // operation mode
56  std::string modeString = regPSet.getParameter<std::string>("mode");
57  if (modeString == "BeamSpotFixed")
59  else if (modeString == "BeamSpotSigma")
61  else if (modeString == "VerticesFixed")
63  else if (modeString == "VerticesSigma")
65  else
66  edm::LogError("PointSeededTrackingRegionsProducer") << "Unknown mode string: " << modeString;
67 
68  // basic inputsi
70  etaPoints = points.getParameter<std::vector<double>>("eta");
71  phiPoints = points.getParameter<std::vector<double>>("phi");
72  if (!(etaPoints.size() == phiPoints.size()))
73  throw edm::Exception(edm::errors::Configuration) << "The parameters 'eta' and 'phi' must have the same size";
74  ;
75  m_maxNRegions = regPSet.getParameter<int>("maxNRegions");
76  token_beamSpot = iC.consumes<reco::BeamSpot>(regPSet.getParameter<edm::InputTag>("beamSpot"));
77  m_maxNVertices = 1;
79  token_vertex = iC.consumes<reco::VertexCollection>(regPSet.getParameter<edm::InputTag>("vertexCollection"));
80  m_maxNVertices = regPSet.getParameter<int>("maxNVertices");
81  }
82 
83  // RectangularEtaPhiTrackingRegion parameters:
84  m_ptMin = regPSet.getParameter<double>("ptMin");
85  m_originRadius = regPSet.getParameter<double>("originRadius");
86  m_zErrorBeamSpot = regPSet.getParameter<double>("zErrorBeamSpot");
87  m_deltaEta = regPSet.getParameter<double>("deltaEta");
88  m_deltaPhi = regPSet.getParameter<double>("deltaPhi");
89  m_precise = regPSet.getParameter<bool>("precise");
91  regPSet.getParameter<std::string>("whereToUseMeasurementTracker"));
92  if (m_whereToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever) {
94  iC.consumes<MeasurementTrackerEvent>(regPSet.getParameter<edm::InputTag>("measurementTrackerName"));
95  }
96  m_searchOpt = false;
97  if (regPSet.exists("searchOpt"))
98  m_searchOpt = regPSet.getParameter<bool>("searchOpt");
99 
100  // mode-dependent z-halflength of tracking regions
101  if (m_mode == VERTICES_SIGMA)
102  m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
103  if (m_mode == VERTICES_FIXED)
104  m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
105  m_nSigmaZBeamSpot = -1.;
106  if (m_mode == BEAM_SPOT_SIGMA) {
107  m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
108  if (m_nSigmaZBeamSpot < 0.)
109  edm::LogError("PointSeededTrackingRegionsProducer")
110  << "nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
111  }
112  }
113 
115 
118 
119  edm::ParameterSetDescription descPoints;
120  descPoints.add<std::vector<double>>("eta", {0.});
121  descPoints.add<std::vector<double>>("phi", {0.});
122  desc.add<edm::ParameterSetDescription>("points", descPoints);
123 
124  desc.add<std::string>("mode", "BeamSpotFixed");
125  desc.add<int>("maxNRegions", 10);
126  desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
127  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
128  desc.add<int>("maxNVertices", 1);
129 
130  desc.add<double>("ptMin", 0.9);
131  desc.add<double>("originRadius", 0.2);
132  desc.add<double>("zErrorBeamSpot", 24.2);
133  desc.add<double>("deltaEta", 0.5);
134  desc.add<double>("deltaPhi", 0.5);
135  desc.add<bool>("precise", true);
136 
137  desc.add<double>("nSigmaZVertex", 3.);
138  desc.add<double>("zErrorVetex", 0.2);
139  desc.add<double>("nSigmaZBeamSpot", 4.);
140 
141  desc.add<std::string>("whereToUseMeasurementTracker", "ForSiStrips");
142  desc.add<edm::InputTag>("measurementTrackerName", edm::InputTag(""));
143 
144  desc.add<bool>("searchOpt", false);
145 
146  // Only for backwards-compatibility
147  edm::ParameterSetDescription descRegion;
148  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
149  //edm::ParameterSetDescription descPoint;
150  //descPoint.add<edm::ParameterSetDescription>("point_input", desc);
151 
152  descriptions.add("pointSeededTrackingRegion", descRegion);
153  }
154 
155  std::vector<std::unique_ptr<TrackingRegion>> regions(const edm::Event& e, const edm::EventSetup& es) const override {
156  std::vector<std::unique_ptr<TrackingRegion>> result;
157 
158  // always need the beam spot (as a fall back strategy for vertex modes)
160  e.getByToken(token_beamSpot, bs);
161  if (!bs.isValid())
162  return result;
163 
164  // this is a default origin for all modes
165  GlobalPoint default_origin(bs->x0(), bs->y0(), bs->z0());
166 
167  // vector of origin & halfLength pairs:
168  std::vector<std::pair<GlobalPoint, float>> origins;
169 
170  // fill the origins and halfLengths depending on the mode
172  origins.push_back(std::make_pair(
173  default_origin, (m_mode == BEAM_SPOT_FIXED) ? m_zErrorBeamSpot : m_nSigmaZBeamSpot * bs->sigmaZ()));
174  } else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA) {
176  e.getByToken(token_vertex, vertices);
177  int n_vert = 0;
178  for (reco::VertexCollection::const_iterator iv = vertices->begin(), ev = vertices->end();
179  iv != ev && n_vert < m_maxNVertices;
180  ++iv) {
181  if (iv->isFake() || !iv->isValid())
182  continue;
183 
184  origins.push_back(std::make_pair(GlobalPoint(iv->x(), iv->y(), iv->z()),
185  (m_mode == VERTICES_FIXED) ? m_zErrorVetex : m_nSigmaZVertex * iv->zError()));
186  ++n_vert;
187  }
188  // no-vertex fall-back case:
189  if (origins.empty()) {
190  origins.push_back(std::make_pair(
191  default_origin, (m_nSigmaZBeamSpot > 0.) ? m_nSigmaZBeamSpot * bs->z0Error() : m_zErrorBeamSpot));
192  }
193  }
194 
199  measurementTracker = hmte.product();
200  }
201 
202  // create tracking regions (maximum MaxNRegions of them) in directions of the
203  // points of interest
204  size_t n_points = etaPoints.size();
205  int n_regions = 0;
206  for (size_t i = 0; i < n_points && n_regions < m_maxNRegions; ++i) {
207  double x = std::cos(phiPoints[i]);
208  double y = std::sin(phiPoints[i]);
209  double theta = 2 * std::atan(std::exp(-etaPoints[i]));
210  double z = 1. / std::tan(theta);
211 
212  GlobalVector direction(x, y, z);
213 
214  for (size_t j = 0; j < origins.size() && n_regions < m_maxNRegions; ++j) {
215  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(direction, // GlobalVector
216  origins[j].first, // GlobalPoint
217  m_ptMin,
219  origins[j].second,
220  m_deltaEta,
221  m_deltaPhi,
223  m_precise,
224  measurementTracker,
225  m_searchOpt));
226  ++n_regions;
227  }
228  }
229  edm::LogInfo("PointSeededTrackingRegionsProducer") << "produced " << n_regions << " regions";
230 
231  return result;
232  }
233 
234 private:
235  Mode m_mode;
236 
241 
242  std::vector<double> etaPoints;
243  std::vector<double> phiPoints;
244 
245  float m_ptMin;
248  float m_deltaEta;
249  float m_deltaPhi;
250  bool m_precise;
254 
258 };
259 
260 #endif
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:65
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
PointSeededTrackingRegionsProducer(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Geom::Theta< T > theta() const
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
U second(std::pair< T, U > const &p)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::VertexCollection > token_vertex
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
double z0Error() const
error on z
Definition: BeamSpot.h:90
T const * product() const
Definition: Handle.h:69
double sigmaZ() const
sigma z
Definition: BeamSpot.h:76
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double y0() const
y coordinate
Definition: BeamSpot.h:63
bool isUninitialized() const
Definition: EDGetToken.h:70
double x0() const
x coordinate
Definition: BeamSpot.h:61