CMS 3D CMS Logo

HcalGeometry.cc
Go to the documentation of this file.
5 #include <algorithm>
6 #include <iostream>
7 
8 #include <Math/Transform3D.h>
9 #include <Math/EulerAngles.h>
10 
15 
16 //#define EDM_ML_DEBUG
17 
19  : m_topology(topology), m_mergePosition(topology.getMergePositionFlag()) {
20  init();
21 }
22 
24 
27  m_mergePosition = false;
28 #ifdef EDM_ML_DEBUG
29  edm::LogVerbatim("HCalGeom") << "HcalGeometry_init(): HBSize " << m_topology.getHBSize() << " HESize "
30  << m_topology.getHESize() << " HOSize " << m_topology.getHOSize() << " HFSize "
31  << m_topology.getHFSize();
32 #endif
33 
38 }
39 
41  // this must test the last record filled to avoid a race condition
42  if (!m_emptyIds.isSet()) {
43  std::unique_ptr<std::vector<DetId>> p_hbIds{new std::vector<DetId>};
44  std::unique_ptr<std::vector<DetId>> p_heIds{new std::vector<DetId>};
45  std::unique_ptr<std::vector<DetId>> p_hoIds{new std::vector<DetId>};
46  std::unique_ptr<std::vector<DetId>> p_hfIds{new std::vector<DetId>};
47  std::unique_ptr<std::vector<DetId>> p_emptyIds{new std::vector<DetId>};
48 
49  const std::vector<DetId>& baseIds(CaloSubdetectorGeometry::getValidDetIds());
50  for (unsigned int i(0); i != baseIds.size(); ++i) {
51  const DetId id(baseIds[i]);
52  if (id.subdetId() == HcalBarrel) {
53  p_hbIds->emplace_back(id);
54  } else if (id.subdetId() == HcalEndcap) {
55  p_heIds->emplace_back(id);
56  } else if (id.subdetId() == HcalOuter) {
57  p_hoIds->emplace_back(id);
58  } else if (id.subdetId() == HcalForward) {
59  p_hfIds->emplace_back(id);
60  }
61  }
62  std::sort(p_hbIds->begin(), p_hbIds->end());
63  std::sort(p_heIds->begin(), p_heIds->end());
64  std::sort(p_hoIds->begin(), p_hoIds->end());
65  std::sort(p_hfIds->begin(), p_hfIds->end());
66  p_emptyIds->resize(0);
67 
68  m_hbIds.set(std::move(p_hbIds));
69  m_heIds.set(std::move(p_heIds));
70  m_hoIds.set(std::move(p_hoIds));
71  m_hfIds.set(std::move(p_hfIds));
72  m_emptyIds.set(std::move(p_emptyIds));
73  }
74 }
75 
76 const std::vector<DetId>& HcalGeometry::getValidDetIds(DetId::Detector det, int subdet) const {
77  if (0 != subdet)
78  fillDetIds();
79  return (0 == subdet
81  : (HcalBarrel == subdet
82  ? *m_hbIds.load()
83  : (HcalEndcap == subdet
84  ? *m_heIds.load()
85  : (HcalOuter == subdet ? *m_hoIds.load()
86  : (HcalForward == subdet ? *m_hfIds.load() : *m_emptyIds.load())))));
87 }
88 
89 std::shared_ptr<const CaloCellGeometry> HcalGeometry::getGeometry(const DetId& id) const {
90 #ifdef EDM_ML_DEBUG
91  std::ostringstream st1;
92  st1 << "HcalGeometry::getGeometry for " << HcalDetId(id) << " " << m_mergePosition << " ";
93 #endif
94  if (!m_mergePosition) {
95 #ifdef EDM_ML_DEBUG
96  st1 << m_topology.detId2denseId(id) << " " << getGeometryBase(id) << "\n";
97  edm::LogVerbatim("HCalGeom") << st1.str();
98 #endif
99  return getGeometryBase(id);
100  } else {
101 #ifdef EDM_ML_DEBUG
103  edm::LogVerbatim("HCalGeom") << st1.str();
104 #endif
105  return getGeometryBase(m_topology.idFront(id));
106  }
107 }
108 
110 
111 DetId HcalGeometry::getClosestCell(const GlobalPoint& r, bool ignoreCorrect) const {
112  // Now find the closest eta_bin, eta value of a bin i is average
113  // of eta[i] and eta[i-1]
114  static const double z_long = 1100.0;
115  double abseta = fabs(r.eta());
116  double absz = fabs(r.z());
117 
118  // figure out subdetector, giving preference to HE in HE/HF overlap region
120  if (abseta <= m_topology.etaMax(HcalBarrel)) {
121  bc = HcalBarrel;
122  } else if (absz >= z_long) {
123  bc = HcalForward;
124  } else if (abseta <= m_topology.etaMax(HcalEndcap)) {
125  bc = HcalEndcap;
126  } else {
127  bc = HcalForward;
128  }
129 
130  // find eta bin
131  int etaring = etaRing(bc, abseta);
132 
133  int phibin = phiBin(bc, etaring, r.phi());
134 
135  // add a sign to the etaring
136  int etabin = (r.z() > 0) ? etaring : -etaring;
137 
138  if (bc == HcalForward) {
139  static const double z_short = 1137.0;
140  // Next line is premature depth 1 and 2 can coexist for large z-extent
141  // HcalDetId bestId(bc,etabin,phibin,((fabs(r.z())>=z_short)?(2):(1)));
142  // above line is no good with finite precision
143  HcalDetId bestId(bc, etabin, phibin, ((fabs(r.z()) - z_short > -0.1) ? (2) : (1)));
144  return correctId(bestId);
145  } else {
146  //Now do depth if required
147  int zside = (r.z() > 0) ? 1 : -1;
148  int dbin = m_topology.dddConstants()->getMinDepth(((int)(bc)-1), etaring, phibin, zside);
149  double pointrz(0), drz(99999.);
150  HcalDetId currentId(bc, etabin, phibin, dbin);
151  if (bc == HcalBarrel)
152  pointrz = r.mag();
153  else
154  pointrz = std::abs(r.z());
155  HcalDetId bestId;
156  for (; currentId != HcalDetId(); m_topology.incrementDepth(currentId)) {
157  std::shared_ptr<const CaloCellGeometry> cell = getGeometry(currentId);
158  if (cell == nullptr) {
159  assert(bestId != HcalDetId());
160  break;
161  } else {
162  double rz;
163  if (bc == HcalEndcap)
164  rz = std::abs(cell->getPosition().z());
165  else
166  rz = cell->getPosition().mag();
167  if (std::abs(pointrz - rz) < drz) {
168  bestId = currentId;
169  drz = std::abs(pointrz - rz);
170  }
171  }
172  }
173 #ifdef EDM_ML_DEBUG
174  edm::LogVerbatim("HCalGeom") << bestId << " Corrected to " << HcalDetId(correctId(bestId));
175 #endif
176 
177  return (ignoreCorrect ? bestId : correctId(bestId));
178  }
179 }
180 
182  if (!m_mergePosition) {
183  return (getGeometryBase(id)->getPosition());
184  } else {
186  }
187 }
188 
190  if (test) {
191  int subdet, z, depth, eta, phi, lay;
193  int sign = (z == 0) ? (-1) : (1);
194  auto id = m_topology.dddConstants()->getHCID(subdet, (sign * eta), phi, lay, depth);
195  auto hcId = ((id.subdet == static_cast<int>(HcalBarrel))
196  ? HcalDetId(HcalBarrel, sign * id.eta, id.phi, id.depth)
197  : ((id.subdet == static_cast<int>(HcalEndcap))
198  ? HcalDetId(HcalEndcap, sign * id.eta, id.phi, id.depth)
199  : ((id.subdet == static_cast<int>(HcalOuter))
200  ? HcalDetId(HcalOuter, sign * id.eta, id.phi, id.depth)
201  : ((id.subdet == static_cast<int>(HcalForward))
202  ? HcalDetId(HcalForward, sign * id.eta, id.phi, id.depth)
203  : HcalDetId()))));
204  return getPosition(DetId(hcId));
205  } else {
206  return getPosition(DetId(idx));
207  }
208 }
209 
211  if (!m_mergePosition) {
212  return (getGeometryBase(id)->getBackPoint());
213  } else {
214  std::vector<HcalDetId> ids;
216  return (getGeometryBase(ids.back())->getBackPoint());
217  }
218 }
219 
221  if (test) {
222  int subdet, z, depth, eta, phi, lay;
224  int sign = (z == 0) ? (-1) : (1);
225  auto id = m_topology.dddConstants()->getHCID(subdet, (sign * eta), phi, lay, depth);
226  auto hcId = ((id.subdet == static_cast<int>(HcalBarrel))
227  ? HcalDetId(HcalBarrel, sign * id.eta, id.phi, id.depth)
228  : ((id.subdet == static_cast<int>(HcalEndcap))
229  ? HcalDetId(HcalEndcap, sign * id.eta, id.phi, id.depth)
230  : ((id.subdet == static_cast<int>(HcalOuter))
231  ? HcalDetId(HcalOuter, sign * id.eta, id.phi, id.depth)
232  : ((id.subdet == static_cast<int>(HcalForward))
233  ? HcalDetId(HcalForward, sign * id.eta, id.phi, id.depth)
234  : HcalDetId()))));
235  return getBackPosition(DetId(hcId));
236  } else {
237  return getBackPosition(DetId(idx));
238  }
239 }
240 
242  if (!m_mergePosition) {
243  return (getGeometryBase(id)->getCorners());
244  } else {
245  std::vector<HcalDetId> ids;
250  for (unsigned int k = 0; k < 4; ++k) {
251  mcorners[k] = mcf[k];
252  mcorners[k + 4] = mcb[k + 4];
253  }
254  return mcorners;
255  }
256 }
257 
258 int HcalGeometry::etaRing(HcalSubdetector bc, double abseta) const { return m_topology.etaRing(bc, abseta); }
259 
260 int HcalGeometry::phiBin(HcalSubdetector bc, int etaring, double phi) const {
261  return m_topology.phiBin(bc, etaring, phi);
262 }
263 
265  CaloSubdetectorGeometry::DetIdSet dis; // this is the return object
266 
267  if (0.000001 < dR) {
268  if (dR > M_PI / 2.) { // this version needs "small" dR
269  dis = CaloSubdetectorGeometry::getCells(r, dR); // base class version
270  } else {
271  const double dR2(dR * dR);
272  const double reta(r.eta());
273  const double rphi(r.phi());
274  const double lowEta(reta - dR);
275  const double highEta(reta + dR);
276  const double lowPhi(rphi - dR);
277  const double highPhi(rphi + dR);
278 
279  const double hfEtaHi(m_topology.etaMax(HcalForward));
280 
281  if (highEta > -hfEtaHi && lowEta < hfEtaHi) { // in hcal
283 
284  for (unsigned int is(0); is != 4; ++is) {
285  const int sign(reta > 0 ? 1 : -1);
286  const int ieta_center(sign * etaRing(hs[is], fabs(reta)));
287  const int ieta_lo((0 < lowEta * sign ? sign : -sign) * etaRing(hs[is], fabs(lowEta)));
288  const int ieta_hi((0 < highEta * sign ? sign : -sign) * etaRing(hs[is], fabs(highEta)));
289  const int iphi_lo(phiBin(hs[is], ieta_center, lowPhi));
290  const int iphi_hi(phiBin(hs[is], ieta_center, highPhi));
291  const int jphi_lo(iphi_lo > iphi_hi ? iphi_lo - 72 : iphi_lo);
292  const int jphi_hi(iphi_hi);
293 
294  const int idep_lo(1 == is ? 4 : 1);
295  const int idep_hi(m_topology.maxDepth(hs[is]));
296  for (int ieta(ieta_lo); ieta <= ieta_hi; ++ieta) { // over eta limits
297  if (ieta != 0) {
298  for (int jphi(jphi_lo); jphi <= jphi_hi; ++jphi) { // over phi limits
299  const int iphi(1 > jphi ? jphi + 72 : jphi);
300  for (int idep(idep_lo); idep <= idep_hi; ++idep) {
301  const HcalDetId did(hs[is], ieta, iphi, idep);
302  if (m_topology.valid(did)) {
303  std::shared_ptr<const CaloCellGeometry> cell(getGeometryBase(did));
304  if (nullptr != cell) {
305  const GlobalPoint& p(cell->getPosition());
306  const double eta(p.eta());
307  const double phi(p.phi());
308  if (reco::deltaR2(eta, phi, reta, rphi) < dR2)
309  dis.insert(did);
310  }
311  }
312  }
313  }
314  }
315  }
316  }
317  }
318  }
319  }
320  return dis;
321 }
322 
323 unsigned int HcalGeometry::getHxSize(const int type) const {
324  unsigned int hxsize(0);
325  if (!m_hbIds.isSet())
326  fillDetIds();
327  if (type == 1)
328  hxsize = (m_hbIds.isSet()) ? m_hbIds->size() : 0;
329  else if (type == 2)
330  hxsize = (m_heIds.isSet()) ? m_heIds->size() : 0;
331  else if (type == 3)
332  hxsize = (m_hoIds.isSet()) ? m_hoIds->size() : 0;
333  else if (type == 4)
334  hxsize = (m_hfIds.isSet()) ? m_hfIds->size() : 0;
335  else if (type == 0)
336  hxsize = (m_emptyIds.isSet()) ? m_emptyIds->size() : 0;
337  return hxsize;
338 }
339 
342  const int ieta(i < numberOfBarrelAlignments() / 2 ? -1 : 1);
343  const int iphi(1 + (4 * i) % 72);
344  return HcalDetId(HcalBarrel, ieta, iphi, 1);
345 }
346 
349  const int ieta(i < numberOfEndcapAlignments() / 2 ? -16 : 16);
350  const int iphi(1 + (4 * i) % 72);
351  return HcalDetId(HcalEndcap, ieta, iphi, 1);
352 }
353 
356  const int ieta(i < numberOfForwardAlignments() / 2 ? -29 : 29);
357  const int iphi(1 + (4 * i) % 72);
358  return HcalDetId(HcalForward, ieta, iphi, 1);
359 }
360 
363  const int ring(i / 12);
364  const int ieta(0 == ring ? -11 : 1 == ring ? -5 : 2 == ring ? 1 : 3 == ring ? 5 : 11);
365  const int iphi(1 + (i - ring * 12) * 6);
366  return HcalDetId(HcalOuter, ieta, iphi, 4);
367 }
368 
371 
372  const unsigned int nB(numberOfBarrelAlignments());
373  const unsigned int nE(numberOfEndcapAlignments());
374  const unsigned int nF(numberOfForwardAlignments());
375  // const unsigned int nO ( numberOfOuterAlignments() ) ;
376 
377  return (i < nB ? detIdFromBarrelAlignmentIndex(i)
378  : i < nB + nE ? detIdFromEndcapAlignmentIndex(i - nB)
379  : i < nB + nE + nF ? detIdFromForwardAlignmentIndex(i - nB - nE)
380  : detIdFromOuterAlignmentIndex(i - nB - nE - nF));
381 }
382 
383 unsigned int HcalGeometry::alignmentBarEndForIndexLocal(const DetId& id, unsigned int nD) {
384  const HcalDetId hid(id);
385  const unsigned int iphi(hid.iphi());
386  const int ieta(hid.ieta());
387  const unsigned int index((0 < ieta ? nD / 2 : 0) + (iphi + 1) % 72 / 4);
388  assert(index < nD);
389  return index;
390 }
391 
394 }
395 
398 }
399 
402 }
403 
405  const HcalDetId hid(id);
406  const int ieta(hid.ieta());
407  const int iphi(hid.iphi());
408  const int ring(ieta < -10 ? 0 : (ieta < -4 ? 1 : (ieta < 5 ? 2 : (ieta < 11 ? 3 : 4))));
409 
410  const unsigned int index(12 * ring + (iphi - 1) / 6);
412  return index;
413 }
414 
416  assert(id.det() == DetId::Hcal);
417 
418  const HcalDetId hid(id);
419  bool isHB = (hid.subdet() == HcalBarrel);
420  bool isHE = (hid.subdet() == HcalEndcap);
421  bool isHF = (hid.subdet() == HcalForward);
422  // bool isHO = (hid.subdet() == HcalOuter);
423 
424  const unsigned int nB(numberOfBarrelAlignments());
425  const unsigned int nE(numberOfEndcapAlignments());
426  const unsigned int nF(numberOfForwardAlignments());
427  // const unsigned int nO ( numberOfOuterAlignments() ) ;
428 
429  const unsigned int index(isHB ? alignmentBarrelIndexLocal(id)
430  : isHE ? alignmentEndcapIndexLocal(id) + nB
431  : isHF ? alignmentForwardIndexLocal(id) + nB + nE
432  : alignmentOuterIndexLocal(id) + nB + nE + nF);
433 
435  return index;
436 }
437 
438 unsigned int HcalGeometry::alignmentTransformIndexGlobal(const DetId& id) { return (unsigned int)DetId::Hcal - 1; }
439 
440 void HcalGeometry::localCorners(Pt3DVec& lc, const CCGFloat* pv, unsigned int i, Pt3D& ref) {
442 
443  if (hid.subdet() == HcalForward) {
444  IdealZPrism::localCorners(lc, pv, ref);
445  } else {
447  }
448 }
449 
451  const GlobalPoint& f1, const GlobalPoint& f2, const GlobalPoint& f3, const CCGFloat* parm, const DetId& detId) {
452  assert(detId.det() == DetId::Hcal);
453 
454  const HcalDetId hid(detId);
455  unsigned int din = m_topology.detId2denseId(detId);
456 
457 #ifdef EDM_ML_DEBUG
458  edm::LogVerbatim("HCalGeom") << "HcalGeometry: newCell subdet " << detId.subdetId() << ", raw ID " << detId.rawId()
459  << ", hid " << hid << ", din " << din << ", index ";
460 #endif
461 
462  if (hid.subdet() == HcalBarrel) {
464  } else if (hid.subdet() == HcalEndcap) {
465  const unsigned int index(din - m_hbCellVec.size());
467  } else if (hid.subdet() == HcalOuter) {
468  const unsigned int index(din - m_hbCellVec.size() - m_heCellVec.size());
470  } else {
471  const unsigned int index(din - m_hbCellVec.size() - m_heCellVec.size() - m_hoCellVec.size());
472  m_hfCellVec.at(index) = IdealZPrism(f1, cornersMgr(), parm, hid.depth() == 1 ? IdealZPrism::EM : IdealZPrism::HADR);
473  }
474 
475  return din;
476 }
477 
479  const GlobalPoint& f1, const GlobalPoint& f2, const GlobalPoint& f3, const CCGFloat* parm, const DetId& detId) {
480  unsigned int din = newCellImpl(f1, f2, f3, parm, detId);
481 
482  addValidID(detId);
483  m_dins.emplace_back(din);
484 }
485 
487  const GlobalPoint& f1, const GlobalPoint& f2, const GlobalPoint& f3, const CCGFloat* parm, const DetId& detId) {
488  unsigned int din = newCellImpl(f1, f2, f3, parm, detId);
489 
490  m_validIds.emplace_back(detId);
491  m_dins.emplace_back(din);
492 }
493 
495  // Modify the RawPtr class
496  const CaloCellGeometry* cell(nullptr);
497  if (m_hbCellVec.size() > din) {
498  cell = (&m_hbCellVec[din]);
499  } else if (m_hbCellVec.size() + m_heCellVec.size() > din) {
500  const unsigned int ind(din - m_hbCellVec.size());
501  cell = (&m_heCellVec[ind]);
502  } else if (m_hbCellVec.size() + m_heCellVec.size() + m_hoCellVec.size() > din) {
503  const unsigned int ind(din - m_hbCellVec.size() - m_heCellVec.size());
504  cell = (&m_hoCellVec[ind]);
505  } else if (m_hbCellVec.size() + m_heCellVec.size() + m_hoCellVec.size() + m_hfCellVec.size() > din) {
506  const unsigned int ind(din - m_hbCellVec.size() - m_heCellVec.size() - m_hoCellVec.size());
507  cell = (&m_hfCellVec[ind]);
508  }
509 
510  return ((nullptr == cell || nullptr == cell->param()) ? nullptr : cell);
511 }
512 
516  CaloSubdetectorGeometry::IVec& dinsVec) const {
517  tVec.reserve(m_topology.ncells() * numberOfTransformParms());
518  iVec.reserve(numberOfShapes() == 1 ? 1 : m_topology.ncells());
519  dVec.reserve(numberOfShapes() * numberOfParametersPerShape());
520  dinsVec.reserve(m_topology.ncells());
521 
522  for (const auto& pv : parVecVec()) {
523  for (float iv : pv) {
524  dVec.emplace_back(iv);
525  }
526  }
527 
528  for (auto i : m_dins) {
529  Tr3D tr;
530  auto ptr = cellGeomPtr(i);
531 
532  if (nullptr != ptr) {
533  dinsVec.emplace_back(i);
534 
535  const CCGFloat* par(ptr->param());
536 
537  unsigned int ishape(9999);
538 
539  for (unsigned int ivv(0); ivv != parVecVec().size(); ++ivv) {
540  bool ok(true);
541  const CCGFloat* pv(&(*parVecVec()[ivv].begin()));
542  for (unsigned int k(0); k != numberOfParametersPerShape(); ++k) {
543  ok = ok && (fabs(par[k] - pv[k]) < 1.e-6);
544  }
545  if (ok) {
546  ishape = ivv;
547  break;
548  }
549  }
550  assert(9999 != ishape);
551 
552  const unsigned int nn((numberOfShapes() == 1) ? (unsigned int)1 : m_dins.size());
553  if (iVec.size() < nn)
554  iVec.emplace_back(ishape);
555 
556  ptr->getTransform(tr, (Pt3DVec*)nullptr);
557 
558  if (Tr3D() == tr) { // for preshower there is no rotation
559  const GlobalPoint& gp(ptr->getPosition());
560  tr = HepGeom::Translate3D(gp.x(), gp.y(), gp.z());
561  }
562 
563  const CLHEP::Hep3Vector tt(tr.getTranslation());
564  tVec.emplace_back(tt.x());
565  tVec.emplace_back(tt.y());
566  tVec.emplace_back(tt.z());
567  if (6 == numberOfTransformParms()) {
568  const CLHEP::HepRotation rr(tr.getRotation());
569  const ROOT::Math::Transform3D rtr(
570  rr.xx(), rr.xy(), rr.xz(), tt.x(), rr.yx(), rr.yy(), rr.yz(), tt.y(), rr.zx(), rr.zy(), rr.zz(), tt.z());
572  rtr.GetRotation(ea);
573  tVec.emplace_back(ea.Phi());
574  tVec.emplace_back(ea.Theta());
575  tVec.emplace_back(ea.Psi());
576  }
577  }
578  }
579 }
580 
582  if (m_mergePosition) {
583  HcalDetId hid(id);
584  return ((DetId)(m_topology.mergedDepthDetId(hid)));
585  } else {
586  return id;
587  }
588 }
589 
590 void HcalGeometry::increaseReserve(unsigned int extra) { m_validIds.reserve(m_validIds.size() + extra); }
591 
static unsigned int numberOfBarrelAlignments()
Definition: HcalGeometry.h:67
int getMinDepth(const int &itype, const int &ieta, const int &iphi, const int &zside) const
Log< level::Info, true > LogVerbatim
bool getMergePositionFlag() const
Definition: HcalTopology.h:167
static unsigned int alignmentTransformIndexLocal(const DetId &id)
DetId correctId(const DetId &id) const
std::vector< CCGFloat > DimVec
CaloCellGeometry::Pt3DVec Pt3DVec
Definition: HcalGeometry.cc:13
GlobalPoint getBackPosition(const DetId &id) const
edm::AtomicPtrCache< std::vector< DetId > > m_heIds
Definition: HcalGeometry.h:150
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HcalGeometry.cc:76
int32_t *__restrict__ iv
static unsigned int alignmentEndcapIndexLocal(const DetId &id)
std::shared_ptr< const CaloCellGeometry > getGeometryBase(const DetId &id) const
Definition: HcalGeometry.h:122
bool set(std::unique_ptr< T > iNewValue) const
static DetId detIdFromOuterAlignmentIndex(unsigned int i)
DetId getClosestCell(const GlobalPoint &r) const override
virtual unsigned int numberOfTransformParms() const
bool valid(const DetId &id) const override
unsigned int getHESize() const
Definition: HcalTopology.h:133
const HcalTopology & m_topology
Definition: HcalGeometry.h:146
std::vector< unsigned int > IVec
edm::AtomicPtrCache< std::vector< DetId > > m_emptyIds
Definition: HcalGeometry.h:153
unsigned int detId2denseId(const DetId &id) const override
return a linear packed id
bool isHE(int etabin, int depth)
std::vector< CCGFloat > TrVec
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.
Definition: HcalGeometry.cc:89
HcalID getHCID(int subdet, int ieta, int iphi, int lay, int idepth) const
HepGeom::Transform3D Tr3D
std::vector< IdealObliquePrism > HECellVec
Definition: HcalGeometry.h:25
int zside(DetId const &)
unsigned int getHBSize() const
Definition: HcalTopology.h:132
std::vector< Pt3D > Pt3DVec
int phiBin(HcalSubdetector bc, int etaring, double phi) const
assert(be >=bs)
bool isHB(int etabin, int depth)
CaloCellGeometry::CCGFloat CCGFloat
bool m_mergePosition
Definition: HcalGeometry.h:147
unsigned int getHOSize() const
Definition: HcalTopology.h:134
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
int etaRing(HcalSubdetector subdet, double eta) const
eta and phi index from eta, phi values
CaloCellGeometry::Tr3D Tr3D
Definition: HcalGeometry.cc:14
T const * load() const
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
void newCellFast(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId)
CaloCellGeometry::Pt3DVec Pt3DVec
Definition: HcalGeometry.h:31
CaloCellGeometry::CornersVec getCorners(const DetId &id) const
Definition: TTTypes.h:54
virtual DetIdSet getCells(const GlobalPoint &r, double dR) const
Get a list of all cells within a dR of the given cell.
HBCellVec m_hbCellVec
Definition: HcalGeometry.h:156
void fillDetIds() const
Definition: HcalGeometry.cc:40
HOCellVec m_hoCellVec
Definition: HcalGeometry.h:158
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
HcalDetId idFront(const HcalDetId &id) const
Definition: HcalTopology.h:170
void getSummary(CaloSubdetectorGeometry::TrVec &trVector, CaloSubdetectorGeometry::IVec &iVector, CaloSubdetectorGeometry::DimVec &dimVector, CaloSubdetectorGeometry::IVec &dinsVector) const override
HFCellVec m_hfCellVec
Definition: HcalGeometry.h:159
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
std::vector< DetId > m_validIds
std::vector< IdealObliquePrism > HOCellVec
Definition: HcalGeometry.h:26
static unsigned int numberOfOuterAlignments()
Definition: HcalGeometry.h:73
void sortValidIds()
static DetId detIdFromForwardAlignmentIndex(unsigned int i)
static unsigned int alignmentTransformIndexGlobal(const DetId &id)
static unsigned int alignmentBarrelIndexLocal(const DetId &id)
unsigned int numberOfParametersPerShape() const override
Definition: HcalGeometry.h:43
int maxDepth(void) const
HcalSubdetector
Definition: HcalAssistant.h:31
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
virtual const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const
Get a list of valid detector ids (for the given subdetector)
static DetId detIdFromBarrelAlignmentIndex(unsigned int i)
DetId denseId2detId(unsigned int) const override
return a linear packed id
CaloCellGeometry::CCGFloat CCGFloat
Definition: HcalGeometry.cc:11
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
unsigned int getHFSize() const
Definition: HcalTopology.h:135
static DetId detIdFromLocalAlignmentIndex(unsigned int i)
CaloSubdetectorGeometry::IVec m_dins
Definition: HcalGeometry.h:154
HcalGeometry(const HcalTopology &topology)
Definition: HcalGeometry.cc:18
unsigned int ncells() const override
return a count of valid cells (for dense indexing use)
void localCorners(Pt3DVec &lc, const CCGFloat *pv, unsigned int i, Pt3D &ref)
#define M_PI
CaloSubdetectorGeometry::DetIdSet getCells(const GlobalPoint &r, double dR) const override
Get a list of all cells within a dR of the given cell.
int etaRing(HcalSubdetector bc, double abseta) const
helper methods for getClosestCell
static unsigned int alignmentOuterIndexLocal(const DetId &id)
std::vector< IdealZPrism > HFCellVec
Definition: HcalGeometry.h:27
~HcalGeometry() override
The HcalGeometry will delete all its cell geometries at destruction time.
Definition: HcalGeometry.cc:23
Definition: DetId.h:17
static unsigned int numberOfEndcapAlignments()
Definition: HcalGeometry.h:69
CaloCellGeometry::Pt3D Pt3D
Definition: HcalGeometry.cc:12
AlgebraicVector EulerAngles
Definition: Definitions.h:34
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void addValidID(const DetId &id)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
unsigned int numberOfShapes() const override
Definition: HcalGeometry.h:42
CaloCellGeometry::CornersMgr * cornersMgr()
Detector
Definition: DetId.h:24
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
Definition: IdealZPrism.cc:73
unsigned int getHxSize(const int type) const
HepGeom::Point3D< CCGFloat > Pt3D
bool isHF(int etabin, int depth)
std::vector< IdealObliquePrism > HBCellVec
Definition: HcalGeometry.h:24
HECellVec m_heCellVec
Definition: HcalGeometry.h:157
CaloCellGeometry::Tr3D Tr3D
const CaloCellGeometry * getGeometryRawPtr(uint32_t index) const override
void increaseReserve(unsigned int extra)
HcalDetId mergedDepthDetId(const HcalDetId &id) const
Definition: HcalTopology.h:166
edm::AtomicPtrCache< std::vector< DetId > > m_hfIds
Definition: HcalGeometry.h:152
GlobalPoint getPosition(const DetId &id) const
void unmergeDepthDetId(const HcalDetId &id, std::vector< HcalDetId > &ids) const
Definition: HcalTopology.h:168
void newCell(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId) override
static unsigned int numberOfForwardAlignments()
Definition: HcalGeometry.h:71
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:164
static unsigned int alignmentForwardIndexLocal(const DetId &id)
bool incrementDepth(HcalDetId &id) const
int phiBin(HcalSubdetector subdet, int etaRing, double phi) const
virtual std::shared_ptr< const CaloCellGeometry > cellGeomPtr(uint32_t index) const
static DetId detIdFromEndcapAlignmentIndex(unsigned int i)
static unsigned int alignmentBarEndForIndexLocal(const DetId &id, unsigned int nD)
edm::AtomicPtrCache< std::vector< DetId > > m_hoIds
Definition: HcalGeometry.h:151
CaloCellGeometry::Pt3D Pt3D
Definition: HcalGeometry.h:30
double etaMax(HcalSubdetector subdet) const
unsigned int newCellImpl(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId)
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
edm::AtomicPtrCache< std::vector< DetId > > m_hbIds
Definition: HcalGeometry.h:149
static unsigned int numberOfAlignments()
Definition: HcalGeometry.h:77
def move(src, dest)
Definition: eostools.py:511
const CCGFloat * param() const