2 #include <alpaka/alpaka.hpp> 17 #ifdef HCAL_MAHI_GPUDEBUG 18 #define DETID_TO_DEBUG 1125647428 28 #if defined(__CUDA_ARCH__) or defined(__HIP_DEVICE_COMPILE__) 29 return __uint_as_float(
val);
31 return edm::bit_cast<
float>(
val);
36 #if defined(__CUDA_ARCH__) or defined(__HIP_DEVICE_COMPILE__) 37 return __float_as_uint(
val);
39 return edm::bit_cast<
unsigned int>(
val);
48 return rawDelay < 0 ? 0 : (rawDelay >
tmax ?
tmax : rawDelay);
53 ALPAKA_STATIC_ACC_MEM_CONSTANT
float const qie8shape[129] = {
54 -1, 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 16,
55 18, 20, 22, 24, 26, 28, 31, 34, 37, 40, 44, 48, 52, 57, 62, 57, 62,
56 67, 72, 77, 82, 87, 92, 97, 102, 107, 112, 117, 122, 127, 132, 142, 152, 162,
57 172, 182, 192, 202, 217, 232, 247, 262, 282, 302, 322, 347, 372, 347, 372, 397, 422,
58 447, 472, 497, 522, 547, 572, 597, 622, 647, 672, 697, 722, 772, 822, 872, 922, 972,
59 1022, 1072, 1147, 1222, 1297, 1372, 1472, 1572, 1672, 1797, 1922, 1797, 1922, 2047, 2172, 2297, 2422,
60 2547, 2672, 2797, 2922, 3047, 3172, 3297, 3422, 3547, 3672, 3922, 4172, 4422, 4672, 4922, 5172, 5422,
61 5797, 6172, 6547, 6922, 7422, 7922, 8422, 9047, 9672, 10297};
63 ALPAKA_STATIC_ACC_MEM_CONSTANT
float const qie11shape[257] = {
64 -0.5, 0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5,
65 11.5, 12.5, 13.5, 14.5, 15.5, 17.5, 19.5, 21.5, 23.5, 25.5, 27.5, 29.5,
66 31.5, 33.5, 35.5, 37.5, 39.5, 41.5, 43.5, 45.5, 47.5, 49.5, 51.5, 53.5,
67 55.5, 59.5, 63.5, 67.5, 71.5, 75.5, 79.5, 83.5, 87.5, 91.5, 95.5, 99.5,
68 103.5, 107.5, 111.5, 115.5, 119.5, 123.5, 127.5, 131.5, 135.5, 139.5, 147.5, 155.5,
69 163.5, 171.5, 179.5, 187.5, 171.5, 179.5, 187.5, 195.5, 203.5, 211.5, 219.5, 227.5,
70 235.5, 243.5, 251.5, 259.5, 267.5, 275.5, 283.5, 291.5, 299.5, 315.5, 331.5, 347.5,
71 363.5, 379.5, 395.5, 411.5, 427.5, 443.5, 459.5, 475.5, 491.5, 507.5, 523.5, 539.5,
72 555.5, 571.5, 587.5, 603.5, 619.5, 651.5, 683.5, 715.5, 747.5, 779.5, 811.5, 843.5,
73 875.5, 907.5, 939.5, 971.5, 1003.5, 1035.5, 1067.5, 1099.5, 1131.5, 1163.5, 1195.5, 1227.5,
74 1259.5, 1291.5, 1355.5, 1419.5, 1483.5, 1547.5, 1611.5, 1675.5, 1547.5, 1611.5, 1675.5, 1739.5,
75 1803.5, 1867.5, 1931.5, 1995.5, 2059.5, 2123.5, 2187.5, 2251.5, 2315.5, 2379.5, 2443.5, 2507.5,
76 2571.5, 2699.5, 2827.5, 2955.5, 3083.5, 3211.5, 3339.5, 3467.5, 3595.5, 3723.5, 3851.5, 3979.5,
77 4107.5, 4235.5, 4363.5, 4491.5, 4619.5, 4747.5, 4875.5, 5003.5, 5131.5, 5387.5, 5643.5, 5899.5,
78 6155.5, 6411.5, 6667.5, 6923.5, 7179.5, 7435.5, 7691.5, 7947.5, 8203.5, 8459.5, 8715.5, 8971.5,
79 9227.5, 9483.5, 9739.5, 9995.5, 10251.5, 10507.5, 11019.5, 11531.5, 12043.5, 12555.5, 13067.5, 13579.5,
80 12555.5, 13067.5, 13579.5, 14091.5, 14603.5, 15115.5, 15627.5, 16139.5, 16651.5, 17163.5, 17675.5, 18187.5,
81 18699.5, 19211.5, 19723.5, 20235.5, 20747.5, 21771.5, 22795.5, 23819.5, 24843.5, 25867.5, 26891.5, 27915.5,
82 28939.5, 29963.5, 30987.5, 32011.5, 33035.5, 34059.5, 35083.5, 36107.5, 37131.5, 38155.5, 39179.5, 40203.5,
83 41227.5, 43275.5, 45323.5, 47371.5, 49419.5, 51467.5, 53515.5, 55563.5, 57611.5, 59659.5, 61707.5, 63755.5,
84 65803.5, 67851.5, 69899.5, 71947.5, 73995.5, 76043.5, 78091.5, 80139.5, 82187.5, 84235.5, 88331.5, 92427.5,
85 96523.5, 100620, 104716, 108812, 112908};
90 uint32_t
const didraw,
int const maxDepthHB,
int const firstHBRing,
int const lastHBRing,
int const nEtaHB) {
92 uint32_t
const value = (did.depth() - 1) +
maxDepthHB * (did.iphi() - 1);
100 int const firstHERing,
101 int const lastHERing,
104 uint32_t
const value = (did.depth() - 1) +
maxDepthHE * (did.iphi() - 1);
105 return did.ieta() > 0 ?
value +
maxDepthHE * maxPhiHE * (did.ieta() - firstHERing)
109 return capid * 4 +
range;
121 int const qieType, uint8_t
const adc, uint8_t
const capid,
float const* qieOffsets,
float const* qieSlopes) {
122 auto const range = qieType == 0 ? (
adc >> 5) & 0x3 : (
adc >> 6) & 0x3;
124 auto const nbins = qieType == 0 ? 32 : 64;
125 auto const center =
adc %
nbins ==
nbins - 1 ? 0.5 * (3 * qieShapeToUse[
adc] - qieShapeToUse[
adc - 1])
126 : 0.5 * (qieShapeToUse[
adc] + qieShapeToUse[
adc + 1]);
128 return (center - qieOffsets[
index]) / qieSlopes[
index];
137 float const* qieOffsets,
138 float const* qieSlopes,
139 bool const isqie11) {
140 constexpr uint32_t mantissaMaskQIE8 = 0x1fu;
141 constexpr uint32_t mantissaMaskQIE11 = 0x3f;
142 auto const mantissaMask = isqie11 ? mantissaMaskQIE11 : mantissaMaskQIE8;
144 auto const mantissa =
adc & mantissaMask;
146 if (mantissa == 0u || mantissa == mantissaMask - 1u)
148 else if (mantissa == 1u || mantissa == mantissaMask)
153 auto const upgain = qup -
q;
154 auto const downgain =
q - qdown;
155 auto const averagegain = (qup - qdown) / 2.
f;
156 if (
std::abs(upgain - downgain) < 0.01f * averagegain)
161 auto const upgain2 = q2up - qup;
162 auto const downgain2 = qdown - q2down;
175 float const pulse_time,
178 auto const& acc25nsVec = pulseShape.acc25nsVec();
179 auto const& diff25nsItvlVec = pulseShape.diff25nsItvlVec();
180 auto const& accVarLenIdxMinusOneVec = pulseShape.accVarLenIdxMinusOneVec();
181 auto const& diffVarItvlIdxMinusOneVec = pulseShape.diffVarItvlIdxMinusOneVec();
182 auto const& accVarLenIdxZeroVec = pulseShape.accVarLenIdxZEROVec();
183 auto const& diffVarItvlIdxZeroVec = pulseShape.diffVarItvlIdxZEROVec();
193 int i_start =
static_cast<int>(i_start_float);
197 if (offset_start == 1.0
f) {
202 int const bin_start =
static_cast<int>(offset_start);
203 float const bin_start_up =
static_cast<float>(bin_start) + 0.5
f;
204 int const bin_0_start = offset_start < bin_start_up ? bin_start - 1 : bin_start;
205 int const its_start = i_start / ns_per_bx;
206 int const distTo25ns_start = ns_per_bx - 1 - i_start % ns_per_bx;
207 auto const factor = offset_start -
static_cast<float>(bin_0_start) - 0.5;
211 if (sample_over10ts == its_start) {
212 value = bin_0_start == -1
213 ? accVarLenIdxMinusOneVec[distTo25ns_start] +
factor * diffVarItvlIdxMinusOneVec[distTo25ns_start]
214 : accVarLenIdxZeroVec[distTo25ns_start] +
factor * diffVarItvlIdxZeroVec[distTo25ns_start];
215 }
else if (sample_over10ts > its_start) {
216 int const bin_idx = distTo25ns_start + 1 + (sample_over10ts - its_start - 1) * ns_per_bx + bin_0_start;
217 value = acc25nsVec[bin_idx] +
factor * diff25nsItvlVec[bin_idx];
232 float const* shrChargeMinusPedestal,
236 int32_t
const nsamplesForCompute,
242 bool const isqie11) {
252 sipmq += shrChargeMinusPedestal[lch * nsamplesForCompute + ts];
253 auto const effectivePixelsFired = sipmq / fcByPE;
257 #ifdef HCAL_MAHI_GPUDEBUG 258 printf(
"first = %d last = %d sipmQ = %f factor = %f rawCharge = %f\n",
first,
last, sipmq,
factor, rawCharge);
266 return a.second >
b.second;
271 template <
typename TAcc,
typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
274 IProductTypef01::ConstView f01HEDigis,
275 IProductTypef5::ConstView f5HBDigis,
276 IProductTypef3::ConstView f3HBDigis,
277 HcalMahiConditionsPortableDevice::ConstView mahi,
278 HcalSiPMCharacteristicsPortableDevice::ConstView sipmCharacteristics,
287 float* electronicNoiseTerms,
289 int const windowSize)
const {
290 auto const nchannelsf015 = f01HEDigis.size() + f5HBDigis.size();
291 auto const startingSample = compute_nsamples<Flavor1>(f01HEDigis.stride()) - windowSize;
292 auto const nchannels = f01HEDigis.size() + f5HBDigis.size() + f3HBDigis.size();
295 auto const nchannels_per_block(alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[0u]);
297 auto const nsamplesForCompute(alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[1u]);
300 float* shrEnergyM0PerTS = alpaka::getDynSharedMem<float>(acc);
301 float* shrChargeMinusPedestal = shrEnergyM0PerTS + nsamplesForCompute * nchannels_per_block;
302 float* shrMethod0EnergyAccum = shrChargeMinusPedestal + nsamplesForCompute * nchannels_per_block;
303 float* shrEnergyM0TotalAccum = shrMethod0EnergyAccum + nchannels_per_block;
304 unsigned long long int* shrMethod0EnergySamplePair =
305 reinterpret_cast<unsigned long long int*
>(shrEnergyM0TotalAccum + nchannels_per_block);
311 auto const gch = channel.global;
312 auto const lch = channel.local;
315 auto const sampleWithinWindow = i_sample;
316 auto const sample = i_sample + startingSample;
317 auto const linearThPerBlock = i_sample + channel.local * nsamplesForCompute;
319 if (sampleWithinWindow == 0) {
320 soiSamples[gch] = -1;
321 shrMethod0EnergyAccum[lch] = 0;
322 shrMethod0EnergySamplePair[lch] = 0;
323 shrEnergyM0TotalAccum[lch] = 0;
327 if (gch < f01HEDigis.size()) {
331 auto const soibit = soibit_for_sample<Flavor1>(&(f01HEDigis.data()[gch][0]),
sample);
333 soiSamples[gch] = sampleWithinWindow;
334 }
else if (gch >= nchannelsf015) {
335 auto const soibit = soibit_for_sample<Flavor3>(&(f3HBDigis.data()[gch - nchannelsf015][0]),
sample);
337 soiSamples[gch] = sampleWithinWindow;
341 auto const id = gch < f01HEDigis.size()
342 ? f01HEDigis.ids()[gch]
343 : (gch < nchannelsf015 ? f5HBDigis.ids()[gch - f01HEDigis.size()]
344 : f3HBDigis.ids()[gch - nchannelsf015]);
348 gch < f01HEDigis.size()
349 ? adc_for_sample<Flavor1>(&(f01HEDigis.data()[gch][0]),
sample)
350 : (gch < nchannelsf015
352 : adc_for_sample<Flavor3>(&(f3HBDigis.data()[gch - nchannelsf015][0]),
sample));
354 gch < f01HEDigis.size()
355 ? capid_for_sample<Flavor1>(&(f01HEDigis.data()[gch][0]),
sample)
356 : (gch < nchannelsf015
357 ? capid_for_sample<Flavor5>(&(f5HBDigis.data()[gch - f01HEDigis.size()][0]),
sample)
362 auto const hashedId =
364 ?
did2linearIndexHB(
id, mahi.maxDepthHB(), mahi.firstHBRing(), mahi.lastHBRing(), mahi.nEtaHB())
371 mahi.offsetForHashes();
374 auto const qieType = mahi.qieTypes_values()[hashedId] > 0 ? 1 : 0;
375 auto const* qieOffsets = mahi.qieCoders_offsets()[hashedId].data();
376 auto const* qieSlopes = mahi.qieCoders_slopes()[hashedId].data();
377 auto const pedestal = mahi.pedestals_value()[hashedId][capid];
385 alpaka::syncBlockThreads(acc);
389 auto const gch = channel.global;
390 auto const lch = channel.local;
392 auto const sampleWithinWindow = i_sample;
393 auto const sample = i_sample + startingSample;
396 if (sampleWithinWindow == 0) {
397 outputGPU.detId()[gch] = 0;
398 outputGPU.energyM0()[gch] = 0;
399 outputGPU.timeM0()[gch] = 0;
400 outputGPU.energy()[gch] = 0;
401 outputGPU.chi2()[gch] = 0;
405 auto* amplitudesForChannel = amplitudes + nsamplesForCompute * gch;
406 auto* noiseTermsForChannel = noiseTerms + nsamplesForCompute * gch;
407 auto* electronicNoiseTermsForChannel = electronicNoiseTerms + nsamplesForCompute * gch;
410 auto const stride = gch < f01HEDigis.size()
411 ? f01HEDigis.stride()
412 : (gch < f5HBDigis.size() ? f5HBDigis.stride() : f3HBDigis.stride());
413 auto const nsamples =
414 gch < f01HEDigis.size()
415 ? compute_nsamples<Flavor1>(
stride)
416 : (gch < nchannelsf015 ? compute_nsamples<Flavor5>(
stride) : compute_nsamples<Flavor3>(
stride));
418 ALPAKA_ASSERT_ACC(nsamples == nsamplesForCompute || nsamples - startingSample == nsamplesForCompute);
420 auto const id = gch < f01HEDigis.size()
421 ? f01HEDigis.ids()[gch]
422 : (gch < nchannelsf015 ? f5HBDigis.ids()[gch - f01HEDigis.size()]
423 : f3HBDigis.ids()[gch - nchannelsf015]);
427 gch < f01HEDigis.size()
428 ? adc_for_sample<Flavor1>(&(f01HEDigis.data()[gch][0]),
sample)
429 : (gch < nchannelsf015
431 : adc_for_sample<Flavor3>(&(f3HBDigis.data()[gch - nchannelsf015][0]),
sample));
433 gch < f01HEDigis.size()
434 ? capid_for_sample<Flavor1>(&(f01HEDigis.data()[gch][0]),
sample)
435 : (gch < nchannelsf015
436 ? capid_for_sample<Flavor5>(&(f5HBDigis.data()[gch - f01HEDigis.size()][0]),
sample)
441 auto const hashedId =
443 ?
did2linearIndexHB(
id, mahi.maxDepthHB(), mahi.firstHBRing(), mahi.lastHBRing(), mahi.nEtaHB())
450 mahi.offsetForHashes();
453 auto const qieType = mahi.qieTypes_values()[hashedId] > 0 ? 1 : 0;
454 auto const* qieOffsets = mahi.qieCoders_offsets()[hashedId].data();
455 auto const* qieSlopes = mahi.qieCoders_slopes()[hashedId].data();
456 auto const* pedestalWidthsForChannel =
458 ? mahi.effectivePedestalWidths()[hashedId].data()
459 : mahi.pedestals_width()[hashedId].data();
461 auto const gain = mahi.gains_value()[hashedId][capid];
462 auto const gain0 = mahi.gains_value()[hashedId][0];
463 auto const respCorrection = mahi.respCorrs_values()[hashedId];
464 auto const pedestal = mahi.pedestals_value()[hashedId][capid];
465 auto const pedestalWidth = pedestalWidthsForChannel[capid];
467 auto const pedestalToUseForMethod0 =
469 ? mahi.effectivePedestals()[hashedId][capid]
471 auto const sipmType = mahi.sipmPar_type()[hashedId];
472 auto const fcByPE = mahi.sipmPar_fcByPE()[hashedId];
473 auto const recoParam1 = recoParamsWithPS.
recoParamView().param1()[hashedId];
474 auto const recoParam2 = recoParamsWithPS.
recoParamView().param2()[hashedId];
476 #ifdef HCAL_MAHI_GPUDEBUG 477 printf(
"qieType = %d qieOffset0 = %f qieOffset1 = %f qieSlope0 = %f qieSlope1 = %f\n",
489 gch < f01HEDigis.size()
491 : (gch < nchannelsf015 ? f5HBDigis.npresamples()[gch - f01HEDigis.size()] : soiSamples[gch]);
493 bool badSOI = (soi < 0 or static_cast<unsigned>(soi) >= nsamplesForCompute);
494 if (badSOI and sampleWithinWindow == 0) {
496 printf(
"Found HBHE channel %d with invalid SOI %d\n", gch, soi);
499 outputGPU.chi2()[gch] = -9999.f;
504 auto const precisionItem = sipmCharacteristics.precisionItem()[sipmType - 1];
505 auto const parLin1 = precisionItem.parLin1_;
506 auto const parLin2 = precisionItem.parLin2_;
507 auto const parLin3 = precisionItem.parLin3_;
513 if (gch >= f01HEDigis.size() && gch < nchannelsf015 && sampleWithinWindow == 0)
514 soiSamples[gch] = f5HBDigis.npresamples()[gch - f01HEDigis.size()];
523 shrChargeMinusPedestal,
533 gch < f01HEDigis.size() || gch >= nchannelsf015);
536 qieType,
adc, capid, qieOffsets, qieSlopes, gch < f01HEDigis.size() || gch >= nchannelsf015);
538 #ifdef COMPUTE_TDC_TIME 540 if (gch >= f01HEDigis.size() && gch < nchannelsf015) {
543 if (gch < f01HEDigis.size())
546 else if (gch >= nchannelsf015)
548 tdc_for_sample<Flavor3>(&(f3HBDigis.data()[gch - nchannelsf015][0]),
sample));
550 #endif // COMPUTE_TDC_TIME 556 auto const nsamplesToAdd = recoParam1 < 10 ? recoParam2 : (recoParam1 >> 14) & 0xF;
558 auto const startSample = startSampleTmp < 0 ? 0 : startSampleTmp;
559 auto const endSample =
560 startSample + nsamplesToAdd < nsamplesForCompute ? startSample + nsamplesToAdd : nsamplesForCompute;
562 auto const energym0_per_ts =
gain * ((rawCharge - pedestalToUseForMethod0) * respCorrection);
563 auto const energym0_per_ts_gain0 = gain0 * ((rawCharge - pedestalToUseForMethod0) * respCorrection);
565 shrEnergyM0PerTS[lch * nsamplesForCompute + sampleWithinWindow] = energym0_per_ts;
568 acc, &shrEnergyM0TotalAccum[lch], energym0_per_ts_gain0, alpaka::hierarchy::Threads{});
570 #ifdef HCAL_MAHI_GPUDEBUG 572 "id = %u sample = %d gch = %d hashedId = %u adc = %u capid = %u\n" 573 " charge = %f rawCharge = %f dfc = %f pedestal = %f\n" 574 " gain = %f respCorrection = %f energym0_per_ts = %f\n",
584 pedestalToUseForMethod0,
588 printf(
"startSample = %d endSample = %d param1 = %u param2 = %u\n",
595 if (sampleWithinWindow >= static_cast<unsigned>(startSample) &&
596 sampleWithinWindow < static_cast<unsigned>(endSample)) {
598 alpaka::atomicAdd(acc, &shrMethod0EnergyAccum[lch], energym0_per_ts, alpaka::hierarchy::Threads{});
602 acc, &shrMethod0EnergySamplePair[lch], {sampleWithinWindow, energym0_per_ts},
TSenergyCompare);
607 if (sampleWithinWindow == static_cast<unsigned>(soi)) {
613 auto const digiStatus_ = mahi.channelQuality_status()[hashedId];
614 const bool taggedBadByDb = (digiStatus_ / 32770);
617 outputGPU.chi2()[gch] = -9999.f;
621 if (not(energym0_per_ts_gain0 >
ts4Thresh)) {
622 outputGPU.chi2()[gch] = -9999.f;
628 auto const amplitude = rawCharge - pedestalToUseForMethod0;
629 auto const noiseADC = (1. /
std::sqrt(12)) * dfc;
631 auto const noiseTerm = noiseADC * noiseADC + noisePhotoSq + pedestalWidth * pedestalWidth;
634 amplitudesForChannel[sampleWithinWindow] =
amplitude;
635 noiseTermsForChannel[sampleWithinWindow] =
noiseTerm;
636 electronicNoiseTermsForChannel[sampleWithinWindow] = pedestalWidth;
640 alpaka::syncBlockThreads(acc);
644 auto const gch = channel.global;
645 auto const lch = channel.local;
647 auto const sampleWithinWindow = i_sample;
650 gch < f01HEDigis.size()
652 : (gch < nchannelsf015 ? f5HBDigis.npresamples()[gch - f01HEDigis.size()] : soiSamples[gch]);
654 auto const id = gch < f01HEDigis.size()
655 ? f01HEDigis.ids()[gch]
656 : (gch < nchannelsf015 ? f5HBDigis.ids()[gch - f01HEDigis.size()]
657 : f3HBDigis.ids()[gch - nchannelsf015]);
659 auto const hashedId =
661 ?
did2linearIndexHB(
id, mahi.maxDepthHB(), mahi.firstHBRing(), mahi.lastHBRing(), mahi.nEtaHB())
668 mahi.offsetForHashes();
670 auto const recoParam1 = recoParamsWithPS.
recoParamView().param1()[hashedId];
671 auto const recoParam2 = recoParamsWithPS.
recoParamView().param2()[hashedId];
673 auto const nsamplesToAdd = recoParam1 < 10 ? recoParam2 : (recoParam1 >> 14) & 0xF;
676 if (sampleWithinWindow == static_cast<unsigned>(soi)) {
677 auto const method0_energy = shrMethod0EnergyAccum[lch];
678 auto const val = shrMethod0EnergySamplePair[lch];
679 int const max_sample = (
val >> 32) & 0xffffffff;
681 float const max_energy =
uint_as_float(static_cast<uint32_t>(
val & 0xffffffff));
683 float const max_energy_1 =
static_cast<unsigned>(max_sample) < nsamplesForCompute - 1
684 ? shrEnergyM0PerTS[lch * nsamplesForCompute + max_sample + 1]
686 float const position = nsamplesToAdd < nsamplesForCompute ? max_sample - soi : max_sample;
687 auto const sum = max_energy + max_energy_1;
693 max_energy > 0.f && max_energy_1 > 0.f ? 25.f * (
position + max_energy_1 / sum) : 25.
f *
position;
696 outputGPU.detId()[gch] =
id;
697 outputGPU.energyM0()[gch] = method0_energy;
698 outputGPU.timeM0()[gch] =
time;
705 if (not(shrEnergyM0TotalAccum[lch] > 0)) {
706 outputGPU.chi2()[gch] = -9999.f;
709 #ifdef HCAL_MAHI_GPUDEBUG 710 printf(
"tsTOT = %f tstrig = %f ts4Thresh = %f\n",
711 shrEnergyM0TotalAccum[lch],
712 energym0_per_ts_gain0,
716 #ifdef HCAL_MAHI_GPUDEBUG 717 printf(
" method0_energy = %f max_sample = %d max_energy = %f time = %f\n",
724 #ifdef HCAL_MAHI_GPUDEBUG 726 "charge(%d) = %f pedestal(%d) = %f dfc(%d) = %f pedestalWidth(%d) = %f noiseADC(%d) = %f " 727 "noisPhoto(%d) =%f outputGPU chi2[gch] = %f \n",
731 pedestalToUseForMethod0,
740 outputGPU.chi2()[gch]);
751 template <
typename TAcc,
typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
753 float* pulseMatrices,
754 float* pulseMatricesM,
755 float* pulseMatricesP,
756 HcalMahiPulseOffsetsPortableDevice::ConstView
pulseOffsets,
757 float const* amplitudes,
758 IProductTypef01::ConstView f01HEDigis,
759 IProductTypef5::ConstView f5HBDigis,
760 IProductTypef3::ConstView f3HBDigis,
762 HcalMahiConditionsPortableDevice::ConstView mahi,
768 float const tzeroTimeSlew,
769 float const slopeTimeSlew,
770 float const tmaxTimeSlew)
const {
772 auto const npulses(alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[1u]);
774 auto const nsamples(alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[2u]);
776 auto const nchannels = f01HEDigis.size() + f5HBDigis.size() + f3HBDigis.size();
777 auto const nchannelsf015 = f01HEDigis.size() + f5HBDigis.size();
786 auto const id = channel < f01HEDigis.size()
787 ? f01HEDigis.ids()[channel]
788 : (channel < nchannelsf015 ? f5HBDigis.ids()[channel - f01HEDigis.size()]
789 : f3HBDigis.ids()[channel - nchannelsf015]);
795 auto const did =
DetId{
id};
796 auto const hashedId =
798 ?
did2linearIndexHB(
id, mahi.maxDepthHB(), mahi.firstHBRing(), mahi.lastHBRing(), mahi.nEtaHB())
805 mahi.offsetForHashes();
806 auto const pulseShape = recoParamsWithPS.
getPulseShape(hashedId);
809 auto* pulseMatrix = pulseMatrices + nsamples * npulses * channel;
810 auto* pulseMatrixM = pulseMatricesM + nsamples * npulses * channel;
811 auto* pulseMatrixP = pulseMatricesP + nsamples * npulses * channel;
814 int const soi = soiSamples[channel];
816 auto const amplitude = amplitudes[channel * nsamples + pulseOffset + soi];
819 pulseMatrix[ipulse * nsamples +
sample] = 0.f;
820 pulseMatrixM[ipulse * nsamples +
sample] = 0.f;
821 pulseMatrixP[ipulse * nsamples +
sample] = 0.f;
824 #ifdef HCAL_MAHI_GPUDEBUG 825 #ifdef HCAL_MAHI_GPUDEBUG_FILTERDETID 826 if (
id != DETID_TO_DEBUG)
829 if (
sample == 0 && ipulse == 0) {
830 for (
int i = 0;
i < 8;
i++)
831 printf(
"amplitude(%d) = %f\n",
i, amplitudes[channel * nsamples +
i]);
832 printf(
"acc25nsVec and diff25nsItvlVec for recoPulseShapeId = %u\n", recoPulseShapeId);
833 for (
int i = 0;
i < 256;
i++) {
834 printf(
"acc25nsVec(%d) = %f diff25nsItvlVec(%d) = %f\n",
i, acc25nsVec[
i],
i, diff25nsItvlVec[
i]);
836 printf(
"accVarLenIdxZEROVec and accVarLenIdxMinusOneVec\n");
837 for (
int i = 0;
i < 25;
i++) {
838 printf(
"accVarLenIdxZEROVec(%d) = %f accVarLenIdxMinusOneVec(%d) = %f\n",
840 accVarLenIdxZeroVec[
i],
842 accVarLenIdxMinusOneVec[
i]);
844 printf(
"diffVarItvlIdxZEROVec and diffVarItvlIdxMinusOneVec\n");
845 for (
int i = 0;
i < 25;
i++) {
846 printf(
"diffVarItvlIdxZEROVec(%d) = %f diffVarItvlIdxMinusOneVec(%d) = %f\n",
848 diffVarItvlIdxZeroVec[
i],
850 diffVarItvlIdxMinusOneVec[
i]);
862 auto const t0m = -deltaT +
t0;
863 auto const t0p = deltaT +
t0;
865 #ifdef HCAL_MAHI_GPUDEBUG 866 if (
sample == 0 && ipulse == 0) {
867 printf(
"time values: %f %f %f\n",
t0, t0m, t0p);
870 if (
sample == 0 && ipulse == 0) {
887 auto const shift = 4 - soi;
889 int32_t
const idx =
sample - pulseOffset;
890 auto const value =
idx >= 0 &&
static_cast<unsigned>(
idx) < nsamples
893 auto const value_t0m =
idx >= 0 &&
static_cast<unsigned>(
idx) < nsamples
896 auto const value_t0p =
idx >= 0 &&
static_cast<unsigned>(
idx) < nsamples
902 pulseMatrixM[ipulse * nsamples +
sample] = value_t0m;
903 pulseMatrixP[ipulse * nsamples +
sample] = value_t0p;
911 template <
int NSAMPLES,
int NPULSES>
919 for (
int ipulse = 0; ipulse < NPULSES; ipulse++) {
920 auto const resultAmplitude = resultAmplitudesVector(ipulse);
921 if (resultAmplitude == 0)
924 #ifdef HCAL_MAHI_GPUDEBUG 925 printf(
"pulse cov array for ibx = %d\n", ipulse);
929 float pmcol[NSAMPLES], pmpcol[NSAMPLES], pmmcol[NSAMPLES];
937 auto const ampl2 = resultAmplitude * resultAmplitude;
939 for (
int col = 0;
col < NSAMPLES;
col++) {
940 auto const valueP_col = pmpcol[
col];
941 auto const valueM_col = pmmcol[
col];
942 auto const value_col = pmcol[
col];
943 auto const tmppcol = valueP_col - value_col;
944 auto const tmpmcol = valueM_col - value_col;
947 auto tmp_value = 0.5 * (tmppcol * tmppcol + tmpmcol * tmpmcol);
948 covarianceMatrix(
col,
col) += ampl2 * tmp_value;
952 for (
int row =
col + 1; row < NSAMPLES; row++) {
953 float const valueP_row = pmpcol[row];
954 float const value_row = pmcol[row];
955 float const valueM_row = pmmcol[row];
957 float tmpprow = valueP_row - value_row;
958 float tmpmrow = valueM_row - value_row;
960 auto const covValue = 0.5 * (tmppcol * tmpprow + tmpmcol * tmpmrow);
962 covarianceMatrix(row,
col) += ampl2 * covValue;
968 template <
int NSAMPLES,
int NPULSES>
971 template <
typename TAcc,
typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
974 float const* amplitudes,
975 float* pulseMatrices,
976 float* pulseMatricesM,
977 float* pulseMatricesP,
978 HcalMahiPulseOffsetsPortableDevice::ConstView pulseOffsetsView,
980 float* electronicNoiseTerms,
982 HcalMahiConditionsPortableDevice::ConstView mahi,
984 IProductTypef01::ConstView f01HEDigis,
985 IProductTypef5::ConstView f5HBDigis,
986 IProductTypef3::ConstView f3HBDigis)
const {
988 static_assert(NPULSES == NSAMPLES);
990 auto const nchannels = f01HEDigis.size() + f5HBDigis.size() + f3HBDigis.size();
992 auto const nchannelsf015 = f01HEDigis.size() + f5HBDigis.size();
994 auto const threadsPerBlock(alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[0u]);
1000 auto const gch = channel.global;
1001 auto const lch = channel.local;
1004 if (outputGPU.chi2()[gch] == -9999.f)
1008 float* shrmem = alpaka::getDynSharedMem<float>(acc);
1013 auto const id = gch < f01HEDigis.size() ? f01HEDigis.ids()[gch]
1014 : (gch < nchannelsf015 ? f5HBDigis.ids()[gch - f01HEDigis.size()]
1015 : f3HBDigis.ids()[gch - nchannelsf015]);
1016 auto const did =
DetId{
id};
1017 auto const hashedId =
1019 ?
did2linearIndexHB(
id, mahi.maxDepthHB(), mahi.firstHBRing(), mahi.lastHBRing(), mahi.nEtaHB())
1026 mahi.offsetForHashes();
1029 auto const* pedestalWidthsForChannel =
1031 ? mahi.effectivePedestalWidths()[hashedId].data()
1032 : mahi.pedestals_width()[hashedId].data();
1033 auto const averagePedestalWidth2 = 0.25 * (pedestalWidthsForChannel[0] * pedestalWidthsForChannel[0] +
1034 pedestalWidthsForChannel[1] * pedestalWidthsForChannel[1] +
1035 pedestalWidthsForChannel[2] * pedestalWidthsForChannel[2] +
1036 pedestalWidthsForChannel[3] * pedestalWidthsForChannel[3]);
1039 auto const gain = mahi.gains_value()[hashedId][0];
1041 auto const respCorrection = mahi.respCorrs_values()[hashedId];
1042 auto const noisecorr = mahi.sipmPar_auxi2()[hashedId];
1044 #ifdef HCAL_MAHI_GPUDEBUG 1045 #ifdef HCAL_MAHI_GPUDEBUG_FILTERDETID 1046 if (
id != DETID_TO_DEBUG)
1053 for (
int i = 0;
i < NPULSES; ++
i)
1061 Eigen::Map<const calo::multifit::ColumnVector<NSAMPLES>> inputAmplitudesView{amplitudes + gch * NSAMPLES};
1062 Eigen::Map<const calo::multifit::ColumnVector<NSAMPLES>> noiseTermsView{noiseTerms + gch * NSAMPLES};
1063 Eigen::Map<const calo::multifit::ColumnVector<NSAMPLES>> noiseElectronicView{electronicNoiseTerms +
1065 Eigen::Map<const calo::multifit::ColMajorMatrix<NSAMPLES, NPULSES>> glbPulseMatrixMView{
1066 pulseMatricesM + gch * NSAMPLES * NPULSES};
1067 Eigen::Map<const calo::multifit::ColMajorMatrix<NSAMPLES, NPULSES>> glbPulseMatrixPView{
1068 pulseMatricesP + gch * NSAMPLES * NPULSES};
1069 Eigen::Map<const calo::multifit::ColMajorMatrix<NSAMPLES, NPULSES>> glbPulseMatrixView{
1070 pulseMatrices + gch * NSAMPLES * NPULSES};
1071 #ifdef HCAL_MAHI_GPUDEBUG 1072 for (
int i = 0;
i < NSAMPLES;
i++)
1073 printf(
"inputValues(%d) = %f noiseTerms(%d) = %f\n",
i, inputAmplitudesView(
i),
i, noiseTermsView(
i));
1074 for (
int i = 0;
i < NSAMPLES;
i++) {
1075 for (
int j = 0;
j < NPULSES;
j++)
1076 printf(
"%f ", glbPulseMatrixView(
i,
j));
1080 for (
int i = 0;
i < NSAMPLES;
i++) {
1081 for (
int j = 0;
j < NPULSES;
j++)
1082 printf(
"%f ", glbPulseMatrixMView(
i,
j));
1086 for (
int i = 0;
i < NSAMPLES;
i++) {
1087 for (
int j = 0;
j < NPULSES;
j++)
1088 printf(
"%f ", glbPulseMatrixPView(
i,
j));
1094 float chi2 = 0, previous_chi2 = 0.f, chi2_2itersback = 0.f;
1100 float* covarianceMatrixStorage = shrMatrixLFnnlsStorage;
1104 covarianceMatrixStorage[
counter] = (noisecorr != 0.
f) ? 0.f : averagePedestalWidth2;
1110 noisecorr * noiseElectronicView.coeffRef(
counter - 1) * noiseElectronicView.coeffRef(
counter);
1117 glbPulseMatrixMView,
1118 glbPulseMatrixPView);
1120 #ifdef HCAL_MAHI_GPUDEBUG 1121 printf(
"covariance matrix\n");
1122 for (
int i = 0;
i < 8;
i++) {
1123 for (
int j = 0;
j < 8;
j++)
1124 printf(
"%f ", covarianceMatrix(
i,
j));
1151 float reg_b[NSAMPLES];
1163 for (
int icol = 0; icol < NPULSES; icol++) {
1164 float reg_ai[NSAMPLES];
1178 AtA(icol, icol) = sum;
1182 for (
int j = icol + 1;
j < NPULSES;
j++) {
1184 float reg_aj[NSAMPLES];
1207 Atb(icol) = sum_atb;
1210 #ifdef HCAL_MAHI_GPUDEBUG 1212 for (
int i = 0;
i < 8;
i++) {
1213 for (
int j = 0;
j < 8;
j++)
1214 printf(
"%f ", AtA(
i,
j));
1218 for (
int i = 0;
i < 8;
i++)
1219 printf(
"%f ", Atb(
i));
1221 printf(
"result Amplitudes before nnls\n");
1222 for (
int i = 0;
i < 8;
i++)
1223 printf(
"%f ", resultAmplitudesVector(
i));
1233 resultAmplitudesVector,
1242 #ifdef HCAL_MAHI_GPUDEBUG 1243 printf(
"result Amplitudes\n");
1244 for (
int i = 0;
i < 8;
i++)
1245 printf(
"resultAmplitudes(%d) = %f\n",
i, resultAmplitudesVector(
i));
1249 matrixL, glbPulseMatrixView, resultAmplitudesVector, inputAmplitudesView,
chi2);
1252 if (
chi2 == chi2_2itersback &&
chi2 < previous_chi2)
1256 chi2_2itersback = previous_chi2;
1257 previous_chi2 =
chi2;
1264 #ifdef HCAL_MAHI_GPUDEBUG 1265 for (
int i = 0;
i < NPULSES;
i++)
1266 printf(
"pulseOffsets(%d) = %d outputAmplitudes(%d) = %f\n",
1270 resultAmplitudesVector(
i));
1271 printf(
"chi2 = %f\n",
chi2);
1274 outputGPU.chi2()[gch] =
chi2;
1275 auto const idx_for_energy =
std::abs(pulseOffsetsView.offsets()[0]);
1276 outputGPU.energy()[gch] = (
gain * resultAmplitudesVector(idx_for_energy)) * respCorrection;
1286 IProductTypef01::ConstView
const& f01HEDigis,
1287 IProductTypef5::ConstView
const& f5HBDigis,
1288 IProductTypef3::ConstView
const& f3HBDigis,
1290 HcalMahiConditionsPortableDevice::ConstView
const& mahi,
1291 HcalSiPMCharacteristicsPortableDevice::ConstView
const& sipmCharacteristics,
1293 HcalMahiPulseOffsetsPortableDevice::ConstView
const& mahiPulseOffsets,
1295 auto const totalChannels =
1296 f01HEDigis.metadata().size() + f5HBDigis.metadata().size() + f3HBDigis.metadata().size();
1311 Vec2D const blocks_2d{blocks_y, 1u};
1312 Vec2D const threads_2d{nchannels_per_block, windowSize};
1313 auto workDivPrep2D = cms::alpakatools::make_workdiv<Acc2D>(blocks_2d, threads_2d);
1316 auto amplitudes = cms::alpakatools::make_device_buffer<float[]>(
queue, totalChannels * windowSize);
1317 auto noiseTerms = cms::alpakatools::make_device_buffer<float[]>(
queue, totalChannels * windowSize);
1318 auto electronicNoiseTerms = cms::alpakatools::make_device_buffer<float[]>(
queue, totalChannels * windowSize);
1319 auto soiSamples = cms::alpakatools::make_device_buffer<int8_t[]>(
queue, totalChannels * windowSize);
1321 alpaka::exec<Acc2D>(
queue,
1329 sipmCharacteristics,
1338 electronicNoiseTerms.data(),
1344 uint32_t
const channelsPerBlock = 1024 / (windowSize * mahiPulseOffsets.metadata().size());
1348 Vec3D const blocks_3d{blocks_z, 1u, 1u};
1349 Vec3D const threads_3d{channelsPerBlock, mahiPulseOffsets.metadata().size(), windowSize};
1351 auto workDivPrep3D = cms::alpakatools::make_workdiv<Acc3D>(blocks_3d, threads_3d);
1354 auto pulseMatrices = cms::alpakatools::make_device_buffer<float[]>(
1355 queue, totalChannels * windowSize * mahiPulseOffsets.metadata().size());
1356 auto pulseMatricesM = cms::alpakatools::make_device_buffer<float[]>(
1357 queue, totalChannels * windowSize * mahiPulseOffsets.metadata().size());
1358 auto pulseMatricesP = cms::alpakatools::make_device_buffer<float[]>(
1359 queue, totalChannels * windowSize * mahiPulseOffsets.metadata().size());
1361 alpaka::exec<Acc3D>(
queue,
1364 pulseMatrices.data(),
1365 pulseMatricesM.data(),
1366 pulseMatricesP.data(),
1386 auto workDivPrep1D = cms::alpakatools::make_workdiv<Acc1D>(blocks_1d, threadsPerBlock);
1388 alpaka::exec<Acc1D>(
queue,
1393 pulseMatrices.data(),
1394 pulseMatricesM.data(),
1395 pulseMatricesP.data(),
1398 electronicNoiseTerms.data(),
1414 template <
typename TAcc>
1417 template <
typename TVec,
typename... TArgs>
1419 TVec
const& threadsPerBlock,
1420 TVec
const& elemsPerThread,
1421 TArgs
const&...) -> std::size_t {
1427 std::size_t bytes = threadsPerBlock[0u] * elemsPerThread[0u] *
1428 ((2 * threadsPerBlock[1u] * elemsPerThread[1u] + 2) *
sizeof(
float) +
sizeof(
uint64_t));
1434 template <
int NSAMPLES,
int NPULSES,
typename TAcc>
1437 template <
typename TVec,
typename... TArgs>
1439 TVec
const& threadsPerBlock,
1440 TVec
const& elemsPerThread,
1441 TArgs
const&...) -> std::size_t {
ALPAKA_FN_ACC float get_raw_charge(double const charge, double const pedestal, float const *shrChargeMinusPedestal, float const parLin1, float const parLin2, float const parLin3, int32_t const nsamplesForCompute, int32_t const soi, int const sipmQTSShift, int const sipmQNTStoSum, float const fcByPE, int32_t const lch, bool const isqie11)
ALPAKA_STATIC_ACC_MEM_CONSTANT float const qie8shape[129]
constexpr float deltaChi2Threashold
ALPAKA_FN_ACC void operator()(TAcc const &acc, OProductType::View outputGPU, float const *amplitudes, float *pulseMatrices, float *pulseMatricesM, float *pulseMatricesP, HcalMahiPulseOffsetsPortableDevice::ConstView pulseOffsetsView, float *noiseTerms, float *electronicNoiseTerms, int8_t *soiSamples, HcalMahiConditionsPortableDevice::ConstView mahi, bool const useEffectivePedestals, IProductTypef01::ConstView f01HEDigis, IProductTypef5::ConstView f5HBDigis, IProductTypef3::ConstView f3HBDigis) const
static const double slope[3]
constexpr bool TSenergyCompare(std::pair< unsigned int, float > a, std::pair< unsigned int, float > b)
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void compute_decomposition_unrolled(MatrixType1 &L, MatrixType2 const &M)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE void atomicMaxPair(const TAcc &acc, unsigned long long int *address, std::pair< unsigned int, float > value, F comparator)
ALPAKA_FN_ACC ALPAKA_FN_INLINE float compute_time_slew_delay(float const fC, float const tzero, float const slope, float const tmax)
ALPAKA_STATIC_ACC_MEM_CONSTANT float const qie11shape[257]
constexpr int nMaxItersMin
ALPAKA_FN_ACC PulseShapeConstElement getPulseShape(uint32_t const hashedId) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE float compute_coder_charge(int const qieType, uint8_t const adc, uint8_t const capid, float const *qieOffsets, float const *qieSlopes)
ALPAKA_FN_ACC ALPAKA_FN_INLINE float compute_diff_charge_gain(int const qieType, uint8_t adc, uint8_t const capid, float const *qieOffsets, float const *qieSlopes, bool const isqie11)
constexpr double nnlsThresh
Eigen::Matrix< T, SIZE, 1 > ColumnVector
void runMahiAsync(Queue &queue, IProductTypef01::ConstView const &f01HEDigis, IProductTypef5::ConstView const &f5HBDigis, IProductTypef3::ConstView const &f3HBDigis, OProductType::View outputGPU, HcalMahiConditionsPortableDevice::ConstView const &mahi, HcalSiPMCharacteristicsPortableDevice::ConstView const &sipmCharacteristics, HcalRecoParamWithPulseShapeDevice::ConstView const &recoParamsWithPS, HcalMahiPulseOffsetsPortableDevice::ConstView const &mahiPulseOffsets, ConfigParameters const &configParameters)
constexpr uint32_t stride
constexpr float getTDCTime(const int tdc)
ALPAKA_FN_ACC ALPAKA_FN_INLINE float compute_pulse_shape_value(PulseShapeConstElement const &pulseShape, float const pulse_time, int const sample, int const shift)
EIGEN_DEVICE_FUNC void fnnls(MatrixType const &AtA, VectorType const &Atb, VectorType &solution, int &npassive, ColumnVector< VectorType::RowsAtCompileTime, int > &pulseOffsets, MapSymM< float, VectorType::RowsAtCompileTime > &matrixL, double eps, const int maxIterations, const int relaxationPeriod, const int relaxationFactor)
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void calculateChiSq(MatrixType1 const &matrixL, MatrixType2 const &pulseMatrixView, MatrixType3 const &resultAmplitudesVector, MatrixType4 const &inputAmplitudesView, float &chi2)
EIGEN_DEVICE_FUNC void solve_forward_subst_matrix(MatrixType1 &A, MatrixType2 const &pulseMatrixView, MatrixType3 const &matrixL)
constexpr int32_t IPHI_MAX
Abs< T >::type abs(const T &t)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum) ...
std::array< uint32_t, 3 > kernelMinimizeThreads
constexpr uint8_t adc_for_sample< Flavor5 >(uint16_t const *const dfstart, uint32_t const sample)
static const double tmax[3]
typename HcalPulseShapeSoA::ConstView::const_element PulseShapeConstElement
EIGEN_DEVICE_FUNC void solve_forward_subst_vector(float reg_b[MatrixType1::RowsAtCompileTime], MatrixType1 inputAmplitudesView, MatrixType2 matrixL)
static ALPAKA_FN_HOST_ACC auto getBlockSharedMemDynSizeBytes(Kernel_minimize< NSAMPLES, NPULSES > const &, TVec const &threadsPerBlock, TVec const &elemsPerThread, TArgs const &...) -> std::size_t
unsigned long long uint64_t
uint32_t kprep1dChannelsPerBlock
ALPAKA_FN_ACC ALPAKA_FN_INLINE float compute_reco_correction_factor(float const par1, float const par2, float const par3, float const x)
constexpr uint8_t capid_for_sample< Flavor3 >(uint16_t const *const dfstart, uint32_t const sample)
constexpr float UNKNOWN_T_NOTDC
static std::atomic< unsigned int > counter
static int position[264][3]
ALPAKA_FN_ACC ALPAKA_FN_INLINE void update_covariance(calo::multifit::ColumnVector< NPULSES > const &resultAmplitudesVector, calo::multifit::MapSymM< float, NSAMPLES > &covarianceMatrix, Eigen::Map< const calo::multifit::ColMajorMatrix< NSAMPLES, NPULSES >> const &pulseMatrix, Eigen::Map< const calo::multifit::ColMajorMatrix< NSAMPLES, NPULSES >> const &pulseMatrixM, Eigen::Map< const calo::multifit::ColMajorMatrix< NSAMPLES, NPULSES >> const &pulseMatrixP)
static unsigned int const shift
static const double tzero[3]
bool useEffectivePedestals
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE uint32_t float_as_uint(float val)
constexpr float iniTimeShift
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t did2linearIndexHE(uint32_t const didraw, int const maxDepthHE, int const maxPhiHE, int const firstHERing, int const lastHERing, int const nEtaHE)
ALPAKA_ASSERT_ACC(offsets)
Eigen::Matrix< float, NROWS, NCOLS, Eigen::ColMajor > ColMajorMatrix
constexpr int nMaxItersNNLS
T1 atomicAdd(T1 *a, T2 b)
constexpr RecoParamCollection::ConstView recoParamView()
static const int IPHI_MAX
uint16_t *__restrict__ uint16_t const *__restrict__ adc
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t did2linearIndexHB(uint32_t const didraw, int const maxDepthHB, int const firstHBRing, int const lastHBRing, int const nEtaHB)
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t get_qiecoder_index(uint32_t const capid, uint32_t const range)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE float uint_as_float(uint32_t val)