CMS 3D CMS Logo

MkFitGeometryESProducer.cc
Go to the documentation of this file.
3 
11 
14 
16 
17 // mkFit includes
23 
24 #include <sstream>
25 
26 // #define DUMP_MKF_GEO
27 
28 //------------------------------------------------------------------------------
29 
31 public:
33 
34  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
35 
36  std::unique_ptr<MkFitGeometry> produce(const TrackerRecoGeometryRecord &iRecord);
37 
38 private:
39  struct GapCollector {
40  struct Interval {
41  float x, y;
42  };
43 
45  void extend_current(float q) {
48  }
50 
51  void add_interval(float x, float y);
52 
53  void sqrt_elements();
54  bool find_gap(Interval &itvl, float eps);
55  void print_gaps(std::ostream &ostr);
56 
57  std::list<Interval> m_coverage;
59  };
60  typedef std::unordered_map<int, GapCollector> layer_gap_map_t;
61 
62  struct MatHistBin {
63  float weight{0}, xi{0}, rl{0};
64  void add(float w, float x, float r) {
65  weight += w;
66  xi += w * x;
67  rl += w * r;
68  }
69  };
71 
72  void considerPoint(const GlobalPoint &gp, mkfit::LayerInfo &lay_info);
73  void fillShapeAndPlacement(const GeomDet *det,
74  mkfit::TrackerInfo &trk_info,
75  MaterialHistogram &material_histogram,
76  layer_gap_map_t *lgc_map = nullptr);
77  void addPixBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
78  void addPixEGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
79  void addTIBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
80  void addTOBGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
81  void addTIDGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
82  void addTECGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
83 
84  void findRZBox(const GlobalPoint &gp, float &rmin, float &rmax, float &zmin, float &zmax);
85  void aggregateMaterialInfo(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram);
86  void fillLayers(mkfit::TrackerInfo &trk_info);
87 
91 
92  const TrackerTopology *trackerTopo_ = nullptr;
93  const TrackerGeometry *trackerGeom_ = nullptr;
95 };
96 
98  auto cc = setWhatProduced(this);
99  geomToken_ = cc.consumes();
100  ttopoToken_ = cc.consumes();
101  trackerToken_ = cc.consumes();
102 }
103 
106  descriptions.addWithDefaultLabel(desc);
107 }
108 
109 //------------------------------------------------------------------------------
110 
112  if (x > y)
113  std::swap(x, y);
114  bool handled = false;
115  for (auto i = m_coverage.begin(); i != m_coverage.end(); ++i) {
116  if (y < i->x) { // fully on 'left'
117  m_coverage.insert(i, {x, y});
118  handled = true;
119  break;
120  } else if (x > i->y) { // fully on 'right'
121  continue;
122  } else if (x < i->x) { // sticking out on 'left'
123  i->x = x;
124  handled = true;
125  break;
126  } else if (y > i->y) { // sticking out on 'right'
127  i->y = y;
128  // check for overlap with the next interval, potentially merge
129  auto j = i;
130  ++j;
131  if (j != m_coverage.end() && i->y >= j->x) {
132  i->y = j->y;
133  m_coverage.erase(j);
134  }
135  handled = true;
136  break;
137  } else { // contained in current interval
138  handled = true;
139  break;
140  }
141  }
142  if (!handled) {
143  m_coverage.push_back({x, y});
144  }
145 }
146 
148  for (auto &itvl : m_coverage) {
149  itvl.x = std::sqrt(itvl.x);
150  itvl.y = std::sqrt(itvl.y);
151  }
152 }
153 
155  auto i = m_coverage.begin();
156  while (i != m_coverage.end()) {
157  auto j = i;
158  ++j;
159  if (j != m_coverage.end()) {
160  if (j->x - i->y > eps) {
161  itvl = {i->y, j->x};
162  return true;
163  }
164  i = j;
165  } else {
166  break;
167  }
168  }
169  return false;
170 }
171 
173  auto i = m_coverage.begin();
174  while (i != m_coverage.end()) {
175  auto j = i;
176  ++j;
177  if (j != m_coverage.end()) {
178  ostr << "(" << i->y << ", " << j->x << ")->" << j->x - i->y;
179  i = j;
180  } else {
181  break;
182  }
183  }
184 }
185 
186 //------------------------------------------------------------------------------
187 
189  // Use radius squared during bounding-region search.
190  float r = gp.perp2(), z = gp.z();
191  li.extend_limits(r, z);
192 }
193 
195  mkfit::TrackerInfo &trk_info,
196  MaterialHistogram &material_histogram,
197  layer_gap_map_t *lgc_map) {
198  DetId detid = det->geographicalId();
199 
200  float xy[4][2];
201  float half_length, dz;
202  const Bounds *b = &((det->surface()).bounds());
203 
204  if (const TrapezoidalPlaneBounds *b2 = dynamic_cast<const TrapezoidalPlaneBounds *>(b)) {
205  // See sec. "TrapezoidalPlaneBounds parameters" in doc/reco-geom-notes.txt
206  std::array<const float, 4> const &par = b2->parameters();
207  xy[0][0] = -par[0];
208  xy[0][1] = -par[3];
209  xy[1][0] = -par[1];
210  xy[1][1] = par[3];
211  xy[2][0] = par[1];
212  xy[2][1] = par[3];
213  xy[3][0] = par[0];
214  xy[3][1] = -par[3];
215  half_length = par[3];
216  dz = par[2];
217 
218 #ifdef DUMP_MKF_GEO
219  printf("TRAP 0x%x %f %f %f %f ", detid.rawId(), par[0], par[1], par[2], par[3]);
220 #endif
221  } else if (const RectangularPlaneBounds *b2 = dynamic_cast<const RectangularPlaneBounds *>(b)) {
222  // Rectangular
223  float dx = b2->width() * 0.5; // half width
224  float dy = b2->length() * 0.5; // half length
225  xy[0][0] = -dx;
226  xy[0][1] = -dy;
227  xy[1][0] = -dx;
228  xy[1][1] = dy;
229  xy[2][0] = dx;
230  xy[2][1] = dy;
231  xy[3][0] = dx;
232  xy[3][1] = -dy;
233  half_length = dy;
234  dz = b2->thickness() * 0.5; // half thickness
235 
236 #ifdef DUMP_MKF_GEO
237  printf("RECT 0x%x %f %f %f ", detid.rawId(), dx, dy, dz);
238 #endif
239  } else {
240  throw cms::Exception("UnimplementedFeature") << "unsupported Bounds class";
241  }
242 
243  const bool useMatched = false;
244  int lay =
246  trackerTopo_->layer(detid),
247  useMatched,
248  trackerTopo_->isStereo(detid),
249  trackerTopo_->side(detid) == static_cast<unsigned>(TrackerDetSide::PosEndcap));
250 #ifdef DUMP_MKF_GEO
251  printf(" subdet=%d layer=%d side=%d is_stereo=%d --> mkflayer=%d\n",
252  detid.subdetId(),
253  trackerTopo_->layer(detid),
254  trackerTopo_->side(detid),
255  trackerTopo_->isStereo(detid),
256  lay);
257 #endif
258 
259  mkfit::LayerInfo &layer_info = trk_info.layer_nc(lay);
260  if (lgc_map) {
261  (*lgc_map)[lay].reset_current();
262  }
263  float zbox_min = 1000, zbox_max = 0, rbox_min = 1000, rbox_max = 0;
264  for (int i = 0; i < 4; ++i) {
265  Local3DPoint lp1(xy[i][0], xy[i][1], -dz);
266  Local3DPoint lp2(xy[i][0], xy[i][1], dz);
267  GlobalPoint gp1 = det->surface().toGlobal(lp1);
268  GlobalPoint gp2 = det->surface().toGlobal(lp2);
269  considerPoint(gp1, layer_info);
270  considerPoint(gp2, layer_info);
271  findRZBox(gp1, rbox_min, rbox_max, zbox_min, zbox_max);
272  findRZBox(gp2, rbox_min, rbox_max, zbox_min, zbox_max);
273  if (lgc_map) {
274  (*lgc_map)[lay].extend_current(gp1.perp2());
275  (*lgc_map)[lay].extend_current(gp2.perp2());
276  }
277  }
278  if (lgc_map) {
279  (*lgc_map)[lay].add_current();
280  }
281  // Module information
282  const auto &p = det->position();
283  auto z = det->rotation().z();
284  auto x = det->rotation().x();
285  layer_info.register_module(
286  {{p.x(), p.y(), p.z()}, {z.x(), z.y(), z.z()}, {x.x(), x.y(), x.z()}, half_length, detid.rawId()});
287  // Set some layer parameters (repeatedly, would require hard-coding otherwise)
288  layer_info.set_subdet(detid.subdetId());
289  layer_info.set_is_pixel(detid.subdetId() <= 2);
290  layer_info.set_is_stereo(trackerTopo_->isStereo(detid));
291 
292  bool doubleSide = false; //double modules have double material
293  if (detid.subdetId() == SiStripSubdetector::TIB)
294  doubleSide = trackerTopo_->tibIsDoubleSide(detid);
295  else if (detid.subdetId() == SiStripSubdetector::TID)
296  doubleSide = trackerTopo_->tidIsDoubleSide(detid);
297  else if (detid.subdetId() == SiStripSubdetector::TOB)
298  doubleSide = trackerTopo_->tobIsDoubleSide(detid);
299  else if (detid.subdetId() == SiStripSubdetector::TEC)
300  doubleSide = trackerTopo_->tecIsDoubleSide(detid);
301 
302  if (!doubleSide) //fill material
303  {
304  //module material
305  const float bbxi = det->surface().mediumProperties().xi();
306  const float radL = det->surface().mediumProperties().radLen();
307  //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
308  const float iBin = trk_info.mat_range_z() / trk_info.mat_nbins_z();
309  const float jBin = trk_info.mat_range_r() / trk_info.mat_nbins_r();
310  for (int i = std::floor(zbox_min / iBin); i < std::ceil(zbox_max / iBin); i++) {
311  for (int j = std::floor(rbox_min / jBin); j < std::ceil(rbox_max / jBin); j++) {
312  const float iF = i * iBin;
313  const float jF = j * jBin;
314  float overlap = std::max(0.f, std::min(jF + jBin, rbox_max) - std::max(jF, rbox_min)) *
315  std::max(0.f, std::min(iF + iBin, zbox_max) - std::max(iF, zbox_min));
316  if (overlap > 0)
317  material_histogram(i, j).add(overlap, bbxi, radL);
318  }
319  }
320  }
321 }
322 
323 //==============================================================================
324 
325 // These functions do the following:
326 // 0. Detect bounding cylinder of each layer.
327 // 1. Setup LayerInfo data.
328 // 2. Establish short module ids.
329 // 3. Store module normal and strip direction vectors.
330 // 4. Extract stereo coverage holes where they exist (TEC, all but last 3 double-layers).
331 //
332 // See python/dumpMkFitGeometry.py and dumpMkFitGeometryPhase2.py
333 
335 #ifdef DUMP_MKF_GEO
336  printf("\n*** addPixBGeometry\n\n");
337 #endif
338  for (auto &det : trackerGeom_->detsPXB()) {
339  fillShapeAndPlacement(det, trk_info, material_histogram);
340  }
341 }
342 
344 #ifdef DUMP_MKF_GEO
345  printf("\n*** addPixEGeometry\n\n");
346 #endif
347  for (auto &det : trackerGeom_->detsPXF()) {
348  fillShapeAndPlacement(det, trk_info, material_histogram);
349  }
350 }
351 
353 #ifdef DUMP_MKF_GEO
354  printf("\n*** addTIBGeometry\n\n");
355 #endif
356  for (auto &det : trackerGeom_->detsTIB()) {
357  fillShapeAndPlacement(det, trk_info, material_histogram);
358  }
359 }
360 
362 #ifdef DUMP_MKF_GEO
363  printf("\n*** addTOBGeometry\n\n");
364 #endif
365  for (auto &det : trackerGeom_->detsTOB()) {
366  fillShapeAndPlacement(det, trk_info, material_histogram);
367  }
368 }
369 
371 #ifdef DUMP_MKF_GEO
372  printf("\n*** addTIDGeometry\n\n");
373 #endif
374  for (auto &det : trackerGeom_->detsTID()) {
375  fillShapeAndPlacement(det, trk_info, material_histogram);
376  }
377 }
378 
380 #ifdef DUMP_MKF_GEO
381  printf("\n*** addTECGeometry\n\n");
382 #endif
383  // For TEC we also need to discover hole in radial extents.
384  layer_gap_map_t lgc_map;
385  for (auto &det : trackerGeom_->detsTEC()) {
386  fillShapeAndPlacement(det, trk_info, material_histogram, &lgc_map);
387  }
388  // Now loop over the GapCollectors and see if there is a coverage gap.
389  std::ostringstream ostr;
390  ostr << "addTECGeometry() gap report:\n";
392  for (auto &[layer, gcol] : lgc_map) {
393  gcol.sqrt_elements();
394  if (gcol.find_gap(itvl, 0.5)) {
395  ostr << " layer: " << layer << ", gap: " << itvl.x << " -> " << itvl.y << " width = " << itvl.y - itvl.x << "\n";
396  ostr << " all gaps: ";
397  gcol.print_gaps(ostr);
398  ostr << "\n";
399  trk_info.layer_nc(layer).set_r_hole_range(itvl.x, itvl.y);
400  }
401  }
402  edm::LogVerbatim("MkFitGeometryESProducer") << ostr.str();
403 }
404 
405 void MkFitGeometryESProducer::findRZBox(const GlobalPoint &gp, float &rmin, float &rmax, float &zmin, float &zmax) {
406  float r = gp.perp(), z = std::abs(gp.z());
407  rmax = std::max(r, rmax);
408  rmin = std::min(r, rmin);
409  zmax = std::max(z, zmax);
410  zmin = std::min(z, zmin);
411 }
412 
414  MaterialHistogram &material_histogram) {
415  //from histogram (vector of tuples) to grid
416  for (int i = 0; i < trk_info.mat_nbins_z(); i++) {
417  for (int j = 0; j < trk_info.mat_nbins_r(); j++) {
418  const MatHistBin &mhb = material_histogram(i, j);
419  if (mhb.weight > 0) {
420  trk_info.material_bbxi(i, j) = mhb.xi / mhb.weight;
421  trk_info.material_radl(i, j) = mhb.rl / mhb.weight;
422  }
423  }
424  }
425 }
426 
428  mkfit::rectvec<int> rneighbor_map(trk_info.mat_nbins_z(), trk_info.mat_nbins_r());
429  mkfit::rectvec<int> zneighbor_map(trk_info.mat_nbins_z(), trk_info.mat_nbins_r());
430 
431  for (int im = 0; im < trk_info.n_layers(); ++im) {
432  const mkfit::LayerInfo &li = trk_info.layer(im);
433  if (!li.is_barrel() && li.zmax() < 0)
434  continue; // neg endcap covered by pos
435  int rin, rout, zmin, zmax;
436  rin = trk_info.mat_bin_r(li.rin());
437  rout = trk_info.mat_bin_r(li.rout()) + 1;
438  if (li.is_barrel()) {
439  zmin = 0;
440  zmax = trk_info.mat_bin_z(std::max(std::abs(li.zmax()), std::abs(li.zmin()))) + 1;
441  } else {
442  zmin = trk_info.mat_bin_z(li.zmin());
443  zmax = trk_info.mat_bin_z(li.zmax()) + 1;
444  }
445  for (int i = zmin; i < zmax; i++) {
446  for (int j = rin; j < rout; j++) {
447  if (trk_info.material_bbxi(i, j) == 0) {
448  float distancesqmin = 100000;
449  for (int i2 = zmin; i2 < zmax; i2++) {
450  for (int j2 = rin; j2 < rout; j2++) {
451  if (j == j2 && i == i2)
452  continue;
453  auto mydistsq = (i - i2) * (i - i2) + (j - j2) * (j - j2);
454  if (mydistsq < distancesqmin && trk_info.material_radl(i2, j2) > 0) {
455  distancesqmin = mydistsq;
456  zneighbor_map(i, j) = i2;
457  rneighbor_map(i, j) = j2;
458  }
459  }
460  } // can work on speedup here
461  }
462  }
463  }
464  for (int i = zmin; i < zmax; i++) {
465  for (int j = rin; j < rout; j++) {
466  if (trk_info.material_bbxi(i, j) == 0) {
467  int iN = zneighbor_map(i, j);
468  int jN = rneighbor_map(i, j);
469  trk_info.material_bbxi(i, j) = trk_info.material_bbxi(iN, jN);
470  trk_info.material_radl(i, j) = trk_info.material_radl(iN, jN);
471  }
472  }
473  }
474  } //module loop
475 }
476 
477 //------------------------------------------------------------------------------
478 // clang-format off
479 namespace {
480  const float phase1QBins[] = {
481  // PIXB, TIB, TOB
482  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,
483  // PIXE+, TID+, TEC+
484  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,
485  // PIXE-, TID-, TEC-
486  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
487  };
488  const float phase2QBins[] = {
489  // TODO: Review these numbers.
490  // PIXB, TOB
491  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,
492  // PIXE+, TEC+
493  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,
494  // PIXE-, TEC-
495  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
496  };
497 }
498 // clang-format on
499 //------------------------------------------------------------------------------
500 
501 std::unique_ptr<MkFitGeometry> MkFitGeometryESProducer::produce(const TrackerRecoGeometryRecord &iRecord) {
502  auto trackerInfo = std::make_unique<mkfit::TrackerInfo>();
503 
504  trackerGeom_ = &iRecord.get(geomToken_);
505  trackerTopo_ = &iRecord.get(ttopoToken_);
506 
507  const float *qBinDefaults = nullptr;
508 
509  // std::string path = "Geometry/TrackerCommonData/data/";
511  edm::LogInfo("MkFitGeometryESProducer") << "Extracting PhaseI geometry";
512  trackerInfo->create_layers(18, 27, 27);
513  qBinDefaults = phase1QBins;
514 
515  trackerInfo->create_material(300, 300.0f, 120, 120.0f);
518  edm::LogInfo("MkFitGeometryESProducer") << "Extracting PhaseII geometry";
520  trackerInfo->create_layers(16, 22, 22);
521  qBinDefaults = phase2QBins;
522  trackerInfo->create_material(300, 300.0f, 120, 120.0f);
523  } else {
524  throw cms::Exception("UnimplementedFeature") << "unsupported / unknowen geometry version";
525  }
526 
527  // Prepare layer boundaries for bounding-box search
528  for (int i = 0; i < trackerInfo->n_layers(); ++i) {
529  auto &li = trackerInfo->layer_nc(i);
530  li.set_limits(
532  li.reserve_modules(256);
533  }
534 
535  MaterialHistogram material_histogram(trackerInfo->mat_nbins_z(), trackerInfo->mat_nbins_r());
536 
537  // This is sort of CMS-phase1 specific ... but fireworks code uses it for PhaseII as well.
538  addPixBGeometry(*trackerInfo, material_histogram);
539  addPixEGeometry(*trackerInfo, material_histogram);
540  addTIBGeometry(*trackerInfo, material_histogram);
541  addTIDGeometry(*trackerInfo, material_histogram);
542  addTOBGeometry(*trackerInfo, material_histogram);
543  addTECGeometry(*trackerInfo, material_histogram);
544 
545  // r_in/out kept as squares until here, root them
546  unsigned int n_mod = 0;
547  for (int i = 0; i < trackerInfo->n_layers(); ++i) {
548  auto &li = trackerInfo->layer_nc(i);
549  li.set_r_in_out(std::sqrt(li.rin()), std::sqrt(li.rout()));
550  li.set_propagate_to(li.is_barrel() ? li.r_mean() : li.z_mean());
551  li.set_q_bin(qBinDefaults[i]);
552  unsigned int maxsid = li.shrink_modules();
553 
554  n_mod += maxsid;
555 
556  // Make sure the short id fits in the 14 bits...
557  assert(maxsid < 1u << 13);
558  assert(n_mod > 0);
559  }
560 
561  // Material grid
562  aggregateMaterialInfo(*trackerInfo, material_histogram);
563  fillLayers(*trackerInfo);
564 
565  // Propagation configuration
566  {
567  using namespace mkfit;
568  PropagationConfig &pconf = trackerInfo->prop_config_nc();
569  pconf.backward_fit_to_pca = false;
574  else
580  pconf.apply_tracker_info(trackerInfo.get());
581  }
582 
583 #ifdef DUMP_MKF_GEO
584  printf("Total number of modules %u, 14-bits fit up to %u modules\n", n_mod, 1u << 13);
585 #endif
586 
587  return std::make_unique<MkFitGeometry>(iRecord.get(geomToken_),
588  iRecord.get(trackerToken_),
589  iRecord.get(ttopoToken_),
590  std::move(trackerInfo),
591  layerNrConv_);
592 }
593 
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:67
void addTIDGeometry(mkfit::TrackerInfo &trk_info, MaterialHistogram &material_histogram)
void set_r_hole_range(float rh1, float rh2)
Definition: TrackerInfo.cc:49
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)
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:204
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:71
unsigned int layer(const DetId &id) const
int mat_nbins_r() const
Definition: TrackerInfo.h:225
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
float zmin() const
Definition: TrackerInfo.h:70
float material_radl(int binZ, int binR) const
Definition: TrackerInfo.h:233
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:19
float radLen() const
int n_layers() const
Definition: TrackerInfo.h:202
int mat_bin_z(float z) const
Definition: TrackerInfo.h:228
float mat_range_z() const
Definition: TrackerInfo.h:226
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]
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
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)
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:227
int mat_bin_r(float r) const
Definition: TrackerInfo.h:229
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:203
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
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
float material_bbxi(int binZ, int binR) const
Definition: TrackerInfo.h:232
double b
Definition: hdecay.h:120
float rout() const
Definition: TrackerInfo.h:68
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
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
constexpr bool usePropToPlane
Definition: Config.h:51
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:224
bool is_barrel() const
Definition: TrackerInfo.h:77
const DetContainer & detsTID() const
std::unordered_map< int, GapCollector > layer_gap_map_t