CMS 3D CMS Logo

HGCalDDDConstants.cc
Go to the documentation of this file.
2 
13 
14 #include <algorithm>
15 #include <functional>
16 #include <numeric>
17 
18 //#define EDM_ML_DEBUG
19 using namespace geant_units::operators;
20 
21 static const int maxType = 2;
22 static const int minType = 0;
23 
25  const std::string& name)
26  : hgpar_(hp), sqrt3_(std::sqrt(3.0)) {
27  mode_ = hgpar_->mode_;
33  std::cos(30._deg));
34  hexside_ = 2.0 * rmax_ * tan30deg_;
35 #ifdef EDM_ML_DEBUG
36  edm::LogVerbatim("HGCalGeom")
37  << "rmax_ " << rmax_ << ":" << hexside_ << " CellSize "
40 #endif
41  }
42  // init maps and constants
43  modHalf_ = 0;
45  for (int simreco = 0; simreco < 2; ++simreco) {
46  tot_layers_[simreco] = layersInit((bool)simreco);
47  max_modules_layer_[simreco].resize(tot_layers_[simreco] + 1);
48  for (unsigned int layer = 1; layer <= tot_layers_[simreco]; ++layer) {
49  max_modules_layer_[simreco][layer] = modulesInit(layer, (bool)simreco);
50  if (simreco == 1) {
51  modHalf_ += max_modules_layer_[simreco][layer];
54 #ifdef EDM_ML_DEBUG
55  edm::LogVerbatim("HGCalGeom") << "Layer " << layer << " with "
56  << max_modules_layer_[simreco][layer]
57  << ":" << modHalf_ << " modules";
58 #endif
59  }
60  }
61  }
62  tot_wafers_ = wafers();
63 
64 #ifdef EDM_ML_DEBUG
65  edm::LogVerbatim("HGCalGeom")
66  << "HGCalDDDConstants initialized for " << name << " with "
67  << layers(false) << ":" << layers(true) << " layers, " << wafers() << ":"
68  << 2 * modHalf_ << " wafers with maximum " << maxWafersPerLayer_
69  << " per layer and "
70  << "maximum of " << maxCells(false) << ":" << maxCells(true) << " cells";
71 #endif
76  int wminT(9999999), wmaxT(-9999999), kount1(0), kount2(0);
77  for (unsigned int i = 0; i < getTrFormN(); ++i) {
78  int lay0 = getTrForm(i).lay;
79  int wmin(9999999), wmax(-9999999), kount(0);
80  for (int wafer = 0; wafer < sectors(); ++wafer) {
81  bool waferIn = waferInLayer(wafer, lay0, true);
84  int kndx = HGCalWaferIndex::waferIndex(
87  waferIn_[kndx] = waferIn;
88  }
89  if (waferIn) {
90  int waferU =
93  ? wafer
95  if (waferU < wmin) wmin = waferU;
96  if (waferU > wmax) wmax = waferU;
97  ++kount;
98  }
99  }
100  if (wminT > wmin) wminT = wmin;
101  if (wmaxT < wmax) wmaxT = wmax;
102  if (kount1 < kount) kount1 = kount;
103  kount2 += kount;
104 #ifdef EDM_ML_DEBUG
105  int lay1 = getIndex(lay0, true).first;
106  edm::LogVerbatim("HGCalGeom")
107  << "Index " << i << " Layer " << lay0 << ":" << lay1 << " Wafer "
108  << wmin << ":" << wmax << ":" << kount;
109 #endif
110  HGCWaferParam a1{{wmin, wmax, kount}};
111  waferLayer_[lay0] = a1;
112  }
113  waferMax_ = std::array<int, 4>{{wminT, wmaxT, kount1, kount2}};
114 #ifdef EDM_ML_DEBUG
115  edm::LogVerbatim("HGCalGeom")
116  << "Overall wafer statistics: " << wminT << ":" << wmaxT << ":"
117  << kount1 << ":" << kount2;
118 #endif
119  }
120 }
121 
123 
124 std::pair<int, int> HGCalDDDConstants::assignCell(float x, float y, int lay,
125  int subSec, bool reco) const {
126  const auto& index = getIndex(lay, reco);
127  if (index.first < 0) return std::make_pair(-1, -1);
130  float xx = (reco) ? x : HGCalParameters::k_ScaleFromDDD * x;
131  float yy = (reco) ? y : HGCalParameters::k_ScaleFromDDD * y;
132 
133  // First the wafer
134  int wafer = cellHex(xx, yy, rmax_, hgpar_->waferPosX_, hgpar_->waferPosY_);
135  if (wafer < 0 || wafer >= (int)(hgpar_->waferTypeT_.size())) {
136  edm::LogWarning("HGCalGeom")
137  << "Wafer no. out of bound for " << wafer << ":"
138  << (hgpar_->waferTypeT_).size() << ":" << (hgpar_->waferPosX_).size()
139  << ":" << (hgpar_->waferPosY_).size() << " ***** ERROR *****";
140  return std::make_pair(-1, -1);
141  } else {
142  // Now the cell
143  xx -= hgpar_->waferPosX_[wafer];
144  yy -= hgpar_->waferPosY_[wafer];
145  if (hgpar_->waferTypeT_[wafer] == 1)
146  return std::make_pair(wafer,
147  cellHex(xx, yy,
149  hgpar_->cellSize_[0],
151  else
152  return std::make_pair(
153  wafer, cellHex(xx, yy,
155  hgpar_->cellSize_[1],
157  }
158  } else {
159  return std::make_pair(-1, -1);
160  }
161 }
162 
163 std::array<int, 5> HGCalDDDConstants::assignCellHex(float x, float y, int lay,
164  bool reco) const {
165  int waferU(0), waferV(0), waferType(-1), cellU(0), cellV(0);
168  double xx = (reco) ? HGCalParameters::k_ScaleToDDD * x : x;
169  double yy = (reco) ? HGCalParameters::k_ScaleToDDD * y : y;
170  double wt(1.0);
171  waferFromPosition(xx, yy, lay, waferU, waferV, cellU, cellV, waferType, wt);
172  }
173  return std::array<int, 5>{{waferU, waferV, waferType, cellU, cellV}};
174 }
175 
176 std::array<int, 3> HGCalDDDConstants::assignCellTrap(float x, float y, float z,
177  int layer,
178  bool reco) const {
179  int irad(-1), iphi(-1), type(-1);
180  const auto& indx = getIndex(layer, reco);
181  if (indx.first < 0) return std::array<int, 3>{{irad, iphi, type}};
182  double xx = (z > 0) ? x : -x;
183  double r =
184  (reco ? std::sqrt(x * x + y * y)
185  : HGCalParameters::k_ScaleFromDDD * std::sqrt(x * x + y * y));
186  double phi = (r == 0. ? 0. : std::atan2(y, xx));
187  if (phi < 0) phi += (2.0 * M_PI);
188  type = hgpar_->scintType(layer);
189  auto ir = std::lower_bound(hgpar_->radiusLayer_[type].begin(),
190  hgpar_->radiusLayer_[type].end(), r);
191  irad = (int)(ir - hgpar_->radiusLayer_[type].begin());
192  irad = std::min(std::max(irad, hgpar_->iradMinBH_[indx.first]),
193  hgpar_->iradMaxBH_[indx.first]);
194  iphi = 1 + (int)(phi / indx.second);
195 #ifdef EDM_ML_DEBUG
196  edm::LogVerbatim("HGCalGeom")
197  << "assignCellTrap Input " << x << ":" << y << ":" << z << ":" << layer
198  << ":" << reco << " x|r " << xx << ":" << r << " phi " << phi << " o/p "
199  << irad << ":" << iphi << ":" << type;
200 #endif
201  return std::array<int, 3>{{irad, iphi, type}};
202 }
203 
204 std::pair<double, double> HGCalDDDConstants::cellEtaPhiTrap(int type,
205  int irad) const {
206  double dr(0), df(0);
208  double r = 0.5*((hgpar_->radiusLayer_[type][irad-1] +
209  hgpar_->radiusLayer_[type][irad]));
210  dr = (hgpar_->radiusLayer_[type][irad] -
211  hgpar_->radiusLayer_[type][irad-1]);
212  df = r * hgpar_->cellSize_[type];
213  }
214  return std::make_pair(dr,df);
215 }
216 
217 bool HGCalDDDConstants::cellInLayer(int waferU, int waferV, int cellU,
218  int cellV, int lay, bool reco) const {
219  const auto& indx = getIndex(lay, true);
220  if (indx.first >= 0) {
225  const auto& xy = (((mode_ == HGCalGeometryMode::Hexagon8) ||
227  ? locateCell(lay, waferU, waferV, cellU, cellV,
228  reco, true, false)
229  : locateCell(cellU, lay, waferU, reco));
230  double rpos = sqrt(xy.first * xy.first + xy.second * xy.second);
231  return ((rpos >= hgpar_->rMinLayHex_[indx.first]) &&
232  (rpos <= hgpar_->rMaxLayHex_[indx.first]));
233  } else {
234  return true;
235  }
236  } else {
237  return false;
238  }
239 }
240 
241 double HGCalDDDConstants::cellThickness(int layer, int waferU,
242  int waferV) const {
243  double thick(-1);
246  auto itr = hgpar_->typesInLayers_.find(
247  HGCalWaferIndex::waferIndex(layer, waferU, waferV));
248  int type = ((itr == hgpar_->typesInLayers_.end()
249  ? maxType
250  : hgpar_->waferTypeL_[itr->second]));
251  thick = 10000.0 * hgpar_->cellThickness_[type]; // cm to micron
252  } else if ((mode_ == HGCalGeometryMode::Hexagon) ||
254  int type = (((waferU >= 0) && (waferU < (int)(hgpar_->waferTypeL_.size())))
255  ? hgpar_->waferTypeL_[waferU]
256  : minType);
257  thick = 100.0 * type; // type = 1,2,3 for 100,200,300 micron
258  }
259  return thick;
260 }
261 
263  int indx = (((mode_ == HGCalGeometryMode::Hexagon8) ||
265  ? ((type >= 1) ? 1 : 0)
266  : ((type == 1) ? 1 : 0));
267  double cell =
269  ? 0.5 * hgpar_->cellSize_[indx]
271  return cell;
272 }
273 
275  int cellV) const {
276  // type=0: in the middle; 1..6: the edges clocwise from bottom left;
277  // =11..16: the corners clockwise from bottom
278  int N = (type == 0) ? hgpar_->nCellsFine_ : hgpar_->nCellsCoarse_;
279  if (cellU == 0) {
280  if (cellV == 0)
282  else if (cellV - cellU == N - 1)
284  else
286  } else if (cellV == 0) {
287  if (cellU - cellV == N)
289  else
291  } else if (cellU - cellV == N) {
292  if (cellU == 2 * N - 1)
294  else
296  } else if (cellU == 2 * N - 1) {
297  if (cellV == 2 * N - 1)
299  else
301  } else if (cellV == 2 * N - 1) {
302  if (cellV - cellU == N - 1)
304  else
306  } else if (cellV - cellU == N - 1) {
308  } else if ((cellU > 2 * N - 1) || (cellV > 2 * N - 1) ||
309  (cellV >= (cellU + N)) || (cellU > (cellV + N))) {
311  } else {
313  }
314 }
315 
316 double HGCalDDDConstants::distFromEdgeHex(double x, double y, double z) const {
317  // Assming the point is within a hexagonal plane of the wafer, calculate
318  // the shortest distance from the edge
319  if (z < 0) x = -x;
320  double dist(0);
321  // Input x, y in Geant4 unit and transformed to CMSSW standard
324  int sizew = (int)(hgpar_->waferPosX_.size());
325  int wafer = sizew;
326  // Transform to the local coordinate frame of the wafer first
327  for (int k = 0; k < sizew; ++k) {
328  double dx = std::abs(xx - hgpar_->waferPosX_[k]);
329  double dy = std::abs(yy - hgpar_->waferPosY_[k]);
330  if ((dx <= rmax_) && (dy <= hexside_) &&
331  ((dy <= 0.5 * hexside_) || (dx * tan30deg_ <= (hexside_ - dy)))) {
332  wafer = k;
333  xx -= hgpar_->waferPosX_[k];
334  yy -= hgpar_->waferPosY_[k];
335  break;
336  }
337  }
338  // Look at only one quarter (both x,y are positive)
339  if (wafer < sizew) {
340  if (std::abs(yy) < 0.5 * hexside_) {
341  dist = rmax_ - std::abs(xx);
342  } else {
343  dist = 0.5 * ((rmax_ - std::abs(xx)) -
344  sqrt3_ * (std::abs(yy) - 0.5 * hexside_));
345  }
346  } else {
347  dist = 0;
348  }
350 #ifdef EDM_ML_DEBUG
351  edm::LogVerbatim("HGCalGeom")
352  << "DistFromEdgeHex: Local " << xx << ":" << yy << " wafer " << wafer
353  << " flag " << (wafer < sizew) << " Distance " << rmax_ << ":"
354  << (rmax_ - std::abs(xx)) << ":" << (std::abs(yy) - 0.5 * hexside_) << ":"
355  << 0.5 * hexside_ << ":" << dist;
356 #endif
357  return dist;
358 }
359 
360 double HGCalDDDConstants::distFromEdgeTrap(double x, double y, double z) const {
361  // Assming the point is within the eta-phi plane of the scintillator tile,
362  // calculate the shortest distance from the edge
363  int lay = getLayer(z, false);
364  double xx = (z < 0) ? -x : x;
365  int indx = layerIndex(lay, false);
366  double r = HGCalParameters::k_ScaleFromDDD * std::sqrt(x * x + y * y);
367  double phi = (r == 0. ? 0. : std::atan2(y, xx));
368  if (phi < 0) phi += (2.0 * M_PI);
369  int type = hgpar_->scintType(lay);
370  double cell = hgpar_->scintCellSize(lay);
371  // Compare with the center of the tile find distances along R and also phi
372  // Take the smaller value
373  auto ir = std::lower_bound(hgpar_->radiusLayer_[type].begin(),
374  hgpar_->radiusLayer_[type].end(), r);
375  int irad = (int)(ir - hgpar_->radiusLayer_[type].begin());
376  if (irad < hgpar_->iradMinBH_[indx])
377  irad = hgpar_->iradMinBH_[indx];
378  else if (irad > hgpar_->iradMaxBH_[indx])
379  irad = hgpar_->iradMaxBH_[indx];
380  int iphi = 1 + (int)(phi / cell);
381  double dphi =
382  std::max(0.0, (0.5 * cell - std::abs(phi - (iphi - 0.5) * cell)));
383  double dist = std::min((r - hgpar_->radiusLayer_[type][irad - 1]),
384  (hgpar_->radiusLayer_[type][irad] - r));
385 #ifdef EDM_ML_DEBUG
386  edm::LogVerbatim("HGCalGeom")
387  << "DistFromEdgeTrap: Global " << x << ":" << y << ":" << z << " Layer "
388  << lay << " Index " << indx << ":" << type << " xx " << xx << " R " << r
389  << ":" << irad << ":" << hgpar_->radiusLayer_[type][irad - 1] << ":"
390  << hgpar_->radiusLayer_[type][irad] << " Phi " << phi << ":" << iphi
391  << ":" << (iphi - 0.5) * cell << " cell " << cell << " Dphi " << dphi
392  << " Dist " << dist << ":" << r * dphi;
393 #endif
394  return HGCalParameters::k_ScaleToDDD * std::min(r * dphi, dist);
395 }
396 
397 int HGCalDDDConstants::getLayer(double z, bool reco) const {
398  // Get the layer # from the gloabl z coordinate
399  unsigned int k = 0;
400  double zz =
402  const auto& zLayerHex = hgpar_->zLayerHex_;
403  auto itr = std::find_if(zLayerHex.begin() + 1, zLayerHex.end(),
404  [&k, &zz, &zLayerHex](double zLayer) {
405  ++k;
406  return zz < 0.5 * (zLayerHex[k-1] + zLayerHex[k]);
407  });
408  int lay = (itr == zLayerHex.end()) ? static_cast<int>(zLayerHex.size()) : k;
409  if (((mode_ == HGCalGeometryMode::Hexagon) ||
411  reco) {
412  int indx = layerIndex(lay, false);
413  if (indx >= 0) lay = hgpar_->layerGroup_[indx];
414  } else {
415  lay += (hgpar_->firstLayer_ - 1);
416  }
417  return lay;
418 }
419 
421  bool hexType,
422  bool reco) const {
424  if (hexType) {
425  if (indx >= hgpar_->waferTypeL_.size())
426  edm::LogWarning("HGCalGeom")
427  << "Wafer no. out bound for index " << indx << ":"
428  << (hgpar_->waferTypeL_).size() << ":" << (hgpar_->waferPosX_).size()
429  << ":" << (hgpar_->waferPosY_).size() << " ***** ERROR *****";
430  unsigned int type =
431  ((indx < hgpar_->waferTypeL_.size()) ? hgpar_->waferTypeL_[indx] : 3);
432  if (type > 0) --type;
433  mytr = hgpar_->getModule(type, reco);
434  } else {
435  mytr = hgpar_->getModule(indx, reco);
436  }
437  return mytr;
438 }
439 
440 std::vector<HGCalParameters::hgtrap> HGCalDDDConstants::getModules() const {
441  std::vector<HGCalParameters::hgtrap> mytrs;
442  for (unsigned int k = 0; k < hgpar_->moduleLayR_.size(); ++k)
443  mytrs.emplace_back(hgpar_->getModule(k, true));
444  return mytrs;
445 }
446 
447 int HGCalDDDConstants::getPhiBins(int lay) const {
449  : 0);
450 }
451 
452 std::pair<int, int> HGCalDDDConstants::getREtaRange(int lay) const {
453  int irmin(0), irmax(0);
455  int indx = layerIndex(lay, false);
456  if ((indx >= 0) && (indx < static_cast<int>(hgpar_->iradMinBH_.size()))) {
457  irmin = hgpar_->iradMinBH_[indx];
458  irmax = hgpar_->iradMaxBH_[indx];
459  }
460  }
461  return std::make_pair(irmin, irmax);
462 }
463 
464 std::vector<HGCalParameters::hgtrform> HGCalDDDConstants::getTrForms() const {
465  std::vector<HGCalParameters::hgtrform> mytrs;
466  for (unsigned int k = 0; k < hgpar_->trformIndex_.size(); ++k)
467  mytrs.emplace_back(hgpar_->getTrForm(k));
468  return mytrs;
469 }
470 
471 int HGCalDDDConstants::getTypeTrap(int layer) const {
472  // Get the module type for scinitllator
474  return hgpar_->scintType(layer);
475  } else {
476  return -1;
477  }
478 }
479 
480 int HGCalDDDConstants::getTypeHex(int layer, int waferU, int waferV) const {
481  // Get the module type for a silicon wafer
484  auto itr = hgpar_->typesInLayers_.find(
485  HGCalWaferIndex::waferIndex(layer, waferU, waferV));
486  return ((itr == hgpar_->typesInLayers_.end()
487  ? 2
488  : hgpar_->waferTypeL_[itr->second]));
489  } else {
490  return -1;
491  }
492 }
493 
494 bool HGCalDDDConstants::isHalfCell(int waferType, int cell) const {
495  if (waferType < 1 || cell < 0) return false;
496  return waferType == 2 ? hgpar_->cellCoarseHalf_[cell]
497  : hgpar_->cellFineHalf_[cell];
498 }
499 
500 bool HGCalDDDConstants::isValidHex(int lay, int mod, int cell,
501  bool reco) const {
502  // Check validity for a layer|wafer|cell of pre-TDR version
503  bool result(false), resultMod(false);
504  int cellmax(0);
507  int32_t copyNumber = hgpar_->waferCopy_[mod];
508  result = ((lay > 0 && lay <= (int)(layers(reco))));
509  if (result) {
510  const int32_t lay_idx = reco ? (lay - 1) * 3 + 1 : lay;
511  const auto& the_modules = hgpar_->copiesInLayers_[lay_idx];
512  auto moditr = the_modules.find(copyNumber);
513  result = resultMod = (moditr != the_modules.end());
514 #ifdef EDM_ML_DEBUG
515  if (!result)
516  edm::LogVerbatim("HGCalGeom")
517  << "HGCalDDDConstants: Layer " << lay << ":" << lay_idx << " Copy "
518  << copyNumber << ":" << mod << " Flag " << result;
519 #endif
520  if (result) {
521  if (moditr->second >= 0) {
522  if (mod >= (int)(hgpar_->waferTypeT_.size()))
523  edm::LogWarning("HGCalGeom")
524  << "Module no. out of bound for " << mod
525  << " to be compared with " << (hgpar_->waferTypeT_).size()
526  << " ***** ERROR *****";
527  cellmax = ((hgpar_->waferTypeT_[mod] == 1)
528  ? (int)(hgpar_->cellFineX_.size())
529  : (int)(hgpar_->cellCoarseX_.size()));
530  result = (cell >= 0 && cell <= cellmax);
531  } else {
532  result = isValidCell(lay_idx, mod, cell);
533  }
534  }
535  }
536  }
537 
538 #ifdef EDM_ML_DEBUG
539  if (!result)
540  edm::LogVerbatim("HGCalGeom")
541  << "HGCalDDDConstants: Layer " << lay << ":"
542  << (lay > 0 && (lay <= (int)(layers(reco)))) << " Module " << mod << ":"
543  << resultMod << " Cell " << cell << ":" << cellmax << ":"
544  << (cell >= 0 && cell <= cellmax) << ":" << maxCells(reco);
545 #endif
546  return result;
547 }
548 
549 bool HGCalDDDConstants::isValidHex8(int layer, int modU, int modV, int cellU,
550  int cellV) const {
551  // Check validity for a layer|wafer|cell of post-TDR version
552  int indx = HGCalWaferIndex::waferIndex(layer, modU, modV);
553  auto itr = hgpar_->typesInLayers_.find(indx);
554 #ifdef EDM_ML_DEBUG
555  edm::LogVerbatim("HGCalGeom")
556  << "HGCalDDDConstants::isValidHex8:WaferType " << layer << ":" << modU
557  << ":" << modV << ":" << indx << " Test "
558  << (itr == hgpar_->typesInLayers_.end());
559 #endif
560  if (itr == hgpar_->typesInLayers_.end()) return false;
561  auto jtr = waferIn_.find(indx);
562 #ifdef EDM_ML_DEBUG
563  edm::LogVerbatim("HGCalGeom") << "HGCalDDDConstants::isValidHex8:WaferIn "
564  << jtr->first << ":" << jtr->second;
565 #endif
566  if (!(jtr->second)) return false;
567  int N = ((hgpar_->waferTypeL_[itr->second] == 0) ? hgpar_->nCellsFine_
568  : hgpar_->nCellsCoarse_);
569 #ifdef EDM_ML_DEBUG
570  edm::LogVerbatim("HGCalGeom")
571  << "HGCalDDDConstants::isValidHex8:Cell " << cellU << ":" << cellV << ":"
572  << N << " Tests " << (cellU >= 0) << ":" << (cellU < 2 * N) << ":"
573  << (cellV >= 0) << ":" << (cellV < 2 * N) << ":" << ((cellV - cellU) < N)
574  << ":" << ((cellU - cellV) <= N);
575 #endif
576  if ((cellU >= 0) && (cellU < 2 * N) && (cellV >= 0) && (cellV < 2 * N)) {
577  return (((cellV - cellU) < N) && ((cellU - cellV) <= N));
578  } else {
579  return false;
580  }
581 }
582 
583 bool HGCalDDDConstants::isValidTrap(int layer, int irad, int iphi) const {
584  // Check validity for a layer|eta|phi of scintillator
585  const auto& indx = getIndex(layer, true);
586  if (indx.first < 0) return false;
587  return ((irad >= hgpar_->iradMinBH_[indx.first]) &&
588  (irad <= hgpar_->iradMaxBH_[indx.first]) && (iphi > 0) &&
589  (iphi <= hgpar_->scintCells(layer)));
590 }
591 
593  return (hgpar_->firstLayer_ + tot_layers_[(int)reco] - 1);
594 }
595 
596 unsigned int HGCalDDDConstants::layers(bool reco) const {
597  return tot_layers_[(int)reco];
598 }
599 
600 int HGCalDDDConstants::layerIndex(int lay, bool reco) const {
601  int ll = lay - hgpar_->firstLayer_;
602  if (ll < 0 || ll >= (int)(hgpar_->layerIndex_.size())) return -1;
605  if (reco && ll >= (int)(hgpar_->depthIndex_.size())) return -1;
606  return (reco ? hgpar_->depthLayerF_[ll] : hgpar_->layerIndex_[ll]);
607  } else {
608  return (hgpar_->layerIndex_[ll]);
609  }
610 }
611 
612 unsigned int HGCalDDDConstants::layersInit(bool reco) const {
613  return (reco ? hgpar_->depthIndex_.size() : hgpar_->layerIndex_.size());
614 }
615 
616 std::pair<float, float> HGCalDDDConstants::locateCell(int cell, int lay,
617  int type,
618  bool reco) const {
619  // type refers to wafer # for hexagon cell
620  float x(999999.), y(999999.);
621  const auto& index = getIndex(lay, reco);
622  int i = index.first;
623  if (i < 0) return std::make_pair(x, y);
626  x = hgpar_->waferPosX_[type];
627  y = hgpar_->waferPosY_[type];
628  if (hgpar_->waferTypeT_[type] == 1) {
629  x += hgpar_->cellFineX_[cell];
630  y += hgpar_->cellFineY_[cell];
631  } else {
632  x += hgpar_->cellCoarseX_[cell];
633  y += hgpar_->cellCoarseY_[cell];
634  }
635  if (!reco) {
638  }
639  }
640  return std::make_pair(x, y);
641 }
642 
643 std::pair<float, float> HGCalDDDConstants::locateCell(int lay, int waferU,
644  int waferV, int cellU,
645  int cellV, bool reco,
646  bool all,
647  bool
648 #ifdef EDM_ML_DEBUG
649  debug
650 #endif
651  ) const {
652 
653  float x(0), y(0);
654  int indx = HGCalWaferIndex::waferIndex(lay, waferU, waferV);
655  auto itr = hgpar_->typesInLayers_.find(indx);
656  int type = ((itr == hgpar_->typesInLayers_.end())
657  ? 2
658  : hgpar_->waferTypeL_[itr->second]);
659 #ifdef EDM_ML_DEBUG
660  if (debug)
661  edm::LogVerbatim("HGCalGeom")
662  << "LocateCell " << lay << ":" << waferU << ":" << waferV << ":" << indx
663  << ":" << (itr == hgpar_->typesInLayers_.end()) << ":" << type;
664 #endif
665  int kndx = cellV * 100 + cellU;
666  if (type == 0) {
667  auto ktr = hgpar_->cellFineIndex_.find(kndx);
668  if (ktr != hgpar_->cellFineIndex_.end()) {
669  x = hgpar_->cellFineX_[ktr->second];
670  y = hgpar_->cellFineY_[ktr->second];
671  }
672 #ifdef EDM_ML_DEBUG
673  if (debug)
674  edm::LogVerbatim("HGCalGeom")
675  << "Fine " << cellU << ":" << cellV << ":" << kndx << ":" << x << ":"
676  << y << ":" << (ktr != hgpar_->cellFineIndex_.end());
677 #endif
678  } else {
679  auto ktr = hgpar_->cellCoarseIndex_.find(kndx);
680  if (ktr != hgpar_->cellCoarseIndex_.end()) {
681  x = hgpar_->cellCoarseX_[ktr->second];
682  y = hgpar_->cellCoarseY_[ktr->second];
683  }
684 #ifdef EDM_ML_DEBUG
685  if (debug)
686  edm::LogVerbatim("HGCalGeom")
687  << "Coarse " << cellU << ":" << cellV << ":" << kndx << ":" << x
688  << ":" << y << ":" << (ktr != hgpar_->cellCoarseIndex_.end());
689 #endif
690  }
691  if (!reco) {
694  }
695  if (all) {
696  const auto& xy = waferPosition(waferU, waferV, reco);
697  x += xy.first;
698  y += xy.second;
699 #ifdef EDM_ML_DEBUG
700  if (debug)
701  edm::LogVerbatim("HGCalGeom") << "With wafer " << x << ":" << y << ":"
702  << xy.first << ":" << xy.second;
703 #endif
704  }
705  return std::make_pair(x, y);
706 }
707 
708 std::pair<float, float> HGCalDDDConstants::locateCellHex(int cell, int wafer,
709  bool reco) const {
710  float x(0), y(0);
711  if (hgpar_->waferTypeT_[wafer] == 1) {
712  x = hgpar_->cellFineX_[cell];
713  y = hgpar_->cellFineY_[cell];
714  } else {
715  x = hgpar_->cellCoarseX_[cell];
716  y = hgpar_->cellCoarseY_[cell];
717  }
718  if (!reco) {
721  }
722  return std::make_pair(x, y);
723 }
724 
725 std::pair<float, float> HGCalDDDConstants::locateCellTrap(int lay, int irad,
726  int iphi,
727  bool reco) const {
728  float x(0), y(0);
729  const auto& indx = getIndex(lay, reco);
730  if (indx.first >= 0) {
731  int ir = std::abs(irad);
732  int type = hgpar_->scintType(lay);
733  double phi = (iphi - 0.5) * indx.second;
734  double z = hgpar_->zLayerHex_[indx.first];
735  double r = 0.5 * (hgpar_->radiusLayer_[type][ir - 1] +
736  hgpar_->radiusLayer_[type][ir]);
737  std::pair<double, double> range = rangeR(z, true);
738  r = std::max(range.first, std::min(r, range.second));
739  x = r * std::cos(phi);
740  y = r * std::sin(phi);
741  if (irad < 0) x = -x;
742  }
743  if (!reco) {
746  }
747  return std::make_pair(x, y);
748 }
749 
750 bool HGCalDDDConstants::maskCell(const DetId& detId, int corners) const {
751  bool mask(false);
752  if (corners > 2 && corners < (int)(HGCalParameters::k_CornerSize)) {
755  int N(0), layer(0), waferU(0), waferV(0), u(0), v(0);
756  if (detId.det() == DetId::Forward) {
757  HFNoseDetId id(detId);
758  N = getUVMax(id.type());
759  layer = id.layer();
760  waferU = id.waferU();
761  waferV = id.waferV();
762  u = id.cellU();
763  v = id.cellV();
764  } else {
765  HGCSiliconDetId id(detId);
766  N = getUVMax(id.type());
767  layer = id.layer();
768  waferU = id.waferU();
769  waferV = id.waferV();
770  u = id.cellU();
771  v = id.cellV();
772  }
773  int wl = HGCalWaferIndex::waferIndex(layer, waferU, waferV);
774  auto itr = hgpar_->waferTypes_.find(wl);
775 #ifdef EDM_ML_DEBUG
776  edm::LogVerbatim("HGCalGeom")
777  << "MaskCell: Layer " << layer << " Wafer " << waferU << ":" << waferV
778  << " Index " << wl << ":" << (itr != hgpar_->waferTypes_.end());
779 #endif
780  if (itr != hgpar_->waferTypes_.end()) {
781  int ncor = (itr->second).first;
782  int fcor = (itr->second).second;
783  if (ncor < corners) {
784  mask = true;
785  } else {
786  if (ncor == 4) {
787  switch (fcor) {
788  case (0): {
789  mask = (v >= N);
790  break;
791  }
792  case (1): {
793  mask = (u >= N);
794  break;
795  }
796  case (2): {
797  mask = (u <= v);
798  break;
799  }
800  case (3): {
801  mask = (v < N);
802  break;
803  }
804  case (4): {
805  mask = (u < N);
806  break;
807  }
808  default: {
809  mask = (u > v);
810  break;
811  }
812  }
813  } else {
814  switch (fcor) {
815  case (0): {
816  if (ncor == 3) {
817  mask = !((u > 2 * v) && (v < N));
818  } else {
819  mask = ((u >= N) && (v >= N) && ((u + v) > (3 * N - 2)));
820  }
821  break;
822  }
823  case (1): {
824  if (ncor == 3) {
825  mask = !((u + v) < N);
826  } else {
827  mask = ((u >= N) && (u > v) && ((2 * u - v) > 2 * N));
828  }
829  break;
830  }
831  case (2): {
832  if (ncor == 3) {
833  mask = !((u < N) && (v > u) && (v > (2 * u - 1)));
834  } else {
835  mask = ((u > 2 * v) && (v < N));
836  }
837  break;
838  }
839  case (3): {
840  if (ncor == 3) {
841  mask = !((v >= u) && ((2 * v - u) > (2 * N - 2)));
842  } else {
843  mask = ((u + v) < N);
844  }
845  break;
846  }
847  case (4): {
848  if (ncor == 3) {
849  mask = !((u >= N) && (v >= N) && ((u + v) > (3 * N - 2)));
850  } else {
851  mask = ((u < N) && (v > u) && (v > (2 * u - 1)));
852  }
853  break;
854  }
855  default: {
856  if (ncor == 3) {
857  mask = !((u >= N) && (u > v) && ((2 * u - v) > 2 * N));
858  } else {
859  mask = ((v >= u) && ((2 * v - u) > (2 * N - 2)));
860  }
861  break;
862  }
863  }
864  }
865 #ifdef EDM_ML_DEBUG
866  edm::LogVerbatim("HGCalGeom")
867  << "Corners: " << ncor << ":" << fcor << " N " << N << " u " << u
868  << " v " << v << " Mask " << mask;
869 #endif
870  }
871  }
872  }
873  }
874  return mask;
875 }
876 
878  int cells(0);
879  for (unsigned int i = 0; i < layers(reco); ++i) {
880  int lay = reco ? hgpar_->depth_[i] : hgpar_->layer_[i];
881  if (cells < maxCells(lay, reco)) cells = maxCells(lay, reco);
882  }
883  return cells;
884 }
885 
886 int HGCalDDDConstants::maxCells(int lay, bool reco) const {
887  const auto& index = getIndex(lay, reco);
888  if (index.first < 0) return 0;
891  unsigned int cells(0);
892  for (unsigned int k = 0; k < hgpar_->waferTypeT_.size(); ++k) {
893  if (waferInLayerTest(k, index.first, hgpar_->defineFull_)) {
894  unsigned int cell = (hgpar_->waferTypeT_[k] == 1)
895  ? (hgpar_->cellFineX_.size())
896  : (hgpar_->cellCoarseX_.size());
897  if (cell > cells) cells = cell;
898  }
899  }
900  return (int)(cells);
901  } else if ((mode_ == HGCalGeometryMode::Hexagon8) ||
903  int cells(0);
904  for (unsigned int k = 0; k < hgpar_->waferCopy_.size(); ++k) {
905  if (waferInLayerTest(k, index.first, hgpar_->defineFull_)) {
909  int type = ((itr == hgpar_->typesInLayers_.end())
910  ? 2
911  : hgpar_->waferTypeL_[itr->second]);
912  int N = (type == 0) ? hgpar_->nCellsFine_ : hgpar_->nCellsCoarse_;
913  cells = std::max(cells, 3 * N * N);
914  }
915  }
916  return cells;
917  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
918  return hgpar_->scintCells(index.first + hgpar_->firstLayer_);
919  } else {
920  return 0;
921  }
922 }
923 
924 int HGCalDDDConstants::maxRows(int lay, bool reco) const {
925  int kymax(0);
926  const auto& index = getIndex(lay, reco);
927  int i = index.first;
928  if (i < 0) return kymax;
931  for (unsigned int k = 0; k < hgpar_->waferCopy_.size(); ++k) {
932  if (waferInLayerTest(k, i, hgpar_->defineFull_)) {
933  int ky = ((hgpar_->waferCopy_[k]) / 100) % 100;
934  if (ky > kymax) kymax = ky;
935  }
936  }
937  } else if ((mode_ == HGCalGeometryMode::Hexagon8) ||
939  kymax = 1 + 2 * hgpar_->waferUVMaxLayer_[index.first];
940  }
941  return kymax;
942 }
943 
944 int HGCalDDDConstants::modifyUV(int uv, int type1, int type2) const {
945  // Modify u/v for transition of type1 to type2
946  return (((type1 == type2) || (type1 * type2 != 0))
947  ? uv
948  : ((type1 == 0) ? (2 * uv + 1) / 3 : (3 * uv) / 2));
949 }
950 
951 int HGCalDDDConstants::modules(int lay, bool reco) const {
952  if (getIndex(lay, reco).first < 0)
953  return 0;
954  else
955  return max_modules_layer_[(int)reco][lay];
956 }
957 
958 int HGCalDDDConstants::modulesInit(int lay, bool reco) const {
959  int nmod(0);
960  const auto& index = getIndex(lay, reco);
961  if (index.first < 0) return nmod;
963  for (unsigned int k = 0; k < hgpar_->waferPosX_.size(); ++k) {
964  if (waferInLayerTest(k, index.first, hgpar_->defineFull_)) ++nmod;
965  }
966  } else {
967  nmod = 1 + hgpar_->lastModule_[index.first] -
968  hgpar_->firstModule_[index.first];
969  }
970  return nmod;
971 }
972 
973 double HGCalDDDConstants::mouseBite(bool reco) const {
974  return (reco ? hgpar_->mouseBite_
976 }
977 
979  int cells(0);
980  unsigned int nlayer = (reco) ? hgpar_->depth_.size() : hgpar_->layer_.size();
981  for (unsigned k = 0; k < nlayer; ++k) {
982  std::vector<int> ncells =
983  numberCells(((reco) ? hgpar_->depth_[k] : hgpar_->layer_[k]), reco);
984  cells = std::accumulate(ncells.begin(), ncells.end(), cells);
985  }
986  return cells;
987 }
988 
989 std::vector<int> HGCalDDDConstants::numberCells(int lay, bool reco) const {
990  const auto& index = getIndex(lay, reco);
991  int i = index.first;
992  std::vector<int> ncell;
993  if (i >= 0) {
996  for (unsigned int k = 0; k < hgpar_->waferTypeT_.size(); ++k) {
997  if (waferInLayerTest(k, i, hgpar_->defineFull_)) {
998  unsigned int cell = (hgpar_->waferTypeT_[k] == 1)
999  ? (hgpar_->cellFineX_.size())
1000  : (hgpar_->cellCoarseX_.size());
1001  ncell.emplace_back((int)(cell));
1002  }
1003  }
1004  } else if (mode_ == HGCalGeometryMode::Trapezoid) {
1005  int nphi = hgpar_->scintCells(lay);
1006  for (int k = hgpar_->firstModule_[i]; k <= hgpar_->lastModule_[i]; ++k)
1007  ncell.emplace_back(nphi);
1008  } else {
1009  for (unsigned int k = 0; k < hgpar_->waferCopy_.size(); ++k) {
1010  if (waferInLayerTest(k, index.first, hgpar_->defineFull_)) {
1011  int cell = numberCellsHexagon(
1014  ncell.emplace_back(cell);
1015  }
1016  }
1017  }
1018  }
1019  return ncell;
1020 }
1021 
1023  if (wafer >= 0 && wafer < (int)(hgpar_->waferTypeT_.size())) {
1024  if (hgpar_->waferTypeT_[wafer] == 1)
1025  return (int)(hgpar_->cellFineX_.size());
1026  else
1027  return (int)(hgpar_->cellCoarseX_.size());
1028  } else {
1029  return 0;
1030  }
1031 }
1032 
1033 int HGCalDDDConstants::numberCellsHexagon(int lay, int waferU, int waferV,
1034  bool flag) const {
1035  auto itr = hgpar_->typesInLayers_.find(
1036  HGCalWaferIndex::waferIndex(lay, waferU, waferV));
1037  int type = ((itr == hgpar_->typesInLayers_.end())
1038  ? 2
1039  : hgpar_->waferTypeL_[itr->second]);
1040  int N = (type == 0) ? hgpar_->nCellsFine_ : hgpar_->nCellsCoarse_;
1041  if (flag)
1042  return (3 * N * N);
1043  else
1044  return N;
1045 }
1046 
1047 std::pair<double, double> HGCalDDDConstants::rangeR(double z, bool reco) const {
1048  double rmin(0), rmax(0), zz(0);
1049  if (hgpar_->detectorType_ > 0) {
1050  zz = (reco ? std::abs(z) : HGCalParameters::k_ScaleFromDDD * std::abs(z));
1051  if (hgpar_->detectorType_ <= 2) {
1053  hgpar_->slopeMin_);
1054  } else {
1057  hgpar_->zLayerHex_,
1059  }
1060  if ((hgpar_->detectorType_ == 2) &&
1061  (zz >= hgpar_->zLayerHex_[hgpar_->firstMixedLayer_ - 1])) {
1064  hgpar_->zLayerHex_,
1066  } else {
1068  hgpar_->rMaxFront_,
1069  hgpar_->slopeTop_);
1070  }
1071  }
1072  if (!reco) {
1075  }
1076 #ifdef EDM_ML_DEBUG
1077  edm::LogVerbatim("HGCalGeom") << "HGCalDDDConstants:rangeR: " << z << ":"
1078  << zz << " R " << rmin << ":" << rmax;
1079 #endif
1080  return std::make_pair(rmin, rmax);
1081 }
1082 
1083 std::pair<double, double> HGCalDDDConstants::rangeZ(bool reco) const {
1084  double zmin = (hgpar_->zLayerHex_[0] - hgpar_->waferThick_);
1085  double zmax =
1086  (hgpar_->zLayerHex_[hgpar_->zLayerHex_.size() - 1] + hgpar_->waferThick_);
1087 #ifdef EDM_ML_DEBUG
1088  edm::LogVerbatim("HGCalGeom") << "HGCalDDDConstants:rangeZ: " << zmin << ":"
1089  << zmax << ":" << hgpar_->waferThick_;
1090 #endif
1091  if (!reco) {
1094  }
1095  return std::make_pair(zmin, zmax);
1096 }
1097 
1098 std::pair<int, int> HGCalDDDConstants::rowColumnWafer(int wafer) const {
1099  int row(0), col(0);
1100  if (wafer < (int)(hgpar_->waferCopy_.size())) {
1101  int copy = hgpar_->waferCopy_[wafer];
1102  col = copy % 100;
1103  if ((copy / 10000) % 10 != 0) col = -col;
1104  row = (copy / 100) % 100;
1105  if ((copy / 100000) % 10 != 0) row = -row;
1106  }
1107  return std::make_pair(row, col);
1108 }
1109 
1110 std::pair<int, int> HGCalDDDConstants::simToReco(int cell, int lay, int mod,
1111  bool half) const {
1112  if ((mode_ != HGCalGeometryMode::Hexagon) &&
1114  return std::make_pair(cell, lay);
1115  } else {
1116  const auto& index = getIndex(lay, false);
1117  int i = index.first;
1118  if (i < 0) {
1119  edm::LogWarning("HGCalGeom")
1120  << "Wrong Layer # " << lay << " not in the list ***** ERROR *****";
1121  return std::make_pair(-1, -1);
1122  }
1123  if (mod >= (int)(hgpar_->waferTypeL_).size()) {
1124  edm::LogWarning("HGCalGeom")
1125  << "Invalid Wafer # " << mod << "should be < "
1126  << (hgpar_->waferTypeL_).size() << " ***** ERROR *****";
1127  return std::make_pair(-1, -1);
1128  }
1129  int depth(-1);
1130  int kx = cell;
1131  int type = hgpar_->waferTypeL_[mod];
1132  if (type == 1) {
1133  depth = hgpar_->layerGroup_[i];
1134  } else if (type == 2) {
1135  depth = hgpar_->layerGroupM_[i];
1136  } else {
1137  depth = hgpar_->layerGroupO_[i];
1138  }
1139  return std::make_pair(kx, depth);
1140  }
1141 }
1142 
1144  const int ncopies = hgpar_->waferCopy_.size();
1145  int wafer(ncopies);
1146  bool result(false);
1147  for (int k = 0; k < ncopies; ++k) {
1148  if (copy == hgpar_->waferCopy_[k]) {
1149  wafer = k;
1150  result = true;
1151  break;
1152  }
1153  }
1154  if (!result) {
1155  wafer = -1;
1156 #ifdef EDM_ML_DEBUG
1157  edm::LogVerbatim("HGCalGeom")
1158  << "Cannot find " << copy << " in a list of " << ncopies << " members";
1159  for (int k = 0; k < ncopies; ++k)
1160  edm::LogVerbatim("HGCalGeom")
1161  << "[" << k << "] " << hgpar_->waferCopy_[k];
1162 #endif
1163  }
1164  return wafer;
1165 }
1166 
1167 void HGCalDDDConstants::waferFromPosition(const double x, const double y,
1168  int& wafer, int& icell,
1169  int& celltyp) const {
1170  // Input x, y in Geant4 unit and transformed to CMSSW standard
1173  int size_ = (int)(hgpar_->waferCopy_.size());
1174  wafer = size_;
1175  for (int k = 0; k < size_; ++k) {
1176  double dx = std::abs(xx - hgpar_->waferPosX_[k]);
1177  double dy = std::abs(yy - hgpar_->waferPosY_[k]);
1178  if (dx <= rmax_ && dy <= hexside_) {
1179  if ((dy <= 0.5 * hexside_) || (dx * tan30deg_ <= (hexside_ - dy))) {
1180  wafer = k;
1181  celltyp = hgpar_->waferTypeT_[k];
1182  xx -= hgpar_->waferPosX_[k];
1183  yy -= hgpar_->waferPosY_[k];
1184  break;
1185  }
1186  }
1187  }
1188  if (wafer < size_) {
1189  if (celltyp == 1)
1190  icell = cellHex(
1191  xx, yy, 0.5 * HGCalParameters::k_ScaleFromDDD * hgpar_->cellSize_[0],
1193  else
1194  icell = cellHex(
1195  xx, yy, 0.5 * HGCalParameters::k_ScaleFromDDD * hgpar_->cellSize_[1],
1197  } else {
1198  wafer = -1;
1199 #ifdef EDM_ML_DEBUG
1200  edm::LogWarning("HGCalGeom") << "Cannot get wafer type corresponding to "
1201  << x << ":" << y << " " << xx << ":" << yy;
1202 #endif
1203  }
1204 #ifdef EDM_ML_DEBUG
1205  edm::LogVerbatim("HGCalGeom")
1206  << "Position " << x << ":" << y << " Wafer " << wafer << ":" << size_
1207  << " XX " << xx << ":" << yy << " Cell " << icell << " Type " << celltyp;
1208 #endif
1209 }
1210 
1211 void HGCalDDDConstants::waferFromPosition(const double x, const double y,
1212  const int layer, int& waferU,
1213  int& waferV, int& cellU, int& cellV,
1214  int& celltype, double& wt,
1215  bool
1216 #ifdef EDM_ML_DEBUG
1217  debug
1218 #endif
1219  ) const {
1220 
1221  double xx(HGCalParameters::k_ScaleFromDDD * x);
1222  double yy(HGCalParameters::k_ScaleFromDDD * y);
1223  waferU = waferV = 1 + hgpar_->waferUVMax_;
1224  for (unsigned int k = 0; k < hgpar_->waferPosX_.size(); ++k) {
1225  double dx = std::abs(xx - hgpar_->waferPosX_[k]);
1226  double dy = std::abs(yy - hgpar_->waferPosY_[k]);
1227  if (dx <= rmax_ && dy <= hexside_) {
1228  if ((dy <= 0.5 * hexside_) || (dx * tan30deg_ <= (hexside_ - dy))) {
1231  auto itr = hgpar_->typesInLayers_.find(
1232  HGCalWaferIndex::waferIndex(layer, waferU, waferV));
1233  celltype = ((itr == hgpar_->typesInLayers_.end())
1234  ? 2
1235  : hgpar_->waferTypeL_[itr->second]);
1236 #ifdef EDM_ML_DEBUG
1237  if (debug)
1238  edm::LogVerbatim("HGCalGeom")
1239  << "WaferFromPosition:: Input " << xx << ":" << yy
1240  << " compared with " << hgpar_->waferPosX_[k] << ":"
1241  << hgpar_->waferPosY_[k] << " difference " << dx << ":" << dy
1242  << ":" << dx * tan30deg_ << ":" << (hexside_ - dy)
1243  << " comparator " << rmax_ << ":" << hexside_ << " wafer "
1244  << waferU << ":" << waferV << ":" << celltype;
1245 #endif
1246  xx -= hgpar_->waferPosX_[k];
1247  yy -= hgpar_->waferPosY_[k];
1248  break;
1249  }
1250  }
1251  }
1252  if (std::abs(waferU) <= hgpar_->waferUVMax_) {
1253  cellHex(xx, yy, celltype, cellU, cellV
1254 #ifdef EDM_ML_DEBUG
1255  ,
1256  debug
1257 #endif
1258  );
1259  wt = ((celltype < 2)
1260  ? (hgpar_->cellThickness_[celltype] / hgpar_->waferThick_)
1261  : 1.0);
1262  } else {
1263  cellU = cellV = 2 * hgpar_->nCellsFine_;
1264  wt = 1.0;
1265  celltype = -1;
1266  }
1267 #ifdef EDM_ML_DEBUG
1268  if (celltype < 0) {
1269  double x1(HGCalParameters::k_ScaleFromDDD * x);
1270  double y1(HGCalParameters::k_ScaleFromDDD * y);
1271  edm::LogVerbatim("HGCalGeom")
1272  << "waferFromPosition: Bad type for X " << x << ":" << x1 << ":" << xx
1273  << " Y " << y << ":" << y1 << ":" << yy << " Wafer " << waferU << ":"
1274  << waferV << " Cell " << cellU << ":" << cellV;
1275  for (unsigned int k = 0; k < hgpar_->waferPosX_.size(); ++k) {
1276  double dx = std::abs(x1 - hgpar_->waferPosX_[k]);
1277  double dy = std::abs(y1 - hgpar_->waferPosY_[k]);
1278  edm::LogVerbatim("HGCalGeom")
1279  << "Wafer [" << k << "] Position (" << hgpar_->waferPosX_[k] << ", "
1280  << hgpar_->waferPosY_[k] << ") difference " << dx << ":" << dy << ":"
1281  << dx * tan30deg_ << ":" << hexside_ - dy << " Paramerers " << rmax_
1282  << ":" << hexside_;
1283  }
1284  }
1285 #endif
1286 }
1287 
1288 bool HGCalDDDConstants::waferInLayer(int wafer, int lay, bool reco) const {
1289  const auto& indx = getIndex(lay, reco);
1290  if (indx.first < 0) return false;
1291  return waferInLayerTest(wafer, indx.first, hgpar_->defineFull_);
1292 }
1293 
1294 bool HGCalDDDConstants::waferFullInLayer(int wafer, int lay, bool reco) const {
1295  const auto& indx = getIndex(lay, reco);
1296  if (indx.first < 0) return false;
1297  return waferInLayerTest(wafer, indx.first, false);
1298 }
1299 
1300 std::pair<double, double> HGCalDDDConstants::waferPosition(int wafer,
1301  bool reco) const {
1302  double xx(0), yy(0);
1303  if (wafer >= 0 && wafer < (int)(hgpar_->waferPosX_.size())) {
1304  xx = hgpar_->waferPosX_[wafer];
1305  yy = hgpar_->waferPosY_[wafer];
1306  }
1307  if (!reco) {
1310  }
1311  return std::make_pair(xx, yy);
1312 }
1313 
1314 std::pair<double, double> HGCalDDDConstants::waferPosition(int waferU,
1315  int waferV,
1316  bool reco) const {
1317  double xx(0), yy(0);
1318  int indx = HGCalWaferIndex::waferIndex(0, waferU, waferV);
1319  auto itr = hgpar_->wafersInLayers_.find(indx);
1320  if (itr != hgpar_->wafersInLayers_.end()) {
1321  xx = hgpar_->waferPosX_[itr->second];
1322  yy = hgpar_->waferPosY_[itr->second];
1323  }
1324  if (!reco) {
1327  }
1328  return std::make_pair(xx, yy);
1329 }
1330 
1331 int HGCalDDDConstants::waferType(DetId const& id) const {
1332  int type(1);
1335  type = ((id.det() != DetId::Forward) ? (1 + HGCSiliconDetId(id).type())
1336  : (1 + HFNoseDetId(id).type()));
1337  } else if ((mode_ == HGCalGeometryMode::Hexagon) ||
1339  type = waferTypeL(HGCalDetId(id).wafer());
1340  }
1341  return type;
1342 }
1343 
1344 bool HGCalDDDConstants::waferVirtual(int layer, int waferU, int waferV) const {
1345  bool type(false);
1348  int wl = HGCalWaferIndex::waferIndex(layer, waferU, waferV, false);
1349  type = (hgpar_->waferTypes_.find(wl) != hgpar_->waferTypes_.end());
1350  } else if ((mode_ == HGCalGeometryMode::Hexagon) ||
1352  int wl = HGCalWaferIndex::waferIndex(layer, waferU, 0, true);
1353  type = (hgpar_->waferTypes_.find(wl) != hgpar_->waferTypes_.end());
1354  }
1355  return type;
1356 }
1357 
1358 double HGCalDDDConstants::waferZ(int lay, bool reco) const {
1359  const auto& index = getIndex(lay, reco);
1360  if (index.first < 0)
1361  return 0;
1362  else
1363  return (reco ? hgpar_->zLayerHex_[index.first]
1365  hgpar_->zLayerHex_[index.first]);
1366 }
1367 
1369  int wafer(0);
1371  for (unsigned int i = 0; i < layers(true); ++i) {
1372  int lay = hgpar_->depth_[i];
1373  wafer += modules(lay, true);
1374  }
1375  } else {
1376  wafer = (int)(hgpar_->moduleLayR_.size());
1377  }
1378  return wafer;
1379 }
1380 
1381 int HGCalDDDConstants::wafers(int layer, int type) const {
1382  int wafer(0);
1384  auto itr = waferLayer_.find(layer);
1385  if (itr != waferLayer_.end()) {
1386  unsigned ity = (type > 0 && type <= 2) ? type : 0;
1387  wafer = (itr->second)[ity];
1388  }
1389  } else {
1390  const auto& index = getIndex(layer, true);
1391  wafer = 1 + hgpar_->lastModule_[index.first] -
1392  hgpar_->firstModule_[index.first];
1393  }
1394  return wafer;
1395 }
1396 
1397 int HGCalDDDConstants::cellHex(double xx, double yy, const double& cellR,
1398  const std::vector<double>& posX,
1399  const std::vector<double>& posY) const {
1400  int num(0);
1401  const double tol(0.00001);
1402  double cellY = 2.0 * cellR * tan30deg_;
1403  for (unsigned int k = 0; k < posX.size(); ++k) {
1404  double dx = std::abs(xx - posX[k]);
1405  double dy = std::abs(yy - posY[k]);
1406  if (dx <= (cellR + tol) && dy <= (cellY + tol)) {
1407  double xmax = (dy <= 0.5 * cellY)
1408  ? cellR
1409  : (cellR - (dy - 0.5 * cellY) / tan30deg_);
1410  if (dx <= (xmax + tol)) {
1411  num = k;
1412  break;
1413  }
1414  }
1415  }
1416  return num;
1417 }
1418 
1419 void HGCalDDDConstants::cellHex(double xloc, double yloc, int cellType,
1420  int& cellU, int& cellV,
1421  bool
1422 #ifdef EDM_ML_DEBUG
1423  debug
1424 #endif
1425  ) const {
1426  int N = (cellType == 0) ? hgpar_->nCellsFine_ : hgpar_->nCellsCoarse_;
1427  double Rc = 2 * rmax_ / (3 * N);
1428  double rc = 0.5 * Rc * sqrt3_;
1429  double v0 = ((xloc / Rc - 1.0) / 1.5);
1430  int cv0 = (v0 > 0) ? (N + (int)(v0 + 0.5)) : (N - (int)(-v0 + 0.5));
1431  double u0 = (0.5 * yloc / rc + 0.5 * cv0);
1432  int cu0 = (u0 > 0) ? (N / 2 + (int)(u0 + 0.5)) : (N / 2 - (int)(-u0 + 0.5));
1433  cu0 = std::max(0, std::min(cu0, 2 * N - 1));
1434  cv0 = std::max(0, std::min(cv0, 2 * N - 1));
1435  if (cv0 - cu0 >= N) cv0 = cu0 + N - 1;
1436 #ifdef EDM_ML_DEBUG
1437  if (debug)
1438  edm::LogVerbatim("HGCalGeom")
1439  << "cellHex: input " << xloc << ":" << yloc << ":" << cellType
1440  << " parameter " << rc << ":" << Rc << " u0 " << u0 << ":" << cu0
1441  << " v0 " << v0 << ":" << cv0;
1442 #endif
1443  bool found(false);
1444  static const int shift[3] = {0, 1, -1};
1445  for (int i1 = 0; i1 < 3; ++i1) {
1446  cellU = cu0 + shift[i1];
1447  for (int i2 = 0; i2 < 3; ++i2) {
1448  cellV = cv0 + shift[i2];
1449  if (((cellV - cellU) < N) && ((cellU - cellV) <= N) && (cellU >= 0) &&
1450  (cellV >= 0) && (cellU < 2 * N) && (cellV < 2 * N)) {
1451  double xc = (1.5 * (cellV - N) + 1.0) * Rc;
1452  double yc = (2 * cellU - cellV - N) * rc;
1453  if ((std::abs(yloc - yc) <= rc) && (std::abs(xloc - xc) <= Rc) &&
1454  ((std::abs(xloc - xc) <= 0.5 * Rc) ||
1455  (std::abs(yloc - yc) <= sqrt3_ * (Rc - std::abs(xloc - xc))))) {
1456 #ifdef EDM_ML_DEBUG
1457  if (debug)
1458  edm::LogVerbatim("HGCalGeom")
1459  << "cellHex: local " << xc << ":" << yc << " difference "
1460  << std::abs(xloc - xc) << ":" << std::abs(yloc - yc) << ":"
1461  << sqrt3_ * (Rc - std::abs(yloc - yc)) << " comparator " << rc
1462  << ":" << Rc << " (u,v) = (" << cellU << "," << cellV << ")";
1463 #endif
1464  found = true;
1465  break;
1466  }
1467  }
1468  }
1469  if (found) break;
1470  }
1471  if (!found) {
1472  cellU = cu0;
1473  cellV = cv0;
1474  }
1475 }
1476 
1477 std::pair<int, float> HGCalDDDConstants::getIndex(int lay, bool reco) const {
1478  int indx = layerIndex(lay, reco);
1479  if (indx < 0) return std::make_pair(-1, 0);
1480  float cell(0);
1481  if ((mode_ == HGCalGeometryMode::Hexagon) ||
1483  cell = (reco ? hgpar_->moduleCellR_[0] : hgpar_->moduleCellS_[0]);
1484  } else {
1487  cell = (reco ? hgpar_->moduleCellR_[0] : hgpar_->moduleCellS_[0]);
1488  } else {
1489  cell = hgpar_->scintCellSize(lay);
1490  }
1491  }
1492  return std::make_pair(indx, cell);
1493 }
1494 
1495 bool HGCalDDDConstants::isValidCell(int lay, int wafer, int cell) const {
1496  // Calculate the position of the cell
1497  // Works for options HGCalHexagon/HGCalHexagonFull
1498  double x = hgpar_->waferPosX_[wafer];
1499  double y = hgpar_->waferPosY_[wafer];
1500  if (hgpar_->waferTypeT_[wafer] == 1) {
1501  x += hgpar_->cellFineX_[cell];
1502  y += hgpar_->cellFineY_[cell];
1503  } else {
1504  x += hgpar_->cellCoarseX_[cell];
1505  y += hgpar_->cellCoarseY_[cell];
1506  }
1507  double rr = sqrt(x * x + y * y);
1508  bool result = ((rr >= hgpar_->rMinLayHex_[lay - 1]) &&
1509  (rr <= hgpar_->rMaxLayHex_[lay - 1]) &&
1510  (wafer < (int)(hgpar_->waferPosX_.size())));
1511 #ifdef EDM_ML_DEBUG
1512  if (!result)
1513  edm::LogVerbatim("HGCalGeom")
1514  << "Input " << lay << ":" << wafer << ":" << cell << " Position " << x
1515  << ":" << y << ":" << rr << " Compare Limits "
1516  << hgpar_->rMinLayHex_[lay - 1] << ":" << hgpar_->rMaxLayHex_[lay - 1]
1517  << " Flag " << result;
1518 #endif
1519  return result;
1520 }
1521 
1522 bool HGCalDDDConstants::waferInLayerTest(int wafer, int lay, bool full) const {
1523  bool flag = ((mode_ == HGCalGeometryMode::Hexagon) ||
1525  ? true
1526  : false;
1527  std::pair<int, int> corner = HGCalGeomTools::waferCorner(
1528  hgpar_->waferPosX_[wafer], hgpar_->waferPosY_[wafer], rmax_, hexside_,
1529  hgpar_->rMinLayHex_[lay], hgpar_->rMaxLayHex_[lay], flag);
1530  bool in = (full ? (corner.first > 0)
1531  : (corner.first == (int)(HGCalParameters::k_CornerSize)));
1532 #ifdef EDM_ML_DEBUG
1533  edm::LogVerbatim("HGCalGeom")
1534  << "WaferInLayerTest: Layer " << lay << " wafer " << wafer << " R-limits "
1535  << hgpar_->rMinLayHex_[lay] << ":" << hgpar_->rMaxLayHex_[lay]
1536  << " Corners " << corner.first << ":" << corner.second << " In " << in;
1537 #endif
1538  return in;
1539 }
1540 
1542 
std::vector< int > iradMaxBH_
size
Write out results.
bool isHalfCell(int waferType, int cell) const
std::vector< double > waferPosY_
type
Definition: HCALResponse.h:21
std::vector< int > layer_
std::vector< int > depthLayerF_
bool isValidTrap(int lay, int ieta, int iphi) const
std::vector< int > depth_
std::vector< double > zFrontMin_
std::vector< double > moduleCellR_
int getTypeTrap(int layer) const
wafer_map cellFineIndex_
layer_map copiesInLayers_
static void radius(double zf, double zb, std::vector< double > const &zFront1, std::vector< double > const &rFront1, std::vector< double > const &slope1, std::vector< double > const &zFront2, std::vector< double > const &rFront2, std::vector< double > const &slope2, int flag, std::vector< double > &zz, std::vector< double > &rin, std::vector< double > &rout)
int getLayer(double z, bool reco) const
std::vector< HGCalParameters::hgtrap > getModules() const
std::vector< bool > cellCoarseHalf_
std::vector< bool > cellFineHalf_
def copy(args, dbName)
int scintType(const int layer) const
void waferFromPosition(const double x, const double y, int &wafer, int &icell, int &celltyp) const
int type() const
get the type
Definition: HFNoseDetId.h:52
std::vector< int > moduleLayR_
double cellSizeHex(int type) const
HGCalParameters::hgtrform getTrForm(unsigned int k) const
Simrecovecs max_modules_layer_
int lastLayer(bool reco) const
double cellThickness(int layer, int waferU, int waferV) const
unsigned int layersInit(bool reco) const
static int32_t waferV(const int32_t index)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool cellInLayer(int waferU, int waferV, int cellU, int cellV, int lay, bool reco) const
int maxRows(int lay, bool reco) const
int getTypeHex(int layer, int waferU, int waferV) const
std::vector< double > cellFineY_
bool waferInLayerTest(int wafer, int lay, bool full) const
std::unordered_map< int32_t, bool > waferIn_
int modulesInit(int lay, bool reco) const
std::pair< float, float > locateCell(int cell, int lay, int type, bool reco) const
std::vector< uint32_t > trformIndex_
bool isValidHex(int lay, int mod, int cell, bool reco) const
std::vector< int > layerGroupM_
HGCalGeometryMode::GeometryMode mode_
HGCalGeometryMode::GeometryMode mode_
double scintCellSize(const int layer) const
int cellHex(double xx, double yy, const double &cellR, const std::vector< double > &posX, const std::vector< double > &posY) const
int layerIndex(int lay, bool reco) const
wafer_map wafersInLayers_
std::pair< float, float > locateCellHex(int cell, int wafer, bool reco) const
U second(std::pair< T, U > const &p)
std::vector< double > cellCoarseX_
std::pair< double, double > rangeR(double z, bool reco) const
std::map< int, HGCWaferParam > waferLayer_
int modules(int lay, bool reco) const
static std::pair< int32_t, int32_t > waferCorner(double xpos, double ypos, double r, double R, double rMin, double rMax, bool oldBug=false)
std::vector< int > firstModule_
std::pair< int, int > simToReco(int cell, int layer, int mod, bool half) const
unsigned int getTrFormN() const
double distFromEdgeTrap(double x, double y, double z) const
std::pair< double, double > rangeZ(bool reco) const
unsigned int layers(bool reco) const
bool isValidHex8(int lay, int modU, int modV, int cellU, int cellV) const
int numberCellsHexagon(int wafer) const
std::vector< double > cellSize_
std::vector< int > waferUVMaxLayer_
std::vector< HGCalParameters::hgtrform > getTrForms() const
std::pair< double, double > cellEtaPhiTrap(int type, int irad) const
std::vector< int > layerIndex_
double mouseBite(bool reco) const
susybsm::HSCParticleRefProd hp
Definition: classes.h:27
hgtrform getTrForm(unsigned int k) const
int type() const
get the type
T sqrt(T t)
Definition: SSEVec.h:18
std::pair< int, float > getIndex(int lay, bool reco) const
static double k_ScaleFromDDD
HGCalDDDConstants(const HGCalParameters *hp, const std::string &name)
std::pair< float, float > locateCellTrap(int lay, int ieta, int iphi, bool reco) const
std::vector< double > rMaxFront_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Definition: GenABIO.cc:168
int scintCells(const int layer) const
hgtrap getModule(unsigned int k, bool reco) const
CellType cellType(int type, int waferU, int waferV) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > slopeTop_
static int32_t waferU(const int32_t index)
bool isValidCell(int layindex, int wafer, int cell) const
int waferTypeL(int wafer) const
T min(T a, T b)
Definition: MathUtil.h:58
int getPhiBins(int lay) const
std::vector< double > rMinLayHex_
static int32_t waferIndex(int32_t layer, int32_t waferU, int32_t waferV, bool old=false)
std::pair< int, int > getREtaRange(int lay) const
std::pair< int, int > rowColumnWafer(const int wafer) const
std::vector< double > zLayerHex_
double waferZ(int layer, bool reco) const
#define M_PI
waferT_map waferTypes_
int k[5][pyjets_maxn]
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:102
int waferType(DetId const &id) const
bool waferVirtual(int layer, int waferU, int waferV) const
static const int maxType
std::vector< double > rMaxLayHex_
static uint32_t k_CornerSize
std::vector< double > slopeMin_
Definition: DetId.h:18
std::vector< int > lastModule_
std::array< int, 4 > waferMax_
int waferFromCopy(int copy) const
std::vector< double > radiusMixBoundary_
#define debug
Definition: HDRShower.cc:19
#define N
Definition: blowfish.cc:9
std::vector< double > cellThickness_
std::vector< int > layerGroup_
std::vector< double > moduleCellS_
static double k_ScaleToDDD
bool maskCell(const DetId &id, int corners) const
int numberCells(bool reco) const
wafer_map cellCoarseIndex_
bool waferFullInLayer(int wafer, int lay, bool reco) const
std::vector< double > rMinFront_
std::vector< int > iradMinBH_
std::vector< double > cellFineX_
wafer_map typesInLayers_
std::vector< int > layerGroupO_
static const int minType
fixed size matrix
std::array< int, 5 > assignCellHex(float x, float y, int lay, bool reco) const
std::array< int, 3 > HGCWaferParam
std::vector< int > waferCopy_
col
Definition: cuy.py:1010
#define EDM_ML_DEBUG
static unsigned int const shift
std::vector< int > depthIndex_
std::pair< double, double > waferPosition(int wafer, bool reco) const
std::pair< int, int > assignCell(float x, float y, int lay, int subSec, bool reco) const
std::array< int, 3 > assignCellTrap(float x, float y, float z, int lay, bool reco) const
std::vector< double > zFrontTop_
std::vector< double > radiusLayer_[2]
std::array< uint32_t, 2 > tot_layers_
std::vector< int > waferTypeT_
std::vector< double > cellCoarseY_
HGCalParameters::hgtrap getModule(unsigned int k, bool hexType, bool reco) const
const HGCalParameters * hgpar_
int maxCells(bool reco) const
std::vector< double > waferPosX_
int modifyUV(int uv, int type1, int type2) const
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
std::vector< int > waferTypeL_
int getUVMax(int type) const
static double tan30deg_
double distFromEdgeHex(double x, double y, double z) const
bool waferInLayer(int wafer, int lay, bool reco) const
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39