CMS 3D CMS Logo

HGCalCell.cc
Go to the documentation of this file.
3 #include <vector>
4 
5 //#define EDM_ML_DEBUG
6 
7 HGCalCell::HGCalCell(double waferSize, int32_t nFine, int32_t nCoarse) {
8  ncell_[0] = nFine;
9  ncell_[1] = nCoarse;
10  for (int k = 0; k < 2; ++k) {
11  cellX_[k] = waferSize / (3 * ncell_[k]);
12  cellY_[k] = sqrt3By2_ * cellX_[k];
13  }
14 #ifdef EDM_ML_DEBUG
15  edm::LogVerbatim("HGCalGeom") << "HGCalCell initialized with waferSize " << waferSize << " number of cells " << nFine
16  << ":" << nCoarse;
17 #endif
18 }
19 
20 std::pair<double, double> HGCalCell::cellUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
21  if (type != 0)
22  type = 1;
23  double x(0), y(0);
24  switch (placementIndex) {
26  x = (1.5 * (v - u) + 0.5) * cellX_[type];
27  y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
28  break;
30  x = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
31  y = (2 * u - v - ncell_[type]) * cellY_[type];
32  break;
34  x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
35  y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
36  break;
38  x = -(1.5 * (v - u) + 0.5) * cellX_[type];
39  y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
40  break;
42  x = -(1.5 * (v - ncell_[type]) + 1) * cellX_[type];
43  y = -(2 * u - v - ncell_[type]) * cellY_[type];
44  break;
46  x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
47  y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
48  break;
50  x = (1.5 * (u - v) - 0.5) * cellX_[type];
51  y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
52  break;
54  x = -(1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
55  y = (2 * u - v - ncell_[type]) * cellY_[type];
56  break;
58  x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
59  y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
60  break;
62  x = -(1.5 * (u - v) - 0.5) * cellX_[type];
63  y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
64  break;
66  x = (1.5 * (v - ncell_[type]) + 1) * cellX_[type];
67  y = -(2 * u - v - ncell_[type]) * cellY_[type];
68  break;
69  default:
70  x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
71  y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
72  break;
73  }
74  return std::make_pair(x, y);
75 }
76 
77 std::pair<double, double> HGCalCell::cellUV2XY2(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
78  if (type != 0)
79  type = 1;
80  double x(0), y(0);
81  if (placementIndex < HGCalCell::cellPlacementExtra) {
82  double x0 = (1.5 * (u - v) - 0.5) * cellX_[type];
83  double y0 = (u + v - 2 * ncell_[type] + 1) * cellY_[type];
84  const std::vector<double> fcos = {1.0, 0.5, -0.5, -1.0, -0.5, 0.5};
85  const std::vector<double> fsin = {0.0, sqrt3By2_, sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_};
86  x = x0 * fcos[placementIndex] - y0 * fsin[placementIndex];
87  y = x0 * fsin[placementIndex] + y0 * fcos[placementIndex];
88  } else {
89  double x0 = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
90  double y0 = (2 * u - v - ncell_[type]) * cellY_[type];
91  const std::vector<double> fcos = {0.5, 1.0, 0.5, -0.5, -1.0, -0.5};
92  const std::vector<double> fsin = {sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_, 0.0, sqrt3By2_};
93  x = x0 * fcos[placementIndex - HGCalCell::cellPlacementExtra] -
94  y0 * fsin[placementIndex - HGCalCell::cellPlacementExtra];
95  y = x0 * fsin[placementIndex - HGCalCell::cellPlacementExtra] +
96  y0 * fcos[placementIndex - HGCalCell::cellPlacementExtra];
97  }
98  return std::make_pair(x, y);
99 }
100 
101 std::pair<int, int> HGCalCell::cellUV2Cell(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
102  if (type != 0)
103  type = 1;
104  int cell(0), cellx(0), cellt(HGCalCell::fullCell);
105  if (placementIndex >= HGCalCell::cellPlacementExtra) {
106  const std::vector<int> itype0 = {0, 7, 8, 9, 10, 11, 6, 3, 4, 5, 4, 5, 3};
107  const std::vector<int> itype1 = {0, 0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 1, 2};
108  const std::vector<int> itype2 = {0, 11, 6, 7, 8, 9, 10, 5, 3, 4, 3, 4, 5};
109  const std::vector<int> itype3 = {0, 4, 5, 0, 1, 2, 3, 2, 0, 1, 2, 0, 1};
110  const std::vector<int> itype4 = {0, 9, 10, 11, 6, 7, 8, 4, 5, 3, 5, 3, 4};
111  const std::vector<int> itype5 = {0, 2, 3, 4, 5, 0, 1, 1, 2, 0, 1, 2, 0};
112  if (u == 0 && v == 0) {
113  cellx = 1;
114  cellt = HGCalCell::cornerCell;
115  } else if (u == 0 && (v - u) == (ncell_[type] - 1)) {
116  cellx = 2;
117  cellt = HGCalCell::cornerCell;
118  } else if ((v - u) == (ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
119  cellx = 3;
120  cellt = HGCalCell::cornerCell;
121  } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
122  cellx = 4;
123  cellt = HGCalCell::cornerCell;
124  } else if (u == (2 * ncell_[type] - 1) && (u - v) == ncell_[type]) {
125  cellx = 5;
126  cellt = HGCalCell::cornerCell;
127  } else if ((u - v) == ncell_[type] && v == 0) {
128  cellx = 6;
129  cellt = HGCalCell::cornerCell;
130  } else if (u == 0) {
131  cellx = 7;
132  cellt = HGCalCell::truncatedCell;
133  } else if ((v - u) == (ncell_[type] - 1)) {
134  cellx = 10;
135  cellt = HGCalCell::extendedCell;
136  } else if (v == (2 * ncell_[type] - 1)) {
137  cellx = 8;
138  cellt = HGCalCell::truncatedCell;
139  } else if (u == (2 * ncell_[type] - 1)) {
140  cellx = 11;
141  cellt = HGCalCell::extendedCell;
142  } else if ((u - v) == ncell_[type]) {
143  cellx = 9;
144  cellt = HGCalCell::truncatedCell;
145  } else if (v == 0) {
146  cellx = 12;
147  cellt = HGCalCell::extendedCell;
148  }
149  switch (placementIndex) {
151  cell = itype0[cellx];
152  break;
154  cell = itype1[cellx];
155  break;
157  cell = itype2[cellx];
158  break;
160  cell = itype3[cellx];
161  break;
163  cell = itype4[cellx];
164  break;
165  default:
166  cell = itype5[cellx];
167  break;
168  }
169  } else {
170  const std::vector<int> itype0 = {0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 0, 1, 2};
171  const std::vector<int> itype1 = {0, 8, 9, 10, 11, 6, 7, 4, 5, 3, 4, 5, 3};
172  const std::vector<int> itype2 = {0, 3, 4, 5, 0, 1, 2, 2, 0, 1, 1, 2, 0};
173  const std::vector<int> itype3 = {0, 10, 11, 6, 7, 8, 9, 5, 3, 4, 5, 3, 4};
174  const std::vector<int> itype4 = {0, 5, 0, 1, 2, 3, 4, 0, 1, 2, 2, 0, 1};
175  const std::vector<int> itype5 = {0, 6, 7, 8, 9, 10, 11, 3, 4, 5, 3, 4, 5};
176  if (u == 0 && v == 0) {
177  cellx = 1;
178  cellt = HGCalCell::cornerCell;
179  } else if (v == 0 && (u - v) == (ncell_[type])) {
180  cellx = 2;
181  cellt = HGCalCell::cornerCell;
182  } else if ((u - v) == (ncell_[type]) && u == (2 * ncell_[type] - 1)) {
183  cellx = 3;
184  cellt = HGCalCell::cornerCell;
185  } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
186  cellx = 4;
187  cellt = HGCalCell::cornerCell;
188  } else if (v == (2 * ncell_[type] - 1) && (v - u) == (ncell_[type] - 1)) {
189  cellx = 5;
190  cellt = HGCalCell::cornerCell;
191  } else if ((v - u) == (ncell_[type] - 1) && u == 0) {
192  cellx = 6;
193  cellt = HGCalCell::cornerCell;
194  } else if (v == 0) {
195  cellx = 10;
196  cellt = HGCalCell::extendedCell;
197  } else if ((u - v) == ncell_[type]) {
198  cellx = 7;
199  cellt = HGCalCell::truncatedCell;
200  } else if (u == (2 * ncell_[type] - 1)) {
201  cellx = 11;
202  cellt = HGCalCell::extendedCell;
203  } else if (v == (2 * ncell_[type] - 1)) {
204  cellx = 8;
205  cellt = HGCalCell::truncatedCell;
206  } else if ((v - u) == (ncell_[type] - 1)) {
207  cellx = 12;
208  cellt = HGCalCell::extendedCell;
209  } else if (u == 0) {
210  cellx = 9;
211  cellt = HGCalCell::truncatedCell;
212  }
213  switch (placementIndex) {
215  cell = itype0[cellx];
216  break;
218  cell = itype1[cellx];
219  break;
221  cell = itype2[cellx];
222  break;
224  cell = itype3[cellx];
225  break;
227  cell = itype4[cellx];
228  break;
229  default:
230  cell = itype5[cellx];
231  break;
232  }
233  }
234  return std::make_pair(cell, cellt);
235 }
236 
237 int HGCalCell::cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient) {
238  int32_t indx = ((iz * frontBack) > 0) ? orient : (orient + HGCalCell::cellPlacementExtra);
239  return indx;
240 }
241 
242 std::pair<int32_t, int32_t> HGCalCell::cellOrient(int32_t placementIndex) {
243  int32_t orient = (placementIndex >= HGCalCell::cellPlacementExtra) ? (placementIndex - HGCalCell::cellPlacementExtra)
244  : placementIndex;
245  int32_t frontBackZside = (placementIndex >= HGCalCell::cellPlacementExtra) ? 1 : -1;
246  return std::make_pair(orient, frontBackZside);
247 }
248 
249 std::pair<int32_t, int32_t> HGCalCell::cellType(int32_t u, int32_t v, int32_t ncell, int32_t placementIndex) {
250  int cell(0), cellx(0), cellt(HGCalCell::fullCell);
251  if (placementIndex >= HGCalCell::cellPlacementExtra) {
252  const std::vector<int> itype0 = {HGCalCell::centralCell,
265  const std::vector<int> itype1 = {HGCalCell::centralCell,
278  const std::vector<int> itype2 = {HGCalCell::centralCell,
291  const std::vector<int> itype3 = {HGCalCell::centralCell,
304  const std::vector<int> itype4 = {HGCalCell::centralCell,
317  const std::vector<int> itype5 = {HGCalCell::centralCell,
330  if (u == 0 && v == 0) {
331  cellx = 1;
332  cellt = HGCalCell::cornerCell;
333  } else if (u == 0 && (v - u) == (ncell - 1)) {
334  cellx = 2;
335  cellt = HGCalCell::cornerCell;
336  } else if ((v - u) == (ncell - 1) && v == (2 * ncell - 1)) {
337  cellx = 3;
338  cellt = HGCalCell::cornerCell;
339  } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
340  cellx = 4;
341  cellt = HGCalCell::cornerCell;
342  } else if (u == (2 * ncell - 1) && (u - v) == ncell) {
343  cellx = 5;
344  cellt = HGCalCell::cornerCell;
345  } else if ((u - v) == ncell && v == 0) {
346  cellx = 6;
347  cellt = HGCalCell::cornerCell;
348  } else if (u == 0) {
349  cellx = 7;
350  cellt = HGCalCell::truncatedCell;
351  } else if ((v - u) == (ncell - 1)) {
352  cellx = 10;
353  cellt = HGCalCell::extendedCell;
354  } else if (v == (2 * ncell - 1)) {
355  cellx = 8;
356  cellt = HGCalCell::truncatedCell;
357  } else if (u == (2 * ncell - 1)) {
358  cellx = 11;
359  cellt = HGCalCell::extendedCell;
360  } else if ((u - v) == ncell) {
361  cellx = 9;
362  cellt = HGCalCell::truncatedCell;
363  } else if (v == 0) {
364  cellx = 12;
365  cellt = HGCalCell::extendedCell;
366  }
367  switch (placementIndex) {
369  cell = itype0[cellx];
370  break;
372  cell = itype1[cellx];
373  break;
375  cell = itype2[cellx];
376  break;
378  cell = itype3[cellx];
379  break;
381  cell = itype4[cellx];
382  break;
383  default:
384  cell = itype5[cellx];
385  break;
386  }
387  } else {
388  const std::vector<int> itype0 = {HGCalCell::centralCell,
401  const std::vector<int> itype1 = {HGCalCell::centralCell,
414  const std::vector<int> itype2 = {HGCalCell::centralCell,
427  const std::vector<int> itype3 = {HGCalCell::centralCell,
440  const std::vector<int> itype4 = {HGCalCell::centralCell,
453  const std::vector<int> itype5 = {HGCalCell::centralCell,
466  if (u == 0 && v == 0) {
467  cellx = 1;
468  cellt = HGCalCell::cornerCell;
469  } else if (v == 0 && (u - v) == (ncell)) {
470  cellx = 2;
471  cellt = HGCalCell::cornerCell;
472  } else if ((u - v) == (ncell) && u == (2 * ncell - 1)) {
473  cellx = 3;
474  cellt = HGCalCell::cornerCell;
475  } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
476  cellx = 4;
477  cellt = HGCalCell::cornerCell;
478  } else if (v == (2 * ncell - 1) && (v - u) == (ncell - 1)) {
479  cellx = 5;
480  cellt = HGCalCell::cornerCell;
481  } else if ((v - u) == (ncell - 1) && u == 0) {
482  cellx = 6;
483  cellt = HGCalCell::cornerCell;
484  } else if (v == 0) {
485  cellx = 10;
486  cellt = HGCalCell::extendedCell;
487  } else if ((u - v) == ncell) {
488  cellx = 7;
489  cellt = HGCalCell::truncatedCell;
490  } else if (u == (2 * ncell - 1)) {
491  cellx = 11;
492  cellt = HGCalCell::extendedCell;
493  } else if (v == (2 * ncell - 1)) {
494  cellx = 8;
495  cellt = HGCalCell::truncatedCell;
496  } else if ((v - u) == (ncell - 1)) {
497  cellx = 12;
498  cellt = HGCalCell::extendedCell;
499  } else if (u == 0) {
500  cellx = 9;
501  cellt = HGCalCell::truncatedCell;
502  }
503  switch (placementIndex) {
505  cell = itype0[cellx];
506  break;
508  cell = itype1[cellx];
509  break;
511  cell = itype2[cellx];
512  break;
514  cell = itype3[cellx];
515  break;
517  cell = itype4[cellx];
518  break;
519  default:
520  cell = itype5[cellx];
521  break;
522  }
523  }
524  return std::make_pair(cell, cellt);
525 }
Log< level::Info, true > LogVerbatim
static int32_t cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient)
Definition: HGCalCell.cc:237
static constexpr int32_t fullCell
Definition: HGCalCell.h:28
const double sqrt3By2_
Definition: HGCalCell.h:60
static constexpr int32_t cellPlacementIndex8
Definition: HGCalCell.h:19
static std::pair< int32_t, int32_t > cellType(int32_t u, int32_t v, int32_t ncell, int32_t placementIndex)
Definition: HGCalCell.cc:249
static constexpr int32_t topCorner
Definition: HGCalCell.h:44
static constexpr int32_t cellPlacementIndex3
Definition: HGCalCell.h:14
static constexpr int32_t leftEdge
Definition: HGCalCell.h:36
int32_t ncell_[2]
Definition: HGCalCell.h:61
static constexpr int32_t cellPlacementIndex10
Definition: HGCalCell.h:21
std::pair< int32_t, int32_t > cellUV2Cell(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
Definition: HGCalCell.cc:101
static constexpr int32_t cellPlacementIndex0
Definition: HGCalCell.h:11
static constexpr int32_t topLeftCorner
Definition: HGCalCell.h:43
static std::pair< int32_t, int32_t > cellOrient(int32_t placementIndex)
Definition: HGCalCell.cc:242
static constexpr int32_t bottomLeftCorner
Definition: HGCalCell.h:42
static constexpr int32_t topRightCorner
Definition: HGCalCell.h:45
static constexpr int32_t bottomRightEdge
Definition: HGCalCell.h:40
static constexpr int32_t cellPlacementIndex9
Definition: HGCalCell.h:20
HGCalCell(double waferSize, int32_t nFine, int32_t nCoarse)
Definition: HGCalCell.cc:7
static constexpr int32_t bottomLeftEdge
Definition: HGCalCell.h:35
static constexpr int32_t cellPlacementIndex7
Definition: HGCalCell.h:18
static constexpr int32_t bottomRightCorner
Definition: HGCalCell.h:46
static constexpr int32_t truncatedCell
Definition: HGCalCell.h:30
static constexpr int32_t centralCell
Definition: HGCalCell.h:34
std::pair< double, double > cellUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
Definition: HGCalCell.cc:20
static constexpr int32_t extendedCell
Definition: HGCalCell.h:31
static constexpr int32_t topRightEdge
Definition: HGCalCell.h:38
static constexpr int32_t cellPlacementIndex11
Definition: HGCalCell.h:22
std::pair< double, double > cellUV2XY2(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
Definition: HGCalCell.cc:77
static constexpr int32_t rightEdge
Definition: HGCalCell.h:39
static constexpr int32_t bottomCorner
Definition: HGCalCell.h:41
static constexpr int32_t cornerCell
Definition: HGCalCell.h:29
static constexpr int32_t cellPlacementIndex4
Definition: HGCalCell.h:15
static constexpr int32_t cellPlacementIndex1
Definition: HGCalCell.h:12
double cellX_[2]
Definition: HGCalCell.h:62
static constexpr int32_t cellPlacementIndex6
Definition: HGCalCell.h:17
double cellY_[2]
Definition: HGCalCell.h:62
static constexpr int32_t cellPlacementIndex2
Definition: HGCalCell.h:13
static constexpr int32_t topLeftEdge
Definition: HGCalCell.h:37
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24