CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  bool good(false);
110  int n2 = n / 2;
111  int n4 = n / 4;
112  int n3 = (n + 1) / 3;
113  switch (type) {
114  case (HGCalTypes::WaferFull): { //WaferFull
115  good = true;
116  break;
117  }
118  case (HGCalTypes::WaferFive): { //WaferFive
119  switch (rotn) {
120  case (HGCalTypes::WaferCorner0): {
121  int u2 = u / 2;
122  good = ((v - u2) < n);
123  break;
124  }
125  case (HGCalTypes::WaferCorner1): {
126  good = ((v + u) < (3 * n - 1));
127  break;
128  }
129  case (HGCalTypes::WaferCorner2): {
130  int v2 = (v + 1) / 2;
131  good = ((u - v2) < n);
132  break;
133  }
134  case (HGCalTypes::WaferCorner3): {
135  int u2 = (u + 1) / 2;
136  good = (u2 <= v);
137  break;
138  }
139  case (HGCalTypes::WaferCorner4): {
140  good = ((v + u) >= n);
141  break;
142  }
143  default: {
144  int v2 = v / 2;
145  good = (u > v2);
146  break;
147  }
148  }
149  break;
150  }
151  case (HGCalTypes::WaferChopTwo): { //WaferChopTwo
152  switch (rotn) {
153  case (HGCalTypes::WaferCorner0): {
154  good = (v < (3 * n2));
155  break;
156  }
157  case (HGCalTypes::WaferCorner1): {
158  good = (u < (3 * n2));
159  break;
160  }
161  case (HGCalTypes::WaferCorner2): {
162  good = ((u - v) <= n2);
163  break;
164  }
165  case (HGCalTypes::WaferCorner3): {
166  good = (v >= n2);
167  break;
168  }
169  case (HGCalTypes::WaferCorner4): {
170  good = (u >= n2);
171  break;
172  }
173  default: {
174  good = ((v - u) < n2);
175  break;
176  }
177  }
178  break;
179  }
180  case (HGCalTypes::WaferChopTwoM): { //WaferChopTwoM
181  switch (rotn) {
182  case (HGCalTypes::WaferCorner0): {
183  good = (v < (5 * n4));
184  break;
185  }
186  case (HGCalTypes::WaferCorner1): {
187  good = (u < (5 * n4));
188  break;
189  }
190  case (HGCalTypes::WaferCorner2): {
191  good = ((u - v) <= n4);
192  break;
193  }
194  case (HGCalTypes::WaferCorner3): {
195  good = (v >= (3 * n4));
196  break;
197  }
198  case (HGCalTypes::WaferCorner4): {
199  good = (u >= (3 * n4));
200  break;
201  }
202  default: {
203  good = ((v - u) < n4);
204  break;
205  }
206  }
207  break;
208  }
209  case (HGCalTypes::WaferHalf): { //WaferHalf
210  switch (rotn) {
211  case (HGCalTypes::WaferCorner0): {
212  good = (v < n);
213  break;
214  }
215  case (HGCalTypes::WaferCorner1): {
216  good = (u < n);
217  break;
218  }
219  case (HGCalTypes::WaferCorner2): {
220  good = (v >= u);
221  break;
222  }
223  case (HGCalTypes::WaferCorner3): {
224  good = (v >= n);
225  break;
226  }
227  case (HGCalTypes::WaferCorner4): {
228  good = (u >= n);
229  break;
230  }
231  default: {
232  good = (u > v);
233  break;
234  }
235  }
236  break;
237  }
238  case (HGCalTypes::WaferSemi): { //WaferSemi
239  switch (rotn) {
240  case (HGCalTypes::WaferCorner0): {
241  good = ((u + v) < (2 * n));
242  break;
243  }
244  case (HGCalTypes::WaferCorner1): {
245  good = ((2 * u - v) < n);
246  break;
247  }
248  case (HGCalTypes::WaferCorner2): {
249  good = ((2 * v - u) >= n);
250  break;
251  }
252  case (HGCalTypes::WaferCorner3): {
253  good = ((u + v) >= (2 * n));
254  break;
255  }
256  case (HGCalTypes::WaferCorner4): {
257  good = ((2 * u - v) > n);
258  break;
259  }
260  default: {
261  good = ((2 * v - u) < n);
262  break;
263  }
264  }
265  break;
266  }
267  case (HGCalTypes::WaferThree): { //WaferThree
268  switch (rotn) {
269  case (HGCalTypes::WaferCorner0): {
270  good = ((v + u) < n);
271  break;
272  }
273  case (HGCalTypes::WaferCorner1): {
274  int v2 = v / 2;
275  good = (u <= v2);
276  break;
277  }
278  case (HGCalTypes::WaferCorner2): {
279  int u2 = (u / 2);
280  good = ((v - u2) >= n);
281  break;
282  }
283  case (HGCalTypes::WaferCorner3): {
284  good = ((v + u) >= (3 * n - 1));
285  break;
286  }
287  case (HGCalTypes::WaferCorner4): {
288  int v2 = ((v + 1) / 2);
289  good = ((u - v2) >= n);
290  break;
291  }
292  default: {
293  int u2 = ((u + 1) / 2);
294  good = (v < u2);
295  break;
296  }
297  }
298  break;
299  }
300  case (HGCalTypes::WaferSemi2): { //WaferSemi2
301  switch (rotn) {
302  case (HGCalTypes::WaferCorner0): {
303  good = ((u + v) < (4 * n3));
304  break;
305  }
306  case (HGCalTypes::WaferCorner1): {
307  good = ((2 * u - v) < n2);
308  break;
309  }
310  case (HGCalTypes::WaferCorner2): {
311  int u2 = ((u + 1) / 2);
312  good = ((v - u2) >= (3 * n4));
313  break;
314  }
315  case (HGCalTypes::WaferCorner3): {
316  good = ((u + v) >= (5 * n2));
317  break;
318  }
319  case (HGCalTypes::WaferCorner4): {
320  good = ((2 * u - v) > (3 * n2));
321  break;
322  }
323  default: {
324  int u2 = ((n == 8) ? ((u + 1) / 2) : (u / 2));
325  good = ((v - u2) < n4);
326  break;
327  }
328  }
329  break;
330  }
331  case (HGCalTypes::WaferFive2): { //WaferFive2
332  switch (rotn) {
333  case (HGCalTypes::WaferCorner0): {
334  good = ((2 * v - u) <= (3 * n2));
335  break;
336  }
337  case (HGCalTypes::WaferCorner1): {
338  good = ((u + v) < (5 * n2));
339  break;
340  }
341  case (HGCalTypes::WaferCorner2): {
342  good = ((2 * u - v) >= (3 * n2));
343  break;
344  }
345  case (HGCalTypes::WaferCorner3): {
346  good = ((2 * v - u) >= n3);
347  break;
348  }
349  case (HGCalTypes::WaferCorner4): {
350  good = ((u + v) > (4 * n3));
351  break;
352  }
353  default: {
354  good = ((2 * u - v) >= n2);
355  break;
356  }
357  }
358  break;
359  }
360  case (HGCalTypes::WaferHalf2): { //WaferHalf2
361  switch (rotn) {
362  case (HGCalTypes::WaferCorner0): {
363  good = (v <= (3 * n4));
364  break;
365  }
366  case (HGCalTypes::WaferCorner1): {
367  good = (u < (3 * n4));
368  break;
369  }
370  case (HGCalTypes::WaferCorner2): {
371  good = ((v - u) >= n4);
372  break;
373  }
374  case (HGCalTypes::WaferCorner3): {
375  good = (v >= (5 * n4));
376  break;
377  }
378  case (HGCalTypes::WaferCorner4): {
379  good = (u >= (5 * n4));
380  break;
381  }
382  default: {
383  good = ((u - v) >= n4);
384  break;
385  }
386  }
387  break;
388  }
389  }
390 #ifdef EDM_ML_DEBUG
391  edm::LogVerbatim("HGCalGeom") << "u|v " << u << ":" << v << " N " << n << " type " << type << " rot " << rotn
392  << " good " << good;
393 #endif
394  return good;
395 }
396 
397 int HGCalWaferMask::getRotation(int zside, int type, int rotn) {
398  if (rotn >= HGCalTypes::WaferCornerMax)
400  int newrotn(rotn);
401  if ((zside < 0) && (type != HGCalTypes::WaferFull)) {
402  if ((type == HGCalTypes::WaferFive) || (type == HGCalTypes::WaferFive2)) { //WaferFive/WaferFive2
403  static constexpr int rot1[HGCalTypes::WaferCornerMax] = {HGCalTypes::WaferCorner4,
409  newrotn = rot1[rotn];
410  } else if ((type == HGCalTypes::WaferThree) || (type == HGCalTypes::WaferSemi) ||
411  (type == HGCalTypes::WaferSemi2)) { //WaferThree/WaferSemi/WaferSemi2
412  static constexpr int rot2[HGCalTypes::WaferCornerMax] = {HGCalTypes::WaferCorner2,
418  newrotn = rot2[rotn];
419  } else { //WaferHalf/WaferChopTwo/WaferChopTwoM/WaferHalf2
420  static constexpr int rot3[HGCalTypes::WaferCornerMax] = {HGCalTypes::WaferCorner3,
426  newrotn = rot3[rotn];
427  }
428  }
429 #ifdef EDM_ML_DEBUG
430  edm::LogVerbatim("HGCalGeom") << "zside " << zside << " type " << type << " rotn " << rotn << ":" << newrotn;
431 #endif
432  return newrotn;
433 }
434 
435 std::pair<int, int> HGCalWaferMask::getTypeMode(const double& xpos,
436  const double& ypos,
437  const double& delX,
438  const double& delY,
439  const double& rin,
440  const double& rout,
441  const int& wType,
442  const int& mode,
443  bool debug) {
444  int ncor(0), iok(0);
446 
447  static constexpr int corners = 6;
448  static constexpr int base = 10;
449  double rin2 = rin * rin;
450  double rout2 = rout * rout;
451  double dx0[corners] = {HGCalTypes::c00 * delX,
452  HGCalTypes::c10 * delX,
453  HGCalTypes::c10 * delX,
454  HGCalTypes::c00 * delX,
455  -HGCalTypes::c10 * delX,
456  -HGCalTypes::c10 * delX};
457  double dy0[corners] = {-HGCalTypes::c10 * delY,
458  -HGCalTypes::c50 * delY,
459  HGCalTypes::c50 * delY,
460  HGCalTypes::c10 * delY,
461  HGCalTypes::c50 * delY,
462  -HGCalTypes::c50 * delY};
463  double xc[corners], yc[corners];
464  for (int k = 0; k < corners; ++k) {
465  xc[k] = xpos + dx0[k];
466  yc[k] = ypos + dy0[k];
467  double rpos2 = (xc[k] * xc[k] + yc[k] * yc[k]);
468  if (rpos2 <= rout2 && rpos2 >= rin2) {
469  ++ncor;
470  iok = iok * base + 1;
471  } else {
472  iok *= base;
473  }
474  }
475  if (debug)
476  edm::LogVerbatim("HGCalGeom") << "I/p " << xpos << ":" << ypos << ":" << delX << ":" << delY << ":" << rin << ":"
477  << rout << ":" << wType << ":" << mode << " Corners " << ncor << " iok " << iok;
478 
479  static constexpr int ipat5[corners] = {101111, 110111, 111011, 111101, 111110, 11111};
480  static constexpr int ipat4[corners] = {100111, 110011, 111001, 111100, 11110, 1111};
481  static constexpr int ipat3[corners] = {100011, 110001, 111000, 11100, 1110, 111};
482  static constexpr int ipat2[corners] = {11, 100001, 110000, 11000, 1100, 110};
483  double dx1[corners] = {HGCalTypes::c50 * delX,
484  HGCalTypes::c10 * delX,
485  HGCalTypes::c50 * delX,
486  -HGCalTypes::c50 * delX,
487  -HGCalTypes::c10 * delX,
488  -HGCalTypes::c50 * delX};
489  double dy1[corners] = {-HGCalTypes::c75 * delY,
490  HGCalTypes::c00 * delY,
491  HGCalTypes::c75 * delY,
492  HGCalTypes::c75 * delY,
493  HGCalTypes::c00 * delY,
494  -HGCalTypes::c75 * delY};
495  double dx2[corners] = {HGCalTypes::c50 * delX,
496  -HGCalTypes::c50 * delX,
497  -HGCalTypes::c10 * delX,
498  -HGCalTypes::c50 * delX,
499  HGCalTypes::c50 * delX,
500  HGCalTypes::c10 * delX};
501  double dy2[corners] = {HGCalTypes::c75 * delY,
502  HGCalTypes::c75 * delY,
503  HGCalTypes::c00 * delY,
504  -HGCalTypes::c75 * delY,
505  -HGCalTypes::c75 * delY,
506  HGCalTypes::c00 * delY};
507  double dx3[corners] = {HGCalTypes::c22 * delX,
508  HGCalTypes::c10 * delX,
509  HGCalTypes::c77 * delX,
510  -HGCalTypes::c22 * delX,
511  -HGCalTypes::c10 * delX,
512  -HGCalTypes::c77 * delX};
513  double dy3[corners] = {-HGCalTypes::c88 * delY,
514  -HGCalTypes::c27 * delY,
515  HGCalTypes::c61 * delY,
516  HGCalTypes::c88 * delY,
517  HGCalTypes::c27 * delY,
518  -HGCalTypes::c61 * delY};
519  double dx4[corners] = {HGCalTypes::c22 * delX,
520  -HGCalTypes::c77 * delX,
521  -HGCalTypes::c10 * delX,
522  -HGCalTypes::c22 * delX,
523  HGCalTypes::c77 * delX,
524  HGCalTypes::c10 * delX};
525  double dy4[corners] = {HGCalTypes::c88 * delY,
526  HGCalTypes::c61 * delY,
527  -HGCalTypes::c27 * delY,
528  -HGCalTypes::c88 * delY,
529  -HGCalTypes::c61 * delY,
530  HGCalTypes::c27 * delY};
531  double dx5[corners] = {-HGCalTypes::c50 * delX,
532  -HGCalTypes::c10 * delX,
533  -HGCalTypes::c50 * delX,
534  HGCalTypes::c50 * delX,
535  HGCalTypes::c10 * delX,
536  HGCalTypes::c50 * delX};
537  double dy5[corners] = {HGCalTypes::c75 * delY,
538  HGCalTypes::c00 * delY,
539  -HGCalTypes::c75 * delY,
540  -HGCalTypes::c75 * delY,
541  HGCalTypes::c00 * delY,
542  HGCalTypes::c75 * delY};
543  double dx6[corners] = {-HGCalTypes::c77 * delX,
544  -HGCalTypes::c10 * delX,
545  -HGCalTypes::c22 * delX,
546  HGCalTypes::c77 * delX,
547  HGCalTypes::c10 * delX,
548  HGCalTypes::c22 * delX};
549  double dy6[corners] = {HGCalTypes::c61 * delY,
550  -HGCalTypes::c27 * delY,
551  -HGCalTypes::c88 * delY,
552  -HGCalTypes::c61 * delY,
553  HGCalTypes::c27 * delY,
554  HGCalTypes::c88 * delY};
555  double dx7[corners] = {-HGCalTypes::c22 * delX,
556  -HGCalTypes::c10 * delX,
557  -HGCalTypes::c77 * delX,
558  HGCalTypes::c22 * delX,
559  HGCalTypes::c10 * delX,
560  HGCalTypes::c77 * delX};
561  double dy7[corners] = {HGCalTypes::c88 * delY,
562  HGCalTypes::c27 * delY,
563  -HGCalTypes::c61 * delY,
564  -HGCalTypes::c88 * delY,
565  -HGCalTypes::c27 * delY,
566  HGCalTypes::c61 * delY};
567  double dx8[corners] = {HGCalTypes::c77 * delX,
568  HGCalTypes::c10 * delX,
569  HGCalTypes::c22 * delX,
570  -HGCalTypes::c77 * delX,
571  -HGCalTypes::c10 * delX,
572  -HGCalTypes::c22 * delX};
573  double dy8[corners] = {-HGCalTypes::c61 * delY,
574  HGCalTypes::c27 * delY,
575  HGCalTypes::c88 * delY,
576  HGCalTypes::c61 * delY,
577  -HGCalTypes::c27 * delY,
578  -HGCalTypes::c88 * delY};
579  double dx9[corners] = {-HGCalTypes::c22 * delX,
580  HGCalTypes::c77 * delX,
581  HGCalTypes::c10 * delX,
582  HGCalTypes::c22 * delX,
583  -HGCalTypes::c77 * delX,
584  -HGCalTypes::c10 * delX};
585  double dy9[corners] = {-HGCalTypes::c88 * delY,
586  -HGCalTypes::c61 * delY,
587  HGCalTypes::c27 * delY,
588  HGCalTypes::c88 * delY,
589  HGCalTypes::c61 * delY,
590  -HGCalTypes::c27 * delY};
591 
592  if (ncor == HGCalGeomTools::k_allCorners) {
593  } else if (ncor == HGCalGeomTools::k_fiveCorners) {
594  rotn = static_cast<int>(std::find(ipat5, ipat5 + 6, iok) - ipat5);
596  } else if (ncor == HGCalGeomTools::k_fourCorners) {
597  rotn = static_cast<int>(std::find(ipat4, ipat4 + 6, iok) - ipat4);
599  double rpos12 = ((xpos + dx1[rotn]) * (xpos + dx1[rotn]) + (ypos + dy1[rotn]) * (ypos + dy1[rotn]));
600  double rpos22(0);
601  if (rpos12 <= rout2 && rpos12 >= rin2) {
602  rpos22 = ((xpos + dx2[rotn]) * (xpos + dx2[rotn]) + (ypos + dy2[rotn]) * (ypos + dy2[rotn]));
603  if (rpos22 <= rout2 && rpos22 >= rin2)
605  }
606  if (debug)
607  edm::LogVerbatim("HGCalGeom") << "Test for Chop2 " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
608  << type;
609  if ((type == HGCalTypes::WaferHalf) && (wType == 0)) {
610  rpos12 = ((xpos + dx3[rotn]) * (xpos + dx3[rotn]) + (ypos + dy3[rotn]) * (ypos + dy3[rotn]));
611  if (rpos12 <= rout2 && rpos12 >= rin2) {
612  rpos22 = ((xpos + dx4[rotn]) * (xpos + dx4[rotn]) + (ypos + dy4[rotn]) * (ypos + dy4[rotn]));
613  if (rpos22 <= rout2 && rpos22 >= rin2)
615  }
616  if (debug)
617  edm::LogVerbatim("HGCalGeom") << "Test for Chop2M " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
618  << type;
619  }
620  } else if (ncor == HGCalGeomTools::k_threeCorners) {
621  rotn = static_cast<int>(std::find(ipat3, ipat3 + 6, iok) - ipat3);
623  double rpos12 = ((xpos + dx7[rotn]) * (xpos + dx7[rotn]) + (ypos + dy7[rotn]) * (ypos + dy7[rotn]));
624  double rpos22(0);
625  if (rpos12 <= rout2 && rpos12 >= rin2) {
626  rpos22 = ((xpos + dx8[rotn]) * (xpos + dx8[rotn]) + (ypos + dy8[rotn]) * (ypos + dy8[rotn]));
627  if (rpos22 <= rout2 && rpos22 >= rin2)
629  }
630  if (debug)
631  edm::LogVerbatim("HGCalGeom") << "Test for Five2 " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
632  << type;
633  if ((type == HGCalTypes::WaferThree) && (wType == 0)) {
634  rpos12 = ((xpos + dx1[rotn]) * (xpos + dx1[rotn]) + (ypos + dy1[rotn]) * (ypos + dy1[rotn]));
635  if (rpos12 <= rout2 && rpos12 >= rin2) {
636  rpos22 = ((xpos + dx5[rotn]) * (xpos + dx5[rotn]) + (ypos + dy5[rotn]) * (ypos + dy5[rotn]));
637  if (rpos22 <= rout2 && rpos22 >= rin2)
639  }
640  if (debug)
641  edm::LogVerbatim("HGCalGeom") << "Test for Semi " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
642  << type;
643  }
644  if ((type == HGCalTypes::WaferThree) && (wType == 0)) {
645  rpos12 = ((xpos + dx3[rotn]) * (xpos + dx3[rotn]) + (ypos + dy3[rotn]) * (ypos + dy3[rotn]));
646  if (rpos12 <= rout2 && rpos12 >= rin2) {
647  rpos22 = ((xpos + dx6[rotn]) * (xpos + dx6[rotn]) + (ypos + dy6[rotn]) * (ypos + dy6[rotn]));
648  if (rpos22 <= rout2 && rpos22 >= rin2)
650  }
651  if (debug)
652  edm::LogVerbatim("HGCalGeom") << "Test for SemiM " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
653  << type;
654  }
655  } else if (ncor == HGCalGeomTools::k_twoCorners) {
656  rotn = static_cast<int>(std::find(ipat2, ipat2 + 6, iok) - ipat2);
658  double rpos12 = ((xpos + dx7[rotn]) * (xpos + dx7[rotn]) + (ypos + dy7[rotn]) * (ypos + dy7[rotn]));
659  double rpos22(0);
660  if (rpos12 <= rout2 && rpos12 >= rin2) {
661  rpos22 = ((xpos + dx9[rotn]) * (xpos + dx9[rotn]) + (ypos + dy9[rotn]) * (ypos + dy9[rotn]));
662  if (rpos22 <= rout2 && rpos22 >= rin2)
664  else
666  }
667  if (debug)
668  edm::LogVerbatim("HGCalGeom") << "Test for Half2 " << std::sqrt(rpos12) << ":" << std::sqrt(rpos22) << " Type "
669  << type;
670  } else {
672  }
673 
674  if (debug)
675  edm::LogVerbatim("HGCalGeom") << "I/p " << xpos << ":" << ypos << ":" << delX << ":" << delY << ":" << rin << ":"
676  << rout << ":" << wType << ":" << mode << " o/p " << iok << ":" << ncor << ":" << type
677  << ":" << rotn;
678  return ((mode == 0) ? std::make_pair(ncor, rotn) : std::make_pair(type, (rotn + HGCalWaferMask::k_OffsetRotation)));
679 }
680 
682  double xpos, double ypos, double delX, double delY, double rin, double rout, int part, int rotn, bool debug) {
684  return false;
685  if (rotn < 0 || rotn > HGCalTypes::WaferCornerMax)
686  return false;
687  double rin2 = rin * rin;
688  double rout2 = rout * rout;
689  double rpos2(0);
690  static constexpr int corners = HGCalTypes::WaferCornerMax;
691  static constexpr int corner2 = 2 * HGCalTypes::WaferCornerMax;
692  static constexpr int base = 10;
693  static constexpr int base2 = 100;
694  double dx0[corners] = {HGCalTypes::c00 * delX,
695  HGCalTypes::c10 * delX,
696  HGCalTypes::c10 * delX,
697  HGCalTypes::c00 * delX,
698  -HGCalTypes::c10 * delX,
699  -HGCalTypes::c10 * delX};
700  double dy0[corners] = {-HGCalTypes::c10 * delY,
701  -HGCalTypes::c50 * delY,
702  HGCalTypes::c50 * delY,
703  HGCalTypes::c10 * delY,
704  HGCalTypes::c50 * delY,
705  -HGCalTypes::c50 * delY};
706  double dx1[corners] = {HGCalTypes::c50 * delX,
707  HGCalTypes::c10 * delX,
708  HGCalTypes::c50 * delX,
709  -HGCalTypes::c50 * delX,
710  -HGCalTypes::c10 * delX,
711  -HGCalTypes::c50 * delX};
712  double dy1[corners] = {-HGCalTypes::c75 * delY,
713  HGCalTypes::c00 * delY,
714  HGCalTypes::c75 * delY,
715  HGCalTypes::c75 * delY,
716  HGCalTypes::c00 * delY,
717  -HGCalTypes::c75 * delY};
718  double dx2[corner2] = {HGCalTypes::c22 * delX,
719  HGCalTypes::c77 * delX,
720  HGCalTypes::c10 * delX,
721  HGCalTypes::c10 * delX,
722  HGCalTypes::c77 * delX,
723  HGCalTypes::c22 * delX,
724  -HGCalTypes::c22 * delX,
725  -HGCalTypes::c77 * delX,
726  -HGCalTypes::c10 * delX,
727  -HGCalTypes::c10 * delX,
728  -HGCalTypes::c77 * delX,
729  -HGCalTypes::c22 * delX};
730  double dy2[corner2] = {-HGCalTypes::c88 * delY,
731  -HGCalTypes::c61 * delY,
732  -HGCalTypes::c27 * delY,
733  HGCalTypes::c27 * delY,
734  HGCalTypes::c61 * delY,
735  HGCalTypes::c88 * delY,
736  HGCalTypes::c88 * delY,
737  HGCalTypes::c61 * delY,
738  HGCalTypes::c27 * delY,
739  -HGCalTypes::c27 * delY,
740  -HGCalTypes::c61 * delY,
741  -HGCalTypes::c88 * delY};
742  bool ok(true);
743  int ncf(-1);
744  switch (part) {
745  case (HGCalTypes::WaferThree): {
746  static constexpr int nc0[corners] = {450, 150, 201, 312, 423, 534};
747  int nc = nc0[rotn];
748  for (int k1 = 0; k1 < 3; ++k1) {
749  int k = nc % base;
750  double xc1 = xpos + dx0[k];
751  double yc1 = ypos + dy0[k];
752  rpos2 = (xc1 * xc1 + yc1 * yc1);
753  if ((rpos2 > rout2) || (rpos2 < rin2)) {
754  ok = false;
755  ncf = k;
756  break;
757  }
758  nc /= base;
759  }
760  break;
761  }
762  case (HGCalTypes::WaferSemi2): {
763  static constexpr int nc10[corners] = {450, 150, 201, 312, 423, 534};
764  static constexpr int nc11[corners] = {700, 902, 1104, 106, 308, 510};
765  int nc = nc10[rotn];
766  for (int k1 = 0; k1 < 3; ++k1) {
767  int k = nc % base;
768  double xc1 = xpos + dx0[k];
769  double yc1 = ypos + dy0[k];
770  rpos2 = (xc1 * xc1 + yc1 * yc1);
771  if ((rpos2 > rout2) || (rpos2 < rin2)) {
772  ok = false;
773  ncf = k;
774  break;
775  }
776  nc /= base;
777  }
778  nc = nc11[rotn];
779  for (int k1 = 0; k1 < 2; ++k1) {
780  int k = nc % base2;
781  double xc1 = xpos + dx2[k];
782  double yc1 = ypos + dy2[k];
783  rpos2 = (xc1 * xc1 + yc1 * yc1);
784  if ((rpos2 > rout2) || (rpos2 < rin2)) {
785  ok = false;
786  ncf = k + base2;
787  break;
788  }
789  nc /= base2;
790  }
791  break;
792  }
793  case (HGCalTypes::WaferSemi): {
794  static constexpr int nc20[corners] = {450, 150, 201, 312, 423, 534};
795  static constexpr int nc21[corners] = {30, 14, 25, 30, 41, 52};
796  int nc = nc20[rotn];
797  for (int k1 = 0; k1 < 3; ++k1) {
798  int k = nc % base;
799  double xc1 = xpos + dx0[k];
800  double yc1 = ypos + dy0[k];
801  rpos2 = (xc1 * xc1 + yc1 * yc1);
802  if ((rpos2 > rout2) || (rpos2 < rin2)) {
803  ok = false;
804  ncf = k;
805  break;
806  }
807  nc /= base;
808  }
809  nc = nc21[rotn];
810  for (int k1 = 0; k1 < 2; ++k1) {
811  int k = nc % base;
812  double xc1 = xpos + dx1[k];
813  double yc1 = ypos + dy1[k];
814  rpos2 = (xc1 * xc1 + yc1 * yc1);
815  if ((rpos2 > rout2) || (rpos2 < rin2)) {
816  ok = false;
817  ncf = k + base2;
818  break;
819  }
820  nc /= base;
821  }
822  break;
823  }
824  case (HGCalTypes::WaferHalf): {
825  static constexpr int nc3[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
826  int nc = nc3[rotn];
827  for (int k1 = 0; k1 < 4; ++k1) {
828  int k = nc % base;
829  double xc1 = xpos + dx0[k];
830  double yc1 = ypos + dy0[k];
831  rpos2 = (xc1 * xc1 + yc1 * yc1);
832  if ((rpos2 > rout2) || (rpos2 < rin2)) {
833  ok = false;
834  ncf = k;
835  break;
836  }
837  nc /= base;
838  }
839  break;
840  }
841  case (HGCalTypes::WaferChopTwoM): {
842  static constexpr int nc40[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
843  static constexpr int nc41[corners] = {500, 702, 904, 1106, 108, 310};
844  int nc = nc40[rotn];
845  for (int k1 = 0; k1 < 4; ++k1) {
846  int k = nc % base;
847  double xc1 = xpos + dx0[k];
848  double yc1 = ypos + dy0[k];
849  rpos2 = (xc1 * xc1 + yc1 * yc1);
850  if ((rpos2 > rout2) || (rpos2 < rin2)) {
851  ok = false;
852  ncf = k;
853  break;
854  }
855  nc /= base;
856  }
857  nc = nc41[rotn];
858  for (int k1 = 0; k1 < 2; ++k1) {
859  int k = nc % base2;
860  double xc1 = xpos + dx2[k];
861  double yc1 = ypos + dy2[k];
862  rpos2 = (xc1 * xc1 + yc1 * yc1);
863  if ((rpos2 > rout2) || (rpos2 < rin2)) {
864  ok = false;
865  ncf = k + base2;
866  break;
867  }
868  nc /= base2;
869  }
870  break;
871  }
872  case (HGCalTypes::WaferChopTwo): {
873  static constexpr int nc50[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
874  static constexpr int nc51[corners] = {20, 13, 24, 35, 40, 51};
875  int nc = nc50[rotn];
876  for (int k1 = 0; k1 < 4; ++k1) {
877  int k = nc % base;
878  double xc1 = xpos + dx0[k];
879  double yc1 = ypos + dy0[k];
880  rpos2 = (xc1 * xc1 + yc1 * yc1);
881  if ((rpos2 > rout2) || (rpos2 < rin2)) {
882  ok = false;
883  ncf = k;
884  break;
885  }
886  nc /= base;
887  }
888  nc = nc51[rotn];
889  for (int k1 = 0; k1 < 2; ++k1) {
890  int k = nc % base;
891  double xc1 = xpos + dx1[k];
892  double yc1 = ypos + dy1[k];
893  rpos2 = (xc1 * xc1 + yc1 * yc1);
894  if ((rpos2 > rout2) || (rpos2 < rin2)) {
895  ok = false;
896  ncf = k + base2;
897  break;
898  }
899  nc /= base;
900  }
901  break;
902  }
903  case (HGCalTypes::WaferFive): {
904  static constexpr int nc6[corners] = {23450, 13450, 24501, 35012, 40123, 51234};
905  int nc = nc6[rotn];
906  for (int k1 = 0; k1 < 5; ++k1) {
907  int k = nc % base;
908  double xc1 = xpos + dx0[k];
909  double yc1 = ypos + dy0[k];
910  rpos2 = (xc1 * xc1 + yc1 * yc1);
911  if ((rpos2 > rout2) || (rpos2 < rin2)) {
912  ok = false;
913  ncf = k;
914  break;
915  }
916  }
917  break;
918  }
919  case (HGCalTypes::WaferFive2): {
920  static constexpr int nc60[corners] = {450, 150, 201, 312, 423, 534};
921  static constexpr int nc61[corners] = {601, 803, 1005, 7, 209, 411};
922  int nc = nc60[rotn];
923  for (int k1 = 0; k1 < 3; ++k1) {
924  int k = nc % base;
925  double xc1 = xpos + dx0[k];
926  double yc1 = ypos + dy0[k];
927  rpos2 = (xc1 * xc1 + yc1 * yc1);
928  if ((rpos2 > rout2) || (rpos2 < rin2)) {
929  ok = false;
930  ncf = k;
931  break;
932  }
933  nc /= base;
934  }
935  nc = nc61[rotn];
936  for (int k1 = 0; k1 < 2; ++k1) {
937  int k = nc % base2;
938  double xc1 = xpos + dx2[k];
939  double yc1 = ypos + dy2[k];
940  rpos2 = (xc1 * xc1 + yc1 * yc1);
941  if ((rpos2 > rout2) || (rpos2 < rin2)) {
942  ok = false;
943  ncf = k + base2;
944  break;
945  }
946  nc /= base2;
947  }
948  break;
949  }
950  case (HGCalTypes::WaferHalf2): {
951  static constexpr int nc70[corners] = {45, 50, 1, 12, 23, 34};
952  static constexpr int nc71[corners] = {611, 801, 1003, 5, 207, 409};
953  int nc = nc70[rotn];
954  for (int k1 = 0; k1 < 2; ++k1) {
955  int k = nc % base;
956  double xc1 = xpos + dx0[k];
957  double yc1 = ypos + dy0[k];
958  rpos2 = (xc1 * xc1 + yc1 * yc1);
959  if ((rpos2 > rout2) || (rpos2 < rin2)) {
960  ok = false;
961  ncf = k;
962  break;
963  }
964  nc /= base;
965  }
966  nc = nc71[rotn];
967  for (int k1 = 0; k1 < 2; ++k1) {
968  int k = nc % base2;
969  double xc1 = xpos + dx2[k];
970  double yc1 = ypos + dy2[k];
971  rpos2 = (xc1 * xc1 + yc1 * yc1);
972  if ((rpos2 > rout2) || (rpos2 < rin2)) {
973  ok = false;
974  ncf = k + base2;
975  break;
976  }
977  nc /= base2;
978  }
979  break;
980  }
981  default: {
982  for (int k = 0; k < corners; ++k) {
983  double xc1 = xpos + dx0[k];
984  double yc1 = ypos + dy0[k];
985  rpos2 = (xc1 * xc1 + yc1 * yc1);
986  if ((rpos2 > rout2) || (rpos2 < rin2)) {
987  ok = false;
988  ncf = k;
989  break;
990  }
991  }
992  break;
993  }
994  }
995  if (debug || (!ok))
996  edm::LogVerbatim("HGCalGeom") << "I/p "
997  << ":" << xpos << ":" << ypos << ":" << delX << ":" << delY << ":" << rin << ":"
998  << rout << ":" << part << ":" << rotn << " Results " << ok << ":" << ncf << " R "
999  << rin2 << ":" << rout2 << ":" << rpos2;
1000  return ok;
1001 }
1002 
1003 std::vector<std::pair<double, double> > HGCalWaferMask::waferXY(
1004  int part, int ori, int zside, double delX, double delY, double xpos, double ypos) {
1005  std::vector<std::pair<double, double> > xy;
1006  int orient = getRotation(-zside, part, ori);
1007 #ifdef EDM_ML_DEBUG
1008  edm::LogVerbatim("HGCalGeom") << "Part " << part << " zSide " << zside << " Orient " << ori << ":" << orient;
1009 #endif
1010  double dx[24] = {HGCalTypes::c00 * delX, HGCalTypes::c10 * delX, HGCalTypes::c10 * delX, HGCalTypes::c00 * delX,
1011  -HGCalTypes::c10 * delX, -HGCalTypes::c10 * delX, HGCalTypes::c50 * delX, HGCalTypes::c10 * delX,
1012  HGCalTypes::c50 * delX, -HGCalTypes::c50 * delX, -HGCalTypes::c10 * delX, -HGCalTypes::c50 * delX,
1013  HGCalTypes::c22 * delX, HGCalTypes::c10 * delX, HGCalTypes::c77 * delX, -HGCalTypes::c22 * delX,
1014  -HGCalTypes::c10 * delX, -HGCalTypes::c77 * delX, HGCalTypes::c22 * delX, -HGCalTypes::c77 * delX,
1015  -HGCalTypes::c10 * delX, -HGCalTypes::c22 * delX, HGCalTypes::c77 * delX, HGCalTypes::c10 * delX};
1016  double dy[24] = {-HGCalTypes::c10 * delY, -HGCalTypes::c50 * delY, HGCalTypes::c50 * delY, HGCalTypes::c10 * delY,
1017  HGCalTypes::c50 * delY, -HGCalTypes::c50 * delY, -HGCalTypes::c75 * delY, HGCalTypes::c00 * delY,
1018  HGCalTypes::c75 * delY, HGCalTypes::c75 * delY, HGCalTypes::c00 * delY, -HGCalTypes::c75 * delY,
1019  -HGCalTypes::c88 * delY, -HGCalTypes::c27 * delY, HGCalTypes::c61 * delY, HGCalTypes::c88 * delY,
1020  HGCalTypes::c27 * delY, -HGCalTypes::c61 * delY, HGCalTypes::c88 * delY, HGCalTypes::c61 * delY,
1021  -HGCalTypes::c27 * delY, -HGCalTypes::c88 * delY, -HGCalTypes::c61 * delY, HGCalTypes::c27 * delY};
1022  if (part == HGCalTypes::WaferFull) {
1023  int np[7] = {0, 1, 2, 3, 4, 5, 0};
1024  for (int k = 0; k < 7; ++k)
1025  xy.push_back(std::make_pair((xpos + dx[np[k]]), (ypos + dy[np[k]])));
1026  } else if (part == HGCalTypes::WaferFive) {
1027  int np[6][6] = {{0, 2, 3, 4, 5, 0},
1028  {1, 3, 4, 5, 0, 1},
1029  {2, 4, 5, 0, 1, 2},
1030  {3, 5, 0, 1, 2, 3},
1031  {4, 0, 1, 2, 3, 4},
1032  {5, 1, 2, 3, 4, 5}};
1033  for (int k = 0; k < 6; ++k) {
1034  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1035 #ifdef EDM_ML_DEBUG
1036  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1037  << dy[np[orient][k]];
1038 #endif
1039  }
1040  } else if (part == HGCalTypes::WaferHalf) {
1041  int np[6][5] = {
1042  {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}};
1043  for (int k = 0; k < 5; ++k) {
1044  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1045 #ifdef EDM_ML_DEBUG
1046  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1047  << dy[np[orient][k]];
1048 #endif
1049  }
1050  } else if (part == HGCalTypes::WaferThree) {
1051  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}};
1052  for (int k = 0; k < 4; ++k) {
1053  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1054 #ifdef EDM_ML_DEBUG
1055  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1056  << dy[np[orient][k]];
1057 #endif
1058  }
1059  } else if (part == HGCalTypes::WaferChopTwo) {
1060  int np[6][7] = {{6, 8, 3, 4, 5, 0, 6},
1061  {7, 9, 4, 5, 0, 1, 7},
1062  {8, 10, 5, 0, 1, 2, 8},
1063  {9, 11, 0, 1, 2, 3, 9},
1064  {10, 6, 1, 2, 3, 4, 10},
1065  {11, 7, 2, 3, 4, 5, 11}};
1066  for (int k = 0; k < 7; ++k) {
1067  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1068 #ifdef EDM_ML_DEBUG
1069  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1070  << dy[np[orient][k]];
1071 #endif
1072  }
1073  } else if (part == HGCalTypes::WaferSemi) {
1074  int np[6][6] = {{6, 9, 4, 5, 0, 6},
1075  {7, 10, 5, 0, 1, 7},
1076  {8, 11, 0, 1, 2, 8},
1077  {9, 6, 1, 2, 3, 9},
1078  {10, 7, 2, 3, 4, 10},
1079  {11, 8, 3, 4, 5, 11}};
1080  for (int k = 0; k < 6; ++k) {
1081  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1082 #ifdef EDM_ML_DEBUG
1083  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1084  << dy[np[orient][k]];
1085 #endif
1086  }
1087  } else if (part == HGCalTypes::WaferChopTwoM) {
1088  int np[6][7] = {{12, 18, 3, 4, 5, 0, 12},
1089  {13, 19, 4, 5, 0, 1, 13},
1090  {14, 20, 5, 0, 1, 2, 14},
1091  {15, 21, 0, 1, 2, 3, 15},
1092  {16, 22, 1, 2, 3, 4, 16},
1093  {17, 23, 2, 3, 4, 5, 17}};
1094  for (int k = 0; k < 7; ++k) {
1095  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1096 #ifdef EDM_ML_DEBUG
1097  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1098  << dy[np[orient][k]];
1099 #endif
1100  }
1101  } else if (part == HGCalTypes::WaferSemi2) {
1102  int np[6][6] = {{12, 19, 4, 5, 0, 12},
1103  {13, 20, 5, 0, 1, 13},
1104  {14, 21, 0, 1, 2, 14},
1105  {15, 22, 1, 2, 3, 15},
1106  {16, 23, 2, 3, 4, 16},
1107  {17, 18, 3, 4, 5, 17}};
1108  for (int k = 0; k < 6; ++k) {
1109  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1110 #ifdef EDM_ML_DEBUG
1111  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1112  << dy[np[orient][k]];
1113 #endif
1114  }
1115  } else if (part == HGCalTypes::WaferFive2) {
1116  int np[6][6] = {{22, 15, 4, 5, 0, 22},
1117  {23, 16, 5, 0, 1, 23},
1118  {18, 17, 0, 1, 2, 18},
1119  {19, 12, 1, 2, 3, 19},
1120  {20, 13, 2, 3, 4, 20},
1121  {21, 14, 3, 4, 5, 21}};
1122  for (int k = 0; k < 6; ++k) {
1123  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1124 #ifdef EDM_ML_DEBUG
1125  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1126  << dy[np[orient][k]];
1127 #endif
1128  }
1129  } else if (part == HGCalTypes::WaferHalf2) {
1130  int np[6][5] = {{21, 15, 4, 5, 21},
1131  {22, 16, 5, 0, 22},
1132  {23, 17, 0, 1, 23},
1133  {18, 12, 1, 2, 18},
1134  {19, 13, 2, 3, 19},
1135  {20, 14, 3, 4, 20}};
1136  for (int k = 0; k < 5; ++k) {
1137  xy.push_back(std::make_pair((xpos + dx[np[orient][k]]), (ypos + dy[np[orient][k]])));
1138 #ifdef EDM_ML_DEBUG
1139  edm::LogVerbatim("HGCalGeom") << k << ":" << np[orient][k] << ":" << dx[np[orient][k]] << ":"
1140  << dy[np[orient][k]];
1141 #endif
1142  }
1143  }
1144 #ifdef EDM_ML_DEBUG
1145  edm::LogVerbatim("HGCalGeom") << "I/p: " << part << ":" << ori << ":" << zside << ":" << delX << ":" << delY << ":"
1146  << xpos << ":" << ypos << " O/p having " << xy.size() << " points:";
1147  std::ostringstream st1;
1148  for (unsigned int i = 0; i < xy.size(); ++i)
1149  st1 << " [" << i << "] " << xy[i].first << ":" << xy[i].second;
1150  edm::LogVerbatim("HGCalGeom") << st1.str();
1151 #endif
1152  return xy;
1153 }
tuple base
Main Program
Definition: newFWLiteAna.py:92
Log< level::Info, true > LogVerbatim
static constexpr int32_t WaferSizeMax
Definition: HGCalTypes.h:92
static constexpr int k_fourCorners
static std::vector< std::pair< double, double > > waferXY(int part, int orient, int zside, double delX, double delY, double xpos, double ypos)
static bool goodCell(int u, int v, int N, int type, int rotn)
static constexpr int k_twoCorners
static constexpr double c50
Definition: HGCalTypes.h:98
static bool goodTypeMode(double xpos, double ypos, double delX, double delY, double rin, double rout, int part, int rotn, bool debug)
static constexpr double c88
Definition: HGCalTypes.h:102
int zside(DetId const &)
static int getRotation(int zside, int type, int rotn)
static constexpr double c27
Definition: HGCalTypes.h:97
static constexpr double c10
Definition: HGCalTypes.h:103
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
U second(std::pair< T, U > const &p)
int np
Definition: AMPTWrapper.h:43
T sqrt(T t)
Definition: SSEVec.h:19
static constexpr int32_t WaferCornerMax
Definition: HGCalTypes.h:91
static constexpr double c00
Definition: HGCalTypes.h:94
Basic2DVector< T > xy() const
#define debug
Definition: HDRShower.cc:19
static constexpr int k_fiveCorners
auto const good
min quality of good
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, bool deug=false)
part
Definition: HCALResponse.h:20
static constexpr double c77
Definition: HGCalTypes.h:101
static constexpr double c61
Definition: HGCalTypes.h:99
static constexpr int k_OffsetRotation
static constexpr int k_threeCorners
static constexpr double c75
Definition: HGCalTypes.h:100
static constexpr double c22
Definition: HGCalTypes.h:95
static constexpr int k_allCorners
static bool maskCell(int u, int v, int N, int ncor, int fcor, int corners)