CbmRoot
Loading...
Searching...
No Matches
Hitfinder.cxx
Go to the documentation of this file.
1/* Copyright (C) 2022 FIAS Frankfurt Institute for Advanced Studies, Frankfurt / Main
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Felix Weiglhofer [committer], Kilian Hunold */
4
5#include "Hitfinder.h"
6
7using namespace cbm::algo;
8
9// Export Constants
11
12
13// Kernel Definitions
15XPU_D void sts::SortDigis::operator()(context& ctx) { ctx.cmem<sts::TheHitfinder>().SortDigisInSpaceAndTime(ctx); }
16
18XPU_D void sts::FindClusters::operator()(context& ctx) { ctx.cmem<sts::TheHitfinder>().FindClustersSingleStep(ctx); }
19
22{
23 ctx.cmem<sts::TheHitfinder>().CalculateOffsetsParallel(ctx);
24}
25
28{
29 ctx.cmem<sts::TheHitfinder>().FindClustersParallel(ctx);
30}
31
34{
35 ctx.cmem<sts::TheHitfinder>().CalculateClustersParallel(ctx);
36}
37
40
42XPU_D void sts::FindHits::operator()(context& ctx) { ctx.cmem<sts::TheHitfinder>().FindHits(ctx); }
43
53{
54 int iModule = ctx.block_idx_x();
55 CbmStsDigi* digis = &digisPerModule[digiOffsetPerModule[iModule]];
56 CbmStsDigi* digisTmp = &digisPerModuleTmp[digiOffsetPerModule[iModule]];
57 int nDigis = GetNDigis(iModule);
58
59 SortDigis::sort_t digiSort(ctx.pos(), ctx.smem());
60
61 CbmStsDigi* digisOut = digiSort.sort(digis, nDigis, digisTmp, [](const CbmStsDigi a) {
62 return ((unsigned long int) a.GetChannel()) << 32 | (unsigned long int) (a.GetTimeU32());
63 });
64
65 if (digisOut == digisTmp) {
66 // Copy back to digis
67 for (int i = ctx.thread_idx_x(); i < nDigis; i += ctx.block_dim_x()) {
68 digis[i] = digisTmp[i];
69 }
70 }
71}
72
74{
75 CalculateOffsetsParallel(ctx);
76 xpu::barrier(ctx.pos());
77 FindClustersParallel(ctx);
78 xpu::barrier(ctx.pos());
79 CalculateClustersParallel(ctx);
80}
81
93 unsigned int* channelOffsets, unsigned int const nDigis) const
94{
95 channelOffsets[0] = 0;
96
97 for (unsigned int pos = ctx.thread_idx_x(); pos < nDigis - 1; pos += (unsigned int) ctx.block_dim_x()) {
98 auto const currChannel = digis[pos].GetChannel();
99 auto const nextChannel = digis[pos + 1].GetChannel();
100 if (currChannel != nextChannel) {
101 for (int i = currChannel + 1; i <= nextChannel; i++) {
102 channelOffsets[i] = pos + 1;
103 }
104 }
105 }
106
107 for (int c = digis[nDigis - 1].GetChannel() + 1; c < nChannels; c++) {
108 channelOffsets[c] = nDigis;
109 }
110}
111
121{
122 int const iModule = ctx.block_idx_x();
123 CbmStsDigi* digis = &digisPerModule[digiOffsetPerModule[iModule]];
124 auto const nDigis = GetNDigis(iModule);
125
126 if (nDigis == 0) return;
127
128 auto* channelOffsets = &channelOffsetPerModule[iModule * nChannels];
129
130 CalculateChannelOffsets(ctx, digis, channelOffsets, nDigis);
131}
132
142{
143 const i32 iGThread = ctx.block_dim_x() * ctx.block_idx_x() + ctx.thread_idx_x();
144 const i32 nDigisTotal = digiOffsetPerModule[2 * nModules];
145 if (iGThread >= nDigisTotal) {
146 return;
147 };
148
149 i32 iModule = 0;
150 while (size_t(iGThread) >= digiOffsetPerModule[iModule + 1]) {
151 iModule++;
152 }
153
154 i32 iLocal = iGThread - digiOffsetPerModule[iModule];
155
156 const CbmStsDigi myDigi = digisPerModule[iGThread];
157
158 u32 nDigis = GetNDigis(iModule);
159 if (nDigis == 0) return;
160
161 xpu::view digis(&digisPerModule[digiOffsetPerModule[iModule]], nDigis);
162 xpu::view digiConnector(&digiConnectorsPerModule[digiOffsetPerModule[iModule]], nDigis);
163 xpu::view channelOffsets(&channelOffsetPerModule[iModule * nChannels], nChannels);
164
165 u16 channel = myDigi.GetChannel();
166 if (channel == nChannels - 1) return;
167
168 //DeltaCalculation
169 const f32 tResol = GetTimeResolution(iModule, channel);
170 const RecoParams::STS& params = ctx.cmem<Params>().sts;
171 const f32 deltaT = (params.timeCutDigiAbs > 0.f ? params.timeCutDigiAbs : params.timeCutDigiSig * 1.4142f * tResol);
172
173 const u32 nextChannelEnd = (channel + 2 < nChannels) ? channelOffsets[channel + 2] : nDigis;
174
175 float firstPossibleTime = float(myDigi.GetTimeU32()) - deltaT;
176
177 // Binary search for the first possible digi
178 u32 start = channelOffsets[channel + 1];
179 u32 end = nextChannelEnd;
180 while (start < end) {
181 u32 mid = (start + end) / 2;
182 if (float(digis[mid].GetTimeU32()) < firstPossibleTime) {
183 start = mid + 1;
184 }
185 else {
186 end = mid;
187 }
188 }
189
190 for (auto nextIter = start; nextIter < nextChannelEnd; nextIter++) {
191
192 const CbmStsDigi otherDigi = digis[nextIter];
193
194 // This should never happen?
195 if (myDigi.GetChannel() >= otherDigi.GetChannel()) continue;
196
197 if (float(myDigi.GetTimeU32()) + deltaT < float(otherDigi.GetTimeU32())) {
198 break;
199 }
200
201 // This check is not necessary? We already found the first possible digi via binary search...
202 if (float(myDigi.GetTimeU32()) - deltaT > float(otherDigi.GetTimeU32())) {
203 continue;
204 }
205
206 // How to handle noise? Digis in same channel within a small time window of the previous digi are discarded by old cluster finder
207
208 // How to handle if multiple digis try to connect to the same digi?
209 digiConnector[iLocal].SetNext(nextIter);
210 digiConnector[nextIter].SetHasPrevious(true);
211 break;
212 }
213}
214
224{
225 const int nModuleSides = 2 * nModules;
226
227 int iModule = 0;
228 int iThread = ctx.block_dim_x() * ctx.block_idx_x() + ctx.thread_idx_x();
229
230 for (; iModule < nModuleSides; iModule++) {
231 i32 nDigis = GetNDigis(iModule);
232 if (iThread < nDigis) {
233 break;
234 }
235 iThread -= nDigis;
236 }
237
238 if (iModule >= nModuleSides) {
239 return;
240 }
241
242 const CbmStsDigi* digis = &digisPerModule[digiOffsetPerModule[iModule]];
243 auto* digiConnector = &digiConnectorsPerModule[digiOffsetPerModule[iModule]];
244
245 // Local index of digi in sensor
246 const int iDigi = iThread;
247
248 if (digiConnector[iDigi].HasPrevious()) {
249 return;
250 }
251
252 if (!digiConnector[iDigi].HasNext()) {
253 // Cluster has 1 element
254 CreateClusterFromConnectors1(iModule, digis, iDigi);
255 }
256 else if (!digiConnector[digiConnector[iDigi].next()].HasNext()) {
257 // Cluster has 2 elements
258 CreateClusterFromConnectors2(iModule, digis, digiConnector, iDigi);
259 }
260 else {
261 // Cluster has >2 elements
262 CreateClusterFromConnectorsN(iModule, digis, digiConnector, iDigi);
263 }
264}
265
288 sts::DigiConnector* digiConnector, unsigned int const nDigis) const
289{
290 int const iModule = ctx.block_idx_x();
291 for (unsigned int currIter = ctx.thread_idx_x(); currIter < nDigis; currIter += (unsigned int) ctx.block_dim_x()) {
292
293 if (digiConnector[currIter].HasPrevious()) continue;
294 }
295}
296
297XPU_D void sts::Hitfinder::CreateClusterFromConnectors1(int const iModule, const CbmStsDigi* digis, int digiIndex) const
298{
299 const CbmStsDigi& digi = digis[digiIndex];
300
301 uint32_t time = digi.GetTimeU32();
302
303 const float timeError = asic.timeResolution;
304
305 sts::Cluster cluster{
306 .fCharge = asic.AdcToCharge(digi.GetChargeU16()),
307 .fSize = 1,
308 .fPosition = float(digi.GetChannel()) + (IsBackside(iModule) ? nChannels : 0.f),
309 .fPositionError = 1.f / xpu::sqrt(24.f),
310 .fTime = time,
311 .fTimeError = timeError,
312 };
313
314 SaveMaxError(timeError, iModule);
315
316 AddCluster(iModule, time, cluster);
317}
318
319XPU_D void sts::Hitfinder::CreateClusterFromConnectors2(int const iModule, const CbmStsDigi* digis,
320 sts::DigiConnector* digiConnector, int const digiIndex) const
321{
322
323 const CbmStsDigi& digi1 = digis[digiIndex];
324 const CbmStsDigi& digi2 = digis[digiConnector[digiIndex].next()];
325
326 float eNoiseSq = asic.noise * asic.noise;
327
328 float chargePerAdc = asic.dynamicRange / float(asic.nAdc);
329 float eDigitSq = chargePerAdc * chargePerAdc / 12.f;
330
331 float x1 = float(digi1.GetChannel());
332 float q1 = asic.AdcToCharge(digi1.GetChargeU16());
333 float q2 = asic.AdcToCharge(digi2.GetChargeU16());
334
335 // Periodic position for clusters round the edge
336 if (digi1.GetChannel() > digi2.GetChannel()) {
337 x1 -= float(nChannels);
338 }
339
340 // Uncertainties of the charge measurements
341 float width1 = LandauWidth(q1);
342 float eq1sq = width1 * width1 + eNoiseSq + eDigitSq;
343 float width2 = LandauWidth(q2);
344 float eq2sq = width2 * width2 + eNoiseSq + eDigitSq;
345
346 // Cluster time
347 float time = 0.5f * (digi1.GetTimeU32() + digi2.GetTimeU32());
348 float timeError = asic.timeResolution * 0.70710678f; // 1/sqrt(2)
349
350 // Cluster position
351 float x = x1 + 0.5f + (q2 - q1) / 3.f / xpu::max(q1, q2);
352
353 // Correct negative position for clusters around the edge
354 if (x < -0.5f) {
355 x += float(nChannels);
356 }
357
358 // Uncertainty on cluster position. See software note.
359 float ex0sq = 0.f;
360 float ex1sq = 0.f;
361 float ex2sq = 0.f;
362 if (q1 < q2) {
363 ex0sq = (q2 - q1) * (q2 - q1) / q2 / q2 / 72.f;
364 ex1sq = eq1sq / q2 / q2 / 9.f;
365 ex2sq = eq2sq * q1 * q1 / q2 / q2 / q2 / q2 / 9.f;
366 }
367 else {
368 ex0sq = (q2 - q1) * (q2 - q1) / q1 / q1 / 72.f;
369 ex1sq = eq1sq * q2 * q2 / q1 / q1 / q1 / q1 / 9.f;
370 ex2sq = eq2sq / q1 / q1 / 9.f;
371 }
372 float xError = xpu::sqrt(ex0sq + ex1sq + ex2sq);
373
374 // Cluster charge
375 float charge = q1 + q2;
376
377 if (IsBackside(iModule)) x += nChannels;
378
379 sts::Cluster cls{
380 .fCharge = charge,
381 .fSize = 2,
382 .fPosition = x,
383 .fPositionError = xError,
384 .fTime = uint32_t(time),
385 .fTimeError = timeError,
386 };
387
388 SaveMaxError(timeError, iModule);
389 AddCluster(iModule, time, cls);
390}
391
392XPU_D void sts::Hitfinder::CreateClusterFromConnectorsN(int iModule, const CbmStsDigi* digis,
393 sts::DigiConnector* digiConnector, int digiIndex) const
394{
396
397 //This Lambda calculates all charges for a single digi inside of a cluster
398 auto calculateClusterCharges = [this, &cProps, &digis, &digiConnector](int index) {
399 CbmStsDigi digi = digis[index];
400 float eNoiseSq = asic.noise * asic.noise;
401 float chargePerAdc = asic.dynamicRange / float(asic.nAdc);
402 float eDigitSq = chargePerAdc * chargePerAdc / 12.f;
403 cProps.tResolSum += asic.timeResolution;
404 cProps.tSum += digi.GetTimeU32();
405 float charge = asic.AdcToCharge(digi.GetChargeU16());
406 float lWidth = LandauWidth(charge);
407 float eChargeSq = lWidth * lWidth + eNoiseSq + eDigitSq;
408
409 if (!digiConnector[index].HasPrevious()) {
410 //begin of Cluster (most left Element)
411 cProps.chanF = digi.GetChannel();
412 cProps.qF = charge;
413 cProps.eqFsq = eChargeSq;
414 }
415 else if (!digiConnector[index].HasNext()) {
416 //end of Cluster (most right Element)
417 cProps.chanL = digi.GetChannel();
418 cProps.qL = charge;
419 cProps.eqLsq = eChargeSq;
420 }
421 else {
422 cProps.qM += charge;
423 cProps.eqMsq += eChargeSq;
424 }
425 cProps.xSum += charge * float(digi.GetChannel());
426 };
427
428 calculateClusterCharges(digiIndex);
429
430 // Calculate ClusterCharges
431 while (digiConnector[digiIndex].HasNext()) {
432 digiIndex = digiConnector[digiIndex].next();
433 calculateClusterCharges(digiIndex);
434 }
435
436 // Periodic channel position for clusters round the edge
437 if (cProps.chanF > cProps.chanL) cProps.chanF -= nChannels;
438
439 // Cluster time and total charge
440 float nDigis = ChanDist(cProps.chanF, cProps.chanL) + 1;
441 cProps.tSum = cProps.tSum / nDigis;
442 float timeError = (cProps.tResolSum / nDigis) / xpu::sqrt(nDigis);
443 float qSum = cProps.qF + cProps.qM + cProps.qL;
444
445 // Average charge in middle strips
446 cProps.qM /= (nDigis - 2.f);
447 cProps.eqMsq /= (nDigis - 2.f);
448
449 // Cluster position
450 float x = 0.5f * (float(cProps.chanF + cProps.chanL) + (cProps.qL - cProps.qF) / cProps.qM);
451
452 // Correct negative cluster position for clusters round the edge
453 if (x < -0.5f) {
454 x += float(nChannels);
455 }
456
457 // Cluster position error
458 float exFsq = cProps.eqFsq / cProps.qM / cProps.qM / 4.f; // error from first charge
459 float exMsq = cProps.eqMsq * (cProps.qL - cProps.qF) * (cProps.qL - cProps.qF) / cProps.qM / cProps.qM / cProps.qM
460 / cProps.qM / 4.f;
461 float exLsq = cProps.eqLsq / cProps.qM / cProps.qM / 4.f;
462 float xError = xpu::sqrt(exFsq + exMsq + exLsq);
463
464 // Correction for corrupt clusters
465 if (x < cProps.chanF || x > cProps.chanL) {
466 x = cProps.xSum / qSum;
467 }
468
469 if (IsBackside(iModule)) {
470 x += nChannels;
471 }
472
473 sts::Cluster cls{
474 .fCharge = qSum,
475 .fSize = int(nDigis),
476 .fPosition = x,
477 .fPositionError = xError,
478 .fTime = uint32_t(cProps.tSum),
479 .fTimeError = timeError,
480 };
481
482 SaveMaxError(timeError, iModule);
483 AddCluster(iModule, cProps.tSum, cls);
484}
485
487{
488 // By default sort all modules in parallel,
489 // but in debug mode sort only one module at a time
490 // to narrow down the error source
491 int firstModule = ctx.block_idx_x();
492 int lastModule = ctx.block_idx_x();
493 if (ctx.grid_dim_x() == 1) {
494 firstModule = 0;
495 lastModule = 2 * nModules - 1;
496 }
497
498 for (int iModule = firstModule; iModule <= lastModule; iModule++) {
499 size_t offset = iModule * maxClustersPerModule;
500 ClusterIdx* clusterIdx = &clusterIdxPerModule[offset];
501 ClusterIdx* clusterIdxTmp = &clusterIdxPerModuleTmp[offset];
502 int nClusters = nClustersPerModule[iModule];
503
504 // if (nClusters == 0 && ctx.thread_idx_x() == 0) {
505 // printf("Module %d: No clusters found, exit early\n", iModule);
506 // }
507 if (nClusters == 0) return;
508
509 // if (ctx.thread_idx_x() == 0) {
510 // printf("Module %d: Start sorting %d clusters\n", iModule, nClusters);
511 // }
512
513 SortClusters::sort_t clsSort(ctx.pos(), ctx.smem());
514 clusterIdx = clsSort.sort(clusterIdx, nClusters, clusterIdxTmp, [](const ClusterIdx a) { return a.fTime; });
515
516 // if (ctx.thread_idx_x() == 0) {
517 // printf("Module %d: Finished sorting %d clusters\n", iModule, nClusters);
518 // }
519
520 if (ctx.thread_idx_x() == 0) clusterIdxSortedPerModule[iModule] = clusterIdx;
521 }
522}
523
525{
526 int iModule = 0;
527 int iThread = ctx.block_dim_x() * ctx.block_idx_x() + ctx.thread_idx_x();
528
529#if XPU_IS_CPU
530 iModule = iThread / kFindHitsChunksPerModule;
531#else
532 for (; iModule < nModules; iModule++) {
533 if (iThread < nClustersPerModule[iModule]) break;
534 iThread -= nClustersPerModule[iModule];
535 }
536#endif
537
538 if (iModule >= nModules) {
539 return;
540 }
541
542 int iModuleF = iModule;
543 int iModuleB = iModuleF + nModules;
544 size_t offsetF = iModuleF * maxClustersPerModule;
545 size_t offsetB = iModuleB * maxClustersPerModule;
546 ClusterIdx* clusterIdxF = clusterIdxSortedPerModule[iModuleF];
547 ClusterIdx* clusterIdxB = clusterIdxSortedPerModule[iModuleB];
548 sts::Cluster* clusterDataF = &clusterDataPerModule[offsetF];
549 sts::Cluster* clusterDataB = &clusterDataPerModule[offsetB];
550 int nClustersF = nClustersPerModule[iModuleF];
551 int nClustersB = nClustersPerModule[iModuleB];
552 size_t nHitsWritten = nHitsPerModule[iModule];
553
554 // global parameters
555 const RecoParams::STS& params = ctx.cmem<Params>().sts;
556
557 const bool doChargeCorrelation = params.doChargeCorrelation;
558 const real ChargeDelta = params.chargeCorrelationDelta; // Allowed charge difference between clusters to create a hit
559
560
561 if (nClustersF == 0 || nClustersB == 0) {
562 return;
563 }
564
565 // Only used on CPU
566 [[maybe_unused]] int nClustersPerChunk = (nClustersF + kFindHitsChunksPerModule - 1) / kFindHitsChunksPerModule;
567 [[maybe_unused]] int iChunk = iThread % kFindHitsChunksPerModule;
568
569 int iClusterF =
570#if XPU_IS_CPU
571 iChunk * nClustersPerChunk;
572#else
573 iThread;
574#endif
575 if (iClusterF >= nClustersF) {
576 return;
577 }
578
579
580 // Stop processing if memory limits are exceed by a large amount
581 // In general we would still want to process all clusters to get an idea
582 // by how much memory estimates are wrong.
583 // However they are currently chosen very generous, so if they are exceeded by a large amount
584 // something must have gone awry. (e.g. monster events in the mCBM data that explode the hit finding combinatorics)
585 if (nHitsWritten > 2 * maxHitsPerModule) {
586 return;
587 }
588
589 HitfinderCache pars;
590 {
591 SensorPar cfg = sensorPars[iModule];
592
593 pars.dY = cfg.dY;
594 pars.pitch = cfg.pitch;
595 pars.stereoF = cfg.stereoF;
596 pars.stereoB = cfg.stereoB;
597 pars.lorentzF = cfg.lorentzF;
598 pars.lorentzB = cfg.lorentzB;
599 pars.tanStereoF = xpu::tan(pars.stereoF * xpu::deg_to_rad());
600 pars.tanStereoB = xpu::tan(pars.stereoB * xpu::deg_to_rad());
601 pars.errorFac = 1.f / (pars.tanStereoB - pars.tanStereoF) / (pars.tanStereoB - pars.tanStereoF);
602 pars.dX = float(nChannels) * pars.pitch;
603 }
604
605 float maxTerrF = maxClusterTimeErrorByModuleSide[iModuleF];
606 float maxTerrB = maxClusterTimeErrorByModuleSide[iModuleB];
607
608 float maxSigmaBoth = 4.f * xpu::sqrt(maxTerrF * maxTerrF + maxTerrB * maxTerrB);
609
610 int startB = 0;
611 float firstPossibleTime = float(clusterIdxF[iClusterF].fTime) - maxSigmaBoth;
612
613 // Use binary search to find the first cluster on back side that can be matched
614 // with the current cluster on front side
615 int start = startB;
616 int end = nClustersB;
617 while (end - start > 1) {
618 int mid = (start + end) / 2;
619 if (float(clusterIdxB[mid].fTime) < firstPossibleTime) {
620 start = mid;
621 }
622 else {
623 end = mid;
624 }
625 }
626 startB = start;
627
628 int endClusterF =
629#if XPU_IS_CPU
630 xpu::min(iClusterF + nClustersPerChunk, nClustersF);
631#else
632 iClusterF + 1;
633#endif
634
635 for (; iClusterF < endClusterF; iClusterF++) {
636 ClusterIdx clsIdxF = clusterIdxF[iClusterF];
637 sts::Cluster clsDataF = clusterDataF[clsIdxF.fIdx];
638
639 float maxSigma = 4.f * xpu::sqrt(clsDataF.fTimeError * clsDataF.fTimeError + maxTerrB * maxTerrB);
640
641 for (int iClusterB = startB; iClusterB < nClustersB; iClusterB++) {
642 ClusterIdx clsIdxB = clusterIdxB[iClusterB];
643 sts::Cluster clsDataB = clusterDataB[clsIdxB.fIdx];
644
645 float timeDiff = float(clsIdxF.fTime) - float(clsIdxB.fTime);
646
647 if (timeDiff > 0 && timeDiff > maxSigmaBoth) {
648 startB++;
649 continue;
650 }
651 else if (timeDiff > 0 && timeDiff > maxSigma) {
652 continue;
653 }
654 else if (timeDiff < 0 && xpu::abs(timeDiff) > maxSigma) {
655 break;
656 }
657
658 // Critical to check for charge correlation AFTER time difference
659 // as time difference advances cluster index on back side
660 // reducing the combinatorics drastically
661 if (doChargeCorrelation) {
662 float chargeDiff = clsDataF.fCharge - clsDataB.fCharge;
663 if (xpu::abs(chargeDiff) > ChargeDelta) {
664 continue;
665 }
666 }
667
668 float timeCut = -1.f;
669 if (params.timeCutClusterAbs > 0.f)
670 timeCut = params.timeCutClusterAbs;
671 else if (params.timeCutClusterSig > 0.f) {
672 float eF = clsDataF.fTimeError;
673 float eB = clsDataB.fTimeError;
674 timeCut = params.timeCutClusterSig * xpu::sqrt(eF * eF + eB * eB);
675 }
676
677 if (xpu::abs(float(clsIdxF.fTime) - float(clsIdxB.fTime)) > timeCut) {
678 continue;
679 }
680
681 IntersectClusters(iModule, pars, clsIdxF, clsDataF, clsIdxB, clsDataB);
682 } // end loop over clusters on back side
683 } // end loop over clusters on front side
684 // clang-format on
685}
686
687XPU_D void sts::Hitfinder::IntersectClusters(int iBlock, const HitfinderCache& pars, const ClusterIdx& idxF,
688 const sts::Cluster& clsF, const ClusterIdx& idxB,
689 const sts::Cluster& clsB) const
690{
691 // --- Calculate cluster centre position on readout edge
692
693 float xF = GetClusterPosition(pars, clsF.fPosition, true);
694 float exF = clsF.fPositionError * pars.pitch;
695 float du = exF * xpu::cos(xpu::deg_to_rad() * pars.stereoF);
696 float xB = GetClusterPosition(pars, clsB.fPosition, false);
697 float exB = clsB.fPositionError * pars.pitch;
698 float dv = exB * xpu::cos(xpu::deg_to_rad() * pars.stereoB);
699
700 // // --- Should be inside active area
701 if (xF < 0.f || xF > pars.dX || xB < 0.f || xB > pars.dX) {
702 return;
703 }
704
705 // // --- Calculate number of line segments due to horizontal
706 // // --- cross-connection. If x(y=0) does not fall on the bottom edge,
707 // // --- the strip is connected to the one corresponding to the line
708 // // --- with top edge coordinate xF' = xF +/- Dx. For odd combinations
709 // // --- of stereo angle and sensor dimensions, this could even happen
710 // // --- multiple times. For each of these lines, the intersection with
711 // // --- the line on the other side is calculated. If inside the active area,
712 // // --- a hit is created.
713 int nF = (xF + pars.dY * pars.tanStereoF) / pars.dX;
714 int nB = (xB + pars.dY * pars.tanStereoB) / pars.dX;
715
716 // --- If n is positive, all lines from 0 to n must be considered,
717 // --- if it is negative (phi negative), all lines from n to 0.
718 int nF1 = xpu::min(0, nF);
719 int nF2 = xpu::max(0, nF);
720 int nB1 = xpu::min(0, nB);
721 int nB2 = xpu::max(0, nB);
722
723 // --- Double loop over possible lines
724 float xC = -1.f;
725 float yC = -1.f;
726 float varX = 0.f;
727 float varY = 0.f;
728 float varXY = 0.f;
729 for (int iF = nF1; iF <= nF2; iF++) {
730 float xFi = xF - float(iF) * pars.dX;
731 for (int iB = nB1; iB <= nB2; iB++) {
732 float xBi = xB - float(iB) * pars.dX;
733
734 // --- Intersect the two lines
735 bool found = Intersect(pars, xFi, exF, xBi, exB, xC, yC, varX, varY, varXY);
736 if (not found) continue;
737
738 // --- Transform into sensor system with origin at sensor centre
739
740 xC -= 0.5f * pars.dX;
741 yC -= 0.5f * pars.dY;
742
743 CreateHit(iBlock, xC, yC, varX, varY, varXY, idxF, clsF, idxB, clsB, du, dv);
744 }
745 }
746}
747
748XPU_D float sts::Hitfinder::GetClusterPosition(const HitfinderCache& pars, float centre, bool isFront) const
749{
750 // Take integer channel
751
752 int iChannel = int(centre);
753 float xDif = centre - float(iChannel);
754
755 // Calculate corresponding strip on sensor
756
757 int iStrip = iChannel - (isFront ? 0 : nChannels);
758
759 // Re-add difference to integer channel. Convert channel to
760 // coordinate
761 float xCluster = (float(iStrip) + xDif + 0.5f) * pars.pitch;
762
763 // // Correct for Lorentz-Shift
764 // // Simplification: The correction uses only the y component of the
765 // // magnetic field. The shift is calculated using the mid-plane of the
766 // // sensor, which is not correct for tracks not traversing the entire
767 // // sensor thickness (i.e., are created or stopped somewhere in the sensor).
768 // // However, this is the best one can do in reconstruction.
769 xCluster -= (isFront ? pars.lorentzF : pars.lorentzB);
770
771 return xCluster;
772}
773
774XPU_D bool sts::Hitfinder::Intersect(const HitfinderCache& pars, float xF, float exF, float xB, float exB, float& x,
775 float& y, float& varX, float& varY, float& varXY) const
776{
777
778 // In the coordinate system with origin at the bottom left corner,
779 // a line along the strips with coordinate x0 at the top edge is
780 // given by the function y(x) = Dy - ( x - x0 ) / tan(phi), if
781 // the stereo angle phi does not vanish. Two lines yF(x), yB(x) with top
782 // edge coordinates xF, xB intersect at
783 // x = ( tan(phiB)*xF - tan(phiF)*xB ) / (tan(phiB) - tan(phiF)
784 // y = Dy + ( xB - xF ) / ( tan(phiB) - tan(phiF) )
785 // For the case that one of the stereo angles vanish (vertical strips),
786 // the calculation of the intersection is straightforward.
787
788 // --- First check whether stereo angles are different. Else there is
789 // --- no intersection.
790 // TODO: if this is true, no hits can be found? So just do this once at the beginning?
791 if (xpu::abs(pars.stereoF - pars.stereoB) < 0.5f) return false;
792
793 // --- Now treat vertical front strips
794 if (xpu::abs(pars.stereoF) < 0.001f) {
795 x = xF;
796 y = pars.dY - (xF - xB) / pars.tanStereoB;
797 varX = exF * exF;
798 varY = (exF * exF + exB * exB) / pars.tanStereoB / pars.tanStereoB;
799 varXY = -1.f * exF * exF / pars.tanStereoB;
800 return IsInside(pars, x - pars.dX / 2.f, y - pars.dY / 2.f);
801 }
802
803 // --- Maybe the back side has vertical strips?
804 if (xpu::abs(pars.stereoB) < 0.001f) {
805 x = xB;
806 y = pars.dY - (xB - xF) / pars.tanStereoF;
807 varX = exB * exB;
808 varY = (exF * exF + exB * exB) / pars.tanStereoF / pars.tanStereoF;
809 varXY = -1.f * exB * exB / pars.tanStereoF;
810 return IsInside(pars, x - pars.dX / 2.f, y - pars.dY / 2.f);
811 }
812
813 // --- OK, both sides have stereo angle
814
815 x = (pars.tanStereoB * xF - pars.tanStereoF * xB) / (pars.tanStereoB - pars.tanStereoF);
816 y = pars.dY + (xB - xF) / (pars.tanStereoB - pars.tanStereoF);
817 varX =
818 pars.errorFac * (exF * exF * pars.tanStereoB * pars.tanStereoB + exB * exB * pars.tanStereoF * pars.tanStereoF);
819 varY = pars.errorFac * (exF * exF + exB * exB);
820 varXY = -1.f * pars.errorFac * (exF * exF * pars.tanStereoB + exB * exB * pars.tanStereoF);
821
822 return IsInside(pars, x - pars.dX / 2.f, y - pars.dY / 2.f);
823}
824
825XPU_D bool sts::Hitfinder::IsInside(const HitfinderCache& pars, float x, float y) const
826{
827 // clang-format off
828 return x >= -pars.dX / 2.f
829 && x <= pars.dX / 2.f
830 && y >= -pars.dY / 2.f
831 && y <= pars.dY / 2.f;
832 // clang-format on
833}
834
835XPU_D void sts::Hitfinder::CreateHit(int iModule, float xLocal, float yLocal, float varX, float varY, float varXY,
836 const ClusterIdx& idxF, const Cluster& clsF, const ClusterIdx& idxB,
837 const Cluster& clsB, float du, float dv) const
838{
839 // --- Transform into global coordinate system
840 float globalX, globalY, globalZ;
841 ToGlobal(iModule, xLocal, yLocal, 0.f, globalX, globalY, globalZ);
842
843 // We assume here that the local-to-global transformations is only translation
844 // plus maybe rotation upside down or front-side back. In that case, the
845 // global covariance matrix is the same as the local one.
846 float errX = xpu::sqrt(varX);
847 float errY = xpu::sqrt(varY);
848
849 // --- Calculate hit time (average of cluster times)
850 float hitTime = 0.5f * (float(idxF.fTime) + float(idxB.fTime));
851 float etF = clsF.fTimeError;
852 float etB = clsB.fTimeError;
853 float hitTimeError = 0.5f * xpu::sqrt(etF * etF + etB * etB);
854
855 sts::Hit hit{
856 .fX = globalX,
857 .fY = globalY,
858 .fZ = globalZ,
859 .fTime = static_cast<u32>(hitTime),
860 .fDx = errX,
861 .fDy = errY,
862 .fDxy = varXY,
863 .fTimeError = hitTimeError,
864 .fDu = du,
865 .fDv = dv,
866
867 .fFrontClusterId = idxF.fIdx,
868 .fBackClusterId = idxB.fIdx,
869 };
870
871 int idx = xpu::atomic_add(&nHitsPerModule[iModule], 1);
872
873 if (size_t(idx) >= maxHitsPerModule) {
874 xpu::atomic_add(&monitor->nHitBucketOverflow, 1);
875 return;
876 }
877
878 hitsPerModule[iModule * hitsAllocatedPerModule + idx] = hit;
879}
880
881XPU_D float sts::Hitfinder::LandauWidth(float charge) const
882{
883 if (charge <= landauStepSize) return landauTable[0];
884 if (charge > landauStepSize * (landauTableSize - 1)) return landauTable[landauTableSize - 1];
885
886 int tableIdx = xpu::ceil(charge / landauStepSize);
887 float e2 = tableIdx * landauStepSize;
888 float v2 = landauTable[tableIdx];
889 tableIdx--;
890 float e1 = tableIdx * landauStepSize;
891 float v1 = landauTable[tableIdx];
892 return v1 + (charge - e1) * (v2 - v1) / (e2 - e1);
893}
894
895XPU_D void sts::Hitfinder::ToGlobal(int iModule, float lx, float ly, float lz, float& gx, float& gy, float& gz) const
896{
897 const float* tr = &localToGlobalTranslationByModule[iModule * 3];
898 const float* rot = &localToGlobalRotationByModule[iModule * 9];
899
900 gx = tr[0] + lx * rot[0] + ly * rot[1] + lz * rot[2];
901 gy = tr[1] + lx * rot[3] + ly * rot[4] + lz * rot[5];
902 gz = tr[2] + lx * rot[6] + ly * rot[7] + lz * rot[8];
903}
XPU_EXPORT(sts::TheHitfinder)
XPU_EXPORT(cbm::algo::Params)
Data class for a single-channel message in the STS.
Definition CbmStsDigi.h:40
XPU_D uint16_t GetChannel() const
Channel number in module @value Channel number.
Definition CbmStsDigi.h:93
XPU_D uint16_t GetChargeU16() const
Definition CbmStsDigi.h:108
XPU_D uint32_t GetTimeU32() const
Definition CbmStsDigi.h:123
XPU_D unsigned int next() const
Definition Hitfinder.h:179
XPU_D void CreateClusterFromConnectors1(int const iModule, const CbmStsDigi *digis, int const digiIndex) const
XPU_D void CalculateChannelOffsets(FindClusters::context &ctx, CbmStsDigi *digis, unsigned int *channelOffsets, unsigned int nDigis) const
Definition Hitfinder.cxx:92
XPU_D void SortClusters(SortClusters::context &) const
XPU_D void CalculateOffsetsParallel(FindClusters::context &) const
XPU_D bool Intersect(const HitfinderCache &pars, float xF, float exF, float xB, float exB, float &x, float &y, float &varX, float &varY, float &varXY) const
XPU_D void FindClustersSingleStep(FindClusters::context &) const
Definition Hitfinder.cxx:73
XPU_D bool IsInside(const HitfinderCache &pars, float x, float y) const
XPU_D float LandauWidth(float charge) const
XPU_D void CreateClusterFromConnectorsN(int const iModule, const CbmStsDigi *digis, DigiConnector *digiConnector, int const digiIndex) const
XPU_D void CalculateClustersDigiWise(FindClusters::context &ctx, CbmStsDigi *digis, DigiConnector *digiConnector, unsigned int const nDigis) const
XPU_D void IntersectClusters(int iBlock, const HitfinderCache &pars, const ClusterIdx &idxF, const sts::Cluster &clsF, const ClusterIdx &idxB, const sts::Cluster &clsB) const
XPU_D void SortDigisInSpaceAndTime(SortDigis::context &) const
Definition Hitfinder.cxx:52
XPU_D void CalculateClustersParallel(FindClusters::context &) const
XPU_D void CreateClusterFromConnectors2(int const iModule, const CbmStsDigi *digis, DigiConnector *digiConnector, int const digiIndex) const
XPU_D void ToGlobal(int iModule, float lx, float ly, float lz, float &gx, float &gy, float &gz) const
XPU_D void CreateHit(int iBlocks, float xLocal, float yLocal, float varX, float varY, float varXY, const ClusterIdx &idxF, const Cluster &clsF, const ClusterIdx &idxB, const sts::Cluster &clsB, float du, float dv) const
XPU_D void FindHits(FindHits::context &) const
XPU_D void FindClustersParallel(FindClusters::context &) const
XPU_D float GetClusterPosition(const HitfinderCache &pars, float centre, bool isFront) const
std::int32_t i32
Definition Definitions.h:20
float f32
Definition Definitions.h:24
std::uint32_t u32
Definition Definitions.h:21
@ kFindHitsChunksPerModule
Definition Hitfinder.h:40
std::uint16_t u16
Definition Definitions.h:19
XPU_D void operator()(context &)
Definition Hitfinder.cxx:21
xpu::kernel_context< shared_memory, constants > context
Definition Hitfinder.h:81
u32 fTime
Cluster time (average of digi times) [ns].
Definition Hitfinder.h:53
real fTimeError
Error of cluster time [ns].
real fPosition
Cluster centre in channel number units.
real fPositionError
Cluster centre error (r.m.s.) in channel number units.
real fCharge
Total charge.
xpu::kernel_context< shared_memory, constants > context
Definition Hitfinder.h:95
XPU_D void operator()(context &)
Definition Hitfinder.cxx:33
xpu::kernel_context< shared_memory, constants > context
Definition Hitfinder.h:88
xpu::kernel_context< shared_memory, constants > context
Definition Hitfinder.h:74
XPU_D void operator()(context &)
Definition Hitfinder.cxx:18
XPU_D void operator()(context &)
Definition Hitfinder.cxx:42
xpu::kernel_context< shared_memory, constants > context
Definition Hitfinder.h:111
XPU_D void operator()(context &)
Definition Hitfinder.cxx:39
xpu::block_sort< u32, ClusterIdx, kSortClustersBlockSize, kSortClustersItemsPerThread > sort_t
Definition Hitfinder.h:101
xpu::kernel_context< shared_memory, constants > context
Definition Hitfinder.h:104
XPU_D void operator()(context &)
Definition Hitfinder.cxx:15
xpu::block_sort< u64, CbmStsDigi, kSortDigisBlockSize, kSortDigisItemsPerThread > sort_t
Definition Hitfinder.h:64
xpu::kernel_context< shared_memory, constants > context
Definition Hitfinder.h:67