forked from cms-patatrack/pixeltrack-standalone
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpixelCPEforGPU.h
344 lines (284 loc) · 12.6 KB
/
pixelCPEforGPU.h
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
#ifndef RecoLocalTracker_SiPixelRecHits_pixelCPEforGPU_h
#define RecoLocalTracker_SiPixelRecHits_pixelCPEforGPU_h
#include <cassert>
#include <cmath>
#include <cstdint>
#include <iterator>
#include "CUDADataFormats/gpuClusteringConstants.h"
#include "DataFormats/SOARotation.h"
#include "Geometry/phase1PixelTopology.h"
#include "CUDACore/cudaCompat.h"
namespace pixelCPEforGPU {
using Frame = SOAFrame<float>;
using Rotation = SOARotation<float>;
// all modules are identical!
struct CommonParams {
float theThicknessB;
float theThicknessE;
float thePitchX;
float thePitchY;
};
struct DetParams {
bool isBarrel;
bool isPosZ;
uint16_t layer;
uint16_t index;
uint32_t rawId;
float shiftX;
float shiftY;
float chargeWidthX;
float chargeWidthY;
// CMSSW 11.2.x adds
//uint16_t pixmx; // max pix charge
// which would break reading the binary dumps
float x0, y0, z0; // the vertex in the local coord of the detector
float sx[3], sy[3]; // the errors...
Frame frame;
};
using phase1PixelTopology::AverageGeometry;
struct LayerGeometry {
uint32_t layerStart[phase1PixelTopology::numberOfLayers + 1];
uint8_t layer[phase1PixelTopology::layerIndexSize];
};
struct ParamsOnGPU {
CommonParams const* m_commonParams;
DetParams const* m_detParams;
LayerGeometry const* m_layerGeometry;
AverageGeometry const* m_averageGeometry;
constexpr CommonParams const& __restrict__ commonParams() const {
CommonParams const* __restrict__ l = m_commonParams;
return *l;
}
constexpr DetParams const& __restrict__ detParams(int i) const {
DetParams const* __restrict__ l = m_detParams;
return l[i];
}
constexpr LayerGeometry const& __restrict__ layerGeometry() const { return *m_layerGeometry; }
constexpr AverageGeometry const& __restrict__ averageGeometry() const { return *m_averageGeometry; }
__device__ uint8_t layer(uint16_t id) const {
return __ldg(m_layerGeometry->layer + id / phase1PixelTopology::maxModuleStride);
};
};
// SOA (on device)
template <uint32_t N>
struct ClusParamsT {
uint32_t minRow[N];
uint32_t maxRow[N];
uint32_t minCol[N];
uint32_t maxCol[N];
int32_t Q_f_X[N];
int32_t Q_l_X[N];
int32_t Q_f_Y[N];
int32_t Q_l_Y[N];
int32_t charge[N];
float xpos[N];
float ypos[N];
float xerr[N];
float yerr[N];
int16_t xsize[N]; // clipped at 127 if negative is edge....
int16_t ysize[N];
};
constexpr int32_t MaxHitsInIter = gpuClustering::maxHitsInIter();
using ClusParams = ClusParamsT<MaxHitsInIter>;
constexpr inline void computeAnglesFromDet(
DetParams const& __restrict__ detParams, float const x, float const y, float& cotalpha, float& cotbeta) {
// x,y local position on det
auto gvx = x - detParams.x0;
auto gvy = y - detParams.y0;
auto gvz = -1.f / detParams.z0;
// normalization not required as only ratio used...
// calculate angles
cotalpha = gvx * gvz;
cotbeta = gvy * gvz;
}
constexpr inline float correction(int sizeM1,
int Q_f, //!< Charge in the first pixel.
int Q_l, //!< Charge in the last pixel.
uint16_t upper_edge_first_pix, //!< As the name says.
uint16_t lower_edge_last_pix, //!< As the name says.
float lorentz_shift, //!< L-shift at half thickness
float theThickness, //detector thickness
float cot_angle, //!< cot of alpha_ or beta_
float pitch, //!< thePitchX or thePitchY
bool first_is_big, //!< true if the first is big
bool last_is_big) //!< true if the last is big
{
if (0 == sizeM1) // size 1
return 0;
float W_eff = 0;
bool simple = true;
if (1 == sizeM1) { // size 2
//--- Width of the clusters minus the edge (first and last) pixels.
//--- In the note, they are denoted x_F and x_L (and y_F and y_L)
// assert(lower_edge_last_pix >= upper_edge_first_pix);
auto W_inner = pitch * float(lower_edge_last_pix - upper_edge_first_pix); // in cm
//--- Predicted charge width from geometry
auto W_pred = theThickness * cot_angle // geometric correction (in cm)
- lorentz_shift; // (in cm) &&& check fpix!
W_eff = std::abs(W_pred) - W_inner;
//--- If the observed charge width is inconsistent with the expectations
//--- based on the track, do *not* use W_pred-W_inner. Instead, replace
//--- it with an *average* effective charge width, which is the average
//--- length of the edge pixels.
simple =
(W_eff < 0.0f) | (W_eff > pitch); // this produces "large" regressions for very small numeric differences...
}
if (simple) {
//--- Total length of the two edge pixels (first+last)
float sum_of_edge = 2.0f;
if (first_is_big)
sum_of_edge += 1.0f;
if (last_is_big)
sum_of_edge += 1.0f;
W_eff = pitch * 0.5f * sum_of_edge; // ave. length of edge pixels (first+last) (cm)
}
//--- Finally, compute the position in this projection
float Qdiff = Q_l - Q_f;
float Qsum = Q_l + Q_f;
//--- Temporary fix for clusters with both first and last pixel with charge = 0
if (Qsum == 0)
Qsum = 1.0f;
return 0.5f * (Qdiff / Qsum) * W_eff;
}
constexpr inline void position(CommonParams const& __restrict__ comParams,
DetParams const& __restrict__ detParams,
ClusParams& cp,
uint32_t ic) {
//--- Upper Right corner of Lower Left pixel -- in measurement frame
uint16_t llx = cp.minRow[ic] + 1;
uint16_t lly = cp.minCol[ic] + 1;
//--- Lower Left corner of Upper Right pixel -- in measurement frame
uint16_t urx = cp.maxRow[ic];
uint16_t ury = cp.maxCol[ic];
auto llxl = phase1PixelTopology::localX(llx);
auto llyl = phase1PixelTopology::localY(lly);
auto urxl = phase1PixelTopology::localX(urx);
auto uryl = phase1PixelTopology::localY(ury);
auto mx = llxl + urxl;
auto my = llyl + uryl;
auto xsize = int(urxl) + 2 - int(llxl);
auto ysize = int(uryl) + 2 - int(llyl);
assert(xsize >= 0); // 0 if bixpix...
assert(ysize >= 0);
if (phase1PixelTopology::isBigPixX(cp.minRow[ic]))
++xsize;
if (phase1PixelTopology::isBigPixX(cp.maxRow[ic]))
++xsize;
if (phase1PixelTopology::isBigPixY(cp.minCol[ic]))
++ysize;
if (phase1PixelTopology::isBigPixY(cp.maxCol[ic]))
++ysize;
int unbalanceX = 8. * std::abs(float(cp.Q_f_X[ic] - cp.Q_l_X[ic])) / float(cp.Q_f_X[ic] + cp.Q_l_X[ic]);
int unbalanceY = 8. * std::abs(float(cp.Q_f_Y[ic] - cp.Q_l_Y[ic])) / float(cp.Q_f_Y[ic] + cp.Q_l_Y[ic]);
xsize = 8 * xsize - unbalanceX;
ysize = 8 * ysize - unbalanceY;
cp.xsize[ic] = std::min(xsize, 1023);
cp.ysize[ic] = std::min(ysize, 1023);
if (cp.minRow[ic] == 0 || cp.maxRow[ic] == phase1PixelTopology::lastRowInModule)
cp.xsize[ic] = -cp.xsize[ic];
if (cp.minCol[ic] == 0 || cp.maxCol[ic] == phase1PixelTopology::lastColInModule)
cp.ysize[ic] = -cp.ysize[ic];
// apply the lorentz offset correction
auto xPos = detParams.shiftX + comParams.thePitchX * (0.5f * float(mx) + float(phase1PixelTopology::xOffset));
auto yPos = detParams.shiftY + comParams.thePitchY * (0.5f * float(my) + float(phase1PixelTopology::yOffset));
float cotalpha = 0, cotbeta = 0;
computeAnglesFromDet(detParams, xPos, yPos, cotalpha, cotbeta);
auto thickness = detParams.isBarrel ? comParams.theThicknessB : comParams.theThicknessE;
auto xcorr = correction(cp.maxRow[ic] - cp.minRow[ic],
cp.Q_f_X[ic],
cp.Q_l_X[ic],
llxl,
urxl,
detParams.chargeWidthX, // lorentz shift in cm
thickness,
cotalpha,
comParams.thePitchX,
phase1PixelTopology::isBigPixX(cp.minRow[ic]),
phase1PixelTopology::isBigPixX(cp.maxRow[ic]));
auto ycorr = correction(cp.maxCol[ic] - cp.minCol[ic],
cp.Q_f_Y[ic],
cp.Q_l_Y[ic],
llyl,
uryl,
detParams.chargeWidthY, // lorentz shift in cm
thickness,
cotbeta,
comParams.thePitchY,
phase1PixelTopology::isBigPixY(cp.minCol[ic]),
phase1PixelTopology::isBigPixY(cp.maxCol[ic]));
cp.xpos[ic] = xPos + xcorr;
cp.ypos[ic] = yPos + ycorr;
}
constexpr inline void errorFromSize(CommonParams const& __restrict__ comParams,
DetParams const& __restrict__ detParams,
ClusParams& cp,
uint32_t ic) {
// Edge cluster errors
cp.xerr[ic] = 0.0050;
cp.yerr[ic] = 0.0085;
// FIXME these are errors form Run1
constexpr float xerr_barrel_l1[] = {0.00115, 0.00120, 0.00088};
constexpr float xerr_barrel_l1_def = 0.00200; // 0.01030;
constexpr float yerr_barrel_l1[] = {
0.00375, 0.00230, 0.00250, 0.00250, 0.00230, 0.00230, 0.00210, 0.00210, 0.00240};
constexpr float yerr_barrel_l1_def = 0.00210;
constexpr float xerr_barrel_ln[] = {0.00115, 0.00120, 0.00088};
constexpr float xerr_barrel_ln_def = 0.00200; // 0.01030;
constexpr float yerr_barrel_ln[] = {
0.00375, 0.00230, 0.00250, 0.00250, 0.00230, 0.00230, 0.00210, 0.00210, 0.00240};
constexpr float yerr_barrel_ln_def = 0.00210;
constexpr float xerr_endcap[] = {0.0020, 0.0020};
constexpr float xerr_endcap_def = 0.0020;
constexpr float yerr_endcap[] = {0.00210};
constexpr float yerr_endcap_def = 0.00210;
auto sx = cp.maxRow[ic] - cp.minRow[ic];
auto sy = cp.maxCol[ic] - cp.minCol[ic];
// is edgy ?
bool isEdgeX = cp.minRow[ic] == 0 or cp.maxRow[ic] == phase1PixelTopology::lastRowInModule;
bool isEdgeY = cp.minCol[ic] == 0 or cp.maxCol[ic] == phase1PixelTopology::lastColInModule;
// is one and big?
bool isBig1X = (0 == sx) && phase1PixelTopology::isBigPixX(cp.minRow[ic]);
bool isBig1Y = (0 == sy) && phase1PixelTopology::isBigPixY(cp.minCol[ic]);
if (!isEdgeX && !isBig1X) {
if (not detParams.isBarrel) {
cp.xerr[ic] = sx < std::size(xerr_endcap) ? xerr_endcap[sx] : xerr_endcap_def;
} else if (detParams.layer == 1) {
cp.xerr[ic] = sx < std::size(xerr_barrel_l1) ? xerr_barrel_l1[sx] : xerr_barrel_l1_def;
} else {
cp.xerr[ic] = sx < std::size(xerr_barrel_ln) ? xerr_barrel_ln[sx] : xerr_barrel_ln_def;
}
}
if (!isEdgeY && !isBig1Y) {
if (not detParams.isBarrel) {
cp.yerr[ic] = sy < std::size(yerr_endcap) ? yerr_endcap[sy] : yerr_endcap_def;
} else if (detParams.layer == 1) {
cp.yerr[ic] = sy < std::size(yerr_barrel_l1) ? yerr_barrel_l1[sy] : yerr_barrel_l1_def;
} else {
cp.yerr[ic] = sy < std::size(yerr_barrel_ln) ? yerr_barrel_ln[sy] : yerr_barrel_ln_def;
}
}
}
constexpr inline void errorFromDB(CommonParams const& __restrict__ comParams,
DetParams const& __restrict__ detParams,
ClusParams& cp,
uint32_t ic) {
// Edge cluster errors
cp.xerr[ic] = 0.0050f;
cp.yerr[ic] = 0.0085f;
auto sx = cp.maxRow[ic] - cp.minRow[ic];
auto sy = cp.maxCol[ic] - cp.minCol[ic];
// is edgy ?
bool isEdgeX = cp.minRow[ic] == 0 or cp.maxRow[ic] == phase1PixelTopology::lastRowInModule;
bool isEdgeY = cp.minCol[ic] == 0 or cp.maxCol[ic] == phase1PixelTopology::lastColInModule;
// is one and big?
uint32_t ix = (0 == sx);
uint32_t iy = (0 == sy);
ix += (0 == sx) && phase1PixelTopology::isBigPixX(cp.minRow[ic]);
iy += (0 == sy) && phase1PixelTopology::isBigPixY(cp.minCol[ic]);
if (not isEdgeX)
cp.xerr[ic] = detParams.sx[ix];
if (not isEdgeY)
cp.yerr[ic] = detParams.sy[iy];
}
} // namespace pixelCPEforGPU
#endif // RecoLocalTracker_SiPixelRecHits_pixelCPEforGPU_h