CMS 3D CMS Logo

AmplitudeComputationCommonKernels.h
Go to the documentation of this file.
1 #ifndef RecoLocalCalo_EcalRecProducers_plugins_alpaka_AmplitudeComputationCommonKernels_h
2 #define RecoLocalCalo_EcalRecProducers_plugins_alpaka_AmplitudeComputationCommonKernels_h
3 
4 #include <cstdlib>
5 #include <limits>
6 #include <alpaka/alpaka.hpp>
7 
19 
20 #include "DeclsForKernels.h"
21 #include "KernelHelpers.h"
22 
24 
31  public:
32  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
33  ALPAKA_FN_ACC void operator()(TAcc const& acc,
34  EcalDigiDeviceCollection::ConstView digisDevEB,
35  EcalDigiDeviceCollection::ConstView digisDevEE,
38  EcalMultifitConditionsDevice::ConstView conditionsDev,
39  ::ecal::multifit::SampleVector* amplitudes,
41  bool* hasSwitchToGain6,
42  bool* hasSwitchToGain1,
43  bool* isSaturated,
44  char* acState,
46  bool const gainSwitchUseMaxSampleEB,
47  bool const gainSwitchUseMaxSampleEE) const {
48  constexpr bool dynamicPedestal = false; //---- default to false, ok
49  constexpr auto nsamples = EcalDataFrame::MAXSAMPLES;
50  constexpr int sample_max = 5;
51  constexpr int full_pulse_max = 9;
52  auto const offsetForHashes = conditionsDev.offsetEE();
53 
54  auto const nchannelsEB = digisDevEB.size();
55  auto const nchannelsEE = digisDevEE.size();
56  auto const nchannels = nchannelsEB + nchannelsEE;
57  auto const totalElements = nchannels * nsamples;
58 
59  auto const elemsPerBlock = alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[0u];
60 
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;
68 
69  for (auto block : cms::alpakatools::blocks_with_stride(acc, totalElements)) {
70  for (auto idx : cms::alpakatools::elements_in_block(acc, block, totalElements)) {
71  // set the output collection size scalars
72  if (idx.global == 0) {
73  uncalibRecHitsEB.size() = nchannelsEB;
74  uncalibRecHitsEE.size() = nchannelsEE;
75  }
76 
77  auto const ch = idx.global / nsamples;
78  // for accessing input arrays
79  int const inputTx = ch >= nchannelsEB ? idx.global - nchannelsEB * nsamples : idx.global;
80  // eb is first and then ee
81  auto const* digis_in = ch >= nchannelsEB ? digisDevEE.data()->data() : digisDevEB.data()->data();
82  auto const gainId = ecalMGPA::gainId(digis_in[inputTx]);
83 
84  // store into shared mem for initialization
85  shr_hasSwitchToGain6[idx.local] = gainId == EcalMgpaBitwiseGain6;
86  shr_hasSwitchToGain1[idx.local] = gainId == EcalMgpaBitwiseGain1;
87  shr_hasSwitchToGain0_tmp[idx.local] = gainId == EcalMgpaBitwiseGain0;
88  shr_hasSwitchToGain0[idx.local] = shr_hasSwitchToGain0_tmp[idx.local];
89  shr_counts[idx.local] = 0;
90  }
91 
92  alpaka::syncBlockThreads(acc);
93 
94  for (auto idx : cms::alpakatools::elements_in_block(acc, block, totalElements)) {
95  auto const sample = idx.local % nsamples;
96 
97  // non-divergent branch (except for the last 4 threads)
98  if (idx.local <= elemsPerBlock - 5) {
100  for (int i = 0; i < 5; ++i)
101  shr_counts[idx.local] += shr_hasSwitchToGain0[idx.local + i];
102  }
103  shr_isSaturated[idx.local] = shr_counts[idx.local] == 5;
104 
105  //
106  // unrolled reductions
107  //
108  if (sample < 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];
111 
112  // duplication of hasSwitchToGain0 in order not to
113  // introduce another syncthreads
114  shr_hasSwitchToGain0_tmp[idx.local] =
115  shr_hasSwitchToGain0_tmp[idx.local] || shr_hasSwitchToGain0_tmp[idx.local + 5];
116  }
117  }
118 
119  alpaka::syncBlockThreads(acc);
120 
121  for (auto idx : cms::alpakatools::elements_in_block(acc, block, totalElements)) {
122  auto const sample = idx.local % nsamples;
123 
124  if (sample < 2) {
125  // note, both threads per channel take value [3] twice to avoid another if
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];
130 
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];
134 
135  // sample < 2 -> first 2 threads of each channel will be used here
136  // => 0 -> will compare 3 and 4 and put into 0
137  // => 1 -> will compare 4 and 5 and put into 1
138  shr_isSaturated[idx.local] = shr_isSaturated[idx.local + 3] || shr_isSaturated[idx.local + 4];
139  }
140  }
141 
142  alpaka::syncBlockThreads(acc);
143 
144  for (auto idx : cms::alpakatools::elements_in_block(acc, block, totalElements)) {
145  auto const ch = idx.global / nsamples;
146  auto const sample = idx.local % nsamples;
147 
148  if (sample == 0) {
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];
153 
154  hasSwitchToGain6[ch] = shr_hasSwitchToGain6[idx.local];
155  hasSwitchToGain1[ch] = shr_hasSwitchToGain1[idx.local];
156 
157  shr_isSaturated[idx.local + 3] = shr_isSaturated[idx.local] || shr_isSaturated[idx.local + 1];
158  isSaturated[ch] = shr_isSaturated[idx.local + 3];
159  }
160  }
161 
162  // TODO: w/o this sync, there is a race
163  // if (idx.local == sample_max) below uses max sample thread, not for 0 sample
164  // check if we can remove it
165  alpaka::syncBlockThreads(acc);
166 
167  for (auto idx : cms::alpakatools::elements_in_block(acc, block, totalElements)) {
168  auto const ch = idx.global / nsamples;
169  auto const sample = idx.local % nsamples;
170 
171  // for accessing input arrays
172  int const inputCh = ch >= nchannelsEB ? ch - nchannelsEB : ch;
173  int const inputTx = ch >= nchannelsEB ? idx.global - nchannelsEB * nsamples : idx.global;
174 
175  auto const* dids = ch >= nchannelsEB ? digisDevEE.id() : digisDevEB.id();
176  auto const did = DetId{dids[inputCh]};
177  auto const isBarrel = did.subdetId() == EcalBarrel;
178  // TODO offset for ee, 0 for eb
179  auto const hashedId = isBarrel ? reconstruction::hashedIndexEB(did.rawId())
180  : offsetForHashes + reconstruction::hashedIndexEE(did.rawId());
181 
182  // eb is first and then ee
183  auto const* digis_in = ch >= nchannelsEB ? digisDevEE.data()->data() : digisDevEB.data()->data();
184 
185  auto* amplitudesForMinimization = reinterpret_cast<::ecal::multifit::SampleVector*>(
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();
193 
194  auto const adc = ecalMGPA::adc(digis_in[inputTx]);
195  auto const gainId = ecalMGPA::gainId(digis_in[inputTx]);
199 
200  // TODO: divergent branch
201  if (gainId == 0 || gainId == 3) {
202  pedestal = conditionsDev.pedestals_mean_x1()[hashedId];
203  gainratio = conditionsDev.gain6Over1()[hashedId] * conditionsDev.gain12Over6()[hashedId];
204  gainsNoise[ch](sample) = 2;
205  } else if (gainId == 1) {
206  pedestal = conditionsDev.pedestals_mean_x12()[hashedId];
207  gainratio = 1.;
208  gainsNoise[ch](sample) = 0;
209  } else if (gainId == 2) {
210  pedestal = conditionsDev.pedestals_mean_x6()[hashedId];
211  gainratio = conditionsDev.gain12Over6()[hashedId];
212  gainsNoise[ch](sample) = 1;
213  }
214 
215  // TODO: compile time constant -> branch should be non-divergent
216  if (dynamicPedestal)
217  amplitude = static_cast<::ecal::multifit::SampleVector::Scalar>(adc) * gainratio;
218  else
219  amplitude = (static_cast<::ecal::multifit::SampleVector::Scalar>(adc) - pedestal) * gainratio;
220  amplitudes[ch][sample] = amplitude;
221 
222 #ifdef ECAL_RECO_ALPAKA_DEBUG
223  printf("%d %d %d %d %f %f %f\n", idx.global, ch, sample, adc, amplitude, pedestal, gainratio);
224  if (adc == 0)
225  printf("adc is zero\n");
226 #endif
227 
228  //
229  // initialization
230  //
231  amplitudesForMinimization[inputCh](sample) = 0;
232  bxs[ch](sample) = sample - 5;
233 
234  // select the thread for the max sample
235  //---> hardcoded above to be 5th sample, ok
236  if (sample == sample_max) {
237  //
238  // initialization
239  //
240  acState[ch] = static_cast<char>(MinimizationState::NotFinished);
241  energies[inputCh] = 0;
242  chi2[inputCh] = 0;
243  g_pedestal[inputCh] = 0;
244  uint32_t flag = 0;
245  dids_out[inputCh] = did.rawId();
246 
247  // start of this channel in shared mem
248  auto const chStart = idx.local - sample_max;
249  // thread for the max sample in shared mem
250  auto const threadMax = idx.local;
251  auto const gainSwitchUseMaxSample = isBarrel ? gainSwitchUseMaxSampleEB : gainSwitchUseMaxSampleEE;
252 
253  // this flag setting is applied to all of the cases
254  if (shr_hasSwitchToGain6[chStart])
256  if (shr_hasSwitchToGain1[chStart])
258 
259  // this corresponds to cpu branching on lastSampleBeforeSaturation
260  // likely false
261  // check only for the idx.local corresponding to sample==0
262  if (sample == 0 && shr_hasSwitchToGain0_tmp[idx.local]) {
263  // assign for the case some sample having gainId == 0
264  //energies[inputCh] = amplitudes[ch][sample_max];
265  energies[inputCh] = amplitude;
266 
267  // check if samples before sample_max have true
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];
272 
273  // if saturation is in the max sample and not in the first 5
274  if (!saturated_before_max && shr_hasSwitchToGain0[threadMax])
275  energies[inputCh] = 49140; // 4095 * 12 (maximum ADC range * MultiGainPreAmplifier (MGPA) gain)
276  // This is the actual maximum range that is set when we saturate.
277  //---- AM FIXME : no pedestal subtraction???
278  //It should be "(4095. - pedestal) * gainratio"
279 
280  // set state flag to terminate further processing of this channel
281  acState[ch] = static_cast<char>(MinimizationState::Precomputed);
283  flags[inputCh] = flag;
284  continue;
285  }
286 
287  // according to cpu version
288  // auto max_amplitude = amplitudes[ch][sample_max];
289  auto const max_amplitude = amplitude;
290  // pulse shape template value
291  auto shape_value = conditionsDev.pulseShapes()[hashedId][full_pulse_max - 7];
292  // note, no syncing as the same thread will be accessing here
293  bool hasGainSwitch =
294  shr_hasSwitchToGain6[chStart] || shr_hasSwitchToGain1[chStart] || shr_isSaturated[chStart + 3];
295 
296  // pedestal is final unconditionally
297  g_pedestal[inputCh] = pedestal;
298  if (hasGainSwitch && gainSwitchUseMaxSample) {
299  // thread for sample=0 will access the right guys
300  energies[inputCh] = max_amplitude / shape_value;
301  acState[ch] = static_cast<char>(MinimizationState::Precomputed);
302  flags[inputCh] = flag;
303  continue;
304  }
305 
306  // will be used in the future for setting state
307  auto const rmsForChecking = conditionsDev.pedestals_rms_x12()[hashedId];
308 
309  // this happens cause sometimes rms_x12 is 0...
310  // needs to be checkec why this is the case
311  // general case here is that noisecov is a Zero matrix
312  if (rmsForChecking == 0) {
313  acState[ch] = static_cast<char>(MinimizationState::Precomputed);
314  flags[inputCh] = flag;
315  continue;
316  }
317 
318  // for the case when no shortcuts were taken
319  flags[inputCh] = flag;
320  }
321  }
322  }
323  }
324  };
325 
331  public:
332  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
333  ALPAKA_FN_ACC void operator()(TAcc const& acc,
334  EcalDigiDeviceCollection::ConstView digisDevEB,
335  EcalDigiDeviceCollection::ConstView digisDevEE,
336  EcalMultifitConditionsDevice::ConstView conditionsDev,
337  ::ecal::multifit::SampleGainVector const* gainsNoise,
338  ::ecal::multifit::SampleMatrix* noisecov,
339  ::ecal::multifit::PulseMatrixType* pulse_matrix,
340  bool const* hasSwitchToGain6,
341  bool const* hasSwitchToGain1,
342  bool const* isSaturated) const {
343  constexpr auto nsamples = EcalDataFrame::MAXSAMPLES;
344  auto const offsetForHashes = conditionsDev.offsetEE();
345  auto const nchannelsEB = digisDevEB.size();
346  constexpr float addPedestalUncertainty = 0.f;
347  constexpr bool dynamicPedestal = false;
348  constexpr bool simplifiedNoiseModelForGainSwitch = true; //---- default is true
349 
350  // pulse matrix
351  auto const* pulse_shapes = reinterpret_cast<const EcalPulseShape*>(conditionsDev.pulseShapes()->data());
352 
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}; // {y, x} coordinates
357 
358  for (auto ndindex : cms::alpakatools::elements_with_stride_nd(acc, size_2d)) {
359  auto const ch = ndindex[1] / nsamples;
360  auto const tx = ndindex[1] % nsamples;
361  auto const ty = ndindex[0];
362 
363  // to access input arrays (ids and digis only)
364  int const inputCh = ch >= nchannelsEB ? ch - nchannelsEB : ch;
365  auto const* dids = ch >= nchannelsEB ? digisDevEE.id() : digisDevEB.id();
366 
367  auto const did = DetId{dids[inputCh]};
368  auto const isBarrel = did.subdetId() == EcalBarrel;
369  auto const hashedId = isBarrel ? ecal::reconstruction::hashedIndexEB(did.rawId())
370  : offsetForHashes + ecal::reconstruction::hashedIndexEE(did.rawId());
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];
378 
379  auto const vidx = std::abs(static_cast<int>(ty) - static_cast<int>(tx));
380 
381  // non-divergent branch for all threads per block
382  if (hasGainSwitch) {
383  // TODO: did not include simplified noise model
384  float noise_value = 0;
385 
386  // non-divergent branch - all threads per block
387  // TODO: all of these constants indicate that
388  // that these parts could be splitted into completely different
389  // kernels and run one of them only depending on the config
391  constexpr int isample_max = 5; // according to cpu defs
392  auto const gainidx = gainsNoise[ch][isample_max];
393 
394  // non-divergent branches
395  if (gainidx == 0) {
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];
408  }
409  if (!dynamicPedestal && addPedestalUncertainty > 0.f)
410  noise_value += addPedestalUncertainty * addPedestalUncertainty;
411  } else {
412  int gainidx = 0;
413  char mask = gainidx;
414  int pedestal = gainsNoise[ch][ty] == mask ? 1 : 0;
415  // NB: gainratio is 1, that is why it does not appear in the formula
416  auto const rms_x12 = conditionsDev.pedestals_rms_x12()[hashedId];
417  noise_value += rms_x12 * rms_x12 * pedestal * G12SamplesCorrelation[vidx];
418  // non-divergent branch
419  if (!dynamicPedestal && addPedestalUncertainty > 0.f) {
420  noise_value += addPedestalUncertainty * addPedestalUncertainty * pedestal; // gainratio is 1
421  }
422 
423  //
424  gainidx = 1;
425  mask = gainidx;
426  pedestal = gainsNoise[ch][ty] == mask ? 1 : 0;
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];
430  // non-divergent branch
431  if (!dynamicPedestal && addPedestalUncertainty > 0.f) {
432  noise_value += gain12Over6 * gain12Over6 * addPedestalUncertainty * addPedestalUncertainty * pedestal;
433  }
434 
435  //
436  gainidx = 2;
437  mask = gainidx;
438  pedestal = gainsNoise[ch][ty] == mask ? 1 : 0;
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];
443  // non-divergent branch
444  if (!dynamicPedestal && addPedestalUncertainty > 0.f) {
445  noise_value += gain12Over1 * gain12Over1 * addPedestalUncertainty * addPedestalUncertainty * pedestal;
446  }
447  }
448 
449  noisecov[ch](ty, tx) = noise_value;
450  } else {
451  auto const rms = conditionsDev.pedestals_rms_x12()[hashedId];
452  float noise_value = rms * rms * G12SamplesCorrelation[vidx];
453  if (!dynamicPedestal && addPedestalUncertainty > 0.f) {
454  //---- add fully correlated component to noise covariance to inflate pedestal uncertainty
455  noise_value += addPedestalUncertainty * addPedestalUncertainty;
456  }
457  noisecov[ch](ty, tx) = noise_value;
458  }
459 
460  auto const posToAccess = 9 - static_cast<int>(tx) + static_cast<int>(ty); // see cpu for reference
461  float const value = posToAccess >= 7 ? pulse_shapes[hashedId].pdfval[posToAccess - 7] : 0;
462  pulse_matrix[ch](ty, tx) = value;
463  }
464  }
465  };
466 
467 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::ecal::multifit
468 
469 namespace alpaka::trait {
471 
473  template <typename TAcc>
474  struct BlockSharedMemDynSizeBytes<Kernel_prep_1d_and_initialize, TAcc> {
476  template <typename TVec, typename... TArgs>
477  ALPAKA_FN_HOST_ACC static auto getBlockSharedMemDynSizeBytes(Kernel_prep_1d_and_initialize const&,
478  TVec const& threadsPerBlock,
479  TVec const& elemsPerThread,
480  TArgs const&...) -> std::size_t {
481  // return the amount of dynamic shared memory needed
482  std::size_t bytes = threadsPerBlock[0u] * elemsPerThread[0u] * (5 * sizeof(bool) + sizeof(char));
483  return bytes;
484  }
485  };
486 } // namespace alpaka::trait
487 
488 #endif // RecoLocalCalo_EcalRecProducers_plugins_AmplitudeComputationCommonKernels_h
ALPAKA_FN_ACC auto blocks_with_stride(TAcc const &acc, TArgs... args)
Definition: workdivision.h:784
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
double Scalar
Definition: Definitions.h:25
Vec< Dim2D > Vec2D
Definition: config.h:26
static ALPAKA_FN_HOST_ACC auto getBlockSharedMemDynSizeBytes(Kernel_prep_1d_and_initialize const &, TVec const &threadsPerBlock, TVec const &elemsPerThread, TArgs const &...) -> std::size_t
#define CMS_UNROLL_LOOP
Definition: CMSUnrollLoop.h:47
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
Definition: EcalDataFrame.h:10
Eigen::Matrix< char, SampleVectorSize, 1 > SampleGainVector
ALPAKA_FN_ACC uint32_t hashedIndexEE(uint32_t id)
ALPAKA_FN_ACC auto elements_with_stride_nd(TAcc const &acc)
Definition: workdivision.h:565
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
Definition: value.py:1
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
ii
Definition: cuy.py:589
Definition: DetId.h:17
constexpr int gainId(sample_type sample)
get the gainId (2 bits)
ALPAKA_FN_ACC auto elements_in_block(TAcc const &acc, TArgs... args)
Definition: workdivision.h:999
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
Definition: EcalDataFrame.h:9
static constexpr int MAXSAMPLES
Definition: EcalDataFrame.h:48
#define EcalMgpaBitwiseGain0
Definition: EcalDataFrame.h:11
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