1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
|
#include "Geometry/HGCalCommonData/interface/HGCalCellUV.h"
#include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
#include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include <iostream>
#include <array>
#include <algorithm>
#include <cassert>
HGCalCellUV::HGCalCellUV(double waferSize, double separation, int32_t nFine, int32_t nCoarse) : waferSize_(waferSize) {
hgcalcell_ = std::make_unique<HGCalCell>(waferSize, nFine, nCoarse);
assert(nFine > 0 && nCoarse > 0);
ncell_[0] = nFine;
ncell_[1] = nCoarse;
for (int k = 0; k < 2; ++k) {
cellX_[k] = waferSize / (3 * ncell_[k]);
cellY_[k] = 0.5 * sqrt3_ * cellX_[k];
cellXTotal_[k] = (waferSize + separation) / (3 * ncell_[k]);
cellYTotal_[k] = 0.5 * sqrt3_ * cellXTotal_[k];
}
// Fill up the placement vectors
for (int placement = 0; placement < HGCalCell::cellPlacementTotal; ++placement) {
// Fine cells
for (int iu = 0; iu < 2 * ncell_[0]; ++iu) {
for (int iv = 0; iv < 2 * ncell_[0]; ++iv) {
int u = (placement < HGCalCell::cellPlacementExtra) ? iv : iu;
int v = (placement < HGCalCell::cellPlacementExtra) ? iu : iv;
if (((v - u) < ncell_[0]) && (u - v) <= ncell_[0]) {
cellPosFine_[placement][std::pair<int, int>(u, v)] = hgcalcell_->cellUV2XY1(u, v, placement, 0);
}
}
}
// Coarse cells
for (int iu = 0; iu < 2 * ncell_[1]; ++iu) {
for (int iv = 0; iv < 2 * ncell_[1]; ++iv) {
int u = (placement < HGCalCell::cellPlacementExtra) ? iv : iu;
int v = (placement < HGCalCell::cellPlacementExtra) ? iu : iv;
if (((v - u) < ncell_[1]) && (u - v) <= ncell_[1]) {
cellPosCoarse_[placement][std::pair<int, int>(u, v)] = hgcalcell_->cellUV2XY1(u, v, placement, 1);
}
}
}
}
}
std::pair<int32_t, int32_t> HGCalCellUV::cellUVFromXY1(
double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const {
if (type != 0)
type = 1;
//--- Reverse transform to placement=0, if placement index ≠6
double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
int rot = placement % HGCalCell::cellPlacementExtra;
static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
double x = xloc1 * fcos[rot] - yloc * fsin[rot];
double y = xloc1 * fsin[rot] + yloc * fcos[rot];
//--- Calculate coordinates in u,v,w system
double u = x * sin60_ + y * cos60_;
double v = -x * sin60_ + y * cos60_;
double w = y;
//--- Rounding in u, v, w coordinates
int iu = std::floor(u / cellY) + 3 * (ncell_[type]) - 1;
int iv = std::floor(v / cellY) + 3 * (ncell_[type]);
int iw = std::floor(w / cellY) + 1;
int isv = (iu + iw) / 3;
int isu = (iv + iw) / 3;
//--- Taking care of extending cells
if ((iu + iw) < 0) {
isu = (iv + iw + 1) / 3;
isv = 0;
} else if (isv - isu > ncell_[type] - 1) {
isu = (iv + iw + 1) / 3;
isv = (iu + iw - 1) / 3;
} else if (isu > 2 * ncell_[type] - 1) {
isu = 2 * ncell_[type] - 1;
isv = (iu + iw - 1) / 3;
}
if (debug)
edm::LogVerbatim("HGCalGeom") << "cellUVFromXY1: Input " << xloc << ":" << yloc << ":" << extend << " Output "
<< isu << ":" << isv;
return std::make_pair(isu, isv);
}
std::pair<int32_t, int32_t> HGCalCellUV::cellUVFromXY2(
double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const {
//--- Using multiple inequalities to find (u, v)
//--- Reverse transform to placement=0, if placement index ≠7
if (type != 0)
type = 1;
double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -1 * xloc;
int rot = placement % HGCalCell::cellPlacementExtra;
static constexpr std::array<double, 6> fcos = {{cos60_, 1.0, cos60_, -cos60_, -1.0, -cos60_}};
static constexpr std::array<double, 6> fsin = {{-sin60_, 0.0, sin60_, sin60_, 0.0, -sin60_}};
double x = xloc1 * fcos[rot] - yloc * fsin[rot];
double y = xloc1 * fsin[rot] + yloc * fcos[rot];
int32_t u(-100), v(-100);
int ncell = ncell_[type];
double r = cellY;
double l1 = (y / r) + ncell - 1.0;
int l2 = std::floor((0.5 * y + 0.5 * x / sqrt3_) / r + ncell - 4.0 / 3.0);
int l3 = std::floor((x / sqrt3_) / r + ncell - 4.0 / 3.0);
double l4 = (y + sqrt3_ * x) / (2 * r) + 2 * ncell - 2;
double l5 = (y - sqrt3_ * x) / (2 * r) - ncell;
double u1 = (y / r) + ncell + 1.0;
int u2 = std::ceil((0.5 * y + 0.5 * x / sqrt3_) / r + ncell + 2.0 / 3.0);
int u3 = std::ceil((x / sqrt3_) / r + ncell);
double u4 = l4 + 2;
double u5 = l5 + 2;
for (int ui = l2 + 1; ui < u2; ui++) {
for (int vi = l3 + 1; vi < u3; vi++) {
int c1 = 2 * ui - vi;
int c4 = ui + vi;
int c5 = ui - 2 * vi;
if ((c1 < u1) && (c1 > l1) && (c4 < u4) && (c4 > l4) && (c5 < u5) && (c5 > l5)) {
u = ui;
v = vi;
}
}
}
//--- Taking care of extending cells
if (v == -1) {
if (y < (2 * u - v - ncell) * r) {
++v;
} else {
++u;
++v;
}
}
if (v - u == ncell) {
if ((y + sqrt3_ * x) < ((u + v - 2 * ncell + 1) * 2 * r)) {
--v;
} else {
++u;
}
}
if (u == 2 * ncell) {
if ((y - sqrt3_ * x) < ((u - 2 * v + ncell - 1) * 2 * r)) {
--u;
} else {
--u;
--v;
}
}
if (debug)
edm::LogVerbatim("HGCalGeom") << "cellUVFromXY2: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
<< ":" << v;
return std::make_pair(u, v);
}
std::pair<int32_t, int32_t> HGCalCellUV::cellUVFromXY3(
double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const {
//--- Using Cube coordinates to find the (u, v)
//--- Reverse transform to placement=0, if placement index ≠6
if (type != 0)
type = 1;
double cellX = (extend) ? cellXTotal_[type] : cellX_[type];
double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
int rot = placement % HGCalCell::cellPlacementExtra;
static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
double xprime = xloc1 * fcos[rot] - yloc * fsin[rot];
double yprime = xloc1 * fsin[rot] + yloc * fcos[rot];
double x = xprime + cellX;
double y = yprime;
x = x / cellX;
y = y / cellY;
double cu = 2 * x / 3;
double cv = -x / 3 + y / 2;
double cw = -x / 3 - y / 2;
int iu = std::round(cu);
int iv = std::round(cv);
int iw = std::round(cw);
if (iu + iv + iw != 0) {
double arr[] = {std::abs(cu - iu), std::abs(cv - iv), std::abs(cw - iw)};
int i = std::distance(arr, std::max_element(arr, arr + 3));
if (i == 1)
iv = (std::round(cv) == std::floor(cv)) ? std::ceil(cv) : std::floor(cv);
else if (i == 2)
iw = (std::round(cw) == std::floor(cw)) ? std::ceil(cw) : std::floor(cw);
}
//--- Taking care of extending cells
int u(ncell_[type] + iv), v(ncell_[type] - 1 - iw);
double xcell = (1.5 * (v - u) + 0.5) * cellX_[type];
double ycell = (v + u - 2 * ncell_[type] + 1) * cellY;
if (v == -1) {
if ((yprime - sqrt3_ * xprime) < (ycell - sqrt3_ * xcell)) {
++v;
} else {
++u;
++v;
}
}
if (v - u == ncell_[type]) {
if (yprime < ycell) {
--v;
} else {
++u;
}
}
if (u == 2 * ncell_[type]) {
if ((yprime + sqrt3_ * xprime) > (ycell + sqrt3_ * xcell)) {
--u;
} else {
--u;
--v;
}
}
if (debug)
edm::LogVerbatim("HGCalGeom") << "cellUVFromXY3: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
<< ":" << v;
return std::make_pair(u, v);
}
std::pair<int, int> HGCalCellUV::cellUVFromXY4(
double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) {
if (type != 0)
return cellUVFromXY4(
xloc, yloc, ncell_[1], cellX_[1], cellY_[1], cellXTotal_[1], cellY_[1], cellPosCoarse_[placement], extend, debug);
else
return cellUVFromXY4(
xloc, yloc, ncell_[0], cellX_[0], cellY_[0], cellXTotal_[0], cellY_[0], cellPosFine_[placement], extend, debug);
}
std::pair<int, int> HGCalCellUV::cellUVFromXY4(double xloc,
double yloc,
int n,
double cellX,
double cellY,
double cellXTotal,
double cellYTotal,
std::map<std::pair<int, int>, std::pair<double, double> >& cellPos,
bool extend,
bool debug) {
std::pair<int, int> uv = std::make_pair(-1, -1);
std::map<std::pair<int, int>, std::pair<double, double> >::const_iterator itr;
for (itr = cellPos.begin(); itr != cellPos.end(); ++itr) {
double delX = std::abs(xloc - (itr->second).first);
double delY = std::abs(yloc - (itr->second).second);
if ((delX < cellX) && (delY < cellY)) {
if ((delX < (0.5 * cellX)) || (delY < (2.0 * cellY - sqrt3_ * delX))) {
uv = itr->first;
break;
}
}
}
if ((uv.first < 0) && extend) {
for (itr = cellPos.begin(); itr != cellPos.end(); ++itr) {
double delX = std::abs(xloc - (itr->second).first);
double delY = std::abs(yloc - (itr->second).second);
if ((delX < cellXTotal) && (delY < cellYTotal)) {
if ((delX < (0.5 * cellXTotal)) || (delY < (2.0 * cellYTotal - sqrt3_ * delX))) {
uv = itr->first;
break;
}
}
}
}
if (debug)
edm::LogVerbatim("HGCalGeom") << "cellUVFromXY4: Input " << xloc << ":" << yloc << ":" << extend << " Output "
<< uv.first << ":" << uv.second;
return uv;
}
std::pair<int32_t, int32_t> HGCalCellUV::cellUVFromXY1( // for v17
double xloc,
double yloc,
int32_t placement,
int32_t type,
int32_t partial,
bool extend,
bool debug) const {
if (type != 0)
type = 1;
double cellX = (extend) ? cellXTotal_[type] : cellX_[type];
double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
std::pair<int, int> uv = HGCalCellUV::cellUVFromXY1(xloc, yloc, placement, type, extend, debug);
int u = uv.first;
int v = uv.second;
if (partial == HGCalTypes::WaferLDTop) {
if (u * HGCalTypes::edgeWaferLDTop[0] + v * HGCalTypes::edgeWaferLDTop[1] == HGCalTypes::edgeWaferLDTop[2] + 1) {
double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
int rot = placement % HGCalCell::cellPlacementExtra;
static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
double xprime = -1 * (xloc1 * fcos[rot] - yloc * fsin[rot]);
double yprime = xloc1 * fsin[rot] + yloc * fcos[rot];
double xcell = -1 * (1.5 * (v - u) + 0.5) * cellX;
double ycell = (v + u - 2 * ncell_[type] + 1) * cellY;
if ((yprime - sqrt3_ * xprime) > (ycell - sqrt3_ * xcell)) {
--u;
if ((v - u) >= ncell_[1])
--v;
} else {
--u;
--v;
v = std::max(v, 0);
}
}
} else if (partial == HGCalTypes::WaferHDBottom) {
if (u * HGCalTypes::edgeWaferHDBottom[0] + v * HGCalTypes::edgeWaferHDBottom[1] ==
HGCalTypes::edgeWaferHDBottom[2] + 1) {
double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
int rot = placement % HGCalCell::cellPlacementExtra;
static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
double xprime = -1 * (xloc1 * fcos[rot] - yloc * fsin[rot]);
double yprime = xloc1 * fsin[rot] + yloc * fcos[rot];
double xcell = -1 * (1.5 * (v - u) + 0.5) * cellX;
double ycell = (v + u - 2 * ncell_[type] + 1) * cellY;
if ((yprime - sqrt3_ * xprime) > (ycell - sqrt3_ * xcell)) {
++u;
++v;
} else {
++u;
}
}
}
if (debug)
edm::LogVerbatim("HGCalGeom") << "cellUVFromXY5: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
<< ":" << v;
return std::make_pair(u, v);
}
std::pair<int32_t, int32_t> HGCalCellUV::cellUVFromXY2( // for v18
double xloc,
double yloc,
int32_t placement,
int32_t type,
int32_t partial,
bool extend,
bool debug) const {
if (type != 0)
type = 1;
std::pair<int, int> uv = HGCalCellUV::cellUVFromXY1(xloc, yloc, placement, type, extend, debug);
int u = uv.first;
int v = uv.second;
if ((partial != HGCalTypes::WaferFull) && (type == 1)) {
if (u == 1 && v == 8) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(HGCalTypes::WaferLDThree, placement, waferSize_, 0.0, false);
if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
++u;
++v;
}
}
if (u == 15 && v == 15) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(HGCalTypes::WaferLDThree, placement, waferSize_, 0.0, false);
if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
--u;
}
}
if (u * HGCalTypes::edgeWaferLDTop[0] + v * HGCalTypes::edgeWaferLDTop[1] == HGCalTypes::edgeWaferLDTop[2] + 1) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(HGCalTypes::WaferLDTop, placement, waferSize_, 0.0, false);
if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
std::pair<double, double> xy1 = hgcalcell_->cellUV2XY1(u, v, placement, 1);
std::array<double, 4> criterion2 =
HGCalWaferMask::maskCut(HGCalTypes::WaferLDThree, placement, waferSize_, 0.0, false);
if (((criterion2[0] * yloc) + (criterion2[1] * xloc) - (criterion2[0] * xy1.second) -
(criterion2[1] * xy1.first)) < 0.0) {
--u;
if ((v - u) >= ncell_[1])
--v;
} else {
--u;
--v;
v = std::max(v, 0);
}
}
}
} else if ((partial != HGCalTypes::WaferFull) && (type == 0)) {
if (u == 10 && v == 0) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(HGCalTypes::WaferHDTop, placement, waferSize_, 0.0, false);
if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
--u;
}
}
if (u == 10 && v == 21) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(HGCalTypes::WaferHDTop, placement, waferSize_, 0.0, false);
if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
--u;
--v;
}
}
if (u * HGCalTypes::edgeWaferHDBottom[0] + v * HGCalTypes::edgeWaferHDBottom[1] ==
HGCalTypes::edgeWaferHDBottom[2] + 1) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(HGCalTypes::WaferHDBottom, placement, waferSize_, 0.0, false);
if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
std::pair<double, double> xy1 = hgcalcell_->cellUV2XY1(u, v, placement, 0);
std::array<double, 4> criterion2 =
HGCalWaferMask::maskCut(HGCalTypes::WaferHDRight, placement, waferSize_, 0.0, false);
if (((criterion2[0] * yloc) + (criterion2[1] * xloc) - (criterion2[0] * xy1.second) -
(criterion2[1] * xy1.first)) < 0.0) {
++u;
++v;
} else {
++u;
}
}
}
}
if (debug)
edm::LogVerbatim("HGCalGeom") << "cellUVFromXY5: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
<< ":" << v;
return std::make_pair(u, v);
}
|