CMS 3D CMS Logo

HGCalWaferMask.cc
Go to the documentation of this file.
5 
6 #include <algorithm>
7 #include <sstream>
8 
9 //#define EDM_ML_DEBUG
10 
11 bool HGCalWaferMask::maskCell(int u, int v, int n, int ncor, int fcor, int corners) {
12  /*
13 Masks each cell (or not) according to its wafer and cell position (detId) and to the user needs (corners).
14 Each wafer has k_CornerSize corners which are defined in anti-clockwise order starting from the corner at the top, which is always #0. 'ncor' denotes the number of corners inside the physical region. 'fcor' is the defined to be the first corner that appears inside the detector's physical volume in anti-clockwise order.
15 The argument 'corners' controls the types of wafers the user wants: for instance, corners=3 masks all wafers that have at least 3 corners inside the physical region.
16  */
17  bool mask(false);
18  if (ncor < corners) {
19  mask = true;
20  } else {
21  if (ncor == HGCalGeomTools::k_fourCorners) {
22  switch (fcor) {
23  case (0): {
24  mask = (v >= n);
25  break;
26  }
27  case (1): {
28  mask = (u >= n);
29  break;
30  }
31  case (2): {
32  mask = (u > v);
33  break;
34  }
35  case (3): {
36  mask = (v < n);
37  break;
38  }
39  case (4): {
40  mask = (u < n);
41  break;
42  }
43  default: {
44  mask = (u <= v);
45  break;
46  }
47  }
48  } else {
49  switch (fcor) {
50  case (0): {
51  if (ncor == HGCalGeomTools::k_threeCorners) {
52  mask = !((u > 2 * v) && (v < n));
53  } else {
54  mask = ((u >= n) && (v >= n) && ((u + v) > (3 * n - 2)));
55  }
56  break;
57  }
58  case (1): {
59  if (ncor == HGCalGeomTools::k_threeCorners) {
60  mask = !((u + v) < n);
61  } else {
62  mask = ((u >= n) && (u > v) && ((2 * u - v) > 2 * n));
63  }
64  break;
65  }
66  case (2): {
67  if (ncor == HGCalGeomTools::k_threeCorners) {
68  mask = !((u < n) && (v > u) && (v > (2 * u - 1)));
69  } else {
70  mask = ((u > 2 * v) && (v < n));
71  }
72  break;
73  }
74  case (3): {
75  if (ncor == HGCalGeomTools::k_threeCorners) {
76  mask = !((v >= u) && ((2 * v - u) > (2 * n - 2)));
77  } else {
78  mask = ((u + v) < n);
79  }
80  break;
81  }
82  case (4): {
83  if (ncor == HGCalGeomTools::k_threeCorners) {
84  mask = !((u >= n) && (v >= n) && ((u + v) > (3 * n - 2)));
85  } else {
86  mask = ((u < n) && (v > u) && (v > (2 * u - 1)));
87  }
88  break;
89  }
90  default: {
91  if (ncor == HGCalGeomTools::k_threeCorners) {
92  mask = !((u >= n) && (u > v) && ((2 * u - v) > 2 * n));
93  } else {
94  mask = ((v >= u) && ((2 * v - u) > (2 * n - 2)));
95  }
96  break;
97  }
98  }
99  }
100  }
101 #ifdef EDM_ML_DEBUG
102  edm::LogVerbatim("HGCalGeom") << "Corners: " << ncor << ":" << fcor << " N " << n << " u " << u << " v " << v
103  << " Mask " << mask;
104 #endif
105  return mask;
106 }
107 
108 bool HGCalWaferMask::goodCell(int u, int v, int n, int type, int rotn) {
109  // for V15 and V16
110  bool good(false);
111  int n2 = n / 2;
112  int n4 = n / 4;
113  int n3 = (n + 1) / 3;
114  switch (type) {
115  case (HGCalTypes::WaferFull): { //WaferFull
116  good = true;
117  break;
118  }
119  case (HGCalTypes::WaferFive): { //WaferFive
120  switch (rotn) {
121  case (HGCalTypes::WaferCorner0): {
122  int u2 = u / 2;
123  good = ((v - u2) <= n);
124  break;
125  }
126  case (HGCalTypes::WaferCorner1): {
127  good = ((v + u) < (3 * n));
128  break;
129  }
130  case (HGCalTypes::WaferCorner2): {
131  int v2 = (v + 1) / 2;
132  good = ((u - v2) <= n);
133  break;
134  }
135  case (HGCalTypes::WaferCorner3): {
136  int u2 = (u - 1) / 2;
137  good = (u2 <= v);
138  break;
139  }
140  case (HGCalTypes::WaferCorner4): {
141  good = ((v + u) >= n - 1);
142  break;
143  }
144  default: {
145  int v2 = v / 2;
146  good = (u >= v2);
147  break;
148  }
149  }
150  break;
151  }
152  case (HGCalTypes::WaferChopTwo): { //WaferChopTwo
153  switch (rotn) {
154  case (HGCalTypes::WaferCorner0): {
155  good = (v < (3 * n2));
156  break;
157  }
158  case (HGCalTypes::WaferCorner1): {
159  good = (u < (3 * n2));
160  break;
161  }
162  case (HGCalTypes::WaferCorner2): {
163  good = ((u - v) <= n2);
164  break;
165  }
166  case (HGCalTypes::WaferCorner3): {
167  good = (v >= n2);
168  break;
169  }
170  case (HGCalTypes::WaferCorner4): {
171  good = (u >= n2 - 1);
172  break;
173  }
174  default: {
175  good = ((v - u) < n2);
176  break;
177  }
178  }
179  break;
180  }
181  case (HGCalTypes::WaferChopTwoM): { //WaferChopTwoM
182  switch (rotn) {
183  case (HGCalTypes::WaferCorner0): {
184  good = (v < (5 * n4));
185  break;
186  }
187  case (HGCalTypes::WaferCorner1): {
188  good = (u <= (5 * n4));
189  break;
190  }
191  case (HGCalTypes::WaferCorner2): {
192  good = ((u - v) <= n4);
193  break;
194  }
195  case (HGCalTypes::WaferCorner3): {
196  good = (v >= (3 * n4 - 1));
197  break;
198  }
199  case (HGCalTypes::WaferCorner4): {
200  good = (u >= (3 * n4));
201  break;
202  }
203  default: {
204  good = ((v - u) <= n4);
205  break;
206  }
207  }
208  break;
209  }
210  case (HGCalTypes::WaferHalf): { //WaferHalf
211  switch (rotn) {
212  case (HGCalTypes::WaferCorner0): {
213  good = (v < n);
214  break;
215  }
216  case (HGCalTypes::WaferCorner1): {
217  good = (u <= n);
218  break;
219  }
220  case (HGCalTypes::WaferCorner2): {
221  good = (v >= u);
222  break;
223  }
224  case (HGCalTypes::WaferCorner3): {
225  good = (v >= n - 1);
226  break;
227  }
228  case (HGCalTypes::WaferCorner4): {
229  good = (u >= n);
230  break;
231  }
232  default: {
233  good = (u >= v);
234  break;
235  }
236  }
237  break;
238  }
239  case (HGCalTypes::WaferSemi): { //WaferSemi
240  switch (rotn) {
241  case (HGCalTypes::WaferCorner0): {
242  good = ((u + v) <= (2 * n));
243  break;
244  }
245  case (HGCalTypes::WaferCorner1): {
246  good = ((2 * u - v) <= (n + 1));
247  break;
248  }
249  case (HGCalTypes::WaferCorner2): {
250  good = ((2 * v - u) >= (n - 2));
251  break;
252  }
253  case (HGCalTypes::WaferCorner3): {
254  good = ((u + v) >= (2 * n - 2));
255  break;
256  }
257  case (HGCalTypes::WaferCorner4): {
258  good = ((2 * u - v) >= (n - 1));
259  break;
260  }
261  default: {
262  good = ((2 * v - u) <= n);
263  break;
264  }
265  }
266  break;
267  }
268  case (HGCalTypes::WaferThree): { //WaferThree
269  switch (rotn) {
270  case (HGCalTypes::WaferCorner0): {
271  good = ((v + u) <= n);
272  break;
273  }
274  case (HGCalTypes::WaferCorner1): {
275  good = ((2 * u - v) <= 1);
276  break;
277  }
278  case (HGCalTypes::WaferCorner2): {
279  int u2 = ((u > 0) ? (u / 2) : 0);
280  int uv = v - u2;
281  good = (uv >= (n - 1));
282  break;
283  }
284  case (HGCalTypes::WaferCorner3): {
285  good = ((v + u) >= (3 * n - 2));
286  break;
287  }
288  case (HGCalTypes::WaferCorner4): {
289  int uv = 2 * u - v;
290  good = (uv >= (2 * n - 1));
291  break;
292  }
293  default: {
294  int uv = u - 2 * v;
295  good = (uv >= 0);
296  break;
297  }
298  }
299  break;
300  }
301  case (HGCalTypes::WaferSemi2): { //WaferSemi2
302  switch (rotn) {
303  case (HGCalTypes::WaferCorner0): {
304  good = ((u + v) <= (4 * n3 + 1));
305  break;
306  }
307  case (HGCalTypes::WaferCorner1): {
308  good = ((2 * u - v) <= n2);
309  break;
310  }
311  case (HGCalTypes::WaferCorner2): {
312  int u2 = ((u + 1) / 2);
313  good = ((v - u2) >= (3 * n4 - 1));
314  break;
315  }
316  case (HGCalTypes::WaferCorner3): {
317  good = ((u + v) >= (5 * n2 - 1));
318  break;
319  }
320  case (HGCalTypes::WaferCorner4): {
321  good = ((2 * u - v) >= (3 * n2));
322  break;
323  }
324  default: {
325  int u2 = (u + 1) / 2;
326  good = ((v - u2) < n4);
327  break;
328  }
329  }
330  break;
331  }
332  case (HGCalTypes::WaferFive2): { //WaferFive2
333  switch (rotn) {
334  case (HGCalTypes::WaferCorner0): {
335  good = ((2 * v - u) <= (3 * n2));
336  break;
337  }
338  case (HGCalTypes::WaferCorner1): {
339  good = ((u + v) < (5 * n2));
340  break;
341  }
342  case (HGCalTypes::WaferCorner2): {
343  good = ((2 * u - v) >= (3 * n2));
344  break;
345  }
346  case (HGCalTypes::WaferCorner3): {
347  good = ((2 * v - u) >= n3);
348  break;
349  }
350  case (HGCalTypes::WaferCorner4): {
351  good = ((u + v) > (4 * n3));
352  break;
353  }
354  default: {
355  good = ((2 * u - v) >= n2);
356  break;
357  }
358  }
359  break;
360  }
361  case (HGCalTypes::WaferHalf2): { //WaferHalf2
362  switch (rotn) {
363  case (HGCalTypes::WaferCorner0): {
364  good = (v <= (3 * n4));
365  break;
366  }
367  case (HGCalTypes::WaferCorner1): {
368  good = (u <= (3 * n4));
369  break;
370  }
371  case (HGCalTypes::WaferCorner2): {
372  good = ((v - u) >= n4 - 1);
373  break;
374  }
375  case (HGCalTypes::WaferCorner3): {
376  good = (v >= (5 * n4 - 1));
377  break;
378  }
379  case (HGCalTypes::WaferCorner4): {
380  good = (u >= (5 * n4 - 1));
381  break;
382  }
383  default: {
384  good = ((u - v) >= n4);
385  break;
386  }
387  }
388  break;
389  }
390  }
391 #ifdef EDM_ML_DEBUG
392  edm::LogVerbatim("HGCalGeom") << "u|v " << u << ":" << v << " N " << n << " type " << type << " rot " << rotn
393  << " good " << good;
394 #endif
395  return good;
396 }
397 
398 bool HGCalWaferMask::goodCell(int u, int v, int waferType) {
399  // for V17
400  bool good(false);
401  switch (waferType) {
402  case (HGCalTypes::WaferFull): { //WaferFull
403  good = true;
404  break;
405  }
406  case (HGCalTypes::WaferLDTop): {
408  break;
409  }
410  case (HGCalTypes::WaferLDBottom): {
413  break;
414  }
415  case (HGCalTypes::WaferLDLeft): {
416  good =
418  break;
419  }
420  case (HGCalTypes::WaferLDRight): {
423  break;
424  }
425  case (HGCalTypes::WaferLDFive): {
426  good =
428  break;
429  }
430  case (HGCalTypes::WaferLDThree): {
433  break;
434  }
435  case (HGCalTypes::WaferHDTop): {
437  break;
438  }
439  case (HGCalTypes::WaferHDBottom): {
442  break;
443  }
444  case (HGCalTypes::WaferHDLeft): {
445  good =
447  break;
448  }
449  case (HGCalTypes::WaferHDRight): {
452  break;
453  }
454  case (HGCalTypes::WaferHDFive): {
455  good =
457  break;
458  }
459  }
460 #ifdef EDM_ML_DEBUG
461  edm::LogVerbatim("HGCalGeom") << "u|v " << u << ":" << v << " WaferType " << waferType << " good " << good;
462 #endif
463  return good;
464 }
465 
466 int HGCalWaferMask::getRotation(int zside, int type, int rotn) {
467  // Needs extension for V17
468  if (rotn >= HGCalTypes::WaferCornerMax)
470  int newrotn(rotn);
471  if ((zside < 0) && (type != HGCalTypes::WaferFull)) {
472  if ((type == HGCalTypes::WaferFive) || (type == HGCalTypes::WaferFive2)) { //WaferFive/WaferFive2
473  static constexpr int rot1[HGCalTypes::WaferCornerMax] = {HGCalTypes::WaferCorner4,
479  newrotn = rot1[rotn];
480  } else if ((type == HGCalTypes::WaferThree) || (type == HGCalTypes::WaferSemi) ||
481  (type == HGCalTypes::WaferSemi2)) { //WaferThree/WaferSemi/WaferSemi2
482  static constexpr int rot2[HGCalTypes::WaferCornerMax] = {HGCalTypes::WaferCorner2,
488  newrotn = rot2[rotn];
489  } else { //WaferHalf/WaferChopTwo/WaferChopTwoM/WaferHalf2
490  static constexpr int rot3[HGCalTypes::WaferCornerMax] = {HGCalTypes::WaferCorner3,
496  newrotn = rot3[rotn];
497  }
498  }
499 #ifdef EDM_ML_DEBUG
500  edm::LogVerbatim("HGCalGeom") << "zside " << zside << " type " << type << " rotn " << rotn << ":" << newrotn;
501 #endif
502  return newrotn;
503 }
504 
505 std::pair<int, int> HGCalWaferMask::getTypeMode(const double& xpos,
506  const double& ypos,
507  const double& delX,
508  const double& delY,
509  const double& rin,
510  const double& rout,
511  const int& wType,
512  const int& mode,
513  const bool& v17OrLess,
514  const bool& debug) {
515  // No need to extend this for V17 -- use flat file information only
516  int ncor(0), iok(0);
518  double c22(HGCalTypes::c22), c27(HGCalTypes::c27), c61(HGCalTypes::c61);
519  double c77(HGCalTypes::c77), c88(HGCalTypes::c88);
520  if (v17OrLess) {
521  c22 = HGCalTypes::c22O;
522  c27 = HGCalTypes::c27O;
523  c61 = HGCalTypes::c61O;
524  c77 = HGCalTypes::c77O;
525  c88 = HGCalTypes::c88O;
526  }
527 
528  static constexpr int corners = 6;
529  static constexpr int base = 10;
530  double rin2 = rin * rin;
531  double rout2 = rout * rout;
532  double dx0[corners] = {HGCalTypes::c00 * delX,
533  HGCalTypes::c10 * delX,
534  HGCalTypes::c10 * delX,
535  HGCalTypes::c00 * delX,
536  -HGCalTypes::c10 * delX,
537  -HGCalTypes::c10 * delX};
538  double dy0[corners] = {-HGCalTypes::c10 * delY,
539  -HGCalTypes::c50 * delY,
540  HGCalTypes::c50 * delY,
541  HGCalTypes::c10 * delY,
542  HGCalTypes::c50 * delY,
543  -HGCalTypes::c50 * delY};
544  double xc[corners], yc[corners];
545  for (int k = 0; k < corners; ++k) {
546  xc[k] = xpos + dx0[k];
547  yc[k] = ypos + dy0[k];
548  double rpos2 = (xc[k] * xc[k] + yc[k] * yc[k]);
549  if (rpos2 <= rout2 && rpos2 >= rin2) {
550  ++ncor;
551  iok = iok * base + 1;
552  } else {
553  iok *= base;
554  }
555  }
556  if (debug)
557  edm::LogVerbatim("HGCalGeom") << "I/p " << xpos << ":" << ypos << ":" << delX << ":" << delY << ":" << rin << ":"
558  << rout << ":" << wType << ":" << mode << " Corners " << ncor << " iok " << iok;
559 
560  static constexpr int ipat5[corners] = {101111, 110111, 111011, 111101, 111110, 11111};
561  static constexpr int ipat4[corners] = {100111, 110011, 111001, 111100, 11110, 1111};
562  static constexpr int ipat3[corners] = {100011, 110001, 111000, 11100, 1110, 111};
563  static constexpr int ipat2[corners] = {11, 100001, 110000, 11000, 1100, 110};
564  double dx1[corners] = {HGCalTypes::c50 * delX,
565  HGCalTypes::c10 * delX,
566  HGCalTypes::c50 * delX,
567  -HGCalTypes::c50 * delX,
568  -HGCalTypes::c10 * delX,
569  -HGCalTypes::c50 * delX};
570  double dy1[corners] = {-HGCalTypes::c75 * delY,
571  HGCalTypes::c00 * delY,
572  HGCalTypes::c75 * delY,
573  HGCalTypes::c75 * delY,
574  HGCalTypes::c00 * delY,
575  -HGCalTypes::c75 * delY};
576  double dx2[corners] = {HGCalTypes::c50 * delX,
577  -HGCalTypes::c50 * delX,
578  -HGCalTypes::c10 * delX,
579  -HGCalTypes::c50 * delX,
580  HGCalTypes::c50 * delX,
581  HGCalTypes::c10 * delX};
582  double dy2[corners] = {HGCalTypes::c75 * delY,
583  HGCalTypes::c75 * delY,
584  HGCalTypes::c00 * delY,
585  -HGCalTypes::c75 * delY,
586  -HGCalTypes::c75 * delY,
587  HGCalTypes::c00 * delY};
588  double dx3[corners] = {
589  c22 * delX, HGCalTypes::c10 * delX, c77 * delX, -c22 * delX, -HGCalTypes::c10 * delX, -c77 * delX};
590  double dy3[corners] = {-c88 * delY, -c27 * delY, c61 * delY, c88 * delY, c27 * delY, -c61 * delY};
591  double dx4[corners] = {
592  c22 * delX, -c77 * delX, -HGCalTypes::c10 * delX, -c22 * delX, c77 * delX, HGCalTypes::c10 * delX};
593  double dy4[corners] = {c88 * delY, c61 * delY, -c27 * delY, -c88 * delY, -c61 * delY, c27 * delY};
594  double dx5[corners] = {-HGCalTypes::c50 * delX,
595  -HGCalTypes::c10 * delX,
596  -HGCalTypes::c50 * delX,
597  HGCalTypes::c50 * delX,
598  HGCalTypes::c10 * delX,
599  HGCalTypes::c50 * delX};
600  double dy5[corners] = {HGCalTypes::c75 * delY,
601  HGCalTypes::c00 * delY,
602  -HGCalTypes::c75 * delY,
603  -HGCalTypes::c75 * delY,
604  HGCalTypes::c00 * delY,
605  HGCalTypes::c75 * delY};
606  double dx6[corners] = {
607  -c77 * delX, -HGCalTypes::c10 * delX, -c22 * delX, c77 * delX, HGCalTypes::c10 * delX, c22 * delX};
608  double dy6[corners] = {c61 * delY, -c27 * delY, -c88 * delY, -c61 * delY, c27 * delY, c88 * delY};
609  double dx7[corners] = {
610  -c22 * delX, -HGCalTypes::c10 * delX, -c77 * delX, c22 * delX, HGCalTypes::c10 * delX, c77 * delX};
611  double dy7[corners] = {c88 * delY, c27 * delY, -c61 * delY, -c88 * delY, -c27 * delY, c61 * delY};
612  double dx8[corners] = {
613  c77 * delX, HGCalTypes::c10 * delX, c22 * delX, -c77 * delX, -HGCalTypes::c10 * delX, -c22 * delX};
614  double dy8[corners] = {-c61 * delY, c27 * delY, c88 * delY, c61 * delY, -c27 * delY, -c88 * delY};
615  double dx9[corners] = {
616  -c22 * delX, c77 * delX, HGCalTypes::c10 * delX, c22 * delX, -c77 * delX, -HGCalTypes::c10 * delX};
617  double dy9[corners] = {-c88 * delY, -c61 * delY, c27 * delY, c88 * delY, c61 * delY, -c27 * delY};
618 
619  if (ncor == HGCalGeomTools::k_allCorners) {
620  } else if (ncor == HGCalGeomTools::k_fiveCorners) {
621  rotn = static_cast<int>(std::find(ipat5, ipat5 + 6, iok) - ipat5);
623  } else if (ncor == HGCalGeomTools::k_fourCorners) {
624  rotn = static_cast<int>(std::find(ipat4, ipat4 + 6, iok) - ipat4);
626  double rpos12 = ((xpos + dx1[rotn]) * (xpos + dx1[rotn]) + (ypos + dy1[rotn]) * (ypos + dy1[rotn]));
627  double rpos22(0);
628  if (rpos12 <= rout2 && rpos12 >= rin2) {
629  rpos22 = ((xpos + dx2[rotn]) * (xpos + dx2[rotn]) + (ypos + dy2[rotn]) * (ypos + dy2[rotn]));
630  if (rpos22 <= rout2 && rpos22 >= rin2)
632  }
633  if (debug)
634  edm::LogVerbatim("HGCalGeom") << "Test for Chop2 " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
635  << type;
636  if ((type == HGCalTypes::WaferHalf) && (wType == 0)) {
637  rpos12 = ((xpos + dx3[rotn]) * (xpos + dx3[rotn]) + (ypos + dy3[rotn]) * (ypos + dy3[rotn]));
638  if (rpos12 <= rout2 && rpos12 >= rin2) {
639  rpos22 = ((xpos + dx4[rotn]) * (xpos + dx4[rotn]) + (ypos + dy4[rotn]) * (ypos + dy4[rotn]));
640  if (rpos22 <= rout2 && rpos22 >= rin2)
642  }
643  if (debug)
644  edm::LogVerbatim("HGCalGeom") << "Test for Chop2M " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
645  << type;
646  }
647  } else if (ncor == HGCalGeomTools::k_threeCorners) {
648  rotn = static_cast<int>(std::find(ipat3, ipat3 + 6, iok) - ipat3);
650  double rpos12 = ((xpos + dx7[rotn]) * (xpos + dx7[rotn]) + (ypos + dy7[rotn]) * (ypos + dy7[rotn]));
651  double rpos22(0);
652  if (rpos12 <= rout2 && rpos12 >= rin2) {
653  rpos22 = ((xpos + dx8[rotn]) * (xpos + dx8[rotn]) + (ypos + dy8[rotn]) * (ypos + dy8[rotn]));
654  if (rpos22 <= rout2 && rpos22 >= rin2)
656  }
657  if (debug)
658  edm::LogVerbatim("HGCalGeom") << "Test for Five2 " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
659  << type;
660  if ((type == HGCalTypes::WaferThree) && (wType == 0)) {
661  rpos12 = ((xpos + dx1[rotn]) * (xpos + dx1[rotn]) + (ypos + dy1[rotn]) * (ypos + dy1[rotn]));
662  if (rpos12 <= rout2 && rpos12 >= rin2) {
663  rpos22 = ((xpos + dx5[rotn]) * (xpos + dx5[rotn]) + (ypos + dy5[rotn]) * (ypos + dy5[rotn]));
664  if (rpos22 <= rout2 && rpos22 >= rin2)
666  }
667  if (debug)
668  edm::LogVerbatim("HGCalGeom") << "Test for Semi " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
669  << type;
670  }
671  if ((type == HGCalTypes::WaferThree) && (wType == 0)) {
672  rpos12 = ((xpos + dx3[rotn]) * (xpos + dx3[rotn]) + (ypos + dy3[rotn]) * (ypos + dy3[rotn]));
673  if (rpos12 <= rout2 && rpos12 >= rin2) {
674  rpos22 = ((xpos + dx6[rotn]) * (xpos + dx6[rotn]) + (ypos + dy6[rotn]) * (ypos + dy6[rotn]));
675  if (rpos22 <= rout2 && rpos22 >= rin2)
677  }
678  if (debug)
679  edm::LogVerbatim("HGCalGeom") << "Test for SemiM " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
680  << type;
681  }
682  } else if (ncor == HGCalGeomTools::k_twoCorners) {
683  rotn = static_cast<int>(std::find(ipat2, ipat2 + 6, iok) - ipat2);
685  double rpos12 = ((xpos + dx7[rotn]) * (xpos + dx7[rotn]) + (ypos + dy7[rotn]) * (ypos + dy7[rotn]));
686  double rpos22(0);
687  if (rpos12 <= rout2 && rpos12 >= rin2) {
688  rpos22 = ((xpos + dx9[rotn]) * (xpos + dx9[rotn]) + (ypos + dy9[rotn]) * (ypos + dy9[rotn]));
689  if (rpos22 <= rout2 && rpos22 >= rin2)
691  else
693  }
694  if (debug)
695  edm::LogVerbatim("HGCalGeom") << "Test for Half2 " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
696  << type;
697  } else {
699  }
700 
701  if (debug)
702  edm::LogVerbatim("HGCalGeom") << "I/p " << xpos << ":" << ypos << ":" << delX << ":" << delY << ":" << rin << ":"
703  << rout << ":" << wType << ":" << mode << " o/p " << iok << ":" << ncor << ":" << type
704  << ":" << rotn;
705  return ((mode == 0) ? std::make_pair(ncor, rotn) : std::make_pair(type, (rotn + HGCalTypes::k_OffsetRotation)));
706 }
707 
708 bool HGCalWaferMask::goodTypeMode(const double& xpos,
709  const double& ypos,
710  const double& delX,
711  const double& delY,
712  const double& rin,
713  const double& rout,
714  const int& part,
715  const int& rotn,
716  const bool& v17OrLess,
717  const bool& debug) {
718  // Needs extension for V17 or above
719  double c22(HGCalTypes::c22), c27(HGCalTypes::c27), c61(HGCalTypes::c61);
720  double c77(HGCalTypes::c77), c88(HGCalTypes::c88);
721  if (v17OrLess) {
722  c22 = HGCalTypes::c22O;
723  c27 = HGCalTypes::c27O;
724  c61 = HGCalTypes::c61O;
725  c77 = HGCalTypes::c77O;
726  c88 = HGCalTypes::c88O;
727  }
728 
730  return false;
731  if (rotn < 0 || rotn > HGCalTypes::WaferCornerMax)
732  return false;
733  double rin2 = rin * rin;
734  double rout2 = rout * rout;
735  double rpos2(0);
736  static constexpr int corners = HGCalTypes::WaferCornerMax;
737  static constexpr int corner2 = 2 * HGCalTypes::WaferCornerMax;
738  static constexpr int base = 10;
739  static constexpr int base2 = 100;
740  double dx0[corners] = {HGCalTypes::c00 * delX,
741  HGCalTypes::c10 * delX,
742  HGCalTypes::c10 * delX,
743  HGCalTypes::c00 * delX,
744  -HGCalTypes::c10 * delX,
745  -HGCalTypes::c10 * delX};
746  double dy0[corners] = {-HGCalTypes::c10 * delY,
747  -HGCalTypes::c50 * delY,
748  HGCalTypes::c50 * delY,
749  HGCalTypes::c10 * delY,
750  HGCalTypes::c50 * delY,
751  -HGCalTypes::c50 * delY};
752  double dx1[corners] = {HGCalTypes::c50 * delX,
753  HGCalTypes::c10 * delX,
754  HGCalTypes::c50 * delX,
755  -HGCalTypes::c50 * delX,
756  -HGCalTypes::c10 * delX,
757  -HGCalTypes::c50 * delX};
758  double dy1[corners] = {-HGCalTypes::c75 * delY,
759  HGCalTypes::c00 * delY,
760  HGCalTypes::c75 * delY,
761  HGCalTypes::c75 * delY,
762  HGCalTypes::c00 * delY,
763  -HGCalTypes::c75 * delY};
764  double dx2[corner2] = {c22 * delX,
765  c77 * delX,
766  HGCalTypes::c10 * delX,
767  HGCalTypes::c10 * delX,
768  c77 * delX,
769  c22 * delX,
770  -c22 * delX,
771  -c77 * delX,
772  -HGCalTypes::c10 * delX,
773  -HGCalTypes::c10 * delX,
774  -c77 * delX,
775  -c22 * delX};
776  double dy2[corner2] = {-c88 * delY,
777  -c61 * delY,
778  -c27 * delY,
779  c27 * delY,
780  c61 * delY,
781  c88 * delY,
782  c88 * delY,
783  c61 * delY,
784  c27 * delY,
785  -c27 * delY,
786  -c61 * delY,
787  -c88 * delY};
788  bool ok(true);
789  int ncf(-1);
790  switch (part) {
791  case (HGCalTypes::WaferThree): {
792  static constexpr int nc0[corners] = {450, 150, 201, 312, 423, 534};
793  int nc = nc0[rotn];
794  for (int k1 = 0; k1 < 3; ++k1) {
795  int k = nc % base;
796  double xc1 = xpos + dx0[k];
797  double yc1 = ypos + dy0[k];
798  rpos2 = (xc1 * xc1 + yc1 * yc1);
799  if ((rpos2 > rout2) || (rpos2 < rin2)) {
800  ok = false;
801  ncf = k;
802  break;
803  }
804  nc /= base;
805  }
806  break;
807  }
808  case (HGCalTypes::WaferSemi2): {
809  static constexpr int nc10[corners] = {450, 150, 201, 312, 423, 534};
810  static constexpr int nc11[corners] = {700, 902, 1104, 106, 308, 510};
811  int nc = nc10[rotn];
812  for (int k1 = 0; k1 < 3; ++k1) {
813  int k = nc % base;
814  double xc1 = xpos + dx0[k];
815  double yc1 = ypos + dy0[k];
816  rpos2 = (xc1 * xc1 + yc1 * yc1);
817  if ((rpos2 > rout2) || (rpos2 < rin2)) {
818  ok = false;
819  ncf = k;
820  break;
821  }
822  nc /= base;
823  }
824  nc = nc11[rotn];
825  for (int k1 = 0; k1 < 2; ++k1) {
826  int k = nc % base2;
827  double xc1 = xpos + dx2[k];
828  double yc1 = ypos + dy2[k];
829  rpos2 = (xc1 * xc1 + yc1 * yc1);
830  if ((rpos2 > rout2) || (rpos2 < rin2)) {
831  ok = false;
832  ncf = k + base2;
833  break;
834  }
835  nc /= base2;
836  }
837  break;
838  }
839  case (HGCalTypes::WaferSemi): {
840  static constexpr int nc20[corners] = {450, 150, 201, 312, 423, 534};
841  static constexpr int nc21[corners] = {30, 14, 25, 30, 41, 52};
842  int nc = nc20[rotn];
843  for (int k1 = 0; k1 < 3; ++k1) {
844  int k = nc % base;
845  double xc1 = xpos + dx0[k];
846  double yc1 = ypos + dy0[k];
847  rpos2 = (xc1 * xc1 + yc1 * yc1);
848  if ((rpos2 > rout2) || (rpos2 < rin2)) {
849  ok = false;
850  ncf = k;
851  break;
852  }
853  nc /= base;
854  }
855  nc = nc21[rotn];
856  for (int k1 = 0; k1 < 2; ++k1) {
857  int k = nc % base;
858  double xc1 = xpos + dx1[k];
859  double yc1 = ypos + dy1[k];
860  rpos2 = (xc1 * xc1 + yc1 * yc1);
861  if ((rpos2 > rout2) || (rpos2 < rin2)) {
862  ok = false;
863  ncf = k + base2;
864  break;
865  }
866  nc /= base;
867  }
868  break;
869  }
870  case (HGCalTypes::WaferHalf): {
871  static constexpr int nc3[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
872  int nc = nc3[rotn];
873  for (int k1 = 0; k1 < 4; ++k1) {
874  int k = nc % base;
875  double xc1 = xpos + dx0[k];
876  double yc1 = ypos + dy0[k];
877  rpos2 = (xc1 * xc1 + yc1 * yc1);
878  if ((rpos2 > rout2) || (rpos2 < rin2)) {
879  ok = false;
880  ncf = k;
881  break;
882  }
883  nc /= base;
884  }
885  break;
886  }
887  case (HGCalTypes::WaferChopTwoM): {
888  static constexpr int nc40[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
889  static constexpr int nc41[corners] = {500, 702, 904, 1106, 108, 310};
890  int nc = nc40[rotn];
891  for (int k1 = 0; k1 < 4; ++k1) {
892  int k = nc % base;
893  double xc1 = xpos + dx0[k];
894  double yc1 = ypos + dy0[k];
895  rpos2 = (xc1 * xc1 + yc1 * yc1);
896  if ((rpos2 > rout2) || (rpos2 < rin2)) {
897  ok = false;
898  ncf = k;
899  break;
900  }
901  nc /= base;
902  }
903  nc = nc41[rotn];
904  for (int k1 = 0; k1 < 2; ++k1) {
905  int k = nc % base2;
906  double xc1 = xpos + dx2[k];
907  double yc1 = ypos + dy2[k];
908  rpos2 = (xc1 * xc1 + yc1 * yc1);
909  if ((rpos2 > rout2) || (rpos2 < rin2)) {
910  ok = false;
911  ncf = k + base2;
912  break;
913  }
914  nc /= base2;
915  }
916  break;
917  }
918  case (HGCalTypes::WaferChopTwo): {
919  static constexpr int nc50[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
920  static constexpr int nc51[corners] = {20, 13, 24, 35, 40, 51};
921  int nc = nc50[rotn];
922  for (int k1 = 0; k1 < 4; ++k1) {
923  int k = nc % base;
924  double xc1 = xpos + dx0[k];
925  double yc1 = ypos + dy0[k];
926  rpos2 = (xc1 * xc1 + yc1 * yc1);
927  if ((rpos2 > rout2) || (rpos2 < rin2)) {
928  ok = false;
929  ncf = k;
930  break;
931  }
932  nc /= base;
933  }
934  nc = nc51[rotn];
935  for (int k1 = 0; k1 < 2; ++k1) {
936  int k = nc % base;
937  double xc1 = xpos + dx1[k];
938  double yc1 = ypos + dy1[k];
939  rpos2 = (xc1 * xc1 + yc1 * yc1);
940  if ((rpos2 > rout2) || (rpos2 < rin2)) {
941  ok = false;
942  ncf = k + base2;
943  break;
944  }
945  nc /= base;
946  }
947  break;
948  }
949  case (HGCalTypes::WaferFive): {
950  static constexpr int nc6[corners] = {23450, 13450, 24501, 35012, 40123, 51234};
951  int nc = nc6[rotn];
952  for (int k1 = 0; k1 < 5; ++k1) {
953  int k = nc % base;
954  double xc1 = xpos + dx0[k];
955  double yc1 = ypos + dy0[k];
956  rpos2 = (xc1 * xc1 + yc1 * yc1);
957  if ((rpos2 > rout2) || (rpos2 < rin2)) {
958  ok = false;
959  ncf = k;
960  break;
961  }
962  }
963  break;
964  }
965  case (HGCalTypes::WaferFive2): {
966  static constexpr int nc60[corners] = {450, 150, 201, 312, 423, 534};
967  static constexpr int nc61[corners] = {601, 803, 1005, 7, 209, 411};
968  int nc = nc60[rotn];
969  for (int k1 = 0; k1 < 3; ++k1) {
970  int k = nc % base;
971  double xc1 = xpos + dx0[k];
972  double yc1 = ypos + dy0[k];
973  rpos2 = (xc1 * xc1 + yc1 * yc1);
974  if ((rpos2 > rout2) || (rpos2 < rin2)) {
975  ok = false;
976  ncf = k;
977  break;
978  }
979  nc /= base;
980  }
981  nc = nc61[rotn];
982  for (int k1 = 0; k1 < 2; ++k1) {
983  int k = nc % base2;
984  double xc1 = xpos + dx2[k];
985  double yc1 = ypos + dy2[k];
986  rpos2 = (xc1 * xc1 + yc1 * yc1);
987  if ((rpos2 > rout2) || (rpos2 < rin2)) {
988  ok = false;
989  ncf = k + base2;
990  break;
991  }
992  nc /= base2;
993  }
994  break;
995  }
996  case (HGCalTypes::WaferHalf2): {
997  static constexpr int nc70[corners] = {45, 50, 1, 12, 23, 34};
998  static constexpr int nc71[corners] = {611, 801, 1003, 5, 207, 409};
999  int nc = nc70[rotn];
1000  for (int k1 = 0; k1 < 2; ++k1) {
1001  int k = nc % base;
1002  double xc1 = xpos + dx0[k];
1003  double yc1 = ypos + dy0[k];
1004  rpos2 = (xc1 * xc1 + yc1 * yc1);
1005  if ((rpos2 > rout2) || (rpos2 < rin2)) {
1006  ok = false;
1007  ncf = k;
1008  break;
1009  }
1010  nc /= base;
1011  }
1012  nc = nc71[rotn];
1013  for (int k1 = 0; k1 < 2; ++k1) {
1014  int k = nc % base2;
1015  double xc1 = xpos + dx2[k];
1016  double yc1 = ypos + dy2[k];
1017  rpos2 = (xc1 * xc1 + yc1 * yc1);
1018  if ((rpos2 > rout2) || (rpos2 < rin2)) {
1019  ok = false;
1020  ncf = k + base2;
1021  break;
1022  }
1023  nc /= base2;
1024  }
1025  break;
1026  }
1027  default: {
1028  for (int k = 0; k < corners; ++k) {
1029  double xc1 = xpos + dx0[k];
1030  double yc1 = ypos + dy0[k];
1031  rpos2 = (xc1 * xc1 + yc1 * yc1);
1032  if ((rpos2 > rout2) || (rpos2 < rin2)) {
1033  ok = false;
1034  ncf = k;
1035  break;
1036  }
1037  }
1038  break;
1039  }
1040  }
1041  if (debug || (!ok))
1042  edm::LogVerbatim("HGCalGeom") << "I/p "
1043  << ":" << xpos << ":" << ypos << ":" << delX << ":" << delY << ":" << rin << ":"
1044  << rout << ":" << part << ":" << rotn << " Results " << ok << ":" << ncf << " R "
1045  << rin2 << ":" << rout2 << ":" << rpos2;
1046  return ok;
1047 }
1048 
1049 std::vector<std::pair<double, double> > HGCalWaferMask::waferXY(const int& part,
1050  const int& ori,
1051  const int& zside,
1052  const double& waferSize,
1053  const double& offset,
1054  const double& xpos,
1055  const double& ypos,
1056  const bool& v17OrLess) {
1057  // Good for V15 and V16 versions
1058  std::vector<std::pair<double, double> > xy;
1059  int orient = getRotation(-zside, part, ori);
1060 #ifdef EDM_ML_DEBUG
1061  edm::LogVerbatim("HGCalGeom") << "Part " << part << " zSide " << zside << " Orient " << ori << ":" << orient;
1062 #endif
1063  double c22(HGCalTypes::c22), c27(HGCalTypes::c27), c61(HGCalTypes::c61);
1064  double c77(HGCalTypes::c77), c88(HGCalTypes::c88);
1065  if (v17OrLess) {
1066  c22 = HGCalTypes::c22O;
1067  c27 = HGCalTypes::c27O;
1068  c61 = HGCalTypes::c61O;
1069  c77 = HGCalTypes::c77O;
1070  c88 = HGCalTypes::c88O;
1071  }
1072  /*
1073  The exact contour of partial wafers are obtained by joining points on
1074  the circumference of a full wafer.
1075  Numbering the points along the edges of a hexagonal wafer, starting from
1076  the bottom corner:
1077 
1078  3
1079  15 18
1080  9 8
1081  19 14
1082  4 2
1083  16 23
1084  10 7
1085  20 13
1086  5 1
1087  17 22
1088  11 6
1089  21 12
1090  0
1091 
1092  Depending on the wafer type and orientation index, the corners
1093  are chosen in the variable *np*
1094  */
1095  double delX = 0.5 * waferSize;
1096  double delY = delX / sin_60_;
1097  double dx[48] = {HGCalTypes::c00 * delX,
1098  HGCalTypes::c10 * delX,
1099  HGCalTypes::c10 * delX,
1100  HGCalTypes::c00 * delX,
1101  -HGCalTypes::c10 * delX,
1102  -HGCalTypes::c10 * delX,
1103  HGCalTypes::c50 * delX,
1104  HGCalTypes::c10 * delX,
1105  HGCalTypes::c50 * delX,
1106  -HGCalTypes::c50 * delX,
1107  -HGCalTypes::c10 * delX,
1108  -HGCalTypes::c50 * delX,
1109  c22 * delX,
1110  HGCalTypes::c10 * delX,
1111  c77 * delX,
1112  -c22 * delX,
1113  -HGCalTypes::c10 * delX,
1114  -c77 * delX,
1115  c22 * delX,
1116  -c77 * delX,
1117  -HGCalTypes::c10 * delX,
1118  -c22 * delX,
1119  c77 * delX,
1120  HGCalTypes::c10 * delX,
1121  HGCalTypes::c50 * delX,
1122  HGCalTypes::c10 * delX,
1123  HGCalTypes::c50 * delX,
1124  -HGCalTypes::c50 * delX,
1125  -HGCalTypes::c10 * delX,
1126  -HGCalTypes::c50 * delX,
1127  HGCalTypes::c50 * delX,
1128  HGCalTypes::c10 * delX,
1129  HGCalTypes::c50 * delX,
1130  -HGCalTypes::c50 * delX,
1131  -HGCalTypes::c10 * delX,
1132  -HGCalTypes::c50 * delX,
1133  c22 * delX,
1134  HGCalTypes::c10 * delX,
1135  c77 * delX,
1136  -c22 * delX,
1137  -HGCalTypes::c10 * delX,
1138  -c77 * delX,
1139  c22 * delX,
1140  -c77 * delX,
1141  -HGCalTypes::c10 * delX,
1142  -c22 * delX,
1143  c77 * delX,
1144  HGCalTypes::c10 * delX};
1145  double dy[48] = {-HGCalTypes::c10 * delY,
1146  -HGCalTypes::c50 * delY,
1147  HGCalTypes::c50 * delY,
1148  HGCalTypes::c10 * delY,
1149  HGCalTypes::c50 * delY,
1150  -HGCalTypes::c50 * delY,
1151  -HGCalTypes::c75 * delY,
1152  HGCalTypes::c00 * delY,
1153  HGCalTypes::c75 * delY,
1154  HGCalTypes::c75 * delY,
1155  HGCalTypes::c00 * delY,
1156  -HGCalTypes::c75 * delY,
1157  -c88 * delY,
1158  -c27 * delY,
1159  c61 * delY,
1160  c88 * delY,
1161  c27 * delY,
1162  -c61 * delY,
1163  c88 * delY,
1164  c61 * delY,
1165  -c27 * delY,
1166  -c88 * delY,
1167  -c61 * delY,
1168  c27 * delY,
1169  -HGCalTypes::c75 * delY,
1170  HGCalTypes::c00 * delY,
1171  -HGCalTypes::c75 * delY,
1172  HGCalTypes::c00 * delY,
1173  HGCalTypes::c75 * delY,
1174  HGCalTypes::c75 * delY,
1175  HGCalTypes::c00 * delY,
1176  -HGCalTypes::c75 * delY,
1177  HGCalTypes::c75 * delY,
1178  HGCalTypes::c75 * delY,
1179  HGCalTypes::c00 * delY,
1180  -HGCalTypes::c75 * delY,
1181  -c88 * delY,
1182  -c27 * delY,
1183  c61 * delY,
1184  c88 * delY,
1185  c27 * delY,
1186  -c61 * delY,
1187  c88 * delY,
1188  c61 * delY,
1189  -c27 * delY,
1190  -c88 * delY,
1191  -c61 * delY,
1192  c27 * delY};
1193 
1194  double offsetx[48] = {0.0,
1195  -offset,
1196  -offset,
1197  0.0,
1198  offset,
1199  offset,
1200  -offset * cos_60_,
1201  -offset,
1202  -offset * cos_60_,
1203  offset * cos_60_,
1204  offset,
1205  offset * cos_60_,
1206  -offset * cos_60_,
1207  -offset,
1208  -offset * cos_60_,
1209  offset * cos_60_,
1210  offset,
1211  offset * cos_60_,
1212  -offset * cos_60_,
1213  offset * cos_60_,
1214  offset,
1215  offset * cos_60_,
1216  -offset * cos_60_,
1217  -offset,
1218  0.0,
1219  -offset,
1220  -offset,
1221  0.0,
1222  offset,
1223  offset,
1224  0.0,
1225  offset,
1226  offset,
1227  0.0,
1228  -offset,
1229  -offset,
1230  0.0,
1231  -offset,
1232  -offset,
1233  0.0,
1234  offset,
1235  offset,
1236  0.0,
1237  offset,
1238  offset,
1239  0.0,
1240  -offset,
1241  -offset};
1242  double offsety[48] = {offset / sin_60_,
1243  offset / tan_60_,
1244  -offset / tan_60_,
1245  -offset / sin_60_,
1246  -offset / tan_60_,
1247  offset / tan_60_,
1248  offset * sin_60_,
1249  0.0,
1250  -offset * sin_60_,
1251  -offset * sin_60_,
1252  0.0,
1253  offset * sin_60_,
1254  offset * sin_60_,
1255  0.0,
1256  -offset * sin_60_,
1257  -offset * sin_60_,
1258  0.0,
1259  offset * sin_60_,
1260  -offset * sin_60_,
1261  -offset * sin_60_,
1262  0.0,
1263  offset * sin_60_,
1264  offset * sin_60_,
1265  0.0,
1266  offset / sin_60_,
1267  offset / tan_60_,
1268  -offset / tan_60_,
1269  -offset / sin_60_,
1270  -offset / tan_60_,
1271  -offset / sin_60_,
1272  -offset / tan_60_,
1273  offset / tan_60_,
1274  offset / sin_60_,
1275  offset / tan_60_,
1276  -offset / tan_60_,
1277  offset / tan_60_,
1278  offset / sin_60_,
1279  offset / tan_60_,
1280  -offset / tan_60_,
1281  -offset / sin_60_,
1282  -offset / tan_60_,
1283  offset / tan_60_,
1284  -offset / sin_60_,
1285  -offset / tan_60_,
1286  offset / tan_60_,
1287  offset / sin_60_,
1288  offset / tan_60_,
1289  -offset / tan_60_};
1290 
1291  if (part == HGCalTypes::WaferFull) {
1292  int np[7] = {0, 1, 2, 3, 4, 5, 0};
1293  for (int k = 0; k < 7; ++k)
1294  xy.push_back(std::make_pair((xpos + dx[np[k]] + offsetx[np[k]]), (ypos + dy[np[k]] + offsety[np[k]])));
1295  } else if (part == HGCalTypes::WaferFive) {
1296  int np[6][6] = {{0, 2, 3, 4, 5, 0},
1297  {1, 3, 4, 5, 0, 1},
1298  {2, 4, 5, 0, 1, 2},
1299  {3, 5, 0, 1, 2, 3},
1300  {4, 0, 1, 2, 3, 4},
1301  {5, 1, 2, 3, 4, 5}};
1302  for (int k = 0; k < 6; ++k) {
1303  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1304  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1305 #ifdef EDM_ML_DEBUG
1306  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1307  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1308 #endif
1309  }
1310  } else if (part == HGCalTypes::WaferHalf) {
1311  int np[6][5] = {
1312  {0, 3, 4, 5, 0}, {1, 4, 5, 0, 1}, {2, 5, 0, 1, 2}, {3, 0, 1, 2, 3}, {4, 1, 2, 3, 4}, {5, 2, 3, 4, 5}};
1313  for (int k = 0; k < 5; ++k) {
1314  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1315  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1316 #ifdef EDM_ML_DEBUG
1317  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1318  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1319 #endif
1320  }
1321  } else if (part == HGCalTypes::WaferThree) {
1322  int np[6][4] = {{0, 4, 5, 0}, {1, 5, 0, 1}, {2, 0, 1, 2}, {3, 1, 2, 3}, {4, 2, 3, 4}, {5, 3, 4, 5}};
1323  for (int k = 0; k < 4; ++k) {
1324  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1325  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1326 #ifdef EDM_ML_DEBUG
1327  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1328  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1329 #endif
1330  }
1331  } else if (part == HGCalTypes::WaferChopTwo) {
1332  int np[6][7] = {{24, 32, 3, 4, 5, 0, 24},
1333  {25, 33, 4, 5, 0, 1, 25},
1334  {26, 34, 5, 0, 1, 2, 26},
1335  {27, 35, 0, 1, 2, 3, 27},
1336  {28, 30, 1, 2, 3, 4, 28},
1337  {29, 31, 2, 3, 4, 5, 29}};
1338  for (int k = 0; k < 7; ++k) {
1339  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1340  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1341 #ifdef EDM_ML_DEBUG
1342  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1343  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1344 #endif
1345  }
1346  } else if (part == HGCalTypes::WaferSemi) {
1347  int np[6][6] = {{6, 9, 4, 5, 0, 6},
1348  {7, 10, 5, 0, 1, 7},
1349  {8, 11, 0, 1, 2, 8},
1350  {9, 6, 1, 2, 3, 9},
1351  {10, 7, 2, 3, 4, 10},
1352  {11, 8, 3, 4, 5, 11}};
1353  for (int k = 0; k < 6; ++k) {
1354  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1355  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1356 #ifdef EDM_ML_DEBUG
1357  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1358  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1359 #endif
1360  }
1361  } else if (part == HGCalTypes::WaferChopTwoM) {
1362  int np[6][7] = {{36, 42, 3, 4, 5, 0, 36},
1363  {37, 43, 4, 5, 0, 1, 37},
1364  {38, 44, 5, 0, 1, 2, 38},
1365  {39, 45, 0, 1, 2, 3, 39},
1366  {40, 46, 1, 2, 3, 4, 40},
1367  {41, 47, 2, 3, 4, 5, 41}};
1368  for (int k = 0; k < 7; ++k) {
1369  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1370  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1371 #ifdef EDM_ML_DEBUG
1372  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1373  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1374 #endif
1375  }
1376  } else if (part == HGCalTypes::WaferSemi2) {
1377  int np[6][6] = {{12, 19, 4, 5, 0, 12},
1378  {13, 20, 5, 0, 1, 13},
1379  {14, 21, 0, 1, 2, 14},
1380  {15, 22, 1, 2, 3, 15},
1381  {16, 23, 2, 3, 4, 16},
1382  {17, 18, 3, 4, 5, 17}};
1383  for (int k = 0; k < 6; ++k) {
1384  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1385  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1386 #ifdef EDM_ML_DEBUG
1387  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1388  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1389 #endif
1390  }
1391  } else if (part == HGCalTypes::WaferFive2) {
1392  int np[6][6] = {{22, 15, 4, 5, 0, 22},
1393  {23, 16, 5, 0, 1, 23},
1394  {18, 17, 0, 1, 2, 18},
1395  {19, 12, 1, 2, 3, 19},
1396  {20, 13, 2, 3, 4, 20},
1397  {21, 14, 3, 4, 5, 21}};
1398  for (int k = 0; k < 6; ++k) {
1399  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1400  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1401 #ifdef EDM_ML_DEBUG
1402  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1403  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1404 #endif
1405  }
1406  } else if (part == HGCalTypes::WaferHalf2) {
1407  int np[6][5] = {{45, 39, 4, 5, 45},
1408  {46, 40, 5, 0, 46},
1409  {47, 41, 0, 1, 47},
1410  {42, 36, 1, 2, 42},
1411  {43, 37, 2, 3, 43},
1412  {44, 38, 3, 4, 44}};
1413  for (int k = 0; k < 5; ++k) {
1414  xy.push_back(std::make_pair((xpos + dx[np[orient][k]] + offsetx[np[orient][k]]),
1415  (ypos + dy[np[orient][k]] + offsety[np[orient][k]])));
1416 #ifdef EDM_ML_DEBUG
1417  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] + offsetx[np[orient][k]]
1418  << ":" << dy[np[orient][k]] + offsety[np[orient][k]];
1419 #endif
1420  }
1421  }
1422 #ifdef EDM_ML_DEBUG
1423  edm::LogVerbatim("HGCalGeom") << "I/p: " << part << ":" << ori << ":" << zside << ":" << delX << ":" << delY << ":"
1424  << xpos << ":" << ypos << " O/p having " << xy.size() << " points:";
1425  std::ostringstream st1;
1426  for (unsigned int i = 0; i < xy.size(); ++i)
1427  st1 << " [" << i << "] " << xy[i].first << ":" << xy[i].second;
1428  edm::LogVerbatim("HGCalGeom") << st1.str();
1429 #endif
1430  return xy;
1431 }
1432 
1433 std::vector<std::pair<double, double> > HGCalWaferMask::waferXY(const int& part,
1434  const int& place,
1435  const double& waferSize,
1436  const double& offset,
1437  const double& xpos,
1438  const double& ypos,
1439  const bool& v17OrLess) {
1440  std::vector<std::pair<double, double> > xy;
1441  // Good for V17 version and uses partial wafer type & placement index
1442 #ifdef EDM_ML_DEBUG
1443  edm::LogVerbatim("HGCalGeom") << "Part " << part << " Placement Index " << place;
1444 #endif
1445  double c22(HGCalTypes::c22), c27(HGCalTypes::c27), c61(HGCalTypes::c61);
1446  double c77(HGCalTypes::c77), c88(HGCalTypes::c88);
1447  if (v17OrLess) {
1448  c22 = HGCalTypes::c22O;
1449  c27 = HGCalTypes::c27O;
1450  c61 = HGCalTypes::c61O;
1451  c77 = HGCalTypes::c77O;
1452  c88 = HGCalTypes::c88O;
1453  }
1454  /*
1455  The exact contour of partial wafers are obtained by joining points on
1456  the circumference of a full wafer.
1457  Numbering the points along the edges of a hexagonal wafer, starting from
1458  the bottom corner:
1459  3
1460  15 18
1461  9 8
1462  19 14
1463  4 2
1464  16 23
1465  10 7
1466  20 13
1467  5 1
1468  17 22
1469  11 6
1470  21 12
1471  0
1472  Depending on the wafer type and placement index, the corners
1473  are chosen in the variable *np*
1474  The points 24-35 are the same as points 12-23 with different offset
1475  */
1476  double delX = 0.5 * waferSize;
1477  double delY = delX / sin_60_;
1478  double dx[60] = {
1479  HGCalTypes::c00 * delX,
1480  HGCalTypes::c10 * delX,
1481  HGCalTypes::c10 * delX,
1482  HGCalTypes::c00 * delX,
1483  -HGCalTypes::c10 * delX,
1484  -HGCalTypes::c10 * delX,
1485  HGCalTypes::c50 * delX,
1486  HGCalTypes::c10 * delX,
1487  HGCalTypes::c50 * delX,
1488  -HGCalTypes::c50 * delX,
1489  -HGCalTypes::c10 * delX,
1490  -HGCalTypes::c50 * delX,
1491  c22 * delX,
1492  HGCalTypes::c10 * delX,
1493  c77 * delX,
1494  -c22 * delX,
1495  -HGCalTypes::c10 * delX,
1496  -c77 * delX,
1497  c22 * delX,
1498  -c77 * delX,
1499  -HGCalTypes::c10 * delX,
1500  -c22 * delX,
1501  c77 * delX,
1502  HGCalTypes::c10 * delX,
1503  c22 * delX,
1504  HGCalTypes::c10 * delX,
1505  c77 * delX,
1506  -c22 * delX,
1507  -HGCalTypes::c10 * delX,
1508  -c77 * delX,
1509  c22 * delX,
1510  -c77 * delX,
1511  -HGCalTypes::c10 * delX,
1512  -c22 * delX,
1513  c77 * delX,
1514  HGCalTypes::c10 * delX,
1515  HGCalTypes::c00 * delX,
1516  HGCalTypes::c10 * delX,
1517  HGCalTypes::c10 * delX,
1518  HGCalTypes::c00 * delX,
1519  -HGCalTypes::c10 * delX,
1520  -HGCalTypes::c10 * delX,
1521  HGCalTypes::c00 * delX,
1522  HGCalTypes::c10 * delX,
1523  HGCalTypes::c10 * delX,
1524  HGCalTypes::c00 * delX,
1525  -HGCalTypes::c10 * delX,
1526  -HGCalTypes::c10 * delX,
1527  HGCalTypes::c00 * delX,
1528  HGCalTypes::c10 * delX,
1529  HGCalTypes::c10 * delX,
1530  HGCalTypes::c00 * delX,
1531  -HGCalTypes::c10 * delX,
1532  -HGCalTypes::c10 * delX,
1533  HGCalTypes::c00 * delX,
1534  HGCalTypes::c10 * delX,
1535  HGCalTypes::c10 * delX,
1536  HGCalTypes::c00 * delX,
1537  -HGCalTypes::c10 * delX,
1538  -HGCalTypes::c10 * delX,
1539  };
1540  double dy[60] = {
1541  -HGCalTypes::c10 * delY,
1542  -HGCalTypes::c50 * delY,
1543  HGCalTypes::c50 * delY,
1544  HGCalTypes::c10 * delY,
1545  HGCalTypes::c50 * delY,
1546  -HGCalTypes::c50 * delY,
1547  -HGCalTypes::c75 * delY,
1548  HGCalTypes::c00 * delY,
1549  HGCalTypes::c75 * delY,
1550  HGCalTypes::c75 * delY,
1551  HGCalTypes::c00 * delY,
1552  -HGCalTypes::c75 * delY,
1553  -c88 * delY,
1554  -c27 * delY,
1555  c61 * delY,
1556  c88 * delY,
1557  c27 * delY,
1558  -c61 * delY,
1559  c88 * delY,
1560  c61 * delY,
1561  -c27 * delY,
1562  -c88 * delY,
1563  -c61 * delY,
1564  c27 * delY,
1565  -c88 * delY,
1566  -c27 * delY,
1567  c61 * delY,
1568  c88 * delY,
1569  c27 * delY,
1570  -c61 * delY,
1571  c88 * delY,
1572  c61 * delY,
1573  -c27 * delY,
1574  -c88 * delY,
1575  -c61 * delY,
1576  c27 * delY,
1577  -HGCalTypes::c10 * delY,
1578  -HGCalTypes::c50 * delY,
1579  HGCalTypes::c50 * delY,
1580  HGCalTypes::c10 * delY,
1581  HGCalTypes::c50 * delY,
1582  -HGCalTypes::c50 * delY,
1583  -HGCalTypes::c10 * delY,
1584  -HGCalTypes::c50 * delY,
1585  HGCalTypes::c50 * delY,
1586  HGCalTypes::c10 * delY,
1587  HGCalTypes::c50 * delY,
1588  -HGCalTypes::c50 * delY,
1589  -HGCalTypes::c10 * delY,
1590  -HGCalTypes::c50 * delY,
1591  HGCalTypes::c50 * delY,
1592  HGCalTypes::c10 * delY,
1593  HGCalTypes::c50 * delY,
1594  -HGCalTypes::c50 * delY,
1595  -HGCalTypes::c10 * delY,
1596  -HGCalTypes::c50 * delY,
1597  HGCalTypes::c50 * delY,
1598  HGCalTypes::c10 * delY,
1599  HGCalTypes::c50 * delY,
1600  -HGCalTypes::c50 * delY,
1601  };
1602 
1603  double offsetx[60] = {0.0,
1604  -offset,
1605  -offset,
1606  0.0,
1607  offset,
1608  offset,
1609  -offset * cos_60_,
1610  -offset,
1611  -offset * cos_60_,
1612  offset * cos_60_,
1613  offset,
1614  offset * cos_60_,
1615  -offset * cos_60_,
1616  -offset,
1617  -offset * cos_60_,
1618  offset * cos_60_,
1619  offset,
1620  offset * cos_60_,
1621  -offset * cos_60_,
1622  offset * cos_60_,
1623  offset,
1624  offset * cos_60_,
1625  -offset * cos_60_,
1626  -offset,
1627  0.0,
1628  -offset,
1629  -offset,
1630  0.0,
1631  offset,
1632  offset,
1633  0.0,
1634  offset,
1635  offset,
1636  0.0,
1637  -offset,
1638  -offset,
1639  -offset,
1640  -offset / cos_60_,
1641  -offset,
1642  offset,
1643  offset / cos_60_,
1644  offset,
1645  offset,
1646  -offset,
1647  -offset / cos_60_,
1648  -offset,
1649  offset,
1650  offset / cos_60_,
1651  -offset * cos_60_,
1652  -offset,
1653  -offset * cos_60_,
1654  offset * cos_60_,
1655  offset,
1656  offset * cos_60_,
1657  offset * cos_60_,
1658  -offset * cos_60_,
1659  -offset,
1660  -offset * cos_60_,
1661  offset * cos_60_,
1662  offset};
1663  double offsety[60] = {offset / sin_60_,
1664  offset / tan_60_,
1665  -offset / tan_60_,
1666  -offset / sin_60_,
1667  -offset / tan_60_,
1668  offset / tan_60_,
1669  offset * sin_60_,
1670  0.0,
1671  -offset * sin_60_,
1672  -offset * sin_60_,
1673  0.0,
1674  offset * sin_60_,
1675  offset * sin_60_,
1676  0.0,
1677  -offset * sin_60_,
1678  -offset * sin_60_,
1679  0.0,
1680  offset * sin_60_,
1681  -offset * sin_60_,
1682  -offset * sin_60_,
1683  0.0,
1684  offset * sin_60_,
1685  offset * sin_60_,
1686  0.0,
1687  offset / sin_60_,
1688  offset / tan_60_,
1689  -offset / tan_60_,
1690  -offset / sin_60_,
1691  -offset / tan_60_,
1692  offset / tan_60_,
1693  -offset / sin_60_,
1694  -offset / tan_60_,
1695  offset / tan_60_,
1696  offset / sin_60_,
1697  offset / tan_60_,
1698  -offset / tan_60_,
1699  offset * tan_60_,
1700  0,
1701  -offset * tan_60_,
1702  -offset * tan_60_,
1703  0,
1704  offset * tan_60_,
1705  offset * tan_60_,
1706  offset * tan_60_,
1707  0,
1708  -offset * tan_60_,
1709  -offset * tan_60_,
1710  0,
1711  offset * sin_60_,
1712  0,
1713  -offset * sin_60_,
1714  -offset * sin_60_,
1715  0,
1716  offset * sin_60_,
1717  offset * sin_60_,
1718  offset * sin_60_,
1719  0,
1720  -offset * sin_60_,
1721  -offset * sin_60_,
1722  0};
1723 
1724  if (part == HGCalTypes::WaferFull) {
1725  int np[7] = {0, 1, 2, 3, 4, 5, 0};
1726  for (int k = 0; k < 7; ++k)
1727  xy.push_back(std::make_pair((xpos + dx[np[k]] + offsetx[np[k]]), (ypos + dy[np[k]] + offsety[np[k]])));
1728  } else if (part == HGCalTypes::WaferLDTop) {
1729  int np[12][5] = {{0, 1, 4, 5, 0},
1730  {1, 2, 5, 0, 1},
1731  {2, 3, 0, 1, 2},
1732  {3, 4, 1, 2, 3},
1733  {4, 5, 2, 3, 4},
1734  {5, 0, 3, 4, 5},
1735  {0, 1, 2, 5, 0},
1736  {5, 0, 1, 4, 5},
1737  {4, 5, 0, 3, 4},
1738  {3, 4, 5, 2, 3},
1739  {2, 3, 4, 1, 2},
1740  {1, 2, 3, 0, 1}};
1741  for (int k = 0; k < 5; ++k) {
1742  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1743  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1744 #ifdef EDM_ML_DEBUG
1745  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1746  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1747 #endif
1748  }
1749  } else if (part == HGCalTypes::WaferLDBottom) {
1750  int np[12][5] = {{1, 2, 3, 4, 1},
1751  {2, 3, 4, 5, 2},
1752  {3, 4, 5, 0, 3},
1753  {4, 5, 0, 1, 4},
1754  {5, 0, 1, 2, 5},
1755  {0, 1, 2, 3, 0},
1756  {5, 2, 3, 4, 5},
1757  {4, 1, 2, 3, 4},
1758  {3, 0, 1, 2, 3},
1759  {2, 5, 0, 1, 2},
1760  {1, 4, 5, 0, 1},
1761  {0, 3, 4, 5, 0}};
1762  for (int k = 0; k < 5; ++k) {
1763  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1764  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1765 #ifdef EDM_ML_DEBUG
1766  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1767  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1768 #endif
1769  }
1770  } else if (part == HGCalTypes::WaferLDLeft) {
1771  int np[12][6] = {{0, 1, 2, 8, 11, 0},
1772  {1, 2, 3, 9, 6, 1},
1773  {2, 3, 4, 10, 7, 2},
1774  {3, 4, 5, 11, 8, 3},
1775  {4, 5, 0, 6, 9, 4},
1776  {5, 0, 1, 7, 10, 5},
1777  {0, 6, 9, 4, 5, 0},
1778  {5, 11, 8, 3, 4, 5},
1779  {4, 10, 7, 2, 3, 4},
1780  {3, 9, 6, 1, 2, 3},
1781  {2, 8, 11, 0, 1, 2},
1782  {1, 7, 10, 5, 0, 1}};
1783  for (int k = 0; k < 6; ++k) {
1784  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1785  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1786 #ifdef EDM_ML_DEBUG
1787  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1788  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1789 #endif
1790  }
1791  } else if (part == HGCalTypes::WaferLDRight) {
1792  int np[12][6] = {{5, 11, 8, 3, 4, 5},
1793  {0, 6, 9, 4, 5, 0},
1794  {1, 7, 10, 5, 0, 1},
1795  {2, 8, 11, 0, 1, 2},
1796  {3, 9, 6, 1, 2, 3},
1797  {4, 10, 7, 2, 3, 4},
1798  {1, 2, 3, 9, 6, 1},
1799  {0, 1, 2, 8, 11, 0},
1800  {5, 0, 1, 7, 10, 5},
1801  {4, 5, 0, 6, 9, 4},
1802  {3, 4, 5, 11, 8, 3},
1803  {2, 3, 4, 10, 7, 2}};
1804  for (int k = 0; k < 6; ++k) {
1805  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1806  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1807 #ifdef EDM_ML_DEBUG
1808  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1809  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1810 #endif
1811  }
1812  } else if (part == HGCalTypes::WaferLDFive) {
1813  int np[12][6] = {{0, 1, 2, 57, 53, 0},
1814  {1, 2, 3, 58, 48, 1},
1815  {2, 3, 4, 59, 49, 2},
1816  {3, 4, 5, 54, 50, 3},
1817  {4, 5, 0, 55, 51, 4},
1818  {5, 0, 1, 56, 52, 5},
1819  {0, 1, 3, 58, 53, 0},
1820  {5, 0, 2, 57, 52, 5},
1821  {4, 5, 1, 56, 51, 4},
1822  {3, 4, 0, 55, 50, 3},
1823  {2, 3, 5, 54, 49, 2},
1824  {1, 2, 4, 59, 48, 1}};
1825  for (int k = 0; k < 6; ++k) {
1826  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1827  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1828 #ifdef EDM_ML_DEBUG
1829  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1830  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1831 #endif
1832  }
1833  } else if (part == HGCalTypes::WaferLDThree) {
1834  int np[12][4] = {{41, 45, 4, 41},
1835  {36, 46, 5, 36},
1836  {37, 47, 0, 37},
1837  {38, 42, 1, 38},
1838  {39, 43, 2, 39},
1839  {40, 44, 3, 40},
1840  {43, 2, 39, 43},
1841  {42, 1, 38, 42},
1842  {47, 0, 37, 47},
1843  {46, 5, 36, 46},
1844  {45, 4, 41, 45},
1845  {44, 3, 40, 44}};
1846  for (int k = 0; k < 4; ++k) {
1847  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1848  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1849 #ifdef EDM_ML_DEBUG
1850  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1851  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1852 #endif
1853  }
1854  } else if (part == HGCalTypes::WaferHDTop) {
1855  int np[12][5] = {{0, 34, 28, 5, 0},
1856  {1, 35, 29, 0, 1},
1857  {2, 30, 24, 1, 2},
1858  {3, 31, 25, 2, 3},
1859  {4, 32, 26, 3, 4},
1860  {5, 33, 27, 4, 5},
1861  {0, 1, 35, 29, 0},
1862  {5, 0, 34, 28, 5},
1863  {4, 5, 33, 27, 4},
1864  {3, 4, 32, 26, 3},
1865  {2, 3, 31, 25, 2},
1866  {1, 2, 30, 24, 1}};
1867  for (int k = 0; k < 5; ++k) {
1868  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1869  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1870 #ifdef EDM_ML_DEBUG
1871  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1872  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1873 #endif
1874  }
1875  } else if (part == HGCalTypes::WaferHDBottom) {
1876  int np[12][7] = {{1, 2, 3, 4, 28, 34, 1},
1877  {2, 3, 4, 5, 29, 35, 2},
1878  {3, 4, 5, 0, 24, 30, 3},
1879  {4, 5, 0, 1, 25, 31, 4},
1880  {5, 0, 1, 2, 26, 32, 5},
1881  {0, 1, 2, 3, 27, 33, 0},
1882  {5, 29, 35, 2, 3, 4, 5},
1883  {4, 28, 34, 1, 2, 3, 4},
1884  {3, 27, 33, 0, 1, 2, 3},
1885  {2, 26, 32, 5, 0, 1, 2},
1886  {1, 25, 31, 4, 5, 0, 1},
1887  {0, 24, 30, 3, 4, 5, 0}};
1888  for (int k = 0; k < 7; ++k) {
1889  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1890  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1891 #ifdef EDM_ML_DEBUG
1892  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1893  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1894 #endif
1895  }
1896  } else if (part == HGCalTypes::WaferHDLeft) {
1897  int np[12][6] = {{0, 1, 2, 14, 21, 0},
1898  {1, 2, 3, 15, 22, 1},
1899  {2, 3, 4, 16, 23, 2},
1900  {3, 4, 5, 17, 18, 3},
1901  {4, 5, 0, 12, 19, 4},
1902  {5, 0, 1, 13, 20, 5},
1903  {0, 12, 19, 4, 5, 0},
1904  {5, 17, 18, 3, 4, 5},
1905  {4, 16, 23, 2, 3, 4},
1906  {3, 15, 22, 1, 2, 3},
1907  {2, 14, 21, 0, 1, 2},
1908  {1, 13, 20, 5, 0, 1}};
1909  for (int k = 0; k < 6; ++k) {
1910  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1911  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1912 #ifdef EDM_ML_DEBUG
1913  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1914  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1915 #endif
1916  }
1917  } else if (part == HGCalTypes::WaferHDRight) {
1918  int np[12][6] = {{5, 17, 18, 3, 4, 5},
1919  {0, 12, 19, 4, 5, 0},
1920  {1, 13, 20, 5, 0, 1},
1921  {2, 14, 21, 0, 1, 2},
1922  {3, 15, 22, 1, 2, 3},
1923  {4, 16, 23, 2, 3, 4},
1924  {1, 2, 3, 15, 22, 1},
1925  {0, 1, 2, 14, 21, 0},
1926  {5, 0, 1, 13, 20, 5},
1927  {4, 5, 0, 12, 19, 4},
1928  {3, 4, 5, 17, 18, 3},
1929  {2, 3, 4, 16, 23, 2}};
1930  for (int k = 0; k < 6; ++k) {
1931  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1932  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1933 #ifdef EDM_ML_DEBUG
1934  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1935  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1936 #endif
1937  }
1938  } else if (part == HGCalTypes::WaferHDFive) {
1939  int np[12][6] = {{0, 1, 2, 18, 17, 0},
1940  {1, 2, 3, 19, 12, 1},
1941  {2, 3, 4, 20, 13, 2},
1942  {3, 4, 5, 21, 14, 3},
1943  {4, 5, 0, 22, 15, 4},
1944  {5, 0, 1, 23, 16, 5},
1945  {0, 22, 15, 4, 5, 0},
1946  {5, 21, 14, 3, 4, 5},
1947  {4, 20, 13, 2, 3, 4},
1948  {3, 19, 12, 1, 2, 3},
1949  {2, 18, 17, 0, 1, 2},
1950  {1, 23, 16, 5, 0, 1}};
1951  for (int k = 0; k < 6; ++k) {
1952  xy.push_back(std::make_pair((xpos + dx[np[place][k]] + offsetx[np[place][k]]),
1953  (ypos + dy[np[place][k]] + offsety[np[place][k]])));
1954 #ifdef EDM_ML_DEBUG
1955  edm::LogVerbatim("HGCalGeom") << k << ":" << np[place][k] << ":" << dx[np[place][k]] + offsetx[np[place][k]]
1956  << ":" << dy[np[place][k]] + offsety[np[place][k]];
1957 #endif
1958  }
1959  }
1960 #ifdef EDM_ML_DEBUG
1961  edm::LogVerbatim("HGCalGeom") << "I/p: " << part << ":" << place << ":" << delX << ":" << delY << ":" << xpos << ":"
1962  << ypos << " O/p having " << xy.size() << " points:";
1963  std::ostringstream st1;
1964  for (unsigned int i = 0; i < xy.size(); ++i)
1965  st1 << " [" << i << "] " << xy[i].first << ":" << xy[i].second;
1966  edm::LogVerbatim("HGCalGeom") << st1.str();
1967 #endif
1968  return xy;
1969 }
static constexpr std::array< int, 3 > edgeWaferHDBottom
Definition: HGCalTypes.h:86
static constexpr std::array< int, 3 > edgeWaferLDFive
Definition: HGCalTypes.h:83
Log< level::Info, true > LogVerbatim
static constexpr double c27O
Definition: HGCalTypes.h:97
static constexpr int32_t WaferHalf2
Definition: HGCalTypes.h:43
static constexpr std::array< int, 3 > edgeWaferLDLeft
Definition: HGCalTypes.h:81
static constexpr int32_t WaferFive2
Definition: HGCalTypes.h:44
static constexpr std::array< int, 3 > edgeWaferLDTop
Definition: HGCalTypes.h:79
static constexpr std::array< int, 3 > edgeWaferHDLeft
Definition: HGCalTypes.h:87
static constexpr double c77O
Definition: HGCalTypes.h:103
static constexpr int k_OffsetRotation
Definition: HGCalTypes.h:91
static constexpr int32_t WaferSizeMax
Definition: HGCalTypes.h:76
static constexpr int32_t WaferLDThree
Definition: HGCalTypes.h:50
static constexpr int k_fourCorners
static constexpr std::array< int, 3 > edgeWaferLDThree
Definition: HGCalTypes.h:84
static bool goodCell(int u, int v, int N, int type, int rotn)
static constexpr int32_t WaferOut
Definition: HGCalTypes.h:56
base
Main Program
Definition: newFWLiteAna.py:92
static constexpr int k_twoCorners
static constexpr int32_t WaferThree
Definition: HGCalTypes.h:42
static constexpr double c50
Definition: HGCalTypes.h:99
static constexpr int32_t WaferHDFive
Definition: HGCalTypes.h:55
static constexpr double c88
Definition: HGCalTypes.h:106
int zside(DetId const &)
static int getRotation(int zside, int type, int rotn)
static constexpr double c27
Definition: HGCalTypes.h:98
static constexpr double c10
Definition: HGCalTypes.h:107
static constexpr double sin_60_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
static constexpr int32_t WaferSemi2
Definition: HGCalTypes.h:41
static std::vector< std::pair< double, double > > waferXY(const int &part, const int &orient, const int &zside, const double &waferSize, const double &offset, const double &xpos, const double &ypos, const bool &v17)
static constexpr int32_t WaferHDLeft
Definition: HGCalTypes.h:53
constexpr uint32_t mask
Definition: gpuClustering.h:26
U second(std::pair< T, U > const &p)
static constexpr int32_t WaferFull
Definition: HGCalTypes.h:35
static constexpr int32_t WaferHalf
Definition: HGCalTypes.h:39
static constexpr int32_t WaferLDBottom
Definition: HGCalTypes.h:46
int np
Definition: AMPTWrapper.h:43
static constexpr double tan_60_
static constexpr int32_t WaferCorner1
Definition: HGCalTypes.h:14
static std::pair< int, int > getTypeMode(const double &xpos, const double &ypos, const double &delX, const double &delY, const double &rin, const double &rout, const int &waferType, const int &mode, const bool &v17, const bool &debug=false)
T sqrt(T t)
Definition: SSEVec.h:19
static constexpr double c22O
Definition: HGCalTypes.h:94
static constexpr int32_t WaferCorner5
Definition: HGCalTypes.h:18
static constexpr int32_t WaferLDRight
Definition: HGCalTypes.h:48
static constexpr int32_t WaferCorner4
Definition: HGCalTypes.h:17
static constexpr int32_t WaferCornerMax
Definition: HGCalTypes.h:75
static constexpr int32_t WaferCorner2
Definition: HGCalTypes.h:15
static constexpr int32_t WaferLDFive
Definition: HGCalTypes.h:49
static constexpr int32_t WaferCorner0
Definition: HGCalTypes.h:13
static constexpr int32_t WaferLDLeft
Definition: HGCalTypes.h:47
static constexpr double c00
Definition: HGCalTypes.h:93
static constexpr int32_t WaferHDBottom
Definition: HGCalTypes.h:52
static constexpr int32_t WaferCorner3
Definition: HGCalTypes.h:16
static constexpr std::array< int, 3 > edgeWaferHDTop
Definition: HGCalTypes.h:85
static constexpr double c88O
Definition: HGCalTypes.h:105
static bool goodTypeMode(const double &xpos, const double &ypos, const double &delX, const double &delY, const double &rin, const double &rout, const int &part, const int &rotn, const bool &v17, const bool &debug=false)
static constexpr int32_t WaferHDRight
Definition: HGCalTypes.h:54
#define debug
Definition: HDRShower.cc:19
static constexpr int k_fiveCorners
static constexpr int32_t WaferChopTwoM
Definition: HGCalTypes.h:38
static constexpr std::array< int, 3 > edgeWaferLDRight
Definition: HGCalTypes.h:82
part
Definition: HCALResponse.h:20
static constexpr double c77
Definition: HGCalTypes.h:104
static constexpr double c61
Definition: HGCalTypes.h:101
static constexpr double cos_60_
static constexpr double c61O
Definition: HGCalTypes.h:100
static constexpr int k_threeCorners
static constexpr int32_t WaferLDTop
Definition: HGCalTypes.h:45
static constexpr double c75
Definition: HGCalTypes.h:102
static constexpr double c22
Definition: HGCalTypes.h:95
static constexpr std::array< int, 3 > edgeWaferHDFive
Definition: HGCalTypes.h:89
static constexpr int32_t WaferHDTop
Definition: HGCalTypes.h:51
static constexpr int32_t WaferFive
Definition: HGCalTypes.h:36
static constexpr int32_t WaferSemi
Definition: HGCalTypes.h:40
static constexpr std::array< int, 3 > edgeWaferHDRight
Definition: HGCalTypes.h:88
static constexpr int k_allCorners
static constexpr std::array< int, 3 > edgeWaferLDBottom
Definition: HGCalTypes.h:80
static constexpr int32_t WaferChopTwo
Definition: HGCalTypes.h:37
static bool maskCell(int u, int v, int N, int ncor, int fcor, int corners)