CMS 3D CMS Logo

HGCalGeometry.cc
Go to the documentation of this file.
1 /* for High Granularity Calorimeter
2  * This geometry is essentially driven by topology,
3  * which is thus encapsulated in this class.
4  * This makes this geometry not suitable to be loaded
5  * by regular CaloGeometryLoader<T>
6  */
15 
16 #include <cmath>
17 
18 #include <Math/Transform3D.h>
19 #include <Math/EulerAngles.h>
20 
22 typedef std::vector<float> ParmVec;
23 
24 //#define EDM_ML_DEBUG
25 
27  : m_topology( topology_ ),
28  m_validGeomIds( topology_.totalGeomModules()),
29  mode_( topology_.geomMode()),
30  m_det( topology_.detector()),
31  m_subdet( topology_.subDetector()),
32  twoBysqrt3_(2.0/std::sqrt(3.0)) {
33 
34  if (m_det == DetId::HGCalHSc) {
35  m_cellVec2 = CellVec2(topology_.totalGeomModules());
36  } else {
37  m_cellVec = CellVec(topology_.totalGeomModules());
38  }
40 #ifdef EDM_ML_DEBUG
41  edm::LogVerbatim("HGCalGeom") << "Expected total # of Geometry Modules "
43 #endif
44 }
45 
47 
49 
51 }
52 
54  const CCGFloat* pv,
55  unsigned int i,
56  Pt3D& ref) {
57  if (m_det == DetId::HGCalHSc) {
58  FlatTrd::localCorners( lc, pv, ref ) ;
59  } else {
60  FlatHexagon::localCorners( lc, pv, ref ) ;
61  }
62 }
63 
65  const GlobalPoint& f2 ,
66  const GlobalPoint& f3 ,
67  const CCGFloat* parm ,
68  const DetId& detId ) {
69 
70  DetId geomId = getGeometryDetId(detId);
71  int cells (0);
75  cells = m_topology.dddConstants().numberCellsHexagon(id.iSec1);
76 #ifdef EDM_ML_DEBUG
77  edm::LogVerbatim("HGCalGeom") << "NewCell " << HGCalDetId(detId)
78  << " GEOM " << HGCalDetId(geomId);
79 #endif
80  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
81  cells = 1;
82 #ifdef EDM_ML_DEBUG
83  edm::LogVerbatim("HGCalGeom") << "NewCell " << HGCScintillatorDetId(detId)
84  << " GEOM " << HGCScintillatorDetId(geomId);
85 #endif
86  } else {
87  cells = m_topology.dddConstants().numberCellsHexagon(id.iLay,id.iSec1,
88  id.iSec2,false);
89 #ifdef EDM_ML_DEBUG
90  edm::LogVerbatim("HGCalGeom") << "NewCell " << HGCSiliconDetId(detId)
91  << " GEOM " << HGCSiliconDetId(geomId);
92 #endif
93  }
94  const uint32_t cellIndex (m_topology.detId2denseGeomId(geomId));
95 
96  if (m_det == DetId::HGCalHSc) {
97  m_cellVec2.at( cellIndex ) = FlatTrd( cornersMgr(), f1, f2, f3, parm ) ;
98  } else {
99  m_cellVec.at( cellIndex ) = FlatHexagon( cornersMgr(), f1, f2, f3, parm ) ;
100  }
101  m_validGeomIds.at( cellIndex ) = geomId ;
102 
103 #ifdef EDM_ML_DEBUG
104  edm::LogVerbatim("HGCalGeom") << "Store for DetId " << std::hex
105  << detId.rawId() << " GeomId "
106  << geomId.rawId() << std::dec << " Index "
107  << cellIndex << " cells " << cells;
108  unsigned int nOld = m_validIds.size();
109 #endif
110  if ((mode_ == HGCalGeometryMode::Hexagon) ||
112  for (int cell = 0; cell < cells; ++cell) {
113  id.iCell1 = cell;
114  DetId idc = m_topology.encode(id);
115  if (m_topology.valid(idc)) {
116  m_validIds.emplace_back(idc);
117 #ifdef EDM_ML_DEBUG
118  edm::LogVerbatim("HGCalGeom") << "Valid Id [" << cell << "] "
119  << HGCalDetId(idc);
120 #endif
121  }
122  }
123  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
124  DetId idc = m_topology.encode(id);
125  if (m_topology.valid(idc)) {
126  m_validIds.emplace_back(idc);
127 #ifdef EDM_ML_DEBUG
128  edm::LogVerbatim("HGCalGeom") << "Valid Id [0] "
129  << HGCScintillatorDetId(idc);
130 #endif
131  } else {
132  edm::LogWarning("HGCalGeom") << "Check " << HGCScintillatorDetId(idc)
133  << " from " << HGCScintillatorDetId(detId)
134  << " ERROR ???";
135  }
136  } else {
137 #ifdef EDM_ML_DEBUG
138  unsigned int cellAll(0), cellSelect(0);
139 #endif
140  for (int u=0; u<2*cells; ++u) {
141  for (int v=0; v<2*cells; ++v) {
142  if (((v-u) < cells) && (u-v) <= cells) {
143  id.iCell1 = u; id.iCell2 = v;
144  DetId idc = m_topology.encode(id);
145 #ifdef EDM_ML_DEBUG
146  ++cellAll;
147 #endif
148  if (m_topology.dddConstants().cellInLayer(id.iSec1,id.iSec2,u,v,
149  id.iLay,true)) {
150  m_validIds.emplace_back(idc);
151 #ifdef EDM_ML_DEBUG
152  ++cellSelect;
153  edm::LogVerbatim("HGCalGeom") << "Valid Id [" << u << ", " << v
154  << "] " << HGCSiliconDetId(idc);
155 #endif
156  }
157  }
158  }
159  }
160 #ifdef EDM_ML_DEBUG
161  edm::LogVerbatim("HGCalGeom") << "HGCalGeometry keeps " << cellSelect
162  << " out of " << cellAll << " for wafer "
163  << id.iSec1 << ":" << id.iSec2 << " in "
164  << " layer " << id.iLay;
165 #endif
166  }
167 #ifdef EDM_ML_DEBUG
168  if (m_det == DetId::HGCalHSc) {
169  edm::LogVerbatim("HGCalGeom") << "HGCalGeometry::newCell-> [" << cellIndex
170  << "]" << " front:" << f1.x() << '/'
171  << f1.y() << '/' << f1.z() << " back:"
172  << f2.x() << '/' << f2.y() << '/' << f2.z()
173  << " eta|phi "
174  << m_cellVec2[cellIndex].etaPos() << ":"
175  << m_cellVec2[cellIndex].phiPos();
176  } else {
177  edm::LogVerbatim("HGCalGeom") << "HGCalGeometry::newCell-> [" << cellIndex
178  << "]" << " front:" << f1.x() << '/'
179  << f1.y() << '/' << f1.z() << " back:"
180  << f2.x() << '/' << f2.y() << '/' << f2.z()
181  << " eta|phi "
182  << m_cellVec[cellIndex].etaPos() << ":"
183  << m_cellVec[cellIndex].phiPos();
184  }
185  unsigned int nNew = m_validIds.size();
186  if ((mode_ == HGCalGeometryMode::Hexagon) ||
188  edm::LogVerbatim("HGCalGeom") << "ID: " << HGCalDetId(detId)
189  << " with valid DetId from " << nOld
190  << " to " << nNew;
191  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
192  edm::LogVerbatim("HGCalGeom") << "ID: " << HGCScintillatorDetId(detId)
193  << " with valid DetId from " << nOld
194  << " to " << nNew;
195  } else if (m_topology.isHFNose()) {
196  edm::LogVerbatim("HGCalGeom") << "ID: " << HFNoseDetId(detId)
197  << " with valid DetId from " << nOld
198  << " to " << nNew;
199  } else {
200  edm::LogVerbatim("HGCalGeom") << "ID: " << HGCSiliconDetId(detId)
201  << " with valid DetId from " << nOld
202  << " to " << nNew;
203  }
204  edm::LogVerbatim("HGCalGeom") << "Cell[" << cellIndex << "] " << std::hex
205  << geomId.rawId() << ":"
206  << m_validGeomIds[cellIndex].rawId()
207  << std::dec;
208 #endif
209 }
210 
211 std::shared_ptr<const CaloCellGeometry> HGCalGeometry::getGeometry(const DetId& detId) const {
212  if (detId == DetId()) return nullptr; // nothing to get
213  DetId geomId = getGeometryDetId(detId);
214  const uint32_t cellIndex (m_topology.detId2denseGeomId(geomId));
215  const GlobalPoint pos = (detId != geomId) ? getPosition(detId) : GlobalPoint();
216  return cellGeomPtr (cellIndex, pos);
217 
218 }
219 
220 bool HGCalGeometry::present(const DetId& detId) const {
221  if (detId == DetId()) return false;
222  DetId geomId = getGeometryDetId(detId);
223  const uint32_t index (m_topology.detId2denseGeomId(geomId));
224  return (nullptr != getGeometryRawPtr(index)) ;
225 }
226 
228 
229  unsigned int cellIndex = indexFor(detid);
230  GlobalPoint glob;
231  unsigned int maxSize = ((mode_ == HGCalGeometryMode::Trapezoid) ?
232  m_cellVec2.size() : m_cellVec.size());
233  if (cellIndex < maxSize) {
235  std::pair<float,float> xy;
238  xy = m_topology.dddConstants().locateCellHex(id.iCell1,id.iSec1,true);
239  const HepGeom::Point3D<float> lcoord(xy.first,xy.second,0);
240  glob = m_cellVec[cellIndex].getPosition(lcoord);
241 #ifdef EDM_ML_DEBUG
242  edm::LogVerbatim("HGCalGeom") << "getPosition:: index " << cellIndex
243  << " Local " << lcoord.x() << ":"
244  << lcoord.y() << " ID " << id.iCell1
245  << ":" << id.iSec1 << " Global " << glob;
246 #endif
247  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
248  const HepGeom::Point3D<float> lcoord(0,0,0);
249  glob = m_cellVec2[cellIndex].getPosition(lcoord);
250 #ifdef EDM_ML_DEBUG
251  edm::LogVerbatim("HGCalGeom") << "getPositionTrap:: index " << cellIndex
252  << " Local " << lcoord.x() << ":"
253  << lcoord.y() << " ID " << id.iLay << ":"
254  << id.iSec1 << ":" << id.iCell1
255  << " Global " << glob;
256 #endif
257  } else {
258  xy = m_topology.dddConstants().locateCell(id.iLay,id.iSec1,id.iSec2,
259  id.iCell1,id.iCell2,true,false);
260  const HepGeom::Point3D<float> lcoord(xy.first,xy.second,0);
261  glob = m_cellVec[cellIndex].getPosition(lcoord);
262 #ifdef EDM_ML_DEBUG
263  edm::LogVerbatim("HGCalGeom") << "getPositionWafer:: index " << cellIndex
264  << " Local " << lcoord.x() << ":"
265  << lcoord.y() << " ID " << id.iLay << ":"
266  << id.iSec1 << ":" << id.iSec2 << ":"
267  << id.iCell1 << ":" << id.iCell2
268  << " Global " << glob;
269 #endif
270  }
271  }
272  return glob;
273 }
274 
276 
277  unsigned int ncorner = ((m_det == DetId::HGCalHSc) ? FlatTrd::ncorner_ :
279  HGCalGeometry::CornersVec co (ncorner, GlobalPoint(0,0,0));
280  unsigned int cellIndex = indexFor(detid);
282  if (cellIndex < m_cellVec2.size() && m_det == DetId::HGCalHSc) {
283  GlobalPoint v = getPosition(detid);
284  std::pair<double,double> rr = m_topology.dddConstants().cellSizeTrap(id.iType,id.iSec1);
285  float dr = 0.5*(rr.second-rr.first);
286  float dfi= m_cellVec2[cellIndex].param()[FlatTrd::k_Cell];
287  float dz = id.zSide*m_cellVec2[cellIndex].param()[FlatTrd::k_dZ];
288  float r = v.perp();
289  float fi = v.phi();
290  static const int signr[] = {1,1,-1,-1,1,1,-1,-1};
291  static const int signf[] = {-1,1,1,-1,-1,1,1,-1};
292  static const int signz[] = {-1,-1,-1,-1,1,1,1,1};
293  for (unsigned int i = 0; i < ncorner; ++i) {
294  co[i] = GlobalPoint((r+signr[i]*dr)*cos(fi+signf[i]*dfi),
295  (r+signr[i]*dr)*sin(fi+signf[i]*dfi),
296  (v.z()+signz[i]*dz));
297  }
298  } else if (cellIndex < m_cellVec.size() && m_det != DetId::HGCalHSc) {
299  std::pair<float,float> xy;
302  xy = m_topology.dddConstants().locateCellHex(id.iCell1,id.iSec1,true);
303  } else {
304  xy = m_topology.dddConstants().locateCell(id.iLay,id.iSec1,id.iSec2,
305  id.iCell1,id.iCell2,true,false);
306  }
307  float dx = m_cellVec[cellIndex].param()[FlatHexagon::k_r];
308  float dy = m_cellVec[cellIndex].param()[FlatHexagon::k_R];
309  float dz = m_cellVec[cellIndex].param()[FlatHexagon::k_dZ];
310  static const int signx[] = {0,-1,-1,0,1,1,0,-1,-1,0,1,1};
311  static const int signy[] = {-2,-1,1,2,1,-1,-2,-1,1,2,1,-1};
312  static const int signz[] = {-1,-1,-1,-1,-1,-1,1,1,1,1,1,1};
313  for (unsigned int i = 0; i < ncorner; ++i) {
314  const HepGeom::Point3D<float> lcoord(xy.first+signx[i]*dx,xy.second+signy[i]*dy,signz[i]*dz);
315  co[i] = m_cellVec[cellIndex].getPosition(lcoord);
316  }
317  }
318  return co;
319 }
320 
322 
323  unsigned int ncorner = FlatTrd::ncorner_;
324  HGCalGeometry::CornersVec co (ncorner, GlobalPoint(0,0,0));
325  unsigned int cellIndex = indexFor(detid);
327  if (cellIndex < m_cellVec2.size() && m_det == DetId::HGCalHSc) {
328  GlobalPoint v = getPosition(detid);
329  std::pair<double,double> rr = m_topology.dddConstants().cellSizeTrap(id.iType,id.iSec1);
330  float dr = 0.5*(rr.second-rr.first);
331  float dfi= m_cellVec2[cellIndex].param()[FlatTrd::k_Cell];
332  float dz = id.zSide*m_cellVec2[cellIndex].param()[FlatTrd::k_dZ];
333  float r = v.perp();
334  float fi = v.phi();
335  static const int signr[] = {1,1,-1,-1,1,1,-1,-1};
336  static const int signf[] = {-1,1,1,-1,-1,1,1,-1};
337  static const int signz[] = {-1,-1,-1,-1,1,1,1,1};
338  for (unsigned int i = 0; i < ncorner; ++i) {
339  co[i] = GlobalPoint((r+signr[i]*dr)*cos(fi+signf[i]*dfi),
340  (r+signr[i]*dr)*sin(fi+signf[i]*dfi),
341  (v.z()+signz[i]*dz));
342  }
343  } else if (cellIndex < m_cellVec.size() && m_det != DetId::HGCalHSc) {
344  std::pair<float,float> xy;
347  xy = m_topology.dddConstants().locateCellHex(id.iCell1,id.iSec1,true);
348  } else {
349  xy = m_topology.dddConstants().locateCell(id.iLay,id.iSec1,id.iSec2,
350  id.iCell1,id.iCell2,true,false);
351  }
352  static const int signx[] = {-1,-1,1,1,-1,-1,1,1};
353  static const int signy[] = {-1,1,1,-1,-1,1,1,-1};
354  static const int signz[] = {-1,-1,-1,-1,1,1,1,1};
355  float dx = m_cellVec[cellIndex].param()[FlatHexagon::k_r];
356  float dz = m_cellVec[cellIndex].param()[FlatHexagon::k_dZ];
357  for (unsigned int i = 0; i < ncorner; ++i) {
358  const HepGeom::Point3D<float> lcoord(xy.first+signx[i]*dx,
359  xy.second+signy[i]*dx,signz[i]*dz);
360  co[i] = m_cellVec[cellIndex].getPosition(lcoord);
361  }
362  }
363  return co;
364 }
365 
367 
368  unsigned int ncorner = (m_det == DetId::HGCalHSc) ? 5 : 7;
369  HGCalGeometry::CornersVec co (ncorner, GlobalPoint(0,0,0));
370  unsigned int cellIndex = indexFor(detid);
372  if (cellIndex < m_cellVec2.size() && m_det == DetId::HGCalHSc) {
373  GlobalPoint v = getPosition(detid);
374  std::pair<double,double> rr = m_topology.dddConstants().cellSizeTrap(id.iType,id.iSec1);
375  float dr = 0.5*(rr.second-rr.first);
376  float dfi= m_cellVec2[cellIndex].param()[FlatTrd::k_Cell];
377  float dz =-id.zSide*m_cellVec2[cellIndex].param()[FlatTrd::k_dZ];
378  float r = v.perp();
379  float fi = v.phi();
380  static const int signr[] = {1,1,-1,-1};
381  static const int signf[] = {-1,1,1,-1};
382  for (unsigned int i = 0; i < ncorner-1; ++i) {
383  co[i] = GlobalPoint((r+signr[i]*dr)*cos(fi+signf[i]*dfi),
384  (r+signr[i]*dr)*sin(fi+signf[i]*dfi),
385  (v.z()+dz));
386  }
387  co[ncorner-1] = GlobalPoint(0,0,-2*dz);
388  } else if (cellIndex < m_cellVec.size() && m_det != DetId::HGCalHSc) {
389  std::pair<float,float> xy;
392  xy = m_topology.dddConstants().locateCellHex(id.iCell1,id.iSec1,true);
393  } else {
394  xy = m_topology.dddConstants().locateCell(id.iLay,id.iSec1,id.iSec2,
395  id.iCell1,id.iCell2,true,false);
396  }
397  float dx = m_cellVec[cellIndex].param()[FlatHexagon::k_r];
398  float dy = m_cellVec[cellIndex].param()[FlatHexagon::k_R];
399  float dz =-id.zSide*m_cellVec[cellIndex].param()[FlatHexagon::k_dZ];
400  static const int signx[] = {0,-1,-1,0,1,1};
401  static const int signy[] = {-2,-1,1,2,1,-1};
402  for (unsigned int i = 0; i < ncorner-1; ++i) {
403  const HepGeom::Point3D<float> lcoord(xy.first+signx[i]*dx,xy.second+signy[i]*dy,dz);
404  co[i] = m_cellVec[cellIndex].getPosition(lcoord);
405  }
406  co[ncorner-1] = GlobalPoint(0,0,-2*dz);
407  }
408  return co;
409 }
410 
412  const GlobalVector& momentum) const {
413  DetId idnew;
415  int lay = ((momentum.z()*id.zSide > 0) ? (id.iLay+1) : (id.iLay-1));
416 #ifdef EDM_ML_DEBUG
417  edm::LogVerbatim("HGCalGeom") << "neighborz1:: ID " << id.iLay << ":"
418  << id.iSec1 << ":" << id.iSec2 << ":"
419  << id.iCell1 << ":" << id.iCell2
420  << " New Layer " << lay << " Range "
421  << m_topology.dddConstants().firstLayer() <<":"
423  << " pz " << momentum.z();
424 #endif
425  if ((lay >= m_topology.dddConstants().firstLayer()) &&
426  (lay <= m_topology.dddConstants().lastLayer(true)) &&
427  (momentum.z() != 0.0)) {
428  GlobalPoint v = getPosition(idin);
429  double z = id.zSide*m_topology.dddConstants().waferZ(lay,true);
430  double grad = (z - v.z())/momentum.z();
431  GlobalPoint r(v.x()+grad*momentum.x(),v.y()+grad*momentum.y(),z);
432  idnew = getClosestCell(r);
433 #ifdef EDM_ML_DEBUG
434  edm::LogVerbatim("HGCalGeom") << "neighborz1:: Position " << v << " New Z "
435  << z << ":" << grad << " new position " << r;
436 #endif
437  }
438  return idnew;
439 }
440 
442  const MagneticField* bField, int charge,
443  const GlobalVector& momentum) const {
444  DetId idnew;
446  int lay = ((momentum.z()*id.zSide > 0) ? (id.iLay+1) : (id.iLay-1));
447 #ifdef EDM_ML_DEBUG
448  edm::LogVerbatim("HGCalGeom") << "neighborz2:: ID " << id.iLay << ":"
449  << id.iSec1 << ":" << id.iSec2 << ":"
450  << id.iCell1 << ":" << id.iCell2
451  << " New Layer " << lay << " Range "
452  << m_topology.dddConstants().firstLayer() <<":"
454  << " pz " << momentum.z();
455 #endif
456  if ((lay >= m_topology.dddConstants().firstLayer()) &&
457  (lay <= m_topology.dddConstants().lastLayer(true)) &&
458  (momentum.z() != 0.0)) {
459  GlobalPoint v = getPosition(idin);
460  double z = id.zSide*m_topology.dddConstants().waferZ(lay,true);
461  FreeTrajectoryState fts (v, momentum, charge, bField);
464  AnalyticalPropagator myAP (bField, alongMomentum, 2*M_PI);
465  TrajectoryStateOnSurface tsos = myAP.propagate(fts, *nPlane);
466  GlobalPoint r;
467  if (tsos.isValid()) {
468  r = tsos.globalPosition();
469  idnew = getClosestCell(r);
470  }
471 #ifdef EDM_ML_DEBUG
472  edm::LogVerbatim("HGCalGeom") << "neighborz2:: Position " << v << " New Z "
473  << z << ":" << charge << ":" << tsos.isValid()
474  << " new position " << r;
475 #endif
476  }
477  return idnew;
478 }
479 
481  unsigned int cellIndex = getClosestCellIndex(r);
482  if ((cellIndex < m_cellVec.size() && m_det != DetId::HGCalHSc) ||
483  (cellIndex < m_cellVec2.size() && m_det == DetId::HGCalHSc)) {
485  HepGeom::Point3D<float> local;
486  if (r.z() > 0) {
487  local = HepGeom::Point3D<float>(r.x(),r.y(),0);
488  id.zSide = 1;
489  } else {
490  local = HepGeom::Point3D<float>(-r.x(),r.y(),0);
491  id.zSide =-1;
492  }
495  const auto & kxy =
496  m_topology.dddConstants().assignCell(local.x(),local.y(),id.iLay,
497  id.iType,true);
498  id.iCell1 = kxy.second;
499  id.iSec1 = kxy.first;
500  id.iType = m_topology.dddConstants().waferTypeT(kxy.first);
501  if (id.iType != 1) id.iType = -1;
502  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
503  id.iLay = m_topology.dddConstants().getLayer(r.z(),true);
504  const auto & kxy =
505  m_topology.dddConstants().assignCellTrap(r.x(),r.y(),r.z(),
506  id.iLay,true);
507  id.iSec1 = kxy[0];
508  id.iCell1 = kxy[1];
509  id.iType = kxy[2];
510  } else {
511  id.iLay = m_topology.dddConstants().getLayer(r.z(),true);
512  const auto & kxy =
513  m_topology.dddConstants().assignCellHex(local.x(),local.y(),id.iLay,
514  true);
515  id.iSec1 = kxy[0]; id.iSec2 = kxy[1]; id.iType = kxy[2];
516  id.iCell1 = kxy[3]; id.iCell2 = kxy[4];
517  }
518 #ifdef EDM_ML_DEBUG
519  edm::LogVerbatim("HGCalGeom") << "getClosestCell: local " << local
520  << " Id " << id.zSide << ":" << id.iLay
521  << ":" << id.iSec1 << ":" << id.iSec2
522  << ":" << id.iType << ":" << id.iCell1
523  << ":" << id.iCell2;
524 #endif
525 
526  //check if returned cell is valid
527  if (id.iCell1>=0) return m_topology.encode(id);
528  }
529 
530  //if not valid or out of bounds return a null DetId
531  return DetId();
532 }
533 
536  return dss;
537 }
538 
540  if (m_subdet == HGCEE || m_det == DetId::HGCalEE) return "HGCalEE";
541  else if (m_subdet == HGCHEF || m_det == DetId::HGCalHSi) return "HGCalHEFront";
542  else if (m_subdet == HGCHEB || m_det == DetId::HGCalHSc) return "HGCalHEBack";
543  else return "Unknown";
544 }
545 
546 unsigned int HGCalGeometry::indexFor(const DetId& detId) const {
547  unsigned int cellIndex = ((m_det == DetId::HGCalHSc) ? m_cellVec2.size() :
548  m_cellVec.size());
549  if (detId != DetId()) {
550  DetId geomId = getGeometryDetId(detId);
551  cellIndex = m_topology.detId2denseGeomId(geomId);
552 #ifdef EDM_ML_DEBUG
553  edm::LogVerbatim("HGCalGeom") << "indexFor " << std::hex << detId.rawId()
554  << ":" << geomId.rawId() << std::dec
555  << " index " << cellIndex;
556 #endif
557  }
558  return cellIndex;
559 }
560 
561 unsigned int HGCalGeometry::sizeForDenseIndex() const {
562  return m_topology.totalGeomModules();
563 }
564 
566  // Modify the RawPtr class
567  if (m_det == DetId::HGCalHSc) {
568  if (m_cellVec2.size() < index) return nullptr;
569  const CaloCellGeometry* cell(&m_cellVec2[index]);
570  return (nullptr == cell->param() ? nullptr : cell);
571  } else {
572  if (m_cellVec2.size() < index) return nullptr;
573  const CaloCellGeometry* cell(&m_cellVec[index]);
574  return (nullptr == cell->param() ? nullptr : cell);
575  }
576 }
577 
578 std::shared_ptr<const CaloCellGeometry> HGCalGeometry::cellGeomPtr(uint32_t index) const {
579  if ((index >= m_cellVec.size() && m_det != DetId::HGCalHSc) ||
580  (index >= m_cellVec2.size() && m_det == DetId::HGCalHSc) ||
581  (m_validGeomIds[index].rawId() == 0)) return nullptr;
582  static const auto do_not_delete = [](const void*){};
583  if (m_det == DetId::HGCalHSc) {
584  auto cell = std::shared_ptr<const CaloCellGeometry>(&m_cellVec2[index],do_not_delete);
585  if (nullptr == cell->param()) return nullptr;
586  return cell;
587  } else {
588  auto cell = std::shared_ptr<const CaloCellGeometry>(&m_cellVec[index],do_not_delete);
589  if (nullptr == cell->param()) return nullptr;
590  return cell;
591  }
592 }
593 
594 std::shared_ptr<const CaloCellGeometry> HGCalGeometry::cellGeomPtr(uint32_t index, const GlobalPoint& pos) const {
595  if ((index >= m_cellVec.size() && m_det != DetId::HGCalHSc) ||
596  (index >= m_cellVec2.size() && m_det == DetId::HGCalHSc) ||
597  (m_validGeomIds[index].rawId() == 0)) return nullptr;
598  if (pos == GlobalPoint()) return cellGeomPtr(index);
599  if (m_det == DetId::HGCalHSc) {
600  auto cell = std::make_shared<FlatTrd>(m_cellVec2[index]);
601  cell->setPosition(pos);
602 #ifdef EDM_ML_DEBUG
603  edm::LogVerbatim("HGCalGeom") << "cellGeomPtr " << index << ":" << cell;
604 #endif
605  if (nullptr == cell->param()) return nullptr;
606  return cell;
607  } else {
608  auto cell = std::make_shared<FlatHexagon>(m_cellVec[index]);
609  cell->setPosition(pos);
610 #ifdef EDM_ML_DEBUG
611  edm::LogVerbatim("HGCalGeom") << "cellGeomPtr " << index << ":" << cell;
612 #endif
613  if (nullptr == cell->param()) return nullptr;
614  return cell;
615  }
616 }
617 
619  edm::LogError("HGCalGeom") << "HGCalGeometry::addValidID is not implemented";
620 }
621 
622 unsigned int HGCalGeometry::getClosestCellIndex (const GlobalPoint& r) const {
623 
625 }
626 
627 template<class T>
629  const std::vector<T>& vec) const {
630 
631  float phip = r.phi();
632  float zp = r.z();
633  float dzmin(9999), dphimin(9999), dphi10(0.175);
634  unsigned int cellIndex = vec.size();
635  for (unsigned int k=0; k<vec.size(); ++k) {
636  float dphi = phip-vec[k].phiPos();
637  while (dphi > M_PI) dphi -= 2*M_PI;
638  while (dphi <= -M_PI) dphi += 2*M_PI;
639  if (std::abs(dphi) < dphi10) {
640  float dz = std::abs(zp - vec[k].getPosition().z());
641  if (dz < (dzmin+0.001)) {
642  dzmin = dz;
643  if (std::abs(dphi) < (dphimin+0.01)) {
644  cellIndex = k;
645  dphimin = std::abs(dphi);
646  } else {
647  if (cellIndex >= vec.size()) cellIndex = k;
648  }
649  }
650  }
651  }
652 #ifdef EDM_ML_DEBUG
653  edm::LogVerbatim("HGCalGeom") << "getClosestCellIndex::Input " << zp << ":"
654  << phip << " Index " << cellIndex;
655  if (cellIndex < vec.size())
656  edm::LogVerbatim("HGCalGeom") << " Cell z "
657  << vec[cellIndex].getPosition().z()
658  << ":" << dzmin << " phi "
659  << vec[cellIndex].phiPos() << ":" << dphimin;
660 #endif
661  return cellIndex;
662 }
663 
664 
665 // FIXME: Change sorting algorithm if needed
666 namespace {
667  struct rawIdSort {
668  bool operator()( const DetId& a, const DetId& b ) {
669  return( a.rawId() < b.rawId());
670  }
671  };
672 }
673 
675  m_validIds.shrink_to_fit();
676  std::sort( m_validIds.begin(), m_validIds.end(), rawIdSort());
677 }
678 
682  CaloSubdetectorGeometry::IVec& dinsVector ) const {
683 
684  unsigned int numberOfCells = m_topology.totalGeomModules(); // total Geom Modules both sides
685  unsigned int numberOfShapes = k_NumberOfShapes;
686  unsigned int numberOfParametersPerShape =
687  ((m_det == DetId::HGCalHSc) ? (unsigned int)(k_NumberOfParametersPerTrd) :
688  (unsigned int)(k_NumberOfParametersPerHex));
689 
690  trVector.reserve( numberOfCells * numberOfTransformParms());
691  iVector.reserve( numberOfCells );
692  dimVector.reserve( numberOfShapes * numberOfParametersPerShape );
693  dinsVector.reserve( numberOfCells );
694 
695  for (unsigned itr=0; itr<m_topology.dddConstants().getTrFormN(); ++itr) {
697  int layer = mytr.lay;
698 
701  for (int wafer=0; wafer<m_topology.dddConstants().sectors(); ++wafer) {
702  if (m_topology.dddConstants().waferInLayer(wafer,layer,true)) {
703  HGCalParameters::hgtrap vol = m_topology.dddConstants().getModule(wafer, true, true);
704  ParmVec params( numberOfParametersPerShape, 0 );
705  params[FlatHexagon::k_dZ] = vol.dz;
706  params[FlatHexagon::k_r] = vol.cellSize;
707  params[FlatHexagon::k_R] = twoBysqrt3_*params[1];
708  dimVector.insert( dimVector.end(), params.begin(), params.end());
709  }
710  }
711  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
712  int indx = m_topology.dddConstants().layerIndex(layer,true);
713  for (int md=m_topology.dddConstants().getParameter()->firstModule_[indx];
715  ++md) {
717  ParmVec params( numberOfParametersPerShape, 0 );
718  params[FlatTrd::k_dZ] = vol.dz;
719  params[FlatTrd::k_Theta] = params[FlatTrd::k_Phi] = 0;
720  params[FlatTrd::k_dY1] = params[FlatTrd::k_dY2] = vol.h;
721  params[FlatTrd::k_dX1] = params[FlatTrd::k_dX3] = vol.bl;
722  params[FlatTrd::k_dX2] = params[FlatTrd::k_dX4] = vol.tl;
723  params[FlatTrd::k_Alp1] = params[FlatTrd::k_Alp2]= vol.alpha;
724  params[FlatTrd::k_Cell] = vol.cellSize;
725  dimVector.insert( dimVector.end(), params.begin(), params.end());
726  }
727  } else {
728  for (int wafer=0; wafer<m_topology.dddConstants().sectors(); ++wafer) {
729  if (m_topology.dddConstants().waferInLayer(wafer,layer,true)) {
730  HGCalParameters::hgtrap vol = m_topology.dddConstants().getModule(wafer, true, true);
731  ParmVec params( numberOfParametersPerShape, 0 );
732  params[FlatHexagon::k_dZ] = vol.dz;
733  params[FlatHexagon::k_r] = vol.cellSize;
734  params[FlatHexagon::k_R] = twoBysqrt3_*params[1];
735  dimVector.insert( dimVector.end(), params.begin(), params.end());
736  }
737  }
738  }
739  }
740 
741  for (unsigned int i( 0 ); i < numberOfCells; ++i) {
742  DetId detId = m_validGeomIds[i];
743  int layer(0);
746  layer = HGCalDetId(detId).layer();
747  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
748  layer = HGCScintillatorDetId(detId).layer();
749  } else if (m_topology.isHFNose()) {
750  layer = HFNoseDetId(detId).layer();
751  } else {
752  layer = HGCSiliconDetId(detId).layer();
753  }
754  dinsVector.emplace_back(m_topology.detId2denseGeomId( detId ));
755  iVector.emplace_back( layer );
756 
757  Tr3D tr;
758  auto ptr = cellGeomPtr( i );
759  if ( nullptr != ptr ) {
760  ptr->getTransform( tr, ( Pt3DVec* ) nullptr );
761 
762  if( Tr3D() == tr ) { // there is no rotation
763  const GlobalPoint& gp( ptr->getPosition());
764  tr = HepGeom::Translate3D( gp.x(), gp.y(), gp.z());
765  }
766 
767  const CLHEP::Hep3Vector tt( tr.getTranslation());
768  trVector.emplace_back( tt.x());
769  trVector.emplace_back( tt.y());
770  trVector.emplace_back( tt.z());
771  if (6 == numberOfTransformParms()) {
772  const CLHEP::HepRotation rr( tr.getRotation());
773  const ROOT::Math::Transform3D rtr( rr.xx(), rr.xy(), rr.xz(), tt.x(),
774  rr.yx(), rr.yy(), rr.yz(), tt.y(),
775  rr.zx(), rr.zy(), rr.zz(), tt.z());
777  rtr.GetRotation( ea );
778  trVector.emplace_back( ea.Phi());
779  trVector.emplace_back( ea.Theta());
780  trVector.emplace_back( ea.Psi());
781  }
782  }
783  }
784 }
785 
787  DetId geomId;
788  if ((mode_ == HGCalGeometryMode::Hexagon) ||
790  geomId = static_cast<DetId>(HGCalDetId(detId).geometryCell());
791  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
792  geomId = static_cast<DetId>(HGCScintillatorDetId(detId).geometryCell());
793  } else if (m_topology.isHFNose()) {
794  geomId = static_cast<DetId>(HFNoseDetId(detId).geometryCell());
795  } else {
796  geomId = static_cast<DetId>(HGCSiliconDetId(detId).geometryCell());
797  }
798  return geomId;
799 }
800 
802 
std::vector< DetId > m_validGeomIds
HGCalGeometryMode::GeometryMode mode_
A base class to handle the particular shape of HGCal volumes.
Definition: FlatTrd.h:19
int getLayer(double z, bool reco) const
std::vector< CCGFloat > DimVec
T perp() const
Definition: PV3DBase.h:72
std::vector< FlatHexagon > CellVec
Definition: HGCalGeometry.h:34
virtual unsigned int numberOfParametersPerShape() const
std::pair< double, double > cellSizeTrap(int type, int irad) const
std::string cellElement() const
std::shared_ptr< const CaloCellGeometry > cellGeomPtr(uint32_t index) const override
unsigned int sizeForDenseIndex() const
HGCalParameters::hgtrform getTrForm(unsigned int k) const
int lastLayer(bool reco) const
unsigned int getClosestCellIndex(const GlobalPoint &r) const
DetId neighborZ(const DetId &idin, const GlobalVector &p) const
const HGCalParameters * getParameter() const
const HGCalTopology & m_topology
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
~HGCalGeometry() override
bool cellInLayer(int waferU, int waferV, int cellU, int cellV, int lay, bool reco) const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
std::vector< unsigned int > IVec
T y() const
Definition: PV3DBase.h:63
static constexpr uint32_t k_Theta
Definition: FlatTrd.h:28
GlobalPoint globalPosition() const
void newCell(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId) override
std::vector< float > ParmVec
std::vector< GlobalPoint > CornersVec
Definition: HGCalGeometry.h:41
std::vector< CCGFloat > TrVec
GlobalPoint getPosition(const DetId &id) const
HepGeom::Transform3D Tr3D
HFNoseDetId geometryCell() const
Definition: HFNoseDetId.h:49
static constexpr uint32_t k_Cell
Definition: FlatTrd.h:46
static constexpr uint32_t k_dY1
Definition: FlatTrd.h:32
std::pair< float, float > locateCell(int cell, int lay, int type, bool reco) const
CaloCellGeometry::CCGFloat CCGFloat
CornersVec getCorners(const DetId &id) const
Returns the corner points of this cell&#39;s volume.
CaloCellGeometry::Tr3D Tr3D
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
Definition: FlatHexagon.cc:175
int layerIndex(int lay, bool reco) const
static constexpr unsigned int ncorner_
Definition: FlatHexagon.h:82
CaloCellGeometry::Pt3DVec Pt3DVec
Definition: HGCalGeometry.h:38
std::pair< float, float > locateCellHex(int cell, int wafer, bool reco) const
const double twoBysqrt3_
unsigned int totalGeomModules() const
HGCalGeometry(const HGCalTopology &topology)
static constexpr unsigned int ncorner_
Definition: FlatTrd.h:98
std::vector< int > firstModule_
static unsigned int k_NumberOfShapes
Definition: HGCalGeometry.h:49
HGCalDetId geometryCell() const
Definition: HGCalDetId.h:36
static constexpr uint32_t k_R
Definition: FlatHexagon.h:30
void initializeParms() override
unsigned int getTrFormN() const
HGCSiliconDetId geometryCell() const
const CCGFloat * param() const
static unsigned int k_NumberOfParametersPerHex
Definition: HGCalGeometry.h:47
virtual uint32_t detId2denseGeomId(const DetId &id) const
int numberCellsHexagon(int wafer) const
static constexpr uint32_t k_r
Definition: FlatHexagon.h:29
int layer() const
get the layer #
Definition: HFNoseDetId.h:58
virtual unsigned int numberOfShapes() const
std::vector< float > ParmVec
DetId encode(const DecodedDetId &id_) const
T sqrt(T t)
Definition: SSEVec.h:18
int layer() const
get the layer #
static constexpr uint32_t k_dX1
Definition: FlatTrd.h:33
std::vector< DetId > m_validIds
static PlanePointer build(Args &&...args)
Definition: Plane.h:33
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
A base class to handle the hexagonal shape of HGCal silicon volumes.
Definition: FlatHexagon.h:20
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::set< DetId > DetIdSet
Definition: HGCalGeometry.h:40
static constexpr uint32_t k_dY2
Definition: FlatTrd.h:39
HGCScintillatorDetId geometryCell() const
CornersVec get8Corners(const DetId &id) const
CellVec2 m_cellVec2
static constexpr uint32_t k_dX4
Definition: FlatTrd.h:42
static unsigned int k_NumberOfParametersPerTrd
Definition: HGCalGeometry.h:46
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
Definition: FlatTrd.cc:175
DetIdSet getCells(const GlobalPoint &r, double dR) const override
Get a list of all cells within a dR of the given cell.
double waferZ(int layer, bool reco) const
void localCorners(Pt3DVec &lc, const CCGFloat *pv, unsigned int i, Pt3D &ref)
#define M_PI
bool isHFNose() const
int k[5][pyjets_maxn]
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:96
DecodedDetId decode(const DetId &id) const
Definition: DetId.h:18
std::vector< int > lastModule_
AlgebraicVector EulerAngles
Definition: Definitions.h:36
ForwardSubdetector m_subdet
static constexpr uint32_t k_dZ
Definition: FlatHexagon.h:28
std::vector< FlatTrd > CellVec2
Definition: HGCalGeometry.h:35
const HGCalDDDConstants & dddConstants() const
DetId getClosestCell(const GlobalPoint &r) const override
static constexpr uint32_t k_Phi
Definition: FlatTrd.h:30
CaloCellGeometry::CornersMgr * cornersMgr()
CornersVec getNewCorners(const DetId &id) const
void getSummary(CaloSubdetectorGeometry::TrVec &trVector, CaloSubdetectorGeometry::IVec &iVector, CaloSubdetectorGeometry::DimVec &dimVector, CaloSubdetectorGeometry::IVec &dinsVector) const override
double b
Definition: hdecay.h:120
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:53
unsigned int totalModules() const
Definition: HGCalTopology.h:99
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const override
Get the cell geometry of a given detector id. Should return false if not found.
CaloCellGeometry::Tr3D Tr3D
int layer() const
get the layer #
CaloCellGeometry::Pt3D Pt3D
Definition: HGCalGeometry.h:37
std::array< int, 5 > assignCellHex(float x, float y, int lay, bool reco) const
double a
Definition: hdecay.h:121
static constexpr uint32_t k_dZ
Definition: FlatTrd.h:27
virtual void fillNamedParams(DDFilteredView fv)
std::pair< int, int > assignCell(float x, float y, int lay, int subSec, bool reco) const
DetId::Detector m_det
std::array< int, 3 > assignCellTrap(float x, float y, float z, int lay, bool reco) const
const CaloCellGeometry * getGeometryRawPtr(uint32_t index) const override
int firstLayer() const
static constexpr uint32_t k_Alp1
Definition: FlatTrd.h:37
void addValidID(const DetId &id)
bool valid(const DetId &id) const override
Is this a valid cell id.
virtual unsigned int numberOfTransformParms() const
static constexpr uint32_t k_dX3
Definition: FlatTrd.h:40
HGCalParameters::hgtrap getModule(unsigned int k, bool hexType, bool reco) const
T x() const
Definition: PV3DBase.h:62
int waferTypeT(int wafer) const
unsigned int indexFor(const DetId &id) const override
bool present(const DetId &id) const override
is this detid present in the geometry?
int layer() const
get the layer #
Definition: HGCalDetId.h:48
static constexpr uint32_t k_dX2
Definition: FlatTrd.h:35
CellVec m_cellVec
bool waferInLayer(int wafer, int lay, bool reco) const
DetId getGeometryDetId(DetId detId) const
static constexpr uint32_t k_Alp2
Definition: FlatTrd.h:44