1 #ifndef RecoLocalCalo_EcalRecProducers_plugins_alpaka_AmplitudeComputationCommonKernels_h 2 #define RecoLocalCalo_EcalRecProducers_plugins_alpaka_AmplitudeComputationCommonKernels_h 6 #include <alpaka/alpaka.hpp> 32 template <
typename TAcc,
typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
34 EcalDigiDeviceCollection::ConstView digisDevEB,
35 EcalDigiDeviceCollection::ConstView digisDevEE,
38 EcalMultifitConditionsDevice::ConstView conditionsDev,
41 bool* hasSwitchToGain6,
42 bool* hasSwitchToGain1,
52 auto const offsetForHashes = conditionsDev.offsetEE();
54 auto const nchannelsEB = digisDevEB.size();
55 auto const nchannelsEE = digisDevEE.size();
56 auto const nchannels = nchannelsEB + nchannelsEE;
57 auto const totalElements = nchannels * nsamples;
59 auto const elemsPerBlock = alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[0u];
61 char* shared_mem = alpaka::getDynSharedMem<char>(acc);
62 auto* shr_hasSwitchToGain6 =
reinterpret_cast<bool*
>(shared_mem);
63 auto* shr_hasSwitchToGain1 = shr_hasSwitchToGain6 + elemsPerBlock;
64 auto* shr_hasSwitchToGain0 = shr_hasSwitchToGain1 + elemsPerBlock;
65 auto* shr_isSaturated = shr_hasSwitchToGain0 + elemsPerBlock;
66 auto* shr_hasSwitchToGain0_tmp = shr_isSaturated + elemsPerBlock;
67 auto* shr_counts =
reinterpret_cast<char*
>(shr_hasSwitchToGain0_tmp) + elemsPerBlock;
72 if (
idx.global == 0) {
73 uncalibRecHitsEB.size() = nchannelsEB;
74 uncalibRecHitsEE.size() = nchannelsEE;
77 auto const ch =
idx.global / nsamples;
79 int const inputTx = ch >= nchannelsEB ?
idx.global - nchannelsEB * nsamples :
idx.global;
81 auto const* digis_in = ch >= nchannelsEB ? digisDevEE.data()->data() : digisDevEB.data()->data();
88 shr_hasSwitchToGain0[
idx.local] = shr_hasSwitchToGain0_tmp[
idx.local];
89 shr_counts[
idx.local] = 0;
92 alpaka::syncBlockThreads(acc);
98 if (
idx.local <= elemsPerBlock - 5) {
100 for (
int i = 0;
i < 5; ++
i)
101 shr_counts[
idx.local] += shr_hasSwitchToGain0[
idx.local +
i];
103 shr_isSaturated[
idx.local] = shr_counts[
idx.local] == 5;
109 shr_hasSwitchToGain6[
idx.local] = shr_hasSwitchToGain6[
idx.local] || shr_hasSwitchToGain6[
idx.local + 5];
110 shr_hasSwitchToGain1[
idx.local] = shr_hasSwitchToGain1[
idx.local] || shr_hasSwitchToGain1[
idx.local + 5];
114 shr_hasSwitchToGain0_tmp[
idx.local] =
115 shr_hasSwitchToGain0_tmp[
idx.local] || shr_hasSwitchToGain0_tmp[
idx.local + 5];
119 alpaka::syncBlockThreads(acc);
122 auto const sample =
idx.local % nsamples;
126 shr_hasSwitchToGain6[
idx.local] = shr_hasSwitchToGain6[
idx.local] || shr_hasSwitchToGain6[
idx.local + 2] ||
127 shr_hasSwitchToGain6[
idx.local + 3];
128 shr_hasSwitchToGain1[
idx.local] = shr_hasSwitchToGain1[
idx.local] || shr_hasSwitchToGain1[
idx.local + 2] ||
129 shr_hasSwitchToGain1[
idx.local + 3];
131 shr_hasSwitchToGain0_tmp[
idx.local] = shr_hasSwitchToGain0_tmp[
idx.local] ||
132 shr_hasSwitchToGain0_tmp[
idx.local + 2] ||
133 shr_hasSwitchToGain0_tmp[
idx.local + 3];
138 shr_isSaturated[
idx.local] = shr_isSaturated[
idx.local + 3] || shr_isSaturated[
idx.local + 4];
142 alpaka::syncBlockThreads(acc);
145 auto const ch =
idx.global / nsamples;
146 auto const sample =
idx.local % nsamples;
149 shr_hasSwitchToGain6[
idx.local] = shr_hasSwitchToGain6[
idx.local] || shr_hasSwitchToGain6[
idx.local + 1];
150 shr_hasSwitchToGain1[
idx.local] = shr_hasSwitchToGain1[
idx.local] || shr_hasSwitchToGain1[
idx.local + 1];
151 shr_hasSwitchToGain0_tmp[
idx.local] =
152 shr_hasSwitchToGain0_tmp[
idx.local] || shr_hasSwitchToGain0_tmp[
idx.local + 1];
154 hasSwitchToGain6[ch] = shr_hasSwitchToGain6[
idx.local];
155 hasSwitchToGain1[ch] = shr_hasSwitchToGain1[
idx.local];
157 shr_isSaturated[
idx.local + 3] = shr_isSaturated[
idx.local] || shr_isSaturated[
idx.local + 1];
165 alpaka::syncBlockThreads(acc);
168 auto const ch =
idx.global / nsamples;
169 auto const sample =
idx.local % nsamples;
172 int const inputCh = ch >= nchannelsEB ? ch - nchannelsEB : ch;
173 int const inputTx = ch >= nchannelsEB ?
idx.global - nchannelsEB * nsamples :
idx.global;
175 auto const* dids = ch >= nchannelsEB ? digisDevEE.id() : digisDevEB.id();
176 auto const did =
DetId{dids[inputCh]};
183 auto const* digis_in = ch >= nchannelsEB ? digisDevEE.data()->data() : digisDevEB.data()->data();
186 ch >= nchannelsEB ? uncalibRecHitsEE.outOfTimeAmplitudes()->data()
187 : uncalibRecHitsEB.outOfTimeAmplitudes()->data());
188 auto* energies = ch >= nchannelsEB ? uncalibRecHitsEE.amplitude() : uncalibRecHitsEB.amplitude();
189 auto*
chi2 = ch >= nchannelsEB ? uncalibRecHitsEE.chi2() : uncalibRecHitsEB.chi2();
190 auto* g_pedestal = ch >= nchannelsEB ? uncalibRecHitsEE.pedestal() : uncalibRecHitsEB.pedestal();
191 auto* dids_out = ch >= nchannelsEB ? uncalibRecHitsEE.id() : uncalibRecHitsEB.id();
192 auto*
flags = ch >= nchannelsEB ? uncalibRecHitsEE.flags() : uncalibRecHitsEB.flags();
202 pedestal = conditionsDev.pedestals_mean_x1()[hashedId];
203 gainratio = conditionsDev.gain6Over1()[hashedId] * conditionsDev.gain12Over6()[hashedId];
204 gainsNoise[ch](
sample) = 2;
206 pedestal = conditionsDev.pedestals_mean_x12()[hashedId];
208 gainsNoise[ch](
sample) = 0;
210 pedestal = conditionsDev.pedestals_mean_x6()[hashedId];
211 gainratio = conditionsDev.gain12Over6()[hashedId];
212 gainsNoise[ch](
sample) = 1;
222 #ifdef ECAL_RECO_ALPAKA_DEBUG 225 printf(
"adc is zero\n");
231 amplitudesForMinimization[inputCh](
sample) = 0;
236 if (
sample == sample_max) {
241 energies[inputCh] = 0;
243 g_pedestal[inputCh] = 0;
245 dids_out[inputCh] = did.rawId();
248 auto const chStart =
idx.local - sample_max;
250 auto const threadMax =
idx.local;
254 if (shr_hasSwitchToGain6[chStart])
256 if (shr_hasSwitchToGain1[chStart])
262 if (
sample == 0 && shr_hasSwitchToGain0_tmp[
idx.local]) {
268 bool saturated_before_max =
false;
270 for (
char ii = 0;
ii < 5; ++
ii)
271 saturated_before_max = saturated_before_max || shr_hasSwitchToGain0[chStart +
ii];
274 if (!saturated_before_max && shr_hasSwitchToGain0[threadMax])
275 energies[inputCh] = 49140;
291 auto shape_value = conditionsDev.pulseShapes()[hashedId][full_pulse_max - 7];
294 shr_hasSwitchToGain6[chStart] || shr_hasSwitchToGain1[chStart] || shr_isSaturated[chStart + 3];
298 if (hasGainSwitch && gainSwitchUseMaxSample) {
300 energies[inputCh] = max_amplitude / shape_value;
307 auto const rmsForChecking = conditionsDev.pedestals_rms_x12()[hashedId];
312 if (rmsForChecking == 0) {
332 template <
typename TAcc,
typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
334 EcalDigiDeviceCollection::ConstView digisDevEB,
335 EcalDigiDeviceCollection::ConstView digisDevEE,
336 EcalMultifitConditionsDevice::ConstView conditionsDev,
340 bool const* hasSwitchToGain6,
341 bool const* hasSwitchToGain1,
344 auto const offsetForHashes = conditionsDev.offsetEE();
345 auto const nchannelsEB = digisDevEB.size();
346 constexpr float addPedestalUncertainty = 0.f;
351 auto const* pulse_shapes =
reinterpret_cast<const EcalPulseShape*
>(conditionsDev.pulseShapes()->data());
353 auto const blockDimX = alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc)[1u];
354 auto const elemsPerBlockX = alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[1u];
355 auto const elemsPerBlockY = alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[0u];
356 Vec2D const size_2d = {elemsPerBlockY, blockDimX * elemsPerBlockX};
359 auto const ch = ndindex[1] / nsamples;
360 auto const tx = ndindex[1] % nsamples;
361 auto const ty = ndindex[0];
364 int const inputCh = ch >= nchannelsEB ? ch - nchannelsEB : ch;
365 auto const* dids = ch >= nchannelsEB ? digisDevEE.id() : digisDevEB.id();
367 auto const did =
DetId{dids[inputCh]};
371 auto const* G12SamplesCorrelation =
isBarrel ? conditionsDev.sampleCorrelation_EB_G12().data()
372 : conditionsDev.sampleCorrelation_EE_G12().data();
373 auto const* G6SamplesCorrelation =
374 isBarrel ? conditionsDev.sampleCorrelation_EB_G6().data() : conditionsDev.sampleCorrelation_EE_G6().data();
375 auto const* G1SamplesCorrelation =
376 isBarrel ? conditionsDev.sampleCorrelation_EB_G1().data() : conditionsDev.sampleCorrelation_EE_G1().data();
377 auto const hasGainSwitch = hasSwitchToGain6[ch] || hasSwitchToGain1[ch] ||
isSaturated[ch];
379 auto const vidx =
std::abs(static_cast<int>(ty) - static_cast<int>(tx));
384 float noise_value = 0;
392 auto const gainidx = gainsNoise[ch][isample_max];
396 auto const rms_x12 = conditionsDev.pedestals_rms_x12()[hashedId];
397 noise_value = rms_x12 * rms_x12 * G12SamplesCorrelation[vidx];
398 }
else if (gainidx == 1) {
399 auto const gain12Over6 = conditionsDev.gain12Over6()[hashedId];
400 auto const rms_x6 = conditionsDev.pedestals_rms_x6()[hashedId];
401 noise_value = gain12Over6 * gain12Over6 * rms_x6 * rms_x6 * G6SamplesCorrelation[vidx];
402 }
else if (gainidx == 2) {
403 auto const gain12Over6 = conditionsDev.gain12Over6()[hashedId];
404 auto const gain6Over1 = conditionsDev.gain6Over1()[hashedId];
405 auto const gain12Over1 = gain12Over6 * gain6Over1;
406 auto const rms_x1 = conditionsDev.pedestals_rms_x1()[hashedId];
407 noise_value = gain12Over1 * gain12Over1 * rms_x1 * rms_x1 * G1SamplesCorrelation[vidx];
409 if (!dynamicPedestal && addPedestalUncertainty > 0.
f)
410 noise_value += addPedestalUncertainty * addPedestalUncertainty;
416 auto const rms_x12 = conditionsDev.pedestals_rms_x12()[hashedId];
417 noise_value += rms_x12 * rms_x12 *
pedestal * G12SamplesCorrelation[vidx];
419 if (!dynamicPedestal && addPedestalUncertainty > 0.
f) {
420 noise_value += addPedestalUncertainty * addPedestalUncertainty *
pedestal;
427 auto const gain12Over6 = conditionsDev.gain12Over6()[hashedId];
428 auto const rms_x6 = conditionsDev.pedestals_rms_x6()[hashedId];
429 noise_value += gain12Over6 * gain12Over6 * rms_x6 * rms_x6 *
pedestal * G6SamplesCorrelation[vidx];
431 if (!dynamicPedestal && addPedestalUncertainty > 0.
f) {
432 noise_value += gain12Over6 * gain12Over6 * addPedestalUncertainty * addPedestalUncertainty *
pedestal;
439 auto const gain6Over1 = conditionsDev.gain6Over1()[hashedId];
440 auto const gain12Over1 = gain12Over6 * gain6Over1;
441 auto const rms_x1 = conditionsDev.pedestals_rms_x1()[hashedId];
442 noise_value += gain12Over1 * gain12Over1 * rms_x1 * rms_x1 *
pedestal * G1SamplesCorrelation[vidx];
444 if (!dynamicPedestal && addPedestalUncertainty > 0.
f) {
445 noise_value += gain12Over1 * gain12Over1 * addPedestalUncertainty * addPedestalUncertainty *
pedestal;
449 noisecov[ch](ty, tx) = noise_value;
451 auto const rms = conditionsDev.pedestals_rms_x12()[hashedId];
452 float noise_value =
rms *
rms * G12SamplesCorrelation[vidx];
453 if (!dynamicPedestal && addPedestalUncertainty > 0.
f) {
455 noise_value += addPedestalUncertainty * addPedestalUncertainty;
457 noisecov[ch](ty, tx) = noise_value;
460 auto const posToAccess = 9 -
static_cast<int>(tx) + static_cast<int>(ty);
461 float const value = posToAccess >= 7 ? pulse_shapes[hashedId].pdfval[posToAccess - 7] : 0;
462 pulse_matrix[ch](ty, tx) =
value;
473 template <
typename TAcc>
476 template <
typename TVec,
typename... TArgs>
478 TVec
const& threadsPerBlock,
479 TVec
const& elemsPerThread,
480 TArgs
const&...) -> std::size_t {
482 std::size_t bytes = threadsPerBlock[0u] * elemsPerThread[0u] * (5 *
sizeof(
bool) +
sizeof(
char));
488 #endif // RecoLocalCalo_EcalRecProducers_plugins_AmplitudeComputationCommonKernels_h
Eigen::Matrix< data_type, SampleVectorSize, SampleVectorSize > SampleMatrix
ALPAKA_FN_ACC void operator()(TAcc const &acc, EcalDigiDeviceCollection::ConstView digisDevEB, EcalDigiDeviceCollection::ConstView digisDevEE, EcalUncalibratedRecHitDeviceCollection::View uncalibRecHitsEB, EcalUncalibratedRecHitDeviceCollection::View uncalibRecHitsEE, EcalMultifitConditionsDevice::ConstView conditionsDev, ::ecal::multifit::SampleVector *amplitudes, ::ecal::multifit::SampleGainVector *gainsNoise, bool *hasSwitchToGain6, bool *hasSwitchToGain1, bool *isSaturated, char *acState, ::ecal::multifit::BXVectorType *bxs, bool const gainSwitchUseMaxSampleEB, bool const gainSwitchUseMaxSampleEE) const
static ALPAKA_FN_HOST_ACC auto getBlockSharedMemDynSizeBytes(Kernel_prep_1d_and_initialize const &, TVec const &threadsPerBlock, TVec const &elemsPerThread, TArgs const &...) -> std::size_t
Eigen::Matrix< data_type, SampleVectorSize, 1 > SampleVector
Eigen::Matrix< data_type, SampleVectorSize, SampleVectorSize > PulseMatrixType
ALPAKA_FN_ACC uint32_t hashedIndexEB(uint32_t id)
#define EcalMgpaBitwiseGain1
simplifiedNoiseModelForGainSwitch
Eigen::Matrix< char, SampleVectorSize, 1 > SampleGainVector
ALPAKA_FN_ACC uint32_t hashedIndexEE(uint32_t id)
Abs< T >::type abs(const T &t)
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
constexpr int gainId(sample_type sample)
get the gainId (2 bits)
constexpr int gainId(sample_type sample)
get the gainId (2 bits)
bool isSaturated(const Digi &digi, const int &maxADCvalue, int ifirst, int n)
Eigen::Matrix< char, SampleVectorSize, 1 > BXVectorType
#define EcalMgpaBitwiseGain6
static constexpr int MAXSAMPLES
#define EcalMgpaBitwiseGain0
ALPAKA_FN_ACC void operator()(TAcc const &acc, EcalDigiDeviceCollection::ConstView digisDevEB, EcalDigiDeviceCollection::ConstView digisDevEE, EcalMultifitConditionsDevice::ConstView conditionsDev, ::ecal::multifit::SampleGainVector const *gainsNoise, ::ecal::multifit::SampleMatrix *noisecov, ::ecal::multifit::PulseMatrixType *pulse_matrix, bool const *hasSwitchToGain6, bool const *hasSwitchToGain1, bool const *isSaturated) const
uint16_t *__restrict__ uint16_t const *__restrict__ adc