CMS 3D CMS Logo

HcalDDDRecConstants.cc
Go to the documentation of this file.
2 
6 #include "CLHEP/Units/GlobalPhysicalConstants.h"
7 #include "CLHEP/Units/GlobalSystemOfUnits.h"
8 #include <algorithm>
9 #include <cmath>
10 
11 //#define EDM_ML_DEBUG
12 using namespace geant_units::operators;
13 
14 enum { kHOSizePreLS1 = 2160, kHFSizePreLS1 = 1728 };
15 
17  : hpar(hp), hcons(hc) {
18 #ifdef EDM_ML_DEBUG
19  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants::HcalDDDRecConstants (const HcalParameters* hp) constructor";
20 #endif
21  initialize();
22 }
23 
25 #ifdef EDM_ML_DEBUG
26  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants::destructed!!!";
27 #endif
28 }
29 
30 std::vector<int> HcalDDDRecConstants::getDepth(const unsigned int& eta, const bool& extra) const {
31  if (!extra) {
32  std::vector<HcalParameters::LayerItem>::const_iterator last = hpar->layerGroupEtaRec.begin();
33  for (std::vector<HcalParameters::LayerItem>::const_iterator it = hpar->layerGroupEtaRec.begin();
34  it != hpar->layerGroupEtaRec.end();
35  ++it) {
36  if (it->layer == eta + 1)
37  return it->layerGroup;
38  if (it->layer > eta + 1)
39  return last->layerGroup;
40  last = it;
41  }
42  return last->layerGroup;
43  } else {
44  std::map<int, int> layers;
46  std::vector<int> depths;
47  for (unsigned int lay = 0; lay < layers.size(); ++lay)
48  depths.emplace_back(layers[lay + 1]);
49  return depths;
50  }
51 }
52 
53 std::vector<int> HcalDDDRecConstants::getDepth(const int& det,
54  const int& phi,
55  const int& zside,
56  const unsigned int& eta) const {
57  std::map<int, int> layers;
58  hcons.ldMap()->getLayerDepth(det, eta + 1, phi, zside, layers);
59  if (layers.empty()) {
60  return getDepth(eta, false);
61  } else {
62  std::vector<int> depths;
63  for (unsigned int lay = 0; lay < layers.size(); ++lay)
64  depths.emplace_back(layers[lay + 1]);
65  return depths;
66  }
67 }
68 
69 std::vector<HcalDDDRecConstants::HcalEtaBin> HcalDDDRecConstants::getEtaBins(const int& itype) const {
70  std::vector<HcalDDDRecConstants::HcalEtaBin> bins;
71  unsigned int type = (itype == 0) ? 0 : 1;
72  HcalSubdetector subdet = HcalSubdetector(type + 1);
73  std::vector<int> phiSp;
74  HcalSubdetector subdetSp = HcalSubdetector(hcons.ldMap()->validDet(phiSp));
75  std::map<int, int> layers;
76  for (int iz = 0; iz < 2; ++iz) {
77  int zside = 2 * iz - 1;
78  for (int ieta = iEtaMin[type]; ieta <= iEtaMax[type]; ++ieta) {
79  std::vector<std::pair<int, double>> phis = getPhis(subdet, ieta);
80  std::vector<std::pair<int, double>> phiUse;
82  if (subdet == subdetSp) {
83  for (auto& phi : phis) {
84  if (std::find(phiSp.begin(), phiSp.end(), (zside * phi.first)) == phiSp.end()) {
85  phiUse.emplace_back(phi);
86  }
87  }
88  } else {
89  phiUse.insert(phiUse.end(), phis.begin(), phis.end());
90  }
91  if (!phiUse.empty())
92  getOneEtaBin(subdet, ieta, zside, phiUse, layers, false, bins);
93  }
94  }
95  if (subdetSp == subdet) {
96  for (int ieta = iEtaMin[type]; ieta <= iEtaMax[type]; ++ieta) {
97  std::vector<std::pair<int, double>> phis = getPhis(subdet, ieta);
98  for (int iz = 0; iz < 2; ++iz) {
99  int zside = 2 * iz - 1;
100  std::vector<std::pair<int, double>> phiUse;
101  for (int i : phiSp) {
102  for (auto& phi : phis) {
103  if (i == zside * phi.first) {
104  phiUse.emplace_back(phi);
105  break;
106  }
107  }
108  }
109  if (!phiUse.empty()) {
110  hcons.ldMap()->getLayerDepth(subdet, ieta, phiUse[0].first, zside, layers);
111  getOneEtaBin(subdet, ieta, zside, phiUse, layers, true, bins);
112  }
113  }
114  }
115  }
116 #ifdef EDM_ML_DEBUG
117  edm::LogVerbatim("HCalGeom") << "Prepares " << bins.size() << " eta bins for type " << type;
118  for (unsigned int i = 0; i < bins.size(); ++i) {
119  edm::LogVerbatim("HCalGeom") << "Bin[" << i << "]: Eta = (" << bins[i].ieta << ":" << bins[i].etaMin << ":"
120  << bins[i].etaMax << "), Zside = " << bins[i].zside << ", phis = ("
121  << bins[i].phis.size() << ":" << bins[i].dphi << ") and " << bins[i].layer.size()
122  << " depths (start) " << bins[i].depthStart;
123  for (unsigned int k = 0; k < bins[i].layer.size(); ++k)
124  edm::LogVerbatim("HCalGeom") << " [" << k << "] " << bins[i].layer[k].first << ":" << bins[i].layer[k].second;
125  edm::LogVerbatim("HCalGeom") << "Phi sets";
126  for (unsigned int k = 0; k < bins[i].phis.size(); ++k)
127  edm::LogVerbatim("HCalGeom") << "[" << k << "] " << bins[i].phis[k].first << ":" << bins[i].phis[k].second;
128  }
129 #endif
130  return bins;
131 }
132 
133 std::pair<double, double> HcalDDDRecConstants::getEtaPhi(const int& subdet, const int& ieta, const int& iphi) const {
134  int ietaAbs = (ieta > 0) ? ieta : -ieta;
135  double eta(0), phi(0);
136  if ((subdet == static_cast<int>(HcalBarrel)) || (subdet == static_cast<int>(HcalEndcap)) ||
137  (subdet == static_cast<int>(HcalOuter))) { // Use Eta Table
138  int unit = hcons.unitPhi(phibin[ietaAbs - 1]);
139  int kphi = (unit == 2) ? ((iphi - 1) / 2 + 1) : iphi;
140  double foff = (ietaAbs <= iEtaMax[0]) ? hpar->phioff[0] : hpar->phioff[1];
141  eta = 0.5 * (etaTable[ietaAbs - 1] + etaTable[ietaAbs]);
142  phi = foff + (kphi - 0.5) * phibin[ietaAbs - 1];
143  } else {
144  ietaAbs -= iEtaMin[2];
145  int unit = hcons.unitPhi(hpar->phitable[ietaAbs - 1]);
146  int kphi = (unit == 4) ? ((iphi - 3) / 4 + 1) : ((iphi - 1) / 2 + 1);
147  double foff = (unit > 2) ? hpar->phioff[4] : hpar->phioff[2];
148  eta = 0.5 * (hpar->etaTableHF[ietaAbs - 1] + hpar->etaTableHF[ietaAbs]);
149  phi = foff + (kphi - 0.5) * hpar->phitable[ietaAbs - 1];
150  }
151  if (ieta < 0)
152  eta = -eta;
153  if (phi > M_PI)
154  phi -= (2 * M_PI);
155 #ifdef EDM_ML_DEBUG
156  edm::LogVerbatim("HCalGeom") << "getEtaPhi: subdet|ieta|iphi " << subdet << "|" << ieta << "|" << iphi << " eta|phi "
157  << eta << "|" << phi;
158 #endif
159  return std::pair<double, double>(eta, phi);
160 }
161 
162 HcalDDDRecConstants::HcalID HcalDDDRecConstants::getHCID(int subdet, int keta, int iphi, int lay, int idepth) const {
163  int ieta = (keta > 0) ? keta : -keta;
164  int zside = (keta > 0) ? 1 : -1;
165  int eta(ieta), phi(iphi), depth(idepth);
166  if ((subdet == static_cast<int>(HcalOuter)) ||
167  ((subdet == static_cast<int>(HcalBarrel)) && (lay > maxLayerHB_ + 1))) {
168  subdet = static_cast<int>(HcalOuter);
169  depth = 4;
170  } else if (subdet == static_cast<int>(HcalBarrel) || subdet == static_cast<int>(HcalEndcap)) {
171  eta = ietaMap[ieta - 1];
172  int unit = phiUnitS[ieta - 1];
173  int phi0 = (iphi - 1) / (hpar->phigroup[eta - 1]);
174  if (unit == 2) {
175  phi0 = (iphi + 1) / 2;
176  phi0 = (phi0 - 1) / (hpar->phigroup[eta - 1]);
177  } else if (unit == 4) {
178  phi0 = (iphi + 1) / 4;
179  phi0 = (phi0 - 1) / (hpar->phigroup[eta - 1]);
180  }
181  ++phi0;
182  unit = hcons.unitPhi(phibin[eta - 1]);
183  phi = hcons.phiNumber(phi0, unit);
184  depth = hcons.findDepth(subdet, eta, phi, zside, lay - 1);
185  if (depth <= 0)
186  depth = layerGroup(eta - 1, lay - 1);
187  if (eta == iEtaMin[1]) {
188  if (subdet == static_cast<int>(HcalBarrel)) {
189  if (depth > hcons.getDepthEta16(subdet, phi, zside))
190  depth = hcons.getDepthEta16(subdet, phi, zside);
191  } else {
192  if (depth < hcons.getDepthEta16(subdet, phi, zside))
193  depth = hcons.getDepthEta16(subdet, phi, zside);
194  }
195  } else if (eta == hpar->noff[0] && lay > 1) {
196  int kphi = phi + int((hpar->phioff[3] + 0.1) / phibin[eta - 1]);
197  kphi = (kphi - 1) % 4 + 1;
198  if (kphi == 2 || kphi == 3)
199  depth = layerGroup(eta - 1, lay - 2);
200  } else if (eta == hpar->noff[1] && depth > hcons.getDepthEta29(phi, zside, 0)) {
201  eta -= hcons.getDepthEta29(phi, zside, 1);
202  }
203  }
204 #ifdef EDM_ML_DEBUG
205  edm::LogVerbatim("HCalGeom") << "getHCID: input " << subdet << ":" << ieta << ":" << iphi << ":" << idepth << ":"
206  << lay << " output " << eta << ":" << phi << ":" << depth;
207 #endif
208  return HcalDDDRecConstants::HcalID(subdet, eta, phi, depth);
209 }
210 
211 std::vector<HcalDDDRecConstants::HFCellParameters> HcalDDDRecConstants::getHFCellParameters() const {
212  std::vector<HcalDDDRecConstants::HFCellParameters> cells;
213  unsigned int nEta = hcons.getPhiTableHF().size();
214  if (maxDepth[2] > 0) {
215  for (unsigned int k = 0; k < nEta; ++k) {
216  int ieta = iEtaMin[2] + k;
217  int dphi = (int)(0.001 + hcons.getPhiTableHF()[k] / (5._deg));
218  int iphi = (dphi == 4) ? 3 : 1;
219  int nphi = 72 / dphi;
220  double rMin = hcons.getRTableHF()[nEta - k - 1] / CLHEP::cm;
221  double rMax = hcons.getRTableHF()[nEta - k] / CLHEP::cm;
223  cells.emplace_back(cell1);
225  cells.emplace_back(cell2);
226  }
227  }
228  if (maxDepth[2] > 2) {
229  if (!hcons.getIdHF2QIE().empty()) {
230  for (unsigned int k = 0; k < hcons.getIdHF2QIE().size(); ++k) {
231  int ieta = hcons.getIdHF2QIE()[k].ieta();
232  int ind = std::abs(ieta) - iEtaMin[2];
233  int dphi = (int)(0.001 + hcons.getPhiTableHF()[ind] / (5._deg));
234  int iphi = hcons.getIdHF2QIE()[k].iphi();
235  double rMin = hcons.getRTableHF()[nEta - ind - 1] / CLHEP::cm;
236  double rMax = hcons.getRTableHF()[nEta - ind] / CLHEP::cm;
238  cells.emplace_back(cell1);
239  }
240  } else {
241  for (unsigned int k = 0; k < nEta; ++k) {
242  int ieta = iEtaMin[2] + k;
243  int dphi = (int)(0.001 + hcons.getPhiTableHF()[k] / (5._deg));
244  int iphi = (dphi == 4) ? 3 : 1;
245  int nphi = 72 / dphi;
246  double rMin = hcons.getRTableHF()[nEta - k - 1] / CLHEP::cm;
247  double rMax = hcons.getRTableHF()[nEta - k] / CLHEP::cm;
249  cells.emplace_back(cell1);
251  cells.emplace_back(cell2);
252  }
253  }
254  }
255 #ifdef EDM_ML_DEBUG
256  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants returns " << cells.size() << " HF cell parameters";
257  for (unsigned int k = 0; k < cells.size(); ++k)
258  edm::LogVerbatim("HCalGeom") << "Cell[" << k << "] : (" << cells[k].ieta << ", " << cells[k].depth << ", "
259  << cells[k].firstPhi << ", " << cells[k].stepPhi << ", " << cells[k].nPhi << ", "
260  << cells[k].rMin << ", " << cells[k].rMax << ")";
261 #endif
262  return cells;
263 }
264 
265 void HcalDDDRecConstants::getLayerDepth(const int& ieta, std::map<int, int>& layers) const {
266  layers.clear();
267  for (unsigned int l = 0; l < layerGroupSize(ieta - 1); ++l) {
268  int lay = l + 1;
269  layers[lay] = layerGroup(ieta - 1, l);
270  }
271 #ifdef EDM_ML_DEBUG
272  edm::LogVerbatim("HCalGeom") << "getLayerDepth::Input " << ieta << " Output " << layers.size() << " entries";
273  for (std::map<int, int>::iterator itr = layers.begin(); itr != layers.end(); ++itr)
274  edm::LogVerbatim("HCalGeom") << " [" << itr->first << "] " << itr->second;
275 #endif
276 }
277 
278 int HcalDDDRecConstants::getLayerBack(const int& idet, const int& ieta, const int& iphi, const int& depth) const {
279  int subdet = (idet == 1) ? 1 : 2;
280  int zside = (ieta > 0) ? 1 : -1;
281  int eta = zside * ieta;
282  int layBack = hcons.ldMap()->getLayerBack(subdet, eta, iphi, zside, depth);
283  int laymax = hcons.getLastLayer(subdet, ieta);
284  if (layBack < 0 && eta <= hpar->etaMax[1]) {
285  for (unsigned int k = 0; k < layerGroupSize(eta - 1); ++k) {
286  if (depth + 1 == (int)layerGroup(eta - 1, k)) {
287  layBack = k - 1;
288  break;
289  }
290  }
291  }
292  if (layBack < 0 || layBack > laymax)
293  layBack = laymax;
294 #ifdef EDM_ML_DEBUG
295  edm::LogVerbatim("HCalGeom") << "getLayerBack::Input " << idet << ":" << ieta << ":" << iphi << ":" << depth
296  << " Output " << layBack;
297 #endif
298  return layBack;
299 }
300 
301 int HcalDDDRecConstants::getLayerFront(const int& idet, const int& ieta, const int& iphi, const int& depth) const {
302  int subdet = (idet == 1) ? 1 : 2;
303  int zside = (ieta > 0) ? 1 : -1;
304  int eta = zside * ieta;
305  int layFront = hcons.ldMap()->getLayerFront(subdet, eta, iphi, zside, depth);
306  int laymin = hcons.getFrontLayer(subdet, ieta);
307  if ((layFront < 0) || ((subdet == static_cast<int>(HcalEndcap)) && (eta == 16))) {
308  if ((subdet == static_cast<int>(HcalEndcap)) && (eta == 16)) {
309  layFront = laymin;
310  } else if (eta <= hpar->etaMax[1]) {
311  for (unsigned int k = 0; k < layerGroupSize(eta - 1); ++k) {
312  if (depth == (int)layerGroup(eta - 1, k)) {
313  if ((int)(k) >= laymin) {
314  layFront = k;
315  break;
316  }
317  }
318  }
319  }
320  } else {
321  if (layFront < laymin)
322  layFront = laymin;
323  }
324 #ifdef EDM_ML_DEBUG
325  edm::LogVerbatim("HCalGeom") << "getLayerFront::Input " << idet << ":" << ieta << ":" << iphi << ":" << depth
326  << " Output " << layFront;
327 #endif
328  return layFront;
329 }
330 
331 int HcalDDDRecConstants::getMaxDepth(const int& itype, const int& ieta, const int& iphi, const int& zside) const {
332  unsigned int type = (itype == 0) ? 0 : 1;
333  int lmax = hcons.getMaxDepth(type + 1, ieta, iphi, zside, true);
334  if (lmax < 0) {
335  unsigned int lymax = (type == 0) ? maxLayerHB_ + 1 : maxLayer_ + 1;
336  lmax = 0;
337  if (layerGroupSize(ieta - 1) > 0) {
338  if (layerGroupSize(ieta - 1) < lymax)
339  lymax = layerGroupSize(ieta - 1);
340  lmax = (int)(layerGroup(ieta - 1, lymax - 1));
341  if (type == 0 && ieta == iEtaMax[type])
342  lmax = hcons.getDepthEta16M(1);
343  if (type == 1 && ieta >= hpar->noff[1])
344  lmax = hcons.getDepthEta29M(0, false);
345  }
346  }
347 #ifdef EDM_ML_DEBUG
348  edm::LogVerbatim("HCalGeom") << "getMaxDepth::Input " << itype << ":" << ieta << ":" << iphi << ":" << zside
349  << " Output " << lmax;
350 #endif
351  return lmax;
352 }
353 
354 int HcalDDDRecConstants::getMinDepth(const int& itype, const int& ieta, const int& iphi, const int& zside) const {
355  int lmin = hcons.getMinDepth(itype + 1, ieta, iphi, zside, true);
356  if (lmin < 0) {
357  if (itype == 2) { // HFn
358  lmin = 1;
359  } else if (itype == 3) { //HO
360  lmin = maxDepth[3];
361  } else {
362  unsigned int type = (itype == 0) ? 0 : 1;
363  if (layerGroupSize(ieta - 1) > 0) {
364  if (type == 1 && ieta == iEtaMin[type])
365  lmin = hcons.getDepthEta16M(2);
366  else
367  lmin = (int)(layerGroup(ieta - 1, 0));
368  }
369  }
370  }
371  return lmin;
372 }
373 
374 std::vector<std::pair<int, double>> HcalDDDRecConstants::getPhis(const int& subdet, const int& ieta) const {
375  std::vector<std::pair<int, double>> phis;
376  int ietaAbs = (ieta > 0) ? ieta : -ieta;
377  int keta = (subdet != HcalForward) ? etaSimValu[ietaAbs - 1].first : ietaAbs;
378  std::pair<double, double> ficons = hcons.getPhiCons(subdet, keta);
379  double fioff = ficons.first;
380  double dphi = (subdet != HcalForward) ? phibin[ietaAbs - 1] : ficons.second;
381  int nphi = int((2._pi + 0.1 * dphi) / dphi);
382  int units = hcons.unitPhi(subdet, keta);
383  for (int ifi = 0; ifi < nphi; ++ifi) {
384  double phi = -fioff + (ifi + 0.5) * dphi;
385  int iphi = hcons.phiNumber(ifi + 1, units);
386  phis.emplace_back(std::pair<int, double>(iphi, phi));
387  }
388 #ifdef EDM_ML_DEBUG
389  edm::LogVerbatim("HCalGeom") << "getEtaPhi: subdet|ieta|iphi " << subdet << "|" << ieta << " with " << phis.size()
390  << " phi bins";
391  for (unsigned int k = 0; k < phis.size(); ++k)
392  edm::LogVerbatim("HCalGeom") << "[" << k << "] iphi " << phis[k].first << " phi "
393  << convertRadToDeg(phis[k].second);
394 #endif
395  return phis;
396 }
397 
398 int HcalDDDRecConstants::getPhiZOne(std::vector<std::pair<int, int>>& phiz) const {
399  phiz.clear();
400  int subdet = hcons.ldMap()->getSubdet();
401  if (subdet > 0) {
402  std::vector<int> phis = hcons.ldMap()->getPhis();
403  for (int k : phis) {
404  int zside = (k > 0) ? 1 : -1;
405  int phi = (k > 0) ? k : -k;
406  phiz.emplace_back(std::pair<int, int>(phi, zside));
407  }
408  }
409 #ifdef EDM_ML_DEBUG
410  edm::LogVerbatim("HCalGeom") << "Special RBX for detector " << subdet << " with " << phiz.size() << " phi/z bins";
411  for (unsigned int k = 0; k < phiz.size(); ++k)
412  edm::LogVerbatim("HCalGeom") << " [" << k << "] " << phiz[k].first << ":" << phiz[k].second;
413 #endif
414  return subdet;
415 }
416 
417 double HcalDDDRecConstants::getRZ(const int& subdet, const int& ieta, const int& depth) const {
418  return getRZ(subdet, ieta, 1, depth);
419 }
420 
421 double HcalDDDRecConstants::getRZ(const int& subdet, const int& ieta, const int& iphi, const int& depth) const {
422  int layf = getLayerFront(subdet, ieta, iphi, depth);
423  double rz =
424  (layf < 0) ? 0.0 : ((subdet == static_cast<int>(HcalBarrel)) ? (gconsHB[layf].first) : (gconsHE[layf].first));
425 #ifdef EDM_ML_DEBUG
426  edm::LogVerbatim("HCalGeom") << "getRZ: subdet|ieta|ipho|depth " << subdet << "|" << ieta << "|" << iphi << "|"
427  << depth << " lay|rz " << layf << "|" << rz;
428 #endif
429  return rz;
430 }
431 
432 double HcalDDDRecConstants::getRZ(const int& subdet, const int& layer) const {
433  double rz(0);
434  if (layer > 0 && layer <= (int)(layerGroupSize(0)))
435  rz = ((subdet == static_cast<int>(HcalBarrel)) ? (gconsHB[layer - 1].first) : (gconsHE[layer - 1].first));
436 #ifdef EDM_ML_DEBUG
437  edm::LogVerbatim("HCalGeom") << "getRZ: subdet|layer " << subdet << "|" << layer << " rz " << rz;
438 #endif
439  return rz;
440 }
441 
442 std::pair<double, double> HcalDDDRecConstants::getRZ(const HcalDetId& id) const {
443  int subdet = id.subdetId();
444  int ieta = id.ieta();
445  int iphi = id.iphi();
446  int depth = id.depth();
447  int zside = (subdet == static_cast<int>(HcalBarrel)) ? 1 : id.zside();
448  int layf = getLayerFront(subdet, ieta, iphi, depth);
449  double rzf = (layf < 0)
450  ? 0.0
451  : ((subdet == static_cast<int>(HcalBarrel)) ? zside * (gconsHB[layf].first - gconsHB[layf].second)
452  : zside * (gconsHE[layf].first - gconsHE[layf].second));
453  int layb = getLayerBack(subdet, ieta, iphi, depth);
454  double rzb = (layb < 0)
455  ? 0.0
456  : ((subdet == static_cast<int>(HcalBarrel)) ? zside * (gconsHB[layb].first + gconsHB[layb].second)
457  : zside * (gconsHE[layb].first + gconsHE[layb].second));
458 #ifdef EDM_ML_DEBUG
459  edm::LogVerbatim("HCalGeom") << "getRZ: subdet|ieta|ipho|depth " << subdet << "|" << ieta << "|" << iphi << "|"
460  << depth << " lay|rz (front) " << layf << "|" << rzf << " lay|rz (back) " << layb << "|"
461  << rzb;
462 #endif
463  return std::pair<double, double>(rzf, rzb);
464 }
465 
466 std::vector<HcalDDDRecConstants::HcalActiveLength> HcalDDDRecConstants::getThickActive(const int& type) const {
467  std::vector<HcalDDDRecConstants::HcalActiveLength> actives;
468  std::vector<HcalDDDRecConstants::HcalEtaBin> bins = getEtaBins(type);
469 #ifdef EDM_ML_DEBUG
470  unsigned int kount(0);
471 #endif
472  for (auto& bin : bins) {
473  int ieta = bin.ieta;
474  int zside = bin.zside;
475  int stype = (bin.phis.size() > 4) ? 0 : 1;
476  int layf = getLayerFront(type + 1, zside * ieta, bin.phis[0].first, bin.depthStart) + 1;
477  int layl = hcons.getLastLayer(type + 1, zside * ieta) + 1;
478  double eta = 0.5 * (bin.etaMin + bin.etaMax);
479  double theta = 2 * atan(exp(-eta));
480  double scale = 1.0 / ((type == 0) ? sin(theta) : cos(theta));
481  int depth = bin.depthStart;
482 #ifdef EDM_ML_DEBUG
483  edm::LogVerbatim("HCalGeom") << "Eta " << ieta << " zside " << zside << " depth " << depth << " Layers " << layf
484  << ":" << layl << ":" << bin.layer.size();
485  for (auto ll : bin.layer)
486  edm::LogVerbatim("HCalGeom") << "Layer " << ll.first << ":" << ll.second;
487  for (auto phi : bin.phis)
488  edm::LogVerbatim("HCalGeom") << "Phi " << phi.first << ":" << convertRadToDeg(phi.second);
489 #endif
490  for (unsigned int i = 0; i < bin.layer.size(); ++i) {
491  double thick(0);
492  int lmin = (type == 1 && ieta == iEtaMin[1]) ? layf : std::max(bin.layer[i].first, layf);
493  int lmax = std::min(bin.layer[i].second, layl);
494  for (int j = lmin; j <= lmax; ++j) {
495  double t = ((type == 0) ? gconsHB[j - 1].second : gconsHE[j - 1].second);
496  if ((type == 1) && (ieta <= 18))
497  t = gconsHE[j].second;
498  if (t > 0)
499  thick += t;
500  }
501 #ifdef EDM_ML_DEBUG
502  edm::LogVerbatim("HCalGeom") << "Type " << type << " L " << lmin << ":" << lmax << " T " << thick;
503 #endif
504  thick *= (2. * scale);
505  HcalDDDRecConstants::HcalActiveLength active(ieta, depth, zside, stype, zside * eta, thick);
506  for (auto phi : bin.phis)
507  active.iphis.emplace_back(phi.first);
508  actives.emplace_back(active);
509  ++depth;
510 #ifdef EDM_ML_DEBUG
511  kount++;
512  edm::LogVerbatim("HCalGeom") << "getThickActive: [" << kount << "] eta:" << active.ieta << ":" << active.eta
513  << " zside " << active.zside << " depth " << active.depth << " type " << active.stype
514  << " thick " << active.thick;
515 #endif
516  }
517  }
518  return actives;
519 }
520 
521 std::vector<HcalCellType> HcalDDDRecConstants::HcalCellTypes(HcalSubdetector subdet) const {
522  if (subdet == HcalBarrel || subdet == HcalEndcap) {
523  std::vector<HcalCellType> cells;
524  int isub = (subdet == HcalBarrel) ? 0 : 1;
525  std::vector<HcalDDDRecConstants::HcalEtaBin> etabins = getEtaBins(isub);
526  std::vector<int> missPhi;
527  for (const auto& etabin : etabins) {
528  std::vector<HcalCellType> temp;
529  std::vector<int> count;
530  std::vector<double> dmin, dmax;
531  for (unsigned int il = 0; il < etabin.layer.size(); ++il) {
532  HcalCellType cell(subdet, etabin.ieta, etabin.zside, 0, HcalCellType::HcalCell());
533  temp.emplace_back(cell);
534  count.emplace_back(0);
535  dmin.emplace_back(0);
536  dmax.emplace_back(0);
537  }
538  int ieta = etabin.ieta;
539  for (int keta = etaSimValu[ieta - 1].first; keta <= etaSimValu[ieta - 1].second; ++keta) {
540  std::vector<HcalCellType> cellsm = hcons.HcalCellTypes(subdet, keta, -1);
541  for (unsigned int il = 0; il < etabin.layer.size(); ++il) {
542  for (auto& ic : cellsm) {
543  if (ic.depthSegment() >= etabin.layer[il].first && ic.depthSegment() <= etabin.layer[il].second &&
544  ic.etaBin() == temp[il].etaBin() && ic.zside() == temp[il].zside()) {
545  if (count[il] == 0) {
546  temp[il] = ic;
547  dmin[il] = ic.depthMin();
548  dmax[il] = ic.depthMax();
549  }
550  ++count[il];
551  if (ic.depthMin() < dmin[il])
552  dmin[il] = ic.depthMin();
553  if (ic.depthMax() > dmax[il])
554  dmax[il] = ic.depthMax();
555  }
556  }
557  }
558  }
559  for (unsigned int il = 0; il < etabin.layer.size(); ++il) {
560  int depth = etabin.depthStart + (int)(il);
561  temp[il].setEta(ieta, etabin.etaMin, etabin.etaMax);
562  temp[il].setDepth(depth, dmin[il], dmax[il]);
563  double foff = (etabin.ieta <= iEtaMax[0]) ? hpar->phioff[0] : hpar->phioff[1];
564  int unit = hcons.unitPhi(etabin.dphi);
565  temp[il].setPhi(etabin.phis, missPhi, foff, etabin.dphi, unit);
566  cells.emplace_back(temp[il]);
567  }
568  }
569 #ifdef EDM_ML_DEBUG
570  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants: found " << cells.size() << " cells for sub-detector type "
571  << isub;
572  for (unsigned int ic = 0; ic < cells.size(); ++ic)
573  edm::LogVerbatim("HCalGeom") << "Cell[" << ic << "] " << cells[ic];
574 #endif
575  return cells;
576  } else {
577  return hcons.HcalCellTypes(subdet, -1, -1);
578  }
579 }
580 
582  int eta = std::abs(ieta);
583  int zside = (ieta > 0) ? 1 : -1;
584  int etamin = iEtaMax[1] - hcons.getDepthEta29(iphi, zside, 1);
585  if ((eta >= etamin) && (eta <= iEtaMax[1])) {
586  int depthMax = getMaxDepth(1, etamin, iphi, zside);
587  int depthMin = hcons.getDepthEta29(iphi, zside, 0) + 1;
588  if (depth >= depthMin && depth <= depthMax)
589  return true;
590  }
591  return false;
592 }
593 
594 std::vector<int> HcalDDDRecConstants::mergedDepthList29(int ieta, int iphi) const {
595  std::vector<int> depths;
596  int eta = std::abs(ieta);
597  int zside = (ieta > 0) ? 1 : -1;
598  int etamin = iEtaMax[1] - hcons.getDepthEta29(iphi, zside, 1);
599  if ((eta >= etamin) && (eta <= iEtaMax[1])) {
600  int depthMax = getMaxDepth(1, etamin, iphi, zside);
601  int depthMin = hcons.getDepthEta29(iphi, zside, 0) + 1;
602  depths.reserve(depthMax - depthMin + 1);
603  for (int depth = depthMin; depth <= depthMax; ++depth)
604  depths.emplace_back(depth);
605  }
606  return depths;
607 }
608 
610  if (subdet == HcalBarrel || subdet == HcalEndcap) {
611  unsigned int num = 0;
612  std::vector<HcalCellType> cellTypes = HcalCellTypes(subdet);
613  for (auto& cellType : cellTypes) {
614  num += (unsigned int)(cellType.nPhiBins());
615  }
616 #ifdef EDM_ML_DEBUG
617  edm::LogInfo("HCalGeom") << "HcalDDDRecConstants:numberOfCells " << cellTypes.size() << " " << num
618  << " for subdetector " << subdet;
619 #endif
620  return num;
621  } else {
622  return hcons.numberOfCells(subdet);
623  }
624 }
625 
626 unsigned int HcalDDDRecConstants::nCells(HcalSubdetector subdet) const {
627  if (subdet == HcalBarrel || subdet == HcalEndcap) {
628  int isub = (subdet == HcalBarrel) ? 0 : 1;
629  std::vector<HcalDDDRecConstants::HcalEtaBin> etabins = getEtaBins(isub);
630  unsigned int ncell(0);
631  for (auto& etabin : etabins) {
632  ncell += ((etabin.phis.size()) * (etabin.layer.size()));
633  }
634  return ncell;
635  } else if (subdet == HcalOuter) {
636  return kHOSizePreLS1;
637  } else if (subdet == HcalForward) {
638  return (unsigned int)(hcons.numberOfCells(subdet));
639  } else {
640  return 0;
641  }
642 }
643 
644 unsigned int HcalDDDRecConstants::nCells() const {
646 }
647 
649  std::map<HcalDetId, HcalDetId>::const_iterator itr = detIdSp_.find(id);
650  if (itr == detIdSp_.end())
651  return id;
652  else
653  return itr->second;
654 }
655 
657  HcalDetId hid(id);
658  std::map<HcalDetId, std::vector<HcalDetId>>::const_iterator itr = detIdSpR_.find(id);
659  if (itr != detIdSpR_.end())
660  hid = HcalDetId(id.subdet(), id.ieta(), id.iphi(), (itr->second)[0].depth());
661  return hid;
662 }
663 
665  HcalDetId hid(id);
666  std::map<HcalDetId, std::vector<HcalDetId>>::const_iterator itr = detIdSpR_.find(id);
667  if (itr != detIdSpR_.end())
668  hid = HcalDetId(id.subdet(), id.ieta(), id.iphi(), (itr->second).back().depth());
669  return hid;
670 }
671 
672 void HcalDDDRecConstants::unmergeDepthDetId(const HcalDetId& id, std::vector<HcalDetId>& ids) const {
673  ids.clear();
674  std::map<HcalDetId, std::vector<HcalDetId>>::const_iterator itr = detIdSpR_.find(id);
675  if (itr == detIdSpR_.end()) {
676  ids.emplace_back(id);
677  } else {
678  for (auto k : itr->second) {
679  HcalDetId hid(id.subdet(), id.ieta(), id.iphi(), k.depth());
680  ids.emplace_back(hid);
681  }
682  }
683 }
684 
685 void HcalDDDRecConstants::specialRBXHBHE(const std::vector<HcalDetId>& idsOld, std::vector<HcalDetId>& idsNew) const {
686  for (auto k : idsOld) {
687  std::map<HcalDetId, HcalDetId>::const_iterator itr = detIdSp_.find(k);
688  if (itr == detIdSp_.end())
689  idsNew.emplace_back(k);
690  else
691  idsNew.emplace_back(itr->second);
692  }
693 }
694 
695 bool HcalDDDRecConstants::specialRBXHBHE(bool tobemerged, std::vector<HcalDetId>& ids) const {
696  if (tobemerged) {
697  std::map<HcalDetId, HcalDetId>::const_iterator itr;
698  for (itr = detIdSp_.begin(); itr != detIdSp_.end(); ++itr)
699  ids.emplace_back(itr->first);
700  } else {
701  std::map<HcalDetId, std::vector<HcalDetId>>::const_iterator itr;
702  for (itr = detIdSpR_.begin(); itr != detIdSpR_.end(); ++itr)
703  ids.emplace_back(itr->first);
704  }
705  return (!ids.empty());
706 }
707 
709  int ieta,
710  int zside,
711  std::vector<std::pair<int, double>>& phis,
712  std::map<int, int>& layers,
713  bool planOne,
714  std::vector<HcalDDDRecConstants::HcalEtaBin>& bins) const {
715  unsigned int lymax = (subdet == HcalBarrel) ? maxLayerHB_ + 1 : maxLayer_ + 1;
716  int type = (subdet == HcalBarrel) ? 0 : 1;
717  double dphi = phibin[ieta - 1];
720  etabin.phis.insert(etabin.phis.end(), phis.begin(), phis.end());
721  int n = (ieta == iEtaMax[type]) ? 0 : 1;
724  etabin0.depthStart = hcons.getDepthEta29(phis[0].first, zside, 0) + 1;
725  int dstart = -1;
726  int lmin(0), lmax(0);
727 
728  std::map<int, int>::iterator itr = layers.begin();
729  if (!layers.empty()) {
730  int dep = itr->second;
731  if (subdet == HcalEndcap && ieta == iEtaMin[type])
732  dep = hcons.getDepthEta16(subdet, phis[0].first, zside);
733  unsigned lymx0 = (layers.size() > lymax) ? lymax : layers.size();
734 #ifdef EDM_ML_DEBUG
735  edm::LogVerbatim("HCalGeom") << "Eta " << ieta << ":" << hpar->noff[1] << " zside " << zside << " lymax " << lymx0
736  << ":" << lymax << " Depth " << dep << ":" << itr->second;
737  unsigned int l(0);
738  for (itr = layers.begin(); itr != layers.end(); ++itr, ++l)
739  edm::LogVerbatim("HCalGeom") << "Layer [" << l << "] " << itr->first << ":" << itr->second;
740  edm::LogVerbatim("HCalGeom") << "With " << phis.size() << " phis";
741  for (unsigned int l = 0; l < phis.size(); ++l)
742  edm::LogVerbatim("HCalGeom") << "[" << l << "] " << phis[l].first << ":" << convertRadToDeg(phis[l].second);
743 #endif
744  for (itr = layers.begin(); itr != layers.end(); ++itr) {
745  if (itr->first <= (int)(lymx0)) {
746  if (itr->second == dep) {
747  if (lmin == 0)
748  lmin = itr->first;
749  lmax = itr->first;
750  } else if (itr->second > dep) {
751  if (dstart < 0)
752  dstart = dep;
753  int lmax0 = (lmax >= lmin) ? lmax : lmin;
754  if (subdet == HcalEndcap && ieta + 1 == hpar->noff[1] && dep > hcons.getDepthEta29(phis[0].first, zside, 0)) {
755  etabin0.layer.emplace_back(std::pair<int, int>(lmin, lmax0));
756  } else {
757  etabin.layer.emplace_back(std::pair<int, int>(lmin, lmax0));
758  }
759  lmin = itr->first;
760  lmax = lmin - 1;
761  dep = itr->second;
762  }
763  if (subdet == HcalBarrel && ieta == iEtaMax[type] && dep > hcons.getDepthEta16M(1))
764  break;
765  if (subdet == HcalEndcap && ieta == hpar->noff[1] && dep > hcons.getDepthEta29M(0, planOne)) {
766  lmax = lymx0;
767  break;
768  }
769  if (itr->first == (int)(lymx0))
770  lmax = lymx0;
771  }
772  }
773  if (lmax >= lmin) {
774  if (ieta + 1 == hpar->noff[1]) {
775  etabin0.layer.emplace_back(std::pair<int, int>(lmin, lmax));
776  etabin0.phis.insert(etabin0.phis.end(), phis.begin(), phis.end());
777  bins.emplace_back(etabin0);
778 #ifdef EDM_ML_DEBUG
779  edm::LogVerbatim("HCalGeom") << "etabin0: dStatrt " << etabin0.depthStart << " layers " << etabin0.layer.size()
780  << ":" << lmin << ":" << lmax << " phis " << phis.size();
781  for (unsigned int k = 0; k < etabin0.layer.size(); ++k)
782  edm::LogVerbatim("HCalGeom") << " [" << k << "] " << etabin0.layer[k].first << ":" << etabin0.layer[k].second;
783 #endif
784  } else if (ieta == hpar->noff[1]) {
785  } else {
786  etabin.layer.emplace_back(std::pair<int, int>(lmin, lmax));
787  if (dstart < 0)
788  dstart = dep;
789  }
790  }
791  }
792  etabin.depthStart = dstart;
793  bins.emplace_back(etabin);
794 #ifdef EDM_ML_DEBUG
795  edm::LogVerbatim("HCalGeom") << "etabin: dStatrt " << etabin.depthStart << " layers " << etabin.layer.size() << ":"
796  << lmin << ":" << lmax << " phis " << etabin.phis.size();
797  for (unsigned int k = 0; k < etabin.layer.size(); ++k)
798  edm::LogVerbatim("HCalGeom") << "[" << k << "] " << etabin.layer[k].first << ":" << etabin.layer[k].second;
799 #endif
800 }
801 
803  //Eta grouping
804  int nEta = (int)(hpar->etagroup.size());
805  if (nEta != (int)(hpar->phigroup.size())) {
806  edm::LogError("HCalGeom") << "HcalDDDRecConstants: sizes of the vectors "
807  << " etaGroup (" << nEta << ") and phiGroup (" << hpar->phigroup.size()
808  << ") do not match";
809  throw cms::Exception("DDException") << "HcalDDDRecConstants: inconsistent array sizes" << nEta << ":"
810  << hpar->phigroup.size();
811  }
812 
813  // First eta table
814  iEtaMin = hpar->etaMin;
815  iEtaMax = hpar->etaMax;
816  etaTable.clear();
817  ietaMap.clear();
818  etaSimValu.clear();
819  int ieta(0), ietaHB(0), ietaHE(0), ietaHEM(0);
820  etaTable.emplace_back(hpar->etaTable[ieta]);
821  for (int i = 0; i < nEta; ++i) {
822  int ef = ieta + 1;
823  ieta += (hpar->etagroup[i]);
824  if (ieta >= (int)(hpar->etaTable.size())) {
825  edm::LogError("HCalGeom") << "HcalDDDRecConstants: Going beyond the array boundary " << hpar->etaTable.size()
826  << " at index " << i << " of etaTable from SimConstant";
827  throw cms::Exception("DDException")
828  << "HcalDDDRecConstants: Going beyond the array boundary " << hpar->etaTable.size() << " at index " << i
829  << " of etaTable from SimConstant";
830  } else {
831  etaTable.emplace_back(hpar->etaTable[ieta]);
832  etaSimValu.emplace_back(std::pair<int, int>(ef, ieta));
833  }
834  for (int k = 0; k < (hpar->etagroup[i]); ++k)
835  ietaMap.emplace_back(i + 1);
836  if (ieta <= hpar->etaMax[0])
837  ietaHB = i + 1;
838  if (ieta <= hpar->etaMin[1])
839  ietaHE = i + 1;
840  if (ieta <= hpar->etaMax[1])
841  ietaHEM = i + 1;
842  }
843  iEtaMin[1] = ietaHE;
844  iEtaMax[0] = ietaHB;
845  iEtaMax[1] = ietaHEM;
846 
847  // Then Phi bins
848  nPhiBins.clear();
849  for (unsigned int k = 0; k < 4; ++k)
850  nPhiBins.emplace_back(0);
851  ieta = 0;
852  phibin.clear();
853  phiUnitS.clear();
854  for (int i = 0; i < nEta; ++i) {
855  double dphi = (hpar->phigroup[i]) * (hpar->phibin[ieta]);
856  phibin.emplace_back(dphi);
857  int nphi = (int)((2._pi + 0.001) / dphi);
858  if (ieta <= iEtaMax[0]) {
859  if (nphi > nPhiBins[0])
860  nPhiBins[3] = nPhiBins[0] = nphi;
861  }
862  if (ieta >= iEtaMin[1]) {
863  if (nphi > nPhiBins[1])
864  nPhiBins[1] = nphi;
865  }
866  ieta += (hpar->etagroup[i]);
867  }
868  for (unsigned int i = 1; i < hpar->etaTable.size(); ++i) {
869  int unit = hcons.unitPhi(hpar->phibin[i - 1]);
870  phiUnitS.emplace_back(unit);
871  }
872  for (double i : hpar->phitable) {
873  int nphi = (int)((2._pi + 0.001) / i);
874  if (nphi > nPhiBins[2])
875  nPhiBins[2] = nphi;
876  }
877 #ifdef EDM_ML_DEBUG
878  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants: Modified eta/deltaphi table for " << nEta << " bins";
879  for (int i = 0; i < nEta; ++i)
880  edm::LogVerbatim("HCalGeom") << "Eta[" << i << "] = " << etaTable[i] << ":" << etaTable[i + 1] << ":"
881  << etaSimValu[i].first << ":" << etaSimValu[i].second << " PhiBin[" << i
882  << "] = " << convertRadToDeg(phibin[i]);
883  for (unsigned int i = 0; i < phiUnitS.size(); ++i)
884  edm::LogVerbatim("HCalGeom") << " PhiUnitS[" << i << "] = " << phiUnitS[i];
885  for (unsigned int i = 0; i < nPhiBins.size(); ++i)
886  edm::LogVerbatim("HCalGeom") << " nPhiBins[" << i << "] = " << nPhiBins[i];
887  for (unsigned int i = 0; i < hpar->etaTableHF.size(); ++i)
888  edm::LogVerbatim("HCalGeom") << " EtaTableHF[" << i << "] = " << hpar->etaTableHF[i];
889  for (unsigned int i = 0; i < hpar->phitable.size(); ++i)
890  edm::LogVerbatim("HCalGeom") << " PhiBinHF[" << i << "] = " << hpar->phitable[i];
891 #endif
892 
893  //Now the depths
895  maxDepth[0] = maxDepth[1] = 0;
896  for (int i = 0; i < nEta; ++i) {
897  unsigned int imx = layerGroupSize(i);
898  int laymax = (imx > 0) ? layerGroup(i, imx - 1) : 0;
899  if (i < iEtaMax[0]) {
900  int laymax0 = (imx > 16) ? layerGroup(i, 16) : laymax;
901  if (i + 1 == iEtaMax[0])
902  laymax0 = hcons.getDepthEta16M(1);
903 #ifdef EDM_ML_DEBUG
904  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:HB " << i << " " << imx << " " << laymax << " " << laymax0;
905 #endif
906  if (maxDepth[0] < laymax0)
907  maxDepth[0] = laymax0;
908  }
909  if (i >= iEtaMin[1] - 1 && i < iEtaMax[1]) {
910 #ifdef EDM_ML_DEBUG
911  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:HE " << i << " " << imx << " " << laymax;
912 #endif
913  if (maxDepth[1] < laymax)
914  maxDepth[1] = laymax;
915  }
916  }
917 #ifdef EDM_ML_DEBUG
918  for (int i = 0; i < 4; ++i)
919  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:Detector Type[" << i << "] iEta " << iEtaMin[i] << ":"
920  << iEtaMax[i] << " MaxDepth " << maxDepth[i];
921 #endif
922 
923  //Now the geometry constants
924  nModule[0] = hpar->modHB[0];
925  nHalves[0] = hpar->modHB[1];
926  for (unsigned int i = 0; i < hpar->rHB.size(); ++i) {
927  gconsHB.emplace_back(std::pair<double, double>(hpar->rHB[i] / CLHEP::cm, hpar->drHB[i] / CLHEP::cm));
928  }
929 #ifdef EDM_ML_DEBUG
930  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:HB with " << nModule[0] << " modules and " << nHalves[0]
931  << " halves and " << gconsHB.size() << " layers";
932  for (unsigned int i = 0; i < gconsHB.size(); ++i)
933  edm::LogVerbatim("HCalGeom") << "rHB[" << i << "] = " << gconsHB[i].first << " +- " << gconsHB[i].second;
934 #endif
935  nModule[1] = hpar->modHE[0];
936  nHalves[1] = hpar->modHE[1];
937  for (unsigned int i = 0; i < hpar->zHE.size(); ++i) {
938  gconsHE.emplace_back(std::pair<double, double>(hpar->zHE[i] / CLHEP::cm, hpar->dzHE[i] / CLHEP::cm));
939  }
940 #ifdef EDM_ML_DEBUG
941  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:HE with " << nModule[1] << " modules and " << nHalves[1]
942  << " halves and " << gconsHE.size() << " layers";
943  for (unsigned int i = 0; i < gconsHE.size(); ++i)
944  edm::LogVerbatim("HCalGeom") << "zHE[" << i << "] = " << gconsHE[i].first << " +- " << gconsHE[i].second;
945 #endif
946 
947  //Special RBX
949  if (depthMaxSp_.first == 0) {
950  depthMaxSp_ = depthMaxDf_ = std::pair<int, int>(2, maxDepth[1]);
951  } else if (depthMaxSp_.first == 1) {
952  depthMaxDf_ = std::pair<int, int>(1, maxDepth[0]);
953  if (depthMaxSp_.second > maxDepth[0])
954  maxDepth[0] = depthMaxSp_.second;
955  } else {
956  depthMaxDf_ = std::pair<int, int>(2, maxDepth[1]);
957  if (depthMaxSp_.second > maxDepth[1])
958  maxDepth[1] = depthMaxSp_.second;
959  }
960 #ifdef EDM_ML_DEBUG
961  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:Detector type and maximum depth for all RBX "
962  << depthMaxDf_.first << ":" << depthMaxDf_.second << " and for special RBX "
963  << depthMaxSp_.first << ":" << depthMaxSp_.second;
964 #endif
965 
966  //Map of special DetId's
967  std::vector<int> phis;
969  detIdSp_.clear();
970  detIdSpR_.clear();
971  if ((subdet == HcalBarrel) || (subdet == HcalEndcap)) {
972  int phi = (phis[0] > 0) ? phis[0] : -phis[0];
973  int zside = (phis[0] > 0) ? 1 : -1;
974  int lymax = (subdet == HcalBarrel) ? maxLayerHB_ + 1 : maxLayer_ + 1;
975  std::pair<int, int> etas = hcons.ldMap()->validEta();
976  for (int eta = etas.first; eta <= etas.second; ++eta) {
977  std::map<int, std::pair<int, int>> oldDep;
978  int depth(0);
979  int lmin = layerGroup(eta - 1, 0);
980  for (int lay = 0; lay < lymax; ++lay) {
981  int depc = layerGroup(eta - 1, lay);
982  if (depth != depc) {
983  if (depth != 0)
984  oldDep[depth] = std::pair<int, int>(lmin, lay - 1);
985  depth = depc;
986  lmin = lay;
987  }
988  }
989  if (depth != 0)
990  oldDep[depth] = std::pair<int, int>(lmin, lymax - 1);
991 #ifdef EDM_ML_DEBUG
992  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:Eta|Phi|Zside " << eta << ":" << phi << ":" << zside
993  << " with " << oldDep.size() << " old Depths";
994  unsigned int kk(0);
995  for (std::map<int, std::pair<int, int>>::const_iterator itr = oldDep.begin(); itr != oldDep.end(); ++itr, ++kk)
996  edm::LogVerbatim("HCalGeom") << "[" << kk << "] " << itr->first << " --> " << itr->second.first << ":"
997  << itr->second.second;
998 #endif
999  std::pair<int, int> depths = hcons.ldMap()->getDepths(eta);
1000  for (int ndepth = depths.first; ndepth <= depths.second; ++ndepth) {
1001  bool flag = ((subdet == HcalBarrel && eta == iEtaMax[0] && ndepth > hcons.getDepthEta16(subdet, phi, zside)) ||
1002  (subdet == HcalEndcap && eta == iEtaMin[1] && ndepth < hcons.getDepthEta16(subdet, phi, zside)));
1003  if (!flag) {
1004  std::vector<int> count(oldDep.size(), 0);
1005  int layFront = hcons.ldMap()->getLayerFront(subdet, eta, phi, zside, ndepth);
1006  int layBack = hcons.ldMap()->getLayerBack(subdet, eta, phi, zside, ndepth);
1007  for (int lay = layFront; lay <= layBack; ++lay) {
1008  unsigned int l(0);
1009  for (std::map<int, std::pair<int, int>>::iterator itr = oldDep.begin(); itr != oldDep.end(); ++itr, ++l) {
1010  if (lay >= (itr->second).first && lay <= (itr->second).second) {
1011  ++count[l];
1012  break;
1013  }
1014  }
1015  }
1016  int odepth(0), maxlay(0);
1017  unsigned int l(0);
1018  for (std::map<int, std::pair<int, int>>::iterator itr = oldDep.begin(); itr != oldDep.end(); ++itr, ++l) {
1019  if (count[l] > maxlay) {
1020  odepth = itr->first;
1021  maxlay = count[l];
1022  }
1023  }
1024 #ifdef EDM_ML_DEBUG
1025  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:New Depth " << ndepth << " old Depth " << odepth
1026  << " max " << maxlay;
1027 #endif
1028  for (int k : phis) {
1029  zside = (k > 0) ? 1 : -1;
1030  phi = (k > 0) ? k : -k;
1031  if (subdet == HcalEndcap && eta == hpar->noff[1] && ndepth > hcons.getDepthEta29M(0, true))
1032  break;
1033  HcalDetId newId(subdet, zside * eta, phi, ndepth);
1034  HcalDetId oldId(subdet, zside * eta, phi, odepth);
1035  detIdSp_[newId] = oldId;
1036  std::vector<HcalDetId> ids;
1037  std::map<HcalDetId, std::vector<HcalDetId>>::iterator itr = detIdSpR_.find(oldId);
1038  if (itr != detIdSpR_.end())
1039  ids = itr->second;
1040  ids.emplace_back(newId);
1041  detIdSpR_[oldId] = ids;
1042  }
1043  }
1044  }
1045  }
1046 #ifdef EDM_ML_DEBUG
1047  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:Map for merging new channels to old channel"
1048  << " IDs with " << detIdSp_.size() << " entries";
1049  int l(0);
1050  for (auto itr : detIdSp_) {
1051  edm::LogVerbatim("HCalGeom") << "[" << l << "] Special " << itr.first << " Standard " << itr.second;
1052  ++l;
1053  }
1054  edm::LogVerbatim("HCalGeom") << "HcalDDDRecConstants:Reverse Map for mapping old to new IDs with "
1055  << detIdSpR_.size() << " entries";
1056  l = 0;
1057  for (auto itr : detIdSpR_) {
1058  edm::LogVerbatim("HCalGeom") << "[" << l << "] Standard " << itr.first << " Special";
1059  for (auto itr1 : itr.second)
1060  edm::LogVerbatim("HCalGeom") << "ID " << (itr1);
1061  ++l;
1062  }
1063 #endif
1064  }
1065 }
1066 
1067 unsigned int HcalDDDRecConstants::layerGroupSize(int eta) const {
1068  unsigned int k = 0;
1069  for (auto const& it : hpar->layerGroupEtaRec) {
1070  if (it.layer == (unsigned int)(eta + 1)) {
1071  return it.layerGroup.size();
1072  }
1073  if (it.layer > (unsigned int)(eta + 1))
1074  break;
1075  k = it.layerGroup.size();
1076  }
1077  return k;
1078 }
1079 
1080 unsigned int HcalDDDRecConstants::layerGroup(int eta, int i) const {
1081  unsigned int k = 0;
1082  for (auto const& it : hpar->layerGroupEtaRec) {
1083  if (it.layer == (unsigned int)(eta + 1)) {
1084  return it.layerGroup.at(i);
1085  }
1086  if (it.layer > (unsigned int)(eta + 1))
1087  break;
1088  k = it.layerGroup.at(i);
1089  }
1090  return k;
1091 }
int getMinDepth(const int &itype, const int &ieta, const int &iphi, const int &zside) const
Log< level::Info, true > LogVerbatim
const int nphi
std::map< HcalDetId, HcalDetId > detIdSp_
double getRZ(const int &subdet, const int &ieta, const int &depth) const
std::vector< double > etaTable
std::vector< int > iEtaMin
std::vector< int > etagroup
int getFrontLayer(const int &det, const int &eta) const
std::pair< int, int > depthMaxSp_
HcalDetId mergedDepthDetId(const HcalDetId &id) const
std::vector< double > rHB
void getOneEtaBin(HcalSubdetector subdet, int ieta, int zside, std::vector< std::pair< int, double >> &phis, std::map< int, int > &layers, bool planOne, std::vector< HcalDDDRecConstants::HcalEtaBin > &bins) const
int getDepthEta29M(const int &i, const bool &planOne) const
std::vector< std::pair< int, double > > getPhis(const int &subdet, const int &ieta) const
constexpr NumType convertRadToDeg(NumType radians)
Definition: angle_units.h:21
bool mergedDepthList29(int ieta, int iphi, int depth) const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
int getLayerFront(const int subdet, const int ieta, const int iphi, const int zside, const int depth) const
int getDepthEta29(const int &phi, const int &zside, const int &i) const
static const int maxLayer_
HcalID getHCID(int subdet, int ieta, int iphi, int lay, int idepth) const
HcalDetId idBack(const HcalDetId &id) const
std::vector< int > maxDepth
int zside(DetId const &)
std::pair< double, double > getPhiCons(const int &det, const int &ieta) const
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
unsigned int nCells() const
std::vector< HcalEtaBin > getEtaBins(const int &itype) const
std::vector< double > etaTableHF
int getDepthEta16M(const int &det) const
unsigned int numberOfCells(const HcalSubdetector &) const
std::vector< int > etaMax
constexpr std::array< uint8_t, layerIndexSize > layer
std::vector< int > phiUnitS
std::vector< int > modHB
int phiNumber(const int &phi, const int &unit) const
std::vector< HFCellParameters > getHFCellParameters() const
U second(std::pair< T, U > const &p)
std::vector< std::pair< double, double > > gconsHE
const HcalParameters * hpar
std::pair< double, double > getEtaPhi(const int &subdet, const int &ieta, const int &iphi) const
void specialRBXHBHE(const std::vector< HcalDetId > &, std::vector< HcalDetId > &) const
const std::vector< double > & getRTableHF() const
std::pair< int, int > getMaxDepthDet(const int &i) const
std::pair< int, int > depthMaxDf_
std::vector< double > zHE
std::vector< double > phibin
std::pair< int, int > getDepths(const int eta) const
std::vector< double > dzHE
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< int > iEtaMax
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< HcalCellType > HcalCellTypes() const
std::vector< int > ietaMap
std::vector< int > nPhiBins
const std::vector< double > & getPhiTableHF() const
const std::vector< int > & getPhis() const
const HcalLayerDepthMap * ldMap() const
std::vector< int > modHE
std::vector< std::pair< int, int > > layer
std::vector< std::pair< int, int > > etaSimValu
int getMaxDepth(const int &type) const
int findDepth(const int &det, const int &eta, const int &phi, const int &zside, const int &lay) const
unsigned int layerGroupSize(int eta) const
Basic3DVector unit() const
int unitPhi(const int &det, const int &etaR) const
std::vector< HcalCellType > HcalCellTypes(HcalSubdetector) const
std::vector< std::pair< double, double > > gconsHB
int getLayerBack(const int subdet, const int ieta, const int iphi, const int zside, const int depth) const
#define M_PI
Log< level::Info, false > LogInfo
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ cells
std::vector< double > phioff
int getPhiZOne(std::vector< std::pair< int, int >> &phiz) const
std::vector< double > etaTable
const std::vector< HcalDetId > & getIdHF2QIE() const
HcalDDDRecConstants(const HcalParameters *hp, const HcalDDDSimConstants &hc)
unsigned int numberOfCells(HcalSubdetector) const
void getLayerDepth(const int &ieta, std::map< int, int > &layers) const
int getMinDepth(const int &det, const int &eta, const int &phi, const int &zside, const bool &partialOnly) const
std::vector< double > phitable
HcalDetId idFront(const HcalDetId &id) const
TString units(TString variable, Char_t axis)
std::vector< double > phibin
std::vector< LayerItem > layerGroupEtaRec
int getLayerFront(const int &det, const int &eta, const int &phi, const int &depth) const
std::vector< double > drHB
const HcalDDDSimConstants & hcons
std::vector< HcalActiveLength > getThickActive(const int &type) const
unsigned int layerGroup(int eta, int i) const
std::vector< int > noff
std::vector< std::pair< int, double > > phis
const int ndepth
void unmergeDepthDetId(const HcalDetId &id, std::vector< HcalDetId > &ids) const
void getLayerDepth(const int subdet, const int ieta, const int iphi, const int zside, std::map< int, int > &layers) const
std::vector< int > maxDepth
std::vector< int > phigroup
std::pair< int, int > validEta() const
std::vector< int > getDepth(const int &det, const int &phi, const int &zside, const unsigned int &eta) const
int getMaxDepth(const int &type) const
int getLastLayer(const int &det, const int &eta) const
Geom::Theta< T > theta() const
static const int maxLayerHB_
int validDet(std::vector< int > &phis) const
int getDepthEta16(const int &det, const int &phi, const int &zside) const
int getLayerBack(const int &det, const int &eta, const int &phi, const int &depth) const
std::map< HcalDetId, std::vector< HcalDetId > > detIdSpR_
std::vector< int > etaMin