CMS 3D CMS Logo

MkFitGeometryESProducer.cc
Go to the documentation of this file.
3 
11 
14 
15 // mkFit includes
21 
22 #include <sstream>
23 
24 //------------------------------------------------------------------------------
25 
27 public:
29 
30  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
31 
32  std::unique_ptr<MkFitGeometry> produce(const TrackerRecoGeometryRecord &iRecord);
33 
34 private:
35  struct GapCollector {
36  struct Interval {
37  float x, y;
38  };
39 
41  void extend_current(float q) {
44  }
46 
47  void add_interval(float x, float y);
48 
49  void sqrt_elements();
50  bool find_gap(Interval &itvl, float eps);
51  void print_gaps(std::ostream &ostr);
52 
53  std::list<Interval> m_coverage;
55  };
56  typedef std::unordered_map<int, GapCollector> layer_gap_map_t;
57 
58  void considerPoint(const GlobalPoint &gp, mkfit::LayerInfo &lay_info);
59  void fillShapeAndPlacement(const GeomDet *det, mkfit::TrackerInfo &trk_info, layer_gap_map_t *lgc_map = nullptr);
60  void addPixBGeometry(mkfit::TrackerInfo &trk_info);
61  void addPixEGeometry(mkfit::TrackerInfo &trk_info);
62  void addTIBGeometry(mkfit::TrackerInfo &trk_info);
63  void addTOBGeometry(mkfit::TrackerInfo &trk_info);
64  void addTIDGeometry(mkfit::TrackerInfo &trk_info);
65  void addTECGeometry(mkfit::TrackerInfo &trk_info);
66 
70 
71  const TrackerTopology *trackerTopo_ = nullptr;
72  const TrackerGeometry *trackerGeom_ = nullptr;
74 };
75 
77  auto cc = setWhatProduced(this);
78  geomToken_ = cc.consumes();
79  ttopoToken_ = cc.consumes();
80  trackerToken_ = cc.consumes();
81 }
82 
85  descriptions.addWithDefaultLabel(desc);
86 }
87 
88 //------------------------------------------------------------------------------
89 
91  if (x > y)
92  std::swap(x, y);
93  bool handled = false;
94  for (auto i = m_coverage.begin(); i != m_coverage.end(); ++i) {
95  if (y < i->x) { // fully on 'left'
96  m_coverage.insert(i, {x, y});
97  handled = true;
98  break;
99  } else if (x > i->y) { // fully on 'right'
100  continue;
101  } else if (x < i->x) { // sticking out on 'left'
102  i->x = x;
103  handled = true;
104  break;
105  } else if (y > i->y) { // sticking out on 'right'
106  i->y = y;
107  // check for overlap with the next interval, potentially merge
108  auto j = i;
109  ++j;
110  if (j != m_coverage.end() && i->y >= j->x) {
111  i->y = j->y;
112  m_coverage.erase(j);
113  }
114  handled = true;
115  break;
116  } else { // contained in current interval
117  handled = true;
118  break;
119  }
120  }
121  if (!handled) {
122  m_coverage.push_back({x, y});
123  }
124 }
125 
127  for (auto &itvl : m_coverage) {
128  itvl.x = std::sqrt(itvl.x);
129  itvl.y = std::sqrt(itvl.y);
130  }
131 }
132 
134  auto i = m_coverage.begin();
135  while (i != m_coverage.end()) {
136  auto j = i;
137  ++j;
138  if (j != m_coverage.end()) {
139  if (j->x - i->y > eps) {
140  itvl = {i->y, j->x};
141  return true;
142  }
143  i = j;
144  } else {
145  break;
146  }
147  }
148  return false;
149 }
150 
152  auto i = m_coverage.begin();
153  while (i != m_coverage.end()) {
154  auto j = i;
155  ++j;
156  if (j != m_coverage.end()) {
157  ostr << "(" << i->y << ", " << j->x << ")->" << j->x - i->y;
158  i = j;
159  } else {
160  break;
161  }
162  }
163 }
164 
165 //------------------------------------------------------------------------------
166 
168  // Use radius squared during bounding-region search.
169  float r = gp.perp2(), z = gp.z();
170  li.extend_limits(r, z);
171 }
172 
174  mkfit::TrackerInfo &trk_info,
175  layer_gap_map_t *lgc_map) {
176  DetId detid = det->geographicalId();
177 
178  float xy[4][2];
179  float dz;
180  const Bounds *b = &((det->surface()).bounds());
181 
182  if (const TrapezoidalPlaneBounds *b2 = dynamic_cast<const TrapezoidalPlaneBounds *>(b)) {
183  // See sec. "TrapezoidalPlaneBounds parameters" in doc/reco-geom-notes.txt
184  std::array<const float, 4> const &par = b2->parameters();
185  xy[0][0] = -par[0];
186  xy[0][1] = -par[3];
187  xy[1][0] = -par[1];
188  xy[1][1] = par[3];
189  xy[2][0] = par[1];
190  xy[2][1] = par[3];
191  xy[3][0] = par[0];
192  xy[3][1] = -par[3];
193  dz = par[2];
194 
195  // printf("TRAP 0x%x %f %f %f %f\n", detid.rawId(), par[0], par[1], par[2], par[3]);
196  } else if (const RectangularPlaneBounds *b2 = dynamic_cast<const RectangularPlaneBounds *>(b)) {
197  // Rectangular
198  float dx = b2->width() * 0.5; // half width
199  float dy = b2->length() * 0.5; // half length
200  xy[0][0] = -dx;
201  xy[0][1] = -dy;
202  xy[1][0] = -dx;
203  xy[1][1] = dy;
204  xy[2][0] = dx;
205  xy[2][1] = dy;
206  xy[3][0] = dx;
207  xy[3][1] = -dy;
208  dz = b2->thickness() * 0.5; // half thickness
209 
210  // printf("RECT 0x%x %f %f %f\n", detid.rawId(), dx, dy, dz);
211  } else {
212  throw cms::Exception("UnimplementedFeature") << "unsupported Bounds class";
213  }
214 
215  const bool useMatched = false;
216  int lay =
218  trackerTopo_->layer(detid),
219  useMatched,
220  trackerTopo_->isStereo(detid),
221  trackerTopo_->side(detid) == static_cast<unsigned>(TrackerDetSide::PosEndcap));
222 
223  mkfit::LayerInfo &layer_info = trk_info.layer_nc(lay);
224  if (lgc_map) {
225  (*lgc_map)[lay].reset_current();
226  }
227  for (int i = 0; i < 4; ++i) {
228  Local3DPoint lp1(xy[i][0], xy[i][1], -dz);
229  Local3DPoint lp2(xy[i][0], xy[i][1], dz);
230  GlobalPoint gp1 = det->surface().toGlobal(lp1);
231  GlobalPoint gp2 = det->surface().toGlobal(lp2);
232  considerPoint(gp1, layer_info);
233  considerPoint(gp2, layer_info);
234  if (lgc_map) {
235  (*lgc_map)[lay].extend_current(gp1.perp2());
236  (*lgc_map)[lay].extend_current(gp2.perp2());
237  }
238  }
239  if (lgc_map) {
240  (*lgc_map)[lay].add_current();
241  }
242  // Module information
243  const auto &p = det->position();
244  auto z = det->rotation().z();
245  auto x = det->rotation().x();
246  layer_info.register_module({{p.x(), p.y(), p.z()}, {z.x(), z.y(), z.z()}, {x.x(), x.y(), x.z()}, detid.rawId()});
247  // Set some layer parameters (repeatedly, would require hard-coding otherwise)
248  layer_info.set_subdet(detid.subdetId());
249  layer_info.set_is_pixel(detid.subdetId() <= 2);
250  layer_info.set_is_stereo(trackerTopo_->isStereo(detid));
251 }
252 
253 //==============================================================================
254 
255 // Ideally these functions would also:
256 // 0. Setup LayerInfo data (which is now done in auto-generated code).
257 // Some data-members are a bit over specific, esp/ bools for CMS sub-detectors.
258 // 1. Establish short module ids (now done in MkFitGeometry constructor).
259 // 2. Store module normal and strip direction vectors
260 // 3. ? Any other information ?
261 // 4. Extract stereo coverage holes where they exist (TEC, all but last 3 double-layers).
262 //
263 // Plugin DumpMkFitGeometry.cc can then be used to export this for stand-alone.
264 // Would also need to be picked up with tk-ntuple converter (to get module ids as
265 // they will now be used as indices into module info vectors).
266 //
267 // An attempt at export cmsRun config is in python/dumpMkFitGeometry.py
268 
270  for (auto &det : trackerGeom_->detsPXB()) {
271  fillShapeAndPlacement(det, trk_info);
272  }
273 }
274 
276  for (auto &det : trackerGeom_->detsPXF()) {
277  fillShapeAndPlacement(det, trk_info);
278  }
279 }
280 
282  for (auto &det : trackerGeom_->detsTIB()) {
283  fillShapeAndPlacement(det, trk_info);
284  }
285 }
286 
288  for (auto &det : trackerGeom_->detsTOB()) {
289  fillShapeAndPlacement(det, trk_info);
290  }
291 }
292 
294  for (auto &det : trackerGeom_->detsTID()) {
295  fillShapeAndPlacement(det, trk_info);
296  }
297 }
298 
300  // For TEC we also need to discover hole in radial extents.
301  layer_gap_map_t lgc_map;
302  for (auto &det : trackerGeom_->detsTEC()) {
303  fillShapeAndPlacement(det, trk_info, &lgc_map);
304  }
305  // Now loop over the GapCollectors and see if there is a coverage gap.
306  std::ostringstream ostr;
307  ostr << "addTECGeometry() gap report:\n";
309  for (auto &[layer, gcol] : lgc_map) {
310  gcol.sqrt_elements();
311  if (gcol.find_gap(itvl, 0.5)) {
312  ostr << " layer: " << layer << ", gap: " << itvl.x << " -> " << itvl.y << " width = " << itvl.y - itvl.x << "\n";
313  ostr << " all gaps: ";
314  gcol.print_gaps(ostr);
315  ostr << "\n";
316  trk_info.layer_nc(layer).set_r_hole_range(itvl.x, itvl.y);
317  }
318  }
319  edm::LogVerbatim("MkFitGeometryESProducer") << ostr.str();
320 }
321 
322 //------------------------------------------------------------------------------
323 // clang-format off
324 namespace {
325  const float phase1QBins[] = {
326  // PIXB, TIB, TOB
327  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,
328  // PIXE+, TID+, TEC+
329  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,
330  // PIXE-, TID-, TEC-
331  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
332  };
333 }
334 // clang-format on
335 //------------------------------------------------------------------------------
336 
337 std::unique_ptr<MkFitGeometry> MkFitGeometryESProducer::produce(const TrackerRecoGeometryRecord &iRecord) {
338  auto trackerInfo = std::make_unique<mkfit::TrackerInfo>();
339 
340  trackerGeom_ = &iRecord.get(geomToken_);
341  trackerTopo_ = &iRecord.get(ttopoToken_);
342 
343  // std::string path = "Geometry/TrackerCommonData/data/";
345  edm::LogInfo("MkFitGeometryESProducer") << "Extracting PhaseI geometry";
346  trackerInfo->create_layers(18, 27, 27);
349  throw cms::Exception("UnimplementedFeature") << "PhaseII geometry extraction";
350  } else {
351  throw cms::Exception("UnimplementedFeature") << "unsupported / unknowen geometry version";
352  }
353 
354  // Prepare layer boundaries for bounding-box search
355  for (int i = 0; i < trackerInfo->n_layers(); ++i) {
356  auto &li = trackerInfo->layer_nc(i);
357  li.set_limits(
359  li.reserve_modules(256);
360  }
361  // This is sort of CMS-2017 specific ... but fireworks code uses it for PhaseII as well.
362  addPixBGeometry(*trackerInfo);
363  addPixEGeometry(*trackerInfo);
364  addTIBGeometry(*trackerInfo);
365  addTIDGeometry(*trackerInfo);
366  addTOBGeometry(*trackerInfo);
367  addTECGeometry(*trackerInfo);
368 
369  // r_in/out kept as squres until here, root them
370  for (int i = 0; i < trackerInfo->n_layers(); ++i) {
371  auto &li = trackerInfo->layer_nc(i);
372  li.set_r_in_out(std::sqrt(li.rin()), std::sqrt(li.rout()));
373  li.set_propagate_to(li.is_barrel() ? li.r_mean() : li.z_mean());
374  li.set_q_bin(phase1QBins[i]);
375  unsigned int maxsid = li.shrink_modules();
376  // Make sure the short id fits in the 12 bits...
377  assert(maxsid < 1u << 11);
378  }
379 
380  return std::make_unique<MkFitGeometry>(
381  iRecord.get(geomToken_), iRecord.get(trackerToken_), iRecord.get(ttopoToken_), std::move(trackerInfo));
382 }
383 
void set_is_stereo(bool s)
Definition: TrackerInfo.h:54
const DetContainer & detsTIB() const
Log< level::Info, true > LogVerbatim
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
void extend_limits(float r, float z)
Definition: TrackerInfo.cc:15
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
bool find_gap(Interval &itvl, float eps)
void addPixBGeometry(mkfit::TrackerInfo &trk_info)
void set_r_hole_range(float rh1, float rh2)
Definition: TrackerInfo.cc:31
const DetContainer & detsPXB() const
void addTOBGeometry(mkfit::TrackerInfo &trk_info)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::ESGetToken< GeometricSearchTracker, TrackerRecoGeometryRecord > trackerToken_
const DetContainer & detsPXF() const
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:163
constexpr bool useMatched
unsigned int layer(const DetId &id) const
constexpr std::array< uint8_t, layerIndexSize > layer
const Surface::RotationType & rotation() const
The rotation defining the local R.F.
Definition: GeomDet.h:46
Basic3DVector< T > x() const
const DetContainer & detsTOB() const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
const TrackerGeometry * trackerGeom_
bool isThere(GeomDetEnumerators::SubDetector subdet) const
void addTIDGeometry(mkfit::TrackerInfo &trk_info)
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int register_module(ModuleInfo &&mi)
Definition: TrackerInfo.h:105
Basic3DVector< T > z() const
void considerPoint(const GlobalPoint &gp, mkfit::LayerInfo &lay_info)
void set_is_pixel(bool p)
Definition: TrackerInfo.h:53
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
Log< level::Info, false > LogInfo
void fillShapeAndPlacement(const GeomDet *det, mkfit::TrackerInfo &trk_info, layer_gap_map_t *lgc_map=nullptr)
Definition: DetId.h:17
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
double b
Definition: hdecay.h:118
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
T perp2() const
Definition: PV3DBase.h:68
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
void addPixEGeometry(mkfit::TrackerInfo &trk_info)
mkfit::LayerNumberConverter layerNrConv_
float x
std::unique_ptr< MkFitGeometry > produce(const TrackerRecoGeometryRecord &iRecord)
const DetContainer & detsTEC() const
Definition: Bounds.h:18
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
static constexpr float b2
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
int convertLayerNumber(int det, int lay, bool useMatched, int isStereo, bool posZ) const
def move(src, dest)
Definition: eostools.py:511
void addTIBGeometry(mkfit::TrackerInfo &trk_info)
const DetContainer & detsTID() const
void addTECGeometry(mkfit::TrackerInfo &trk_info)
std::unordered_map< int, GapCollector > layer_gap_map_t
void set_subdet(int sd)
Definition: TrackerInfo.h:52