CMS 3D CMS Logo

MkFitGeometryESProducer.cc
Go to the documentation of this file.
3 
11 
14 
16 
17 // mkFit includes
23 
24 #include <list>
25 #include <vector>
26 #include <unordered_map>
27 #include <sstream>
28 
29 // #define DUMP_MKF_GEO
30 
31 //------------------------------------------------------------------------------
32 
34 public:
36 
37  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
38 
39  std::unique_ptr<MkFitGeometry> produce(const TrackerRecoGeometryRecord &iRecord);
40 
41 private:
42  struct GapCollector {
43  struct Interval {
44  float x, y;
45  };
46 
48  void extend_current(float q) {
51  }
53 
54  void add_interval(float x, float y);
55 
56  void sqrt_elements();
57  bool find_gap(Interval &itvl, float eps);
58  void print_gaps(std::ostream &ostr);
59 
60  std::list<Interval> m_coverage;
62  };
63  using layer_gap_map_t = std::unordered_map<int, GapCollector>;
64 
66  std::size_t operator()(const mkfit::ModuleShape &s) const noexcept {
67  return std::hash<float>{}(s.dx1 + s.dx2 + s.dy + s.dz);
68  }
69  };
70  using module_shape_hmap_t = std::unordered_map<mkfit::ModuleShape, unsigned short, ModuleShape_hash>;
71  using layer_module_shape_vec_t = std::vector<module_shape_hmap_t>;
72 
73  struct MatHistBin {
74  float weight{0}, xi{0}, rl{0};
75  void add(float w, float x, float r) {
76  weight += w;
77  xi += w * x;
78  rl += w * r;
79  }
80  };
82 
83  void considerPoint(const GlobalPoint &gp, mkfit::LayerInfo &lay_info);
84  void fillShapeAndPlacement(const GeomDet *det,
85  mkfit::TrackerInfo &trk_info,
86  MaterialHistogram &material_histogram,
87  layer_gap_map_t *lgc_map = nullptr);
88  void addPixBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
89  void addPixEGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
90  void addTIBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
91  void addTOBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
92  void addTIDGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
93  void addTECGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
94 
95  void findRZBox(const GlobalPoint &gp, float &rmin, float &rmax, float &zmin, float &zmax);
96  void aggregateMaterialInfo(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
97  void fillLayers(mkfit::TrackerInfo &trk_info);
98 
102 
103  const TrackerTopology *trackerTopo_ = nullptr;
104  const TrackerGeometry *trackerGeom_ = nullptr;
107 };
108 
110  auto cc = setWhatProduced(this);
111  geomToken_ = cc.consumes();
112  ttopoToken_ = cc.consumes();
113  trackerToken_ = cc.consumes();
114 }
115 
118  descriptions.addWithDefaultLabel(desc);
119 }
120 
121 //------------------------------------------------------------------------------
122 
124  if (x > y)
125  std::swap(x, y);
126  bool handled = false;
127  for (auto i = m_coverage.begin(); i != m_coverage.end(); ++i) {
128  if (y < i->x) { // fully on 'left'
129  m_coverage.insert(i, {x, y});
130  handled = true;
131  break;
132  } else if (x > i->y) { // fully on 'right'
133  continue;
134  } else if (x < i->x) { // sticking out on 'left'
135  i->x = x;
136  handled = true;
137  break;
138  } else if (y > i->y) { // sticking out on 'right'
139  i->y = y;
140  // check for overlap with the next interval, potentially merge
141  auto j = i;
142  ++j;
143  if (j != m_coverage.end() && i->y >= j->x) {
144  i->y = j->y;
145  m_coverage.erase(j);
146  }
147  handled = true;
148  break;
149  } else { // contained in current interval
150  handled = true;
151  break;
152  }
153  }
154  if (!handled) {
155  m_coverage.push_back({x, y});
156  }
157 }
158 
160  for (auto &itvl : m_coverage) {
161  itvl.x = std::sqrt(itvl.x);
162  itvl.y = std::sqrt(itvl.y);
163  }
164 }
165 
167  auto i = m_coverage.begin();
168  while (i != m_coverage.end()) {
169  auto j = i;
170  ++j;
171  if (j != m_coverage.end()) {
172  if (j->x - i->y > eps) {
173  itvl = {i->y, j->x};
174  return true;
175  }
176  i = j;
177  } else {
178  break;
179  }
180  }
181  return false;
182 }
183 
185  auto i = m_coverage.begin();
186  while (i != m_coverage.end()) {
187  auto j = i;
188  ++j;
189  if (j != m_coverage.end()) {
190  ostr << "(" << i->y << ", " << j->x << ")->" << j->x - i->y;
191  i = j;
192  } else {
193  break;
194  }
195  }
196 }
197 
198 //------------------------------------------------------------------------------
199 
201  // Use radius squared during bounding-region search.
202  float r = gp.perp2(), z = gp.z();
203  li.extend_limits(r, z);
204 }
205 
207  mkfit::TrackerInfo &trk_info,
208  MaterialHistogram &material_histogram,
209  layer_gap_map_t *lgc_map) {
210  const DetId detid = det->geographicalId();
211 
212  bool doubleSide = false; //double modules have double material
213  if (detid.subdetId() == SiStripSubdetector::TIB)
214  doubleSide = trackerTopo_->tibIsDoubleSide(detid);
215  else if (detid.subdetId() == SiStripSubdetector::TID)
216  doubleSide = trackerTopo_->tidIsDoubleSide(detid);
217  else if (detid.subdetId() == SiStripSubdetector::TOB)
218  doubleSide = trackerTopo_->tobIsDoubleSide(detid);
219  else if (detid.subdetId() == SiStripSubdetector::TEC)
220  doubleSide = trackerTopo_->tecIsDoubleSide(detid);
221 
222  float xy[4][2];
223  float dz;
224  const Bounds *b = &((det->surface()).bounds());
226 
227  if (const TrapezoidalPlaneBounds *b2 = dynamic_cast<const TrapezoidalPlaneBounds *>(b)) {
228  // See sec. "TrapezoidalPlaneBounds parameters" in doc/reco-geom-notes.txt
229  std::array<const float, 4> const &par = b2->parameters();
230  xy[0][0] = -par[0];
231  xy[0][1] = -par[3];
232  xy[1][0] = -par[1];
233  xy[1][1] = par[3];
234  xy[2][0] = par[1];
235  xy[2][1] = par[3];
236  xy[3][0] = par[0];
237  xy[3][1] = -par[3];
238  dz = par[2];
239  ms.round_assign(par[0], par[1], par[3], par[2]);
240 
241 #ifdef DUMP_MKF_GEO
242  printf("TRAP 0x%x %f %f %f %f ", detid.rawId(), par[0], par[1], par[2], par[3]);
243 #endif
244  } else if (const RectangularPlaneBounds *b2 = dynamic_cast<const RectangularPlaneBounds *>(b)) {
245  // Rectangular
246  float dx = b2->width() * 0.5; // half width
247  float dy = b2->length() * 0.5; // half length
248  xy[0][0] = -dx;
249  xy[0][1] = -dy;
250  xy[1][0] = -dx;
251  xy[1][1] = dy;
252  xy[2][0] = dx;
253  xy[2][1] = dy;
254  xy[3][0] = dx;
255  xy[3][1] = -dy;
256  dz = b2->thickness() * 0.5; // half thickness
257  ms.round_assign(dx, 0.0f, dy, dz);
258 
259 #ifdef DUMP_MKF_GEO
260  printf("RECT 0x%x %f %f %f ", detid.rawId(), dx, dy, dz);
261 #endif
262  } else {
263  throw cms::Exception("UnimplementedFeature") << "unsupported Bounds class";
264  }
265 
266  const bool useMatched = false;
267  int lay =
270  useMatched,
272  trackerTopo_->side(detid) == static_cast<unsigned>(TrackerDetSide::PosEndcap));
273 
274  unsigned short shape_id = 9999;
275  if (!doubleSide) {
277  auto bhmi = bhm.find(ms);
278  if (bhmi == bhm.end()) {
279  bhmi = bhm.insert({ms, (unsigned short)bhm.size()}).first;
280  }
281  shape_id = bhmi->second;
282  }
283 
284 #ifdef DUMP_MKF_GEO
285  printf(" subdet=%d layer=%d side=%d is_stereo=%d is_double_side=%d --> mkflayer=%d; unique shape id=%hu\n",
286  detid.subdetId(),
290  doubleSide,
291  lay,
292  shape_id);
293 #endif
294 
295  mkfit::LayerInfo &layer_info = trk_info.layer_nc(lay);
296  if (lgc_map) {
297  (*lgc_map)[lay].reset_current();
298  }
299  float zbox_min = 1000, zbox_max = 0, rbox_min = 1000, rbox_max = 0;
300  for (int i = 0; i < 4; ++i) {
301  Local3DPoint lp1(xy[i][0], xy[i][1], -dz);
302  Local3DPoint lp2(xy[i][0], xy[i][1], dz);
303  GlobalPoint gp1 = det->surface().toGlobal(lp1);
304  GlobalPoint gp2 = det->surface().toGlobal(lp2);
305  considerPoint(gp1, layer_info);
306  considerPoint(gp2, layer_info);
307  findRZBox(gp1, rbox_min, rbox_max, zbox_min, zbox_max);
308  findRZBox(gp2, rbox_min, rbox_max, zbox_min, zbox_max);
309  if (lgc_map) {
310  (*lgc_map)[lay].extend_current(gp1.perp2());
311  (*lgc_map)[lay].extend_current(gp2.perp2());
312  }
313  }
314  if (lgc_map) {
315  (*lgc_map)[lay].add_current();
316  }
317 
318  // Double-sided module (join of two modules) information is not used in mkFit and
319  // also not needed for the material calculation.
320  // NOTE: This check should actually be performed even before the bounding box calculation
321  // as double-side module bounding box is a box-shaped union of bounding shapes of the
322  // constituent modules -- thus extending the layer bounding box unnecessarily.
323  // The 'if' is here as mkFit layer hit/miss logic has been tuned with this extended
324  // bounding boxes. To be moved higher up once this detection is improved through usage of
325  // the thin-thick layer concept (expected towards the end of 2024).
326  if (doubleSide)
327  return;
328 
329  // Module information
330  const auto &p = det->position();
331  auto z = det->rotation().z();
332  auto x = det->rotation().x();
333  layer_info.register_module(
334  {{p.x(), p.y(), p.z()}, {z.x(), z.y(), z.z()}, {x.x(), x.y(), x.z()}, detid.rawId(), shape_id});
335  // Set some layer parameters (repeatedly, would require hard-coding otherwise)
336  layer_info.set_subdet(detid.subdetId());
337  layer_info.set_is_pixel(detid.subdetId() <= 2);
338  layer_info.set_is_stereo(trackerTopo_->isStereo(detid));
339  if (layerNrConv_.isPhase2() && !layer_info.is_pixel())
340  layer_info.set_has_charge(false);
341 
342  // Fill material
343  {
344  // module material
345  const float bbxi = det->surface().mediumProperties().xi();
346  const float radL = det->surface().mediumProperties().radLen();
347  // loop over bins to fill histogram with bbxi, radL and their weight, which the overlap surface in r-z with the cmsquare of a bin
348  const float iBin = trk_info.mat_range_z() / trk_info.mat_nbins_z();
349  const float jBin = trk_info.mat_range_r() / trk_info.mat_nbins_r();
350  for (int i = std::floor(zbox_min / iBin); i < std::ceil(zbox_max / iBin); i++) {
351  for (int j = std::floor(rbox_min / jBin); j < std::ceil(rbox_max / jBin); j++) {
352  const float iF = i * iBin;
353  const float jF = j * jBin;
354  float overlap = std::max(0.f, std::min(jF + jBin, rbox_max) - std::max(jF, rbox_min)) *
355  std::max(0.f, std::min(iF + iBin, zbox_max) - std::max(iF, zbox_min));
356  if (overlap > 0)
357  material_histogram(i, j).add(overlap, bbxi, radL);
358  }
359  }
360  }
361 }
362 
363 //==============================================================================
364 
365 // These functions do the following:
366 // 0. Detect bounding cylinder of each layer.
367 // 1. Setup LayerInfo data.
368 // 2. Establish short module ids.
369 // 3. Store module normal and strip direction vectors.
370 // 4. Extract stereo coverage holes where they exist (TEC, all but last 3 double-layers).
371 //
372 // See python/dumpMkFitGeometry.py and dumpMkFitGeometryPhase2.py
373 
375 #ifdef DUMP_MKF_GEO
376  printf("\n*** addPixBGeometry\n\n");
377 #endif
378  for (auto &det : trackerGeom_->detsPXB()) {
379  fillShapeAndPlacement(det, trk_info, material_histogram);
380  }
381 }
382 
384 #ifdef DUMP_MKF_GEO
385  printf("\n*** addPixEGeometry\n\n");
386 #endif
387  for (auto &det : trackerGeom_->detsPXF()) {
388  fillShapeAndPlacement(det, trk_info, material_histogram);
389  }
390 }
391 
393 #ifdef DUMP_MKF_GEO
394  printf("\n*** addTIBGeometry\n\n");
395 #endif
396  for (auto &det : trackerGeom_->detsTIB()) {
397  fillShapeAndPlacement(det, trk_info, material_histogram);
398  }
399 }
400 
402 #ifdef DUMP_MKF_GEO
403  printf("\n*** addTOBGeometry\n\n");
404 #endif
405  for (auto &det : trackerGeom_->detsTOB()) {
406  fillShapeAndPlacement(det, trk_info, material_histogram);
407  }
408 }
409 
411 #ifdef DUMP_MKF_GEO
412  printf("\n*** addTIDGeometry\n\n");
413 #endif
414  for (auto &det : trackerGeom_->detsTID()) {
415  fillShapeAndPlacement(det, trk_info, material_histogram);
416  }
417 }
418 
420 #ifdef DUMP_MKF_GEO
421  printf("\n*** addTECGeometry\n\n");
422 #endif
423  // For TEC we also need to discover hole in radial extents.
424  layer_gap_map_t lgc_map;
425  for (auto &det : trackerGeom_->detsTEC()) {
426  fillShapeAndPlacement(det, trk_info, material_histogram, &lgc_map);
427  }
428  // Now loop over the GapCollectors and see if there is a coverage gap.
429  std::ostringstream ostr;
430  ostr << "addTECGeometry() gap report:\n";
432  for (auto &[layer, gcol] : lgc_map) {
433  gcol.sqrt_elements();
434  if (gcol.find_gap(itvl, 0.5)) {
435  ostr << " layer: " << layer << ", gap: " << itvl.x << " -> " << itvl.y << " width = " << itvl.y - itvl.x << "\n";
436  ostr << " all gaps: ";
437  gcol.print_gaps(ostr);
438  ostr << "\n";
439  trk_info.layer_nc(layer).set_r_hole_range(itvl.x, itvl.y);
440  }
441  }
442  edm::LogVerbatim("MkFitGeometryESProducer") << ostr.str();
443 }
444 
445 void MkFitGeometryESProducer::findRZBox(const GlobalPoint &gp, float &rmin, float &rmax, float &zmin, float &zmax) {
446  float r = gp.perp(), z = std::abs(gp.z());
447  rmax = std::max(r, rmax);
448  rmin = std::min(r, rmin);
449  zmax = std::max(z, zmax);
450  zmin = std::min(z, zmin);
451 }
452 
454  MaterialHistogram &material_histogram) {
455  //from histogram (vector of tuples) to grid
456  for (int i = 0; i < trk_info.mat_nbins_z(); i++) {
457  for (int j = 0; j < trk_info.mat_nbins_r(); j++) {
458  const MatHistBin &mhb = material_histogram(i, j);
459  if (mhb.weight > 0) {
460  trk_info.material_bbxi(i, j) = mhb.xi / mhb.weight;
461  trk_info.material_radl(i, j) = mhb.rl / mhb.weight;
462  }
463  }
464  }
465 }
466 
468  mkfit::rectvec<int> rneighbor_map(trk_info.mat_nbins_z(), trk_info.mat_nbins_r());
469  mkfit::rectvec<int> zneighbor_map(trk_info.mat_nbins_z(), trk_info.mat_nbins_r());
470 
471  for (int im = 0; im < trk_info.n_layers(); ++im) {
472  const mkfit::LayerInfo &li = trk_info.layer(im);
473  if (!li.is_barrel() && li.zmax() < 0)
474  continue; // neg endcap covered by pos
475  int rin, rout, zmin, zmax;
476  rin = trk_info.mat_bin_r(li.rin());
477  rout = trk_info.mat_bin_r(li.rout()) + 1;
478  if (li.is_barrel()) {
479  zmin = 0;
480  zmax = trk_info.mat_bin_z(std::max(std::abs(li.zmax()), std::abs(li.zmin()))) + 1;
481  } else {
482  zmin = trk_info.mat_bin_z(li.zmin());
483  zmax = trk_info.mat_bin_z(li.zmax()) + 1;
484  }
485  for (int i = zmin; i < zmax; i++) {
486  for (int j = rin; j < rout; j++) {
487  if (trk_info.material_bbxi(i, j) == 0) {
488  float distancesqmin = 100000;
489  for (int i2 = zmin; i2 < zmax; i2++) {
490  for (int j2 = rin; j2 < rout; j2++) {
491  if (j == j2 && i == i2)
492  continue;
493  auto mydistsq = (i - i2) * (i - i2) + (j - j2) * (j - j2);
494  if (mydistsq < distancesqmin && trk_info.material_radl(i2, j2) > 0) {
495  distancesqmin = mydistsq;
496  zneighbor_map(i, j) = i2;
497  rneighbor_map(i, j) = j2;
498  }
499  }
500  } // can work on speedup here
501  }
502  }
503  }
504  for (int i = zmin; i < zmax; i++) {
505  for (int j = rin; j < rout; j++) {
506  if (trk_info.material_bbxi(i, j) == 0) {
507  int iN = zneighbor_map(i, j);
508  int jN = rneighbor_map(i, j);
509  trk_info.material_bbxi(i, j) = trk_info.material_bbxi(iN, jN);
510  trk_info.material_radl(i, j) = trk_info.material_radl(iN, jN);
511  }
512  }
513  }
514  } //module loop
515 }
516 
517 //------------------------------------------------------------------------------
518 // clang-format off
519 namespace {
520  const float phase1QBins[] = {
521  // PIXB, TIB, TOB
522  2.0, 2.0, 2.0, 2.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0, 9.5, 9.5, 9.5, 9.5, 9.5, 9.5, 9.5, 9.5,
523  // PIXE+, TID+, TEC+
524  1.0, 1.0, 1.0, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5,
525  // PIXE-, TID-, TEC-
526  1.0, 1.0, 1.0, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5, 10.25, 7.5
527  };
528  const float phase2QBins[] = {
529  // TODO: Review these numbers.
530  // PIXB, TOB
531  2.0, 2.0, 2.0, 2.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0, 6.0,
532  // PIXE+, TEC+
533  1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6,
534  // PIXE-, TEC-
535  1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6, 5.6
536  };
537 }
538 // clang-format on
539 //------------------------------------------------------------------------------
540 
541 std::unique_ptr<MkFitGeometry> MkFitGeometryESProducer::produce(const TrackerRecoGeometryRecord &iRecord) {
542  auto trackerInfo = std::make_unique<mkfit::TrackerInfo>();
543 
544  trackerGeom_ = &iRecord.get(geomToken_);
545  trackerTopo_ = &iRecord.get(ttopoToken_);
546 
547  const float *qBinDefaults = nullptr;
548 
549  // std::string path = "Geometry/TrackerCommonData/data/";
551  edm::LogInfo("MkFitGeometryESProducer") << "Extracting PhaseI geometry";
552  trackerInfo->create_layers(18, 27, 27);
553  qBinDefaults = phase1QBins;
554  trackerInfo->create_material(300, 300.0f, 120, 120.0f);
557  edm::LogInfo("MkFitGeometryESProducer") << "Extracting PhaseII geometry";
558 #if !defined(MKFIT_PHASE2CUSTOMFLAGS)
559  // In Phase2, by default use prop-to-plane and pT-dependent MS.
560  // Option is kept to use custom flags, for R&D/test purposes.
561  using namespace mkfit;
562  Config::usePropToPlane = true;
563  Config::usePtMultScat = true;
564 #endif
566  trackerInfo->create_layers(16, 22, 22);
567  qBinDefaults = phase2QBins;
568  trackerInfo->create_material(300, 300.0f, 120, 120.0f);
569  } else {
570  throw cms::Exception("UnimplementedFeature") << "unsupported / unknowen geometry version";
571  }
572 
573  // Prepare layer boundaries for bounding-rz-box search
574  for (int i = 0; i < trackerInfo->n_layers(); ++i) {
575  auto &li = trackerInfo->layer_nc(i);
576  li.set_limits(
578  li.reserve_modules(256);
579  }
580  // Resize the module bounds collector vector
581  layerModuleShapeVec_.resize(trackerInfo->n_layers());
582 
583  MaterialHistogram material_histogram(trackerInfo->mat_nbins_z(), trackerInfo->mat_nbins_r());
584 
585  // This works for both Phase1 and Phase2.
586  // Phase2 TrackerGeometry returns empty det-vectors for TIB and TEC.
587  addPixBGeometry(*trackerInfo, material_histogram);
588  addPixEGeometry(*trackerInfo, material_histogram);
589  addTIBGeometry(*trackerInfo, material_histogram);
590  addTIDGeometry(*trackerInfo, material_histogram);
591  addTOBGeometry(*trackerInfo, material_histogram);
592  addTECGeometry(*trackerInfo, material_histogram);
593 
594  // r_in/out kept as squares until here, root them; fill ModuleShape vectors of layers
595  unsigned int n_mod = 0;
596  for (int i = 0; i < trackerInfo->n_layers(); ++i) {
597  auto &li = trackerInfo->layer_nc(i);
598  li.set_r_in_out(std::sqrt(li.rin()), std::sqrt(li.rout()));
599  li.set_propagate_to(li.is_barrel() ? li.r_mean() : li.z_mean());
600  li.set_q_bin(qBinDefaults[i]);
601  unsigned int maxsid = li.shrink_modules();
602 
603  n_mod += maxsid;
604 
605  // Make sure the short id fits in the 14 bits...
606  assert(maxsid < 1u << 13);
607  assert(n_mod > 0);
608 
609  // Fill ModuleShape vectors
610  int n_shapes = layerModuleShapeVec_[i].size();
611  li.resize_shapes(n_shapes);
612  for (auto &[shape, id] : layerModuleShapeVec_[i]) {
613  li.register_shape(shape, id);
614  }
615  }
616 
617  // Material grid
618  aggregateMaterialInfo(*trackerInfo, material_histogram);
619  fillLayers(*trackerInfo);
620 
621  // Propagation configuration
622  {
623  using namespace mkfit;
624  PropagationConfig &pconf = trackerInfo->prop_config_nc();
625  pconf.backward_fit_to_pca = false;
630  else
636  pconf.apply_tracker_info(trackerInfo.get());
637  }
638 
639 #ifdef DUMP_MKF_GEO
640  printf("Total number of modules %u, 14-bits fit up to %u modules\n", n_mod, 1u << 13);
641 #endif
642 
643  return std::make_unique<MkFitGeometry>(iRecord.get(geomToken_),
644  iRecord.get(trackerToken_),
645  iRecord.get(ttopoToken_),
646  std::move(trackerInfo),
647  layerNrConv_);
648 }
649 
const DetContainer & detsTIB() const
Log< level::Info, true > LogVerbatim
constexpr int32_t ceil(float num)
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:166
bool tibIsDoubleSide(const DetId &id) const
void extend_limits(float r, float z)
Definition: TrackerInfo.cc:33
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
bool tecIsDoubleSide(const DetId &id) const
bool tidIsDoubleSide(const DetId &id) const
bool find_gap(Interval &itvl, float eps)
float rin() const
Definition: TrackerInfo.h:95
void addTIDGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
void set_r_hole_range(float rh1, float rh2)
Definition: TrackerInfo.cc:49
layer_module_shape_vec_t layerModuleShapeVec_
void addTOBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
const DetContainer & detsPXB() const
void addTECGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
T w() const
PropagationFlags forward_fit_pflags
void fillLayers(mkfit::TrackerInfo &trk_info)
void addPixBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< module_shape_hmap_t > layer_module_shape_vec_t
edm::ESGetToken< GeometricSearchTracker, TrackerRecoGeometryRecord > trackerToken_
PropagationFlags backward_fit_pflags
const DetContainer & detsPXF() const
Definition: weight.py:1
bool isStereo(const DetId &id) const
const TrackerTopology * trackerTopo_
unsigned int side(const DetId &id) const
assert(be >=bs)
MkFitGeometryESProducer(const edm::ParameterSet &iConfig)
LayerInfo & layer_nc(int l)
Definition: TrackerInfo.h:242
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:112
bool tobIsDoubleSide(const DetId &id) const
constexpr bool useMatched
float zmax() const
Definition: TrackerInfo.h:99
unsigned int layer(const DetId &id) const
int mat_nbins_r() const
Definition: TrackerInfo.h:263
const Surface::RotationType & rotation() const
The rotation defining the local R.F.
Definition: GeomDet.h:46
PropagationFlags finding_inter_layer_pflags
Basic3DVector< T > x() const
std::unordered_map< mkfit::ModuleShape, unsigned short, ModuleShape_hash > module_shape_hmap_t
float zmin() const
Definition: TrackerInfo.h:98
float material_radl(int binZ, int binR) const
Definition: TrackerInfo.h:271
const DetContainer & detsTOB() const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
const TrackerGeometry * trackerGeom_
bool isThere(GeomDetEnumerators::SubDetector subdet) const
T sqrt(T t)
Definition: SSEVec.h:23
float radLen() const
int n_layers() const
Definition: TrackerInfo.h:240
int mat_bin_z(float z) const
Definition: TrackerInfo.h:266
float mat_range_z() const
Definition: TrackerInfo.h:264
Basic3DVector< T > z() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void considerPoint(const GlobalPoint &gp, mkfit::LayerInfo &lay_info)
PropagationFlags seed_fit_pflags
double f[11][100]
void fillShapeAndPlacement(const GeomDet *det, mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram, layer_gap_map_t *lgc_map=nullptr)
void addTIBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
std::unordered_map< int, GapCollector > layer_gap_map_t
PropagationFlags pca_prop_pflags
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:61
float mat_range_r() const
Definition: TrackerInfo.h:265
int mat_bin_r(float r) const
Definition: TrackerInfo.h:267
void add(float w, float x, float r)
bias2_t b2[25]
Definition: b2.h:9
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
const LayerInfo & layer(int l) const
Definition: TrackerInfo.h:241
bool usePtMultScat
Definition: Config.cc:8
Log< level::Info, false > LogInfo
Definition: DetId.h:17
void aggregateMaterialInfo(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
float material_bbxi(int binZ, int binR) const
Definition: TrackerInfo.h:270
double b
Definition: hdecay.h:120
float rout() const
Definition: TrackerInfo.h:96
void findRZBox(const GlobalPoint &gp, float &rmin, float &rmax, float &zmin, float &zmax)
T perp2() const
Definition: PV3DBase.h:68
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
mkfit::LayerNumberConverter layerNrConv_
void apply_tracker_info(const TrackerInfo *ti)
Definition: TrackerInfo.cc:13
std::size_t operator()(const mkfit::ModuleShape &s) const noexcept
void addPixEGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
float x
PropagationFlags finding_intra_layer_pflags
std::unique_ptr< MkFitGeometry > produce(const TrackerRecoGeometryRecord &iRecord)
const DetContainer & detsTEC() const
Definition: Bounds.h:18
const MediumProperties & mediumProperties() const
Definition: Surface.h:83
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
int convertLayerNumber(int det, int lay, bool useMatched, int isStereo, bool posZ) const
float xi() const
def move(src, dest)
Definition: eostools.py:511
int mat_nbins_z() const
Definition: TrackerInfo.h:262
bool is_barrel() const
Definition: TrackerInfo.h:105
bool usePropToPlane
Definition: Config.cc:7
const DetContainer & detsTID() const
void round_assign(float x1, float x2, float y, float z)
Definition: TrackerInfo.h:39