CMS 3D CMS Logo

HGCalCell.cc
Go to the documentation of this file.
4 #include <vector>
5 #include <iostream>
6 
7 //#define EDM_ML_DEBUG
8 
9 HGCalCell::HGCalCell(double waferSize, int32_t nFine, int32_t nCoarse) {
10  ncell_[0] = nFine;
11  ncell_[1] = nCoarse;
12  for (int k = 0; k < 2; ++k) {
13  cellX_[k] = waferSize / (3 * ncell_[k]);
14  cellY_[k] = sqrt3By2_ * cellX_[k];
15  }
16 #ifdef EDM_ML_DEBUG
17  edm::LogVerbatim("HGCalGeom") << "HGCalCell initialized with waferSize " << waferSize << " number of cells " << nFine
18  << ":" << nCoarse;
19 #endif
20 }
21 
22 std::pair<double, double> HGCalCell::cellUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
23  if (type != 0)
24  type = 1;
25  double x(0), y(0);
26  switch (placementIndex) {
28  x = (1.5 * (v - u) + 0.5) * cellX_[type];
29  y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
30  break;
32  x = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
33  y = (2 * u - v - ncell_[type]) * cellY_[type];
34  break;
36  x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
37  y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
38  break;
40  x = -(1.5 * (v - u) + 0.5) * cellX_[type];
41  y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
42  break;
44  x = -(1.5 * (v - ncell_[type]) + 1) * cellX_[type];
45  y = -(2 * u - v - ncell_[type]) * cellY_[type];
46  break;
48  x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
49  y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
50  break;
52  x = (1.5 * (u - v) - 0.5) * cellX_[type];
53  y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
54  break;
56  x = -(1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
57  y = (2 * u - v - ncell_[type]) * cellY_[type];
58  break;
60  x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
61  y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
62  break;
64  x = -(1.5 * (u - v) - 0.5) * cellX_[type];
65  y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
66  break;
68  x = (1.5 * (v - ncell_[type]) + 1) * cellX_[type];
69  y = -(2 * u - v - ncell_[type]) * cellY_[type];
70  break;
71  default:
72  x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
73  y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
74  break;
75  }
76  return std::make_pair(x, y);
77 }
78 
79 std::pair<double, double> HGCalCell::cellUV2XY2(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
80  if (type != 0)
81  type = 1;
82  double x(0), y(0);
83  if (placementIndex < HGCalCell::cellPlacementExtra) {
84  double x0 = (1.5 * (u - v) - 0.5) * cellX_[type];
85  double y0 = (u + v - 2 * ncell_[type] + 1) * cellY_[type];
86  const std::vector<double> fcos = {1.0, 0.5, -0.5, -1.0, -0.5, 0.5};
87  const std::vector<double> fsin = {0.0, sqrt3By2_, sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_};
88  x = x0 * fcos[placementIndex] - y0 * fsin[placementIndex];
89  y = x0 * fsin[placementIndex] + y0 * fcos[placementIndex];
90  } else {
91  double x0 = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
92  double y0 = (2 * u - v - ncell_[type]) * cellY_[type];
93  const std::vector<double> fcos = {0.5, 1.0, 0.5, -0.5, -1.0, -0.5};
94  const std::vector<double> fsin = {sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_, 0.0, sqrt3By2_};
95  x = x0 * fcos[placementIndex - HGCalCell::cellPlacementExtra] -
96  y0 * fsin[placementIndex - HGCalCell::cellPlacementExtra];
97  y = x0 * fsin[placementIndex - HGCalCell::cellPlacementExtra] +
98  y0 * fcos[placementIndex - HGCalCell::cellPlacementExtra];
99  }
100  return std::make_pair(x, y);
101 }
102 
103 std::pair<int, int> HGCalCell::cellUV2Cell(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
104  if (type != 0)
105  type = 1;
106  int cell(0), cellx(0), cellt(HGCalCell::fullCell);
107  if (placementIndex >= HGCalCell::cellPlacementExtra) {
108  const std::vector<int> itype0 = {0, 7, 8, 9, 10, 11, 6, 3, 4, 5, 4, 5, 3};
109  const std::vector<int> itype1 = {0, 0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 1, 2};
110  const std::vector<int> itype2 = {0, 11, 6, 7, 8, 9, 10, 5, 3, 4, 3, 4, 5};
111  const std::vector<int> itype3 = {0, 4, 5, 0, 1, 2, 3, 2, 0, 1, 2, 0, 1};
112  const std::vector<int> itype4 = {0, 9, 10, 11, 6, 7, 8, 4, 5, 3, 5, 3, 4};
113  const std::vector<int> itype5 = {0, 2, 3, 4, 5, 0, 1, 1, 2, 0, 1, 2, 0};
114  if (u == 0 && v == 0) {
115  cellx = 1;
116  cellt = HGCalCell::cornerCell;
117  } else if (u == 0 && (v - u) == (ncell_[type] - 1)) {
118  cellx = 2;
119  cellt = HGCalCell::cornerCell;
120  } else if ((v - u) == (ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
121  cellx = 3;
122  cellt = HGCalCell::cornerCell;
123  } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
124  cellx = 4;
125  cellt = HGCalCell::cornerCell;
126  } else if (u == (2 * ncell_[type] - 1) && (u - v) == ncell_[type]) {
127  cellx = 5;
128  cellt = HGCalCell::cornerCell;
129  } else if ((u - v) == ncell_[type] && v == 0) {
130  cellx = 6;
131  cellt = HGCalCell::cornerCell;
132  } else if (u == 0) {
133  cellx = 7;
134  cellt = HGCalCell::truncatedCell;
135  } else if ((v - u) == (ncell_[type] - 1)) {
136  cellx = 10;
137  cellt = HGCalCell::extendedCell;
138  } else if (v == (2 * ncell_[type] - 1)) {
139  cellx = 8;
140  cellt = HGCalCell::truncatedCell;
141  } else if (u == (2 * ncell_[type] - 1)) {
142  cellx = 11;
143  cellt = HGCalCell::extendedCell;
144  } else if ((u - v) == ncell_[type]) {
145  cellx = 9;
146  cellt = HGCalCell::truncatedCell;
147  } else if (v == 0) {
148  cellx = 12;
149  cellt = HGCalCell::extendedCell;
150  }
151  switch (placementIndex) {
153  cell = itype0[cellx];
154  break;
156  cell = itype1[cellx];
157  break;
159  cell = itype2[cellx];
160  break;
162  cell = itype3[cellx];
163  break;
165  cell = itype4[cellx];
166  break;
167  default:
168  cell = itype5[cellx];
169  break;
170  }
171  } else {
172  const std::vector<int> itype0 = {0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 0, 1, 2};
173  const std::vector<int> itype1 = {0, 8, 9, 10, 11, 6, 7, 4, 5, 3, 4, 5, 3};
174  const std::vector<int> itype2 = {0, 3, 4, 5, 0, 1, 2, 2, 0, 1, 1, 2, 0};
175  const std::vector<int> itype3 = {0, 10, 11, 6, 7, 8, 9, 5, 3, 4, 5, 3, 4};
176  const std::vector<int> itype4 = {0, 5, 0, 1, 2, 3, 4, 0, 1, 2, 2, 0, 1};
177  const std::vector<int> itype5 = {0, 6, 7, 8, 9, 10, 11, 3, 4, 5, 3, 4, 5};
178  if (u == 0 && v == 0) {
179  cellx = 1;
180  cellt = HGCalCell::cornerCell;
181  } else if (v == 0 && (u - v) == (ncell_[type])) {
182  cellx = 2;
183  cellt = HGCalCell::cornerCell;
184  } else if ((u - v) == (ncell_[type]) && u == (2 * ncell_[type] - 1)) {
185  cellx = 3;
186  cellt = HGCalCell::cornerCell;
187  } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
188  cellx = 4;
189  cellt = HGCalCell::cornerCell;
190  } else if (v == (2 * ncell_[type] - 1) && (v - u) == (ncell_[type] - 1)) {
191  cellx = 5;
192  cellt = HGCalCell::cornerCell;
193  } else if ((v - u) == (ncell_[type] - 1) && u == 0) {
194  cellx = 6;
195  cellt = HGCalCell::cornerCell;
196  } else if (v == 0) {
197  cellx = 10;
198  cellt = HGCalCell::extendedCell;
199  } else if ((u - v) == ncell_[type]) {
200  cellx = 7;
201  cellt = HGCalCell::truncatedCell;
202  } else if (u == (2 * ncell_[type] - 1)) {
203  cellx = 11;
204  cellt = HGCalCell::extendedCell;
205  } else if (v == (2 * ncell_[type] - 1)) {
206  cellx = 8;
207  cellt = HGCalCell::truncatedCell;
208  } else if ((v - u) == (ncell_[type] - 1)) {
209  cellx = 12;
210  cellt = HGCalCell::extendedCell;
211  } else if (u == 0) {
212  cellx = 9;
213  cellt = HGCalCell::truncatedCell;
214  }
215  switch (placementIndex) {
217  cell = itype0[cellx];
218  break;
220  cell = itype1[cellx];
221  break;
223  cell = itype2[cellx];
224  break;
226  cell = itype3[cellx];
227  break;
229  cell = itype4[cellx];
230  break;
231  default:
232  cell = itype5[cellx];
233  break;
234  }
235  }
236  return std::make_pair(cell, cellt);
237 }
238 
239 int HGCalCell::cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient) {
240  int32_t indx = ((iz * frontBack) > 0) ? orient : (orient + HGCalCell::cellPlacementExtra);
241  return indx;
242 }
243 
244 std::pair<int32_t, int32_t> HGCalCell::cellOrient(int32_t placementIndex) {
245  int32_t orient = (placementIndex >= HGCalCell::cellPlacementExtra) ? (placementIndex - HGCalCell::cellPlacementExtra)
246  : placementIndex;
247  int32_t frontBackZside = (placementIndex >= HGCalCell::cellPlacementExtra) ? 1 : -1;
248  return std::make_pair(orient, frontBackZside);
249 }
250 
251 std::pair<int32_t, int32_t> HGCalCell::cellType(int32_t u, int32_t v, int32_t ncell, int32_t placementIndex) {
252  int cell(0), cellx(0), cellt(HGCalCell::fullCell);
253  if (placementIndex >= HGCalCell::cellPlacementExtra) {
254  const std::vector<int> itype0 = {HGCalCell::centralCell,
267  const std::vector<int> itype1 = {HGCalCell::centralCell,
280  const std::vector<int> itype2 = {HGCalCell::centralCell,
293  const std::vector<int> itype3 = {HGCalCell::centralCell,
306  const std::vector<int> itype4 = {HGCalCell::centralCell,
319  const std::vector<int> itype5 = {HGCalCell::centralCell,
332  if (u == 0 && v == 0) {
333  cellx = 1;
334  cellt = HGCalCell::cornerCell;
335  } else if (u == 0 && (v - u) == (ncell - 1)) {
336  cellx = 2;
337  cellt = HGCalCell::cornerCell;
338  } else if ((v - u) == (ncell - 1) && v == (2 * ncell - 1)) {
339  cellx = 3;
340  cellt = HGCalCell::cornerCell;
341  } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
342  cellx = 4;
343  cellt = HGCalCell::cornerCell;
344  } else if (u == (2 * ncell - 1) && (u - v) == ncell) {
345  cellx = 5;
346  cellt = HGCalCell::cornerCell;
347  } else if ((u - v) == ncell && v == 0) {
348  cellx = 6;
349  cellt = HGCalCell::cornerCell;
350  } else if (u == 0) {
351  cellx = 7;
352  cellt = HGCalCell::truncatedCell;
353  if (v == 1) {
354  cellx = 1;
356  } else if (v == ncell - 2) {
357  cellx = 2;
359  }
360  } else if ((v - u) == (ncell - 1)) {
361  cellx = 8;
362  cellt = HGCalCell::extendedCell;
363  if (v == ncell) {
364  cellx = 2;
366  } else if (v == 2 * ncell - 2) {
367  cellx = 3;
369  }
370  } else if (v == (2 * ncell - 1)) {
371  cellx = 9;
372  cellt = HGCalCell::truncatedCell;
373  if (u == ncell + 1) {
374  cellx = 3;
376  } else if (u == 2 * ncell - 2) {
377  cellx = 4;
379  }
380  } else if (u == (2 * ncell - 1)) {
381  cellx = 10;
382  cellt = HGCalCell::extendedCell;
383  if (v == 2 * ncell - 2) {
384  cellx = 4;
386  } else if (v == ncell + 1) {
387  cellx = 5;
389  }
390  } else if ((u - v) == ncell) {
391  cellx = 11;
392  cellt = HGCalCell::truncatedCell;
393  if (u == 2 * ncell - 2) {
394  cellx = 5;
396  } else if (u == ncell + 1) {
397  cellx = 6;
399  }
400  } else if (v == 0) {
401  cellx = 12;
402  cellt = HGCalCell::extendedCell;
403  if (u == ncell - 1) {
404  cellx = 6;
406  } else if (u == 1) {
407  cellx = 1;
409  }
410  }
411  switch (placementIndex) {
413  cell = itype0[cellx];
414  break;
416  cell = itype1[cellx];
417  break;
419  cell = itype2[cellx];
420  break;
422  cell = itype3[cellx];
423  break;
425  cell = itype4[cellx];
426  break;
427  default:
428  cell = itype5[cellx];
429  break;
430  }
431  } else {
432  const std::vector<int> itype0 = {HGCalCell::centralCell,
445  const std::vector<int> itype1 = {HGCalCell::centralCell,
458  const std::vector<int> itype2 = {HGCalCell::centralCell,
471  const std::vector<int> itype3 = {HGCalCell::centralCell,
484  const std::vector<int> itype4 = {HGCalCell::centralCell,
497  const std::vector<int> itype5 = {HGCalCell::centralCell,
510  if (u == 0 && v == 0) {
511  cellx = 1;
512  cellt = HGCalCell::cornerCell;
513  } else if (v == 0 && (u - v) == (ncell)) {
514  cellx = 2;
515  cellt = HGCalCell::cornerCell;
516  } else if ((u - v) == (ncell) && u == (2 * ncell - 1)) {
517  cellx = 3;
518  cellt = HGCalCell::cornerCell;
519  } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
520  cellx = 4;
521  cellt = HGCalCell::cornerCell;
522  } else if (v == (2 * ncell - 1) && (v - u) == (ncell - 1)) {
523  cellx = 5;
524  cellt = HGCalCell::cornerCell;
525  } else if ((v - u) == (ncell - 1) && u == 0) {
526  cellx = 6;
527  cellt = HGCalCell::cornerCell;
528  } else if (v == 0) {
529  cellx = 7;
530  cellt = HGCalCell::extendedCell;
531  if (u == 1) {
532  cellx = 1;
534  } else if (u == ncell - 1) {
535  cellx = 2;
537  }
538  } else if ((u - v) == ncell) {
539  cellx = 8;
540  cellt = HGCalCell::truncatedCell;
541  if (u == 2 * ncell - 2) {
542  cellx = 3;
544  } else if (u == ncell + 1) {
545  cellx = 2;
547  }
548  } else if (u == (2 * ncell - 1)) {
549  cellx = 9;
550  cellt = HGCalCell::extendedCell;
551  if (v == 2 * ncell - 2) {
552  cellx = 4;
554  } else if (v == ncell + 1) {
555  cellx = 3;
557  }
558  } else if (v == (2 * ncell - 1)) {
559  cellx = 10;
560  cellt = HGCalCell::truncatedCell;
561  if (u == ncell + 1) {
562  cellx = 5;
564  } else if (u == 2 * ncell - 2) {
565  cellx = 4;
567  }
568  } else if ((v - u) == (ncell - 1)) {
569  cellx = 11;
570  cellt = HGCalCell::extendedCell;
571  if (v == ncell) {
572  cellx = 6;
574  } else if (v == 2 * ncell - 2) {
575  cellx = 5;
577  }
578  } else if (u == 0) {
579  cellx = 12;
580  cellt = HGCalCell::truncatedCell;
581  if (v == 1) {
582  cellx = 1;
584  } else if (v == ncell - 2) {
585  cellx = 6;
587  }
588  }
589  switch (placementIndex) {
591  cell = itype0[cellx];
592  break;
594  cell = itype1[cellx];
595  break;
597  cell = itype2[cellx];
598  break;
600  cell = itype3[cellx];
601  break;
603  cell = itype4[cellx];
604  break;
605  default:
606  cell = itype5[cellx];
607  break;
608  }
609  }
610  return std::make_pair(cell, cellt);
611 }
612 
613 std::pair<int32_t, int32_t> HGCalCell::cellType(
614  int32_t u, int32_t v, int32_t ncell, int32_t placementIndex, int32_t partialType) {
615  std::pair<int, int> cell = HGCalCell::cellType(u, v, ncell, placementIndex);
616  int cellx = cell.first;
617  int cellt = cell.second;
618  if ((partialType >= HGCalTypes::WaferPartLDOffset) &&
620  if ((u == 7 && v == 14) || (u == 7 && v == 0)) {
622  if (u == 7 && v == 0) {
623  cellx = HGCalCell::leftCell;
624  } else {
625  cellx = HGCalCell::rightCell;
626  }
627  } else if ((u == 8 && v == 15) || (u == 8 && v == 0)) {
629  if (u == 8 && v == 0) {
630  cellx = HGCalCell::leftCell;
631  } else {
632  cellx = HGCalCell::rightCell;
633  }
634  } else if (u == 2 && v == 9) {
636  } else if (u == 0 && v == 7) {
638  } else if (u == 14 && v == 15) {
640  } else if (u == 15 && v == 15) {
642  } else if (u == 7) {
643  cellt = HGCalCell::extendedCell;
644  cellx = HGCalCell::topCell;
645  } else if (u == 8) {
646  cellt = HGCalCell::truncatedCell;
647  cellx = HGCalCell::bottomCell;
648  } else if ((partialType == HGCalTypes::WaferLDLeft) || (partialType == HGCalTypes::WaferLDRight) ||
649  (partialType == HGCalTypes::WaferLDFive) || (partialType == HGCalTypes::WaferLDThree)) {
650  if (((u == 15 && v == 11) || (u == 7 && v == 7)) &&
651  ((partialType == HGCalTypes::WaferLDLeft) || (partialType == HGCalTypes::WaferLDRight))) {
652  cellt = HGCalCell::halfExtCell;
653  if (partialType == HGCalTypes::WaferLDLeft) {
654  cellx = HGCalCell::leftCell;
655  } else {
656  cellx = HGCalCell::rightCell;
657  }
658  } else if ((u == 7 && v == 11) &&
659  ((partialType == HGCalTypes::WaferLDFive) || (partialType == HGCalTypes::WaferLDThree))) {
660  cellt = HGCalCell::halfExtCell;
661  if (partialType == HGCalTypes::WaferLDFive) {
662  cellx = HGCalCell::leftCell;
663  } else {
664  cellx = HGCalCell::rightCell;
665  }
666  } else if (2 * v - u == 7) {
667  cellt = HGCalCell::halfCell;
668  if (partialType == HGCalTypes::WaferLDLeft) {
669  cellx = HGCalCell::leftCell;
670  } else if (partialType == HGCalTypes::WaferLDRight) {
671  cellx = HGCalCell::rightCell;
672  }
673  } else if (2 * v - u == 15) {
674  if (partialType == HGCalTypes::WaferLDFive) {
675  cellx = HGCalCell::leftCell;
676  } else if (partialType == HGCalTypes::WaferLDThree) {
677  cellx = HGCalCell::rightCell;
678  }
679  cellt = HGCalCell::halfCell;
680  std::cout << u << ":" << v << " 1" << std::endl;
681  }
682  }
683  } else if ((partialType >= HGCalTypes::WaferPartHDOffset) &&
685  if ((u == 9 && v == 20) || (u == 9 && v == 0)) {
687  if (u == 9 && v == 0) {
688  cellx = HGCalCell::leftCell;
689  } else {
690  cellx = HGCalCell::rightCell;
691  }
692  } else if ((u == 10 && v == 21) || (u == 10 && v == 0)) {
694  if (u == 10 && v == 0) {
695  cellx = HGCalCell::leftCell;
696  } else {
697  cellx = HGCalCell::rightCell;
698  }
699  } else if (u == 9) {
700  cellt = HGCalCell::truncatedCell;
701  cellx = HGCalCell::topCell;
702  } else if (u == 10) {
703  cellt = HGCalCell::extendedCell;
704  cellx = HGCalCell::bottomCell;
705  } else if ((partialType == HGCalTypes::WaferHDLeft) || (partialType == HGCalTypes::WaferHDRight) ||
706  (partialType == HGCalTypes::WaferHDFive)) {
707  if ((u == 10 && v == 7) || (u == 10 && v == 14)) {
708  cellt = HGCalCell::halfExtCell;
709  if ((partialType == HGCalTypes::WaferHDLeft) || (partialType == HGCalTypes::WaferHDFive)) {
710  cellx = HGCalCell::leftCell;
711  } else {
712  cellx = HGCalCell::rightCell;
713  }
714  } else if ((u == 0 && v == 2) || (u == 0 && v == 9)) {
715  cellt = HGCalCell::halfTrunCell;
716  if ((partialType == HGCalTypes::WaferHDLeft) || (partialType == HGCalTypes::WaferHDFive)) {
717  cellx = HGCalCell::leftCell;
718  } else {
719  cellx = HGCalCell::rightCell;
720  }
721  } else if (2 * v - u == 4) {
722  cellt = HGCalCell::halfCell;
723  if (partialType == HGCalTypes::WaferHDLeft) {
724  cellx = HGCalCell::leftCell;
725  }
726  } else if (2 * v - u == 18) {
727  cellt = HGCalCell::halfCell;
728  if (partialType == HGCalTypes::WaferHDFive) {
729  cellx = HGCalCell::leftCell;
730  } else if (partialType == HGCalTypes::WaferHDRight) {
731  cellx = HGCalCell::rightCell;
732  }
733  }
734  }
735  }
736  return std::make_pair(cellx, cellt);
737 }
static constexpr int32_t LDPartial0007Cell
Definition: HGCalCell.h:43
Log< level::Info, true > LogVerbatim
static int32_t cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient)
Definition: HGCalCell.cc:239
static constexpr int32_t fullCell
Definition: HGCalCell.h:28
static constexpr int32_t WaferPartLDOffset
Definition: HGCalTypes.h:57
const double sqrt3By2_
Definition: HGCalCell.h:86
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:251
static constexpr int32_t WaferLDThree
Definition: HGCalTypes.h:50
static constexpr int32_t topCorner
Definition: HGCalCell.h:62
static constexpr int32_t HDPartial1021Cell
Definition: HGCalCell.h:49
static constexpr int32_t extendedMBCell
Definition: HGCalCell.h:33
static constexpr int32_t cellPlacementIndex3
Definition: HGCalCell.h:14
static constexpr int32_t leftEdge
Definition: HGCalCell.h:54
int32_t ncell_[2]
Definition: HGCalCell.h:87
static constexpr int32_t HDPartial0920Cell
Definition: HGCalCell.h:48
static constexpr int32_t cellPlacementIndex10
Definition: HGCalCell.h:21
static constexpr int32_t WaferHDFive
Definition: HGCalTypes.h:55
static constexpr int32_t truncatedMBCell
Definition: HGCalCell.h:32
std::pair< int32_t, int32_t > cellUV2Cell(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
Definition: HGCalCell.cc:103
static constexpr int32_t leftCell
Definition: HGCalCell.h:66
static constexpr int32_t LDPartial0209Cell
Definition: HGCalCell.h:42
static constexpr int32_t WaferHDLeft
Definition: HGCalTypes.h:53
static constexpr int32_t cellPlacementIndex0
Definition: HGCalCell.h:11
static constexpr int32_t WaferPartLDCount
Definition: HGCalTypes.h:59
static constexpr int32_t topLeftCorner
Definition: HGCalCell.h:61
static constexpr int32_t topCell
Definition: HGCalCell.h:68
static std::pair< int32_t, int32_t > cellOrient(int32_t placementIndex)
Definition: HGCalCell.cc:244
static constexpr int32_t bottomLeftCorner
Definition: HGCalCell.h:60
static constexpr int32_t WaferLDRight
Definition: HGCalTypes.h:48
static constexpr int32_t topRightCorner
Definition: HGCalCell.h:63
static constexpr int32_t bottomRightEdge
Definition: HGCalCell.h:58
static constexpr int32_t cellPlacementIndex9
Definition: HGCalCell.h:20
static constexpr int32_t WaferLDFive
Definition: HGCalTypes.h:49
static constexpr int32_t WaferLDLeft
Definition: HGCalTypes.h:47
static constexpr int32_t halfCell
Definition: HGCalCell.h:36
HGCalCell(double waferSize, int32_t nFine, int32_t nCoarse)
Definition: HGCalCell.cc:9
static constexpr int32_t LDPartial0714Cell
Definition: HGCalCell.h:41
static constexpr int32_t LDPartial1515Cell
Definition: HGCalCell.h:46
static constexpr int32_t bottomLeftEdge
Definition: HGCalCell.h:53
static constexpr int32_t cellPlacementIndex7
Definition: HGCalCell.h:18
static constexpr int32_t WaferPartHDOffset
Definition: HGCalTypes.h:58
static constexpr int32_t WaferHDRight
Definition: HGCalTypes.h:54
static constexpr int32_t bottomRightCorner
Definition: HGCalCell.h:64
static constexpr int32_t truncatedCell
Definition: HGCalCell.h:30
static constexpr int32_t halfTrunCell
Definition: HGCalCell.h:37
static constexpr int32_t centralCell
Definition: HGCalCell.h:52
std::pair< double, double > cellUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
Definition: HGCalCell.cc:22
static constexpr int32_t extendedCell
Definition: HGCalCell.h:31
static constexpr int32_t topRightEdge
Definition: HGCalCell.h:56
static constexpr int32_t cellPlacementIndex11
Definition: HGCalCell.h:22
static constexpr int32_t LDPartial0815Cell
Definition: HGCalCell.h:44
std::pair< double, double > cellUV2XY2(int32_t u, int32_t v, int32_t placementIndex, int32_t type)
Definition: HGCalCell.cc:79
static constexpr int32_t LDPartial1415Cell
Definition: HGCalCell.h:45
static constexpr int32_t halfExtCell
Definition: HGCalCell.h:38
static constexpr int32_t rightEdge
Definition: HGCalCell.h:57
static constexpr int32_t bottomCorner
Definition: HGCalCell.h:59
static constexpr int32_t cornerCell
Definition: HGCalCell.h:29
static constexpr int32_t rightCell
Definition: HGCalCell.h:67
static constexpr int32_t cellPlacementIndex4
Definition: HGCalCell.h:15
static constexpr int32_t cellPlacementIndex1
Definition: HGCalCell.h:12
double cellX_[2]
Definition: HGCalCell.h:88
static constexpr int32_t cellPlacementIndex6
Definition: HGCalCell.h:17
double cellY_[2]
Definition: HGCalCell.h:88
static constexpr int32_t cellPlacementIndex2
Definition: HGCalCell.h:13
static constexpr int32_t topLeftEdge
Definition: HGCalCell.h:55
static constexpr int32_t bottomCell
Definition: HGCalCell.h:69
static constexpr int32_t WaferPartHDCount
Definition: HGCalTypes.h:60
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24