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