Project
Loading...
Searching...
No Matches
TPCFastTransformPOD.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
16
17#ifndef ALICEO2_GPU_TPCFastTransformPOD_H
18#define ALICEO2_GPU_TPCFastTransformPOD_H
19
20#include "GPUCommonRtypes.h"
21#include "TPCFastTransform.h"
24#ifndef GPUCA_GPUCODE
25#include <memory>
26#include <cstdlib>
28#endif
29
30/*
31Binary buffer should be cast to TPCFastTransformPOD class using static TPCFastTransformPOD& t = get(buffer); method,
32so that its head becomes `this` pointer of the object.
33
34First we have all the fixed size data members mentioned explicitly. Part of them is duplicating fixed size
35data members of TPCFastSpaceChargeCorrection but those starting with mOffs... provide the offset in bytes
36(wrt this) for dynamic data which cannot be declared as data member explicitly (since we cannot have any
37pointer except `this`) but obtained via getters using stored offsets wrt `this`.
38This is followed dynamic part itself.
39
40dynamic part layout:
411) size_t[ mNumberOfScenarios ] array starting at offset mOffsScenariosOffsets, each element is the offset
42of distict spline object (scenario in TPCFastSpaceChargeCorrection)
432) size_t[ mNSplineIDs ] array starting at offset mOffsSplineDataOffsets, each element is the offset of the
44beginning of splines data for give splineID
45
46*/
47
48namespace o2::gpu
49{
51{
52 public:
56
61
63 GPUd() static const TPCFastTransformPOD& get(const char* head) { return *reinterpret_cast<const TPCFastTransformPOD*>(head); }
64
67 // Methods taking extra reference transform are legacy compound transforms used to scale corrections.
68 GPUd() void Transform(int32_t sector, int32_t row, float pad, float time, float& x, float& y, float& z, float vertexTime = 0) const;
69 GPUd() void TransformXYZ(int32_t sector, int32_t row, float& x, float& y, float& z) const;
70
72 GPUd() void TransformInTimeFrame(int32_t sector, int32_t row, float pad, float time, float& x, float& y, float& z, float maxTimeBin) const;
73 GPUd() void TransformInTimeFrame(int32_t sector, float time, float& z, float maxTimeBin) const;
74
76 GPUd() void InverseTransformInTimeFrame(int32_t sector, int32_t row, float /*x*/, float y, float z, float& pad, float& time, float maxTimeBin) const;
77 GPUd() float InverseTransformInTimeFrame(int32_t sector, float z, float maxTimeBin) const;
78
80 GPUd() void InverseTransformYZtoX(int32_t sector, int32_t row, float y, float z, float& x) const;
81
83 GPUd() void InverseTransformYZtoNominalYZ(int32_t sector, int32_t row, float y, float z, float& ny, float& nz) const;
84
86 GPUd() void InverseTransformXYZtoNominalXYZ(int32_t sector, int32_t row, float x, float y, float z, float& nx, float& ny, float& nz) const;
87
89 GPUd() void TransformIdeal(int32_t sector, int32_t row, float pad, float time, float& x, float& y, float& z, float vertexTime) const;
90 GPUd() void TransformIdealZ(int32_t sector, float time, float& z, float vertexTime) const;
91
92 GPUd() void convPadTimeToLocal(int32_t sector, int32_t row, float pad, float time, float& y, float& z, float vertexTime) const;
93 GPUd() void convPadTimeToLocalInTimeFrame(int32_t sector, int32_t row, float pad, float time, float& y, float& z, float maxTimeBin) const;
94
95 GPUd() void convLocalToPadTime(int32_t sector, int32_t row, float y, float z, float& pad, float& time, float vertexTime) const;
96 GPUd() void convLocalToPadTimeInTimeFrame(int32_t sector, int32_t row, float y, float z, float& pad, float& time, float maxTimeBin) const;
97
98 GPUd() float convTimeToZinTimeFrame(int32_t sector, float time, float maxTimeBin) const;
99 GPUd() float convZtoTimeInTimeFrame(int32_t sector, float z, float maxTimeBin) const;
100 GPUd() float convDeltaTimeToDeltaZinTimeFrame(int32_t sector, float deltaTime) const;
101 GPUd() float convDeltaZtoDeltaTimeInTimeFrame(int32_t sector, float deltaZ) const;
102 GPUd() float convDeltaZtoDeltaTimeInTimeFrameAbs(float deltaZ) const;
103 GPUd() float convZOffsetToVertexTime(int32_t sector, float zOffset, float maxTimeBin) const;
104 GPUd() float convVertexTimeToZOffset(int32_t sector, float vertexTime, float maxTimeBin) const;
105
107 void setApplyCorrectionOn() { mApplyCorrection = 1; }
108 void setApplyCorrectionOff() { mApplyCorrection = 0; }
109 bool isCorrectionApplied() { return mApplyCorrection; }
110
112 GPUd() const TPCFastTransformGeoPOD& getGeometry() const { return mGeo; }
113
115 GPUd() const RowInfo& getRowInfo(int32_t row) const { return mRowInfos[row]; }
116
118 GPUd() RowInfo& getRowInfo(int32_t row) { return mRowInfos[row]; }
119
121 GPUd() size_t size() const { return mTotalSize; }
122
124 GPUd() long int getTimeStamp() const { return mTimeStamp; }
125
127 GPUd() float getVDrift() const { return mVdrift; }
128
130 GPUd() float getT0() const { return mT0; }
131
133 GPUd() float getIDC() const { return mIDC; }
134
136 GPUd() float getLumi() const { return mLumi; }
137
139 GPUd() float getMaxDriftTime(int32_t sector, int32_t row, float pad) const;
140
142 GPUd() float getMaxDriftTime(int32_t sector, int32_t row) const;
143
145 GPUd() float getMaxDriftTime(int32_t sector) const;
146
148 GPUd() void setTimeStamp(long int v) { mTimeStamp = v; }
149
151 GPUd() void setVDrift(float v) { mVdrift = v; }
152
154 GPUd() void setT0(float v) { mT0 = v; }
155
157 GPUd() void setIDC(float v) { mIDC = v; }
158
160 GPUd() void setLumi(float v) { mLumi = v; }
161
162 GPUd() void setCalibration(int64_t timeStamp, float t0, float vDrift);
163
165 GPUd() const SplineType& getSplineForRow(int32_t row) const { return *reinterpret_cast<const SplineType*>(getThis() + getScenarioOffset(getRowInfo(row).splineScenarioID)); }
166
168 GPUd() const float* getCorrectionData(int32_t sector, int32_t row, int32_t iSpline = 0) const { return reinterpret_cast<const float*>(getThis() + mSplineDataOffsets[sector][iSpline] + getRowInfo(row).dataOffsetBytes[iSpline]); }
169
171 GPUd() const SplineTypeInvX& getSplineInvXforRow(int32_t row) const { return reinterpret_cast<const SplineTypeInvX&>(getSplineForRow(row)); }
172
174 GPUd() const float* getCorrectionDataInvX(int32_t sector, int32_t row) const { return getCorrectionData(sector, row, 1); }
175
177 GPUd() const SplineTypeInvYZ& getSplineInvYZforRow(int32_t row) const { return reinterpret_cast<const SplineTypeInvYZ&>(getSplineForRow(row)); }
178
180 GPUd() const float* getCorrectionDataInvYZ(int32_t sector, int32_t row) const { return getCorrectionData(sector, row, 2); }
181
183 GPUdi() void getCorrectionLocal(int32_t sector, int32_t row, float y, float z, float& dx, float& dy, float& dz) const;
184
186 GPUd() float getCorrectionXatRealYZ(int32_t sector, int32_t row, float realY, float realZ) const;
187
189 GPUd() void getCorrectionYZatRealYZ(int32_t sector, int32_t row, float realY, float realZ, float& measuredY, float& measuredZ) const;
190
192 GPUd() void TransformLocal(int32_t sector, int32_t row, float& x, float& y, float& z) const;
193
195
198 GPUd() void convLocalToGrid(int32_t sector, int32_t row, float y, float z, float& u, float& v, float& s) const;
199
202 GPUd() void convGridToLocal(int32_t sector, int32_t row, float u, float v, float& y, float& z) const;
203
206 GPUd() void convRealLocalToGrid(int32_t sector, int32_t row, float y, float z, float& u, float& v, float& s) const;
207
210 GPUd() void convGridToRealLocal(int32_t sector, int32_t row, float u, float v, float& y, float& z) const;
211
212 GPUd() bool isLocalInsideGrid(int32_t sector, int32_t row, float y, float z) const;
213 GPUd() bool isRealLocalInsideGrid(int32_t sector, int32_t row, float y, float z) const;
214
215#if !defined(GPUCA_GPUCODE)
218
221
223 {
224 destVector.alloc(src.size());
225 std::memcpy(destVector.get(), &src, src.size());
226 return destVector.get();
227 }
228
229 static TPCFastTransformPOD* create(char* buff, size_t buffSize, const TPCFastTransform& src);
230 static TPCFastTransformPOD* create(char* buff, size_t buffSize, const TPCFastSpaceChargeCorrection& src);
231 static size_t estimateSize(const TPCFastTransform& src) { return estimateSize(src.getCorrection()); }
232 static size_t estimateSize(const TPCFastSpaceChargeCorrection& origCorr);
233
234 bool test(const TPCFastTransform& src, int32_t npoints = 100000) const { return test(src.getCorrection(), npoints); }
235 bool test(const TPCFastSpaceChargeCorrection& origCorr, int32_t npoints = 100000) const;
236#endif
237
239 void print() const;
240
241 GPUd() float convDriftLengthToTime(float driftLength, float vertexTime) const;
242
243 static constexpr int NROWS = o2::tpc::constants::MAXGLOBALPADROW;
244 static constexpr int NSECTORS = o2::tpc::constants::MAXSECTOR;
245 static constexpr int NSECTORSA = o2::tpc::constants::MAXSECTOR / 2;
246 static constexpr int NSplineIDs = 3;
247
248 private:
249#if !defined(GPUCA_GPUCODE)
250 static constexpr size_t AlignmentBytes = 8;
251 static size_t alignOffset(size_t offs)
252 {
253 auto res = offs % AlignmentBytes;
254 return res ? offs + (AlignmentBytes - res) : offs;
255 }
256 GPUd() static TPCFastTransformPOD& getNonConst(char* head) { return *reinterpret_cast<TPCFastTransformPOD*>(head); }
257#endif
258
260 GPUd() const char* getThis() const { return reinterpret_cast<const char*>(this); }
261
263 GPUd() size_t getScenarioOffset(int s) const { return (reinterpret_cast<const size_t*>(getThis() + mOffsScenariosOffsets))[s]; }
264
265 GPUd() size_t getFlatBufferOffset(int s) const { return (reinterpret_cast<const size_t*>(getThis() + mOffsFlatBufferOffsets))[s]; }
266
267 // Returns a pointer to the flat buffer of scenario isc, using only the
268 // stored offset array (mOffsFlatBufferOffsets). No stale pointer involved.
269 GPUd() const char* getSplineFlatBuffer(int32_t isc) const
270 {
271 const size_t* offs = reinterpret_cast<const size_t*>(getThis() + mOffsFlatBufferOffsets);
272 return getThis() + offs[isc];
273 }
274
275 // Returns a pointer to mGridX2's flat buffer inside the spline flat buffer.
276 // Reproduces the layout from Spline2DContainer::setActualBufferAddress using
277 // only safe values: getFlatBufferSize() reads mNumberOfKnots/mUmax (plain ints).
278 template <typename SplineT>
279 GPUd() const char* getGridX2FlatBuffer(const SplineT& spline, int32_t isc) const
280 {
281 const size_t g1sz = spline.getGridX1().getFlatBufferSize();
282 const size_t g2align = spline.getGridX2().getBufferAlignmentBytes();
283 return getSplineFlatBuffer(isc) + FlatObject::alignSize(g1sz, g2align);
284 }
285
286 bool mApplyCorrection{};
287 int mNumberOfScenarios{};
288 size_t mTotalSize{};
289 size_t mOffsScenariosOffsets{};
290 size_t mOffsFlatBufferOffsets{};
291 size_t mSplineDataOffsets[TPCFastTransformGeo::getNumberOfSectors()][NSplineIDs];
292 long int mTimeStamp{};
293 float mT0;
294 float mVdrift;
295 float mLumi;
296 float mIDC;
297
298 TPCFastTransformGeoPOD mGeo;
299 RowInfo mRowInfos[NROWS];
300
301 ClassDefNV(TPCFastTransformPOD, 0);
302};
303
304GPUdi() void TPCFastTransformPOD::getCorrectionLocal(int32_t sector, int32_t row, float y, float z, float& dx, float& dy, float& dz) const
305{
306 const auto& info = getRowInfo(row);
307 const int32_t isc = info.splineScenarioID;
308 const SplineType& spline = getSplineForRow(row);
309 const float* splineData = getCorrectionData(sector, row);
310
311 float u, v, s;
312 convLocalToGrid(sector, row, y, z, u, v, s);
313
314 const char* g1buf = getSplineFlatBuffer(isc);
315 const char* g2buf = getGridX2FlatBuffer(spline, isc);
316
317 float dxyz[3];
318 spline.interpolateAtUZeroCopy(g1buf, g2buf, splineData, u, v, dxyz);
319
320 if (CAMath::Abs(dxyz[0]) > TPCFastSpaceChargeCorrection::kMaxCorrection || CAMath::Abs(dxyz[1]) > TPCFastSpaceChargeCorrection::kMaxCorrection || CAMath::Abs(dxyz[2]) > TPCFastSpaceChargeCorrection::kMaxCorrection) {
321 s = 0.f; // TODO: DR: Protect from FPEs, fix upstream and remove once guaranteed that it is fixed
322 }
323
324 dx = s * dxyz[0];
325 dy = s * dxyz[1];
326 dz = s * dxyz[2];
327}
328
329GPUdi() float TPCFastTransformPOD::getCorrectionXatRealYZ(int32_t sector, int32_t row, float realY, float realZ) const
330{
331 const auto& info = getRowInfo(row);
332 float u, v, s;
333 convRealLocalToGrid(sector, row, realY, realZ, u, v, s);
334
335 const int32_t isc = info.splineScenarioID;
336 const auto& spline = getSplineInvXforRow(row);
337 const char* g1buf = getSplineFlatBuffer(isc);
338 const char* g2buf = getGridX2FlatBuffer(spline, isc);
339
340 float dx = 0;
341 spline.interpolateAtUZeroCopy(g1buf, g2buf, getCorrectionDataInvX(sector, row), u, v, &dx);
342 if (CAMath::Abs(dx) > TPCFastSpaceChargeCorrection::kMaxCorrection) {
343 s = 0.f; // TODO: DR: Protect from FPEs, fix upstream and remove once guaranteed that it is fixed
344 }
345 dx = s * dx;
346 return dx;
347}
348
349GPUdi() void TPCFastTransformPOD::getCorrectionYZatRealYZ(int32_t sector, int32_t row, float realY, float realZ, float& y, float& z) const
350{
351 float u, v, s;
352 convRealLocalToGrid(sector, row, realY, realZ, u, v, s);
353 const auto& info = getRowInfo(row);
354 const int32_t isc = info.splineScenarioID;
355 const auto& spline = getSplineInvYZforRow(row);
356 const char* g1buf = getSplineFlatBuffer(isc);
357 const char* g2buf = getGridX2FlatBuffer(spline, isc);
358
359 float dyz[2];
360 spline.interpolateAtUZeroCopy(g1buf, g2buf, getCorrectionDataInvYZ(sector, row), u, v, dyz);
361 if (CAMath::Abs(dyz[0]) > TPCFastSpaceChargeCorrection::kMaxCorrection || CAMath::Abs(dyz[1]) > TPCFastSpaceChargeCorrection::kMaxCorrection) {
362 s = 0.f; // TODO: DR: Protect from FPEs, fix upstream and remove once guaranteed that it is fixed
363 }
364 y = s * dyz[0];
365 z = s * dyz[1];
366}
367
368GPUdi() void TPCFastTransformPOD::convLocalToGrid(int32_t sector, int32_t row, float y, float z, float& u, float& v, float& s) const
369{
372 const SplineType& spline = getSplineForRow(row);
373 getRowInfo(row).gridMeasured.convLocalToGridUntruncated(sector, y, z, u, v, s);
374 // shrink to the grid
375 u = GPUCommonMath::Clamp(u, 0.f, (float)spline.getGridX1().getUmax());
376 v = GPUCommonMath::Clamp(v, 0.f, (float)spline.getGridX2().getUmax());
377}
378
379GPUdi() void TPCFastTransformPOD::convGridToLocal(int32_t sector, int32_t row, float gridU, float gridV, float& y, float& z) const
380{
382 getRowInfo(row).gridMeasured.convGridToLocal(sector, gridU, gridV, y, z);
383}
384
385GPUdi() void TPCFastTransformPOD::convRealLocalToGrid(int32_t sector, int32_t row, float y, float z, float& u, float& v, float& s) const
386{
388 const SplineType& spline = getSplineForRow(row);
389 getRowInfo(row).gridReal.convLocalToGridUntruncated(sector, y, z, u, v, s);
390 // shrink to the grid
391 u = GPUCommonMath::Clamp(u, 0.f, (float)spline.getGridX1().getUmax());
392 v = GPUCommonMath::Clamp(v, 0.f, (float)spline.getGridX2().getUmax());
393}
394
395GPUdi() void TPCFastTransformPOD::convGridToRealLocal(int32_t sector, int32_t row, float gridU, float gridV, float& y, float& z) const
396{
398 getRowInfo(row).gridReal.convGridToLocal(sector, gridU, gridV, y, z);
399}
400
401GPUdi() bool TPCFastTransformPOD::isLocalInsideGrid(int32_t sector, int32_t row, float y, float z) const
402{
404 float u, v, s;
405 getRowInfo(row).gridMeasured.convLocalToGridUntruncated(sector, y, z, u, v, s);
406 const auto& spline = getSplineForRow(row);
407 // shrink to the grid
408 if (u < 0.f || u > (float)spline.getGridX1().getUmax() || //
409 v < 0.f || v > (float)spline.getGridX2().getUmax()) {
410 return false;
411 }
412 return true;
413}
414
415GPUdi() bool TPCFastTransformPOD::isRealLocalInsideGrid(int32_t sector, int32_t row, float y, float z) const
416{
418 float u, v, s;
419 getRowInfo(row).gridReal.convLocalToGridUntruncated(sector, y, z, u, v, s);
420 const auto& spline = getSplineForRow(row);
421 // shrink to the grid
422 if (u < 0.f || u > (float)spline.getGridX1().getUmax() || //
423 v < 0.f || v > (float)spline.getGridX2().getUmax()) {
424 return false;
425 }
426 return true;
427}
428
429GPUdi() void TPCFastTransformPOD::TransformLocal(int32_t sector, int32_t row, float& x, float& y, float& z) const
430{
431 if (!mApplyCorrection) {
432 return;
433 }
434 float dx, dy, dz;
435 getCorrectionLocal(sector, row, y, z, dx, dy, dz);
436
437 GPUCA_DEBUG_STREAMER_CHECK(if (o2::utils::DebugStreamer::checkStream(o2::utils::StreamFlags::streamFastTransform)) {
438 float lx = x, ly = y, lz = z;
439 float gx, gy, gz;
440 getGeometry().convLocalToGlobal(sector, lx, ly, lz, gx, gy, gz);
441 float lxT = lx + dx;
442 float lyT = ly + dy;
443 float lzT = lz + dz;
444 float invYZtoX;
445 InverseTransformYZtoX(sector, row, lyT, lzT, invYZtoX);
446
447 float YZtoNominalY;
448 float YZtoNominalZ;
449 InverseTransformYZtoNominalYZ(sector, row, lyT, lzT, YZtoNominalY, YZtoNominalZ);
450
451 o2::utils::DebugStreamer::instance()->getStreamer("debug_fasttransform", "UPDATE") << o2::utils::DebugStreamer::instance()->getUniqueTreeName("tree_Transform").data()
452 // corrections in x, u, v
453 << "dx=" << dx
454 << "dy=" << dy
455 << "dz=" << dz
456 << "row=" << row
457 << "sector=" << sector
458 // original local coordinates
459 << "ly=" << ly
460 << "lz=" << lz
461 << "lx=" << lx
462 // corrected local coordinated
463 << "lxT=" << lxT
464 << "lyT=" << lyT
465 << "lzT=" << lzT
466 // global uncorrected coordinates
467 << "gx=" << gx
468 << "gy=" << gy
469 << "gz=" << gz
470 // some transformations which are applied
471 << "invYZtoX=" << invYZtoX
472 << "YZtoNominalY=" << YZtoNominalY
473 << "YZtoNominalZ=" << YZtoNominalZ
474 << "\n";
475 })
476
477 x += dx;
478 y += dy;
479 z += dz;
480}
481
482GPUdi() void TPCFastTransformPOD::Transform(int32_t sector, int32_t row, float pad, float time, float& x, float& y, float& z, float vertexTime) const
483{
489
490 x = getGeometry().getRowInfoX(row);
491 convPadTimeToLocal(sector, row, pad, time, y, z, vertexTime);
492 TransformLocal(sector, row, x, y, z);
493}
494
495GPUdi() void TPCFastTransformPOD::TransformXYZ(int32_t sector, int32_t row, float& x, float& y, float& z) const
496{
497
498 TransformLocal(sector, row, x, y, z);
499}
500
501GPUdi() void TPCFastTransformPOD::TransformInTimeFrame(int32_t sector, float time, float& z, float maxTimeBin) const
502{
503 float l = (time - mT0 - maxTimeBin) * mVdrift; // drift length cm
504 z = getGeometry().convDriftLengthToZ(sector, l);
505}
506
507GPUdi() void TPCFastTransformPOD::TransformInTimeFrame(int32_t sector, int32_t row, float pad, float time, float& x, float& y, float& z, float maxTimeBin) const
508{
514
515 x = getGeometry().getRowInfoX(row);
516 convPadTimeToLocalInTimeFrame(sector, row, pad, time, y, z, maxTimeBin);
517}
518
519GPUdi() void TPCFastTransformPOD::InverseTransformInTimeFrame(int32_t sector, int32_t row, float /*x*/, float y, float z, float& pad, float& time, float maxTimeBin) const
520{
522 convLocalToPadTimeInTimeFrame(sector, row, y, z, pad, time, maxTimeBin);
523}
524
525GPUdi() float TPCFastTransformPOD::InverseTransformInTimeFrame(int32_t sector, float z, float maxTimeBin) const
526{
527 float pad, time;
528 InverseTransformInTimeFrame(sector, 0, 0, 0, z, pad, time, maxTimeBin);
529 return time;
530}
531
532GPUdi() void TPCFastTransformPOD::TransformIdealZ(int32_t sector, float time, float& z, float vertexTime) const
533{
540
541 float l = (time - mT0 - vertexTime) * mVdrift; // drift length cm
542 z = getGeometry().convDriftLengthToZ(sector, l);
543}
544
545GPUdi() void TPCFastTransformPOD::TransformIdeal(int32_t sector, int32_t row, float pad, float time, float& x, float& y, float& z, float vertexTime) const
546{
553
554 x = getGeometry().getRowInfoX(row);
555 float driftLength = (time - mT0 - vertexTime) * mVdrift; // drift length cm
556 getGeometry().convPadDriftLengthToLocal(sector, row, pad, driftLength, y, z);
557}
558
559GPUdi() float TPCFastTransformPOD::convTimeToZinTimeFrame(int32_t sector, float time, float maxTimeBin) const
560{
567
568 float v = (time - mT0 - maxTimeBin) * mVdrift; // drift length cm
569 float z = (sector < getGeometry().getNumberOfSectorsA()) ? -v : v;
570 return z;
571}
572
573GPUdi() float TPCFastTransformPOD::convZtoTimeInTimeFrame(int32_t sector, float z, float maxTimeBin) const
574{
576 float v = (sector < getGeometry().getNumberOfSectorsA()) ? -z : z;
577 return mT0 + maxTimeBin + v / mVdrift;
578}
579
580GPUdi() float TPCFastTransformPOD::convDeltaTimeToDeltaZinTimeFrame(int32_t sector, float deltaTime) const
581{
582 float deltaZ = deltaTime * mVdrift;
583 return sector < getGeometry().getNumberOfSectorsA() ? -deltaZ : deltaZ;
584}
585
586GPUdi() float TPCFastTransformPOD::convDeltaZtoDeltaTimeInTimeFrameAbs(float deltaZ) const
587{
588 return deltaZ / mVdrift;
589}
590
591GPUdi() float TPCFastTransformPOD::convDeltaZtoDeltaTimeInTimeFrame(int32_t sector, float deltaZ) const
592{
593 float deltaT = deltaZ / mVdrift;
594 return sector < getGeometry().getNumberOfSectorsA() ? -deltaT : deltaT;
595}
596
597GPUdi() float TPCFastTransformPOD::getMaxDriftTime(int32_t sector, int32_t row, float pad) const
598{
600 return convDriftLengthToTime(getGeometry().getTPCzLength(), 0.f);
601}
602
603GPUdi() float TPCFastTransformPOD::getMaxDriftTime(int32_t sector, int32_t row) const
604{
606 return convDriftLengthToTime(getGeometry().getTPCzLength(), 0.f);
607}
608
609GPUdi() float TPCFastTransformPOD::getMaxDriftTime(int32_t sector) const
610{
612 return convDriftLengthToTime(getGeometry().getTPCzLength(), 0.f);
613}
614
615GPUdi() void TPCFastTransformPOD::InverseTransformYZtoX(int32_t sector, int32_t row, float realY, float realZ, float& realX) const
616{
618 float dx = 0.f;
619 dx = getCorrectionXatRealYZ(sector, row, realY, realZ);
620 realX = getGeometry().getRowInfoX(row) + dx;
621
622 GPUCA_DEBUG_STREAMER_CHECK(if (o2::utils::DebugStreamer::checkStream(o2::utils::StreamFlags::streamFastTransform)) {
623 o2::utils::DebugStreamer::instance()->getStreamer("debug_fasttransform", "UPDATE") << o2::utils::DebugStreamer::instance()->getUniqueTreeName("tree_InverseTransformYZtoX").data()
624 << "sector=" << sector
625 << "row=" << row
626 << "y=" << realY
627 << "z=" << realZ
628 << "x=" << realX
629 << "\n";
630 })
631}
632
633GPUdi() void TPCFastTransformPOD::InverseTransformYZtoNominalYZ(int32_t sector, int32_t row, float realY, float realZ, float& measuredY, float& measuredZ) const
634{
636 float dy, dz;
637 getCorrectionYZatRealYZ(sector, row, realY, realZ, dy, dz);
638 measuredY = realY - dy;
639 measuredZ = realZ - dz;
640
641 GPUCA_DEBUG_STREAMER_CHECK(if (o2::utils::DebugStreamer::checkStream(o2::utils::StreamFlags::streamFastTransform)) {
642 o2::utils::DebugStreamer::instance()->getStreamer("debug_fasttransform", "UPDATE") << o2::utils::DebugStreamer::instance()->getUniqueTreeName("tree_InverseTransformYZtoNominalYZ").data()
643 << "sector=" << sector
644 << "row=" << row
645 << "real y=" << realY
646 << "real z=" << realZ
647 << "measured y=" << measuredY
648 << "measured z=" << measuredZ
649 << "\n";
650 })
651}
652
653GPUdi() void TPCFastTransformPOD::convPadTimeToLocal(int32_t sector, int32_t row, float pad, float time, float& y, float& z, float vertexTime) const
654{
655 float l = (time - mT0 - vertexTime) * mVdrift;
656 getGeometry().convPadDriftLengthToLocal(sector, row, pad, l, y, z);
657}
658
659GPUdi() void TPCFastTransformPOD::convPadTimeToLocalInTimeFrame(int32_t sector, int32_t row, float pad, float time, float& y, float& z, float maxTimeBin) const
660{
661 float l = getGeometry().getTPCzLength() + (time - mT0 - maxTimeBin) * mVdrift;
662 getGeometry().convPadDriftLengthToLocal(sector, row, pad, l, y, z);
663}
664
665GPUdi() void TPCFastTransformPOD::convLocalToPadTimeInTimeFrame(int32_t sector, int32_t row, float y, float z, float& pad, float& time, float maxTimeBin) const
666{
667 float length = 0;
668 getGeometry().convLocalToPadDriftLength(sector, row, y, z, pad, length);
669 time = convDriftLengthToTime(length, maxTimeBin);
670}
671
672GPUdi() float TPCFastTransformPOD::convDriftLengthToTime(float driftLength, float vertexTime) const
673{
674 return (mT0 + vertexTime + driftLength / mVdrift);
675}
676
677GPUdi() float TPCFastTransformPOD::convZOffsetToVertexTime(int32_t sector, float zOffset, float maxTimeBin) const
678{
679 if (sector < getGeometry().getNumberOfSectorsA()) {
680 return maxTimeBin - (getGeometry().getTPCzLength() + zOffset) / mVdrift;
681 } else {
682 return maxTimeBin - (getGeometry().getTPCzLength() - zOffset) / mVdrift;
683 }
684}
685
686GPUdi() float TPCFastTransformPOD::convVertexTimeToZOffset(int32_t sector, float vertexTime, float maxTimeBin) const
687{
688 if (sector < getGeometry().getNumberOfSectorsA()) {
689 return (maxTimeBin - vertexTime) * mVdrift - getGeometry().getTPCzLength();
690 } else {
691 return -((maxTimeBin - vertexTime) * mVdrift - getGeometry().getTPCzLength());
692 }
693}
694
695#ifndef GPUCA_GPUCODE_DEVICE // Functions not needed during GPU processing
696GPUdi() void TPCFastTransformPOD::setCalibration(int64_t timeStamp, float t0, float vDrift)
697{
698 mTimeStamp = timeStamp;
699 mT0 = t0;
700 mVdrift = vDrift;
701}
702
703GPUdi() void TPCFastTransformPOD::InverseTransformXYZtoNominalXYZ(int32_t sector, int32_t row, float x, float y, float z, float& nx, float& ny, float& nz) const
704{
706 int32_t row2 = row + 1;
707 if (row2 >= getGeometry().getNumberOfRows()) {
708 row2 = row - 1;
709 }
710 float nx1, ny1, nz1; // nominal coordinates for row
711 float nx2, ny2, nz2; // nominal coordinates for row2
712 nx1 = getGeometry().getRowInfoX(row);
713 nx2 = getGeometry().getRowInfoX(row2);
714 InverseTransformYZtoNominalYZ(sector, row, y, z, ny1, nz1);
715 InverseTransformYZtoNominalYZ(sector, row2, y, z, ny2, nz2);
716 float c1 = (nx2 - nx) / (nx2 - nx1);
717 float c2 = (nx - nx1) / (nx2 - nx1);
718 nx = x;
719 ny = (ny1 * c1 + ny2 * c2);
720 nz = (nz1 * c1 + nz2 * c2);
721}
722#endif // GPUCA_GPUCODE_DEVICE
723
724} // namespace o2::gpu
725
726#endif
int16_t time
Definition RawEventData.h:4
#define GPUCA_DEBUG_STREAMER_CHECK(...)
uint32_t res
Definition RawData.h:0
Version using constexpr GPUTPCGeometry to be used for TPCFastTransformationPOD.
Definition of TPCFastTransform class.
static constexpr size_t alignSize(size_t sizeBytes, size_t alignmentBytes)
_______________ Generic utilities _______________________________________________
Definition FlatObject.h:275
Forward declaration — specializations below select ClassDefNV based on FlatBase.
Definition Spline2D.h:104
Spline2D< float, 2, NoFlatObject > SlimSplineTypeInvYZ
Spline2D< float, 3, NoFlatObject > SlimSplineTypeXYZ
Slim variants (NoFlatObject base) for use in TPCFastTransformPOD.
Spline2D< float, 1, NoFlatObject > SlimSplineTypeInvX
static constexpr int32_t getNumberOfSectors()
_______________ Getters _________________________________
TPCFastSpaceChargeCorrection::SlimSplineTypeInvYZ SplineTypeInvYZ
static TPCFastTransformPOD * create(aligned_unique_buffer_ptr< TPCFastTransformPOD > &destVector, const TPCFastTransformPOD &src)
bool test(const TPCFastSpaceChargeCorrection &origCorr, int32_t npoints=100000) const
int32_t float float float float float float vertexTime
GPUd() RowInfo &getRowInfo(int32_t row)
Gives TPC sector & row info.
bool test(const TPCFastTransform &src, int32_t npoints=100000) const
GPUd() float getIDC() const
Return IDC estimator.
static TPCFastTransformPOD * create(aligned_unique_buffer_ptr< TPCFastTransformPOD > &destVector, const TPCFastTransform &src)
Create POD transform from old flat-buffer one. Provided vector will serve as a buffer.
GPUd() const RowInfo &getRowInfo(int32_t row) const
Gives TPC sector & row info.
GPUd() const TPCFastTransformGeoPOD &getGeometry() const
TPC geometry information.
GPUd() long int getTimeStamp() const
Gives the time stamp of the current calibaration parameters.
GPUd() float getVDrift() const
Return mVDrift in cm / time bin.
GPUd() void setLumi(float v)
Sets CTP Lumi estimator.
int32_t float float float float & v
GPUd() void setIDC(float v)
Sets IDC estimator.
GPUd() float getLumi() const
Return Lumi estimator.
void print() const
Print method.
int32_t float float float & u
int32_t float float float & ny
GPUd() const float *getCorrectionData(int32_t sector
Gives pointer to spline data.
GPUd() void setVDrift(float v)
Sets current vdrift.
GPUd() static const TPCFastTransformPOD &get(const char *head)
convert prefilled buffer to TPCFastTransformPOD
void setApplyCorrectionOn()
_______________ methods a la TPCFastSpaceChargeCorrection: cluster correction _______________________
int32_t float float float &z const
GPUd() void Transform(int32_t sector
GPUd() float getT0() const
Return T0 in time bin units.
GPUd() void setCalibration(int64_t timeStamp
GPUd() void setT0(float v)
Sets current T0.
static size_t estimateSize(const TPCFastTransform &src)
int32_t float float float float & nx
GPUd() size_t size() const
Gives its own size including dynamic part.
TPCFastSpaceChargeCorrection::SlimSplineTypeInvX SplineTypeInvX
static constexpr int NSplineIDs
number of spline data sets for each sector/row
TPCFastSpaceChargeCorrection::RowInfo RowInfo
int32_t float float float & measuredY
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum src
Definition glcorearb.h:1767
GLsizeiptr size
Definition glcorearb.h:659
const GLdouble * v
Definition glcorearb.h:832
GLint y
Definition glcorearb.h:270
GLuint GLsizei GLsizei * length
Definition glcorearb.h:790
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat t0
Definition glcorearb.h:5034
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
GPUdi() o2
Definition TrackTRD.h:39
@ streamFastTransform
stream tpc fast transform
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::vector< int > row