CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
EcalUncalibRecHitMultiFitAlgoPortable.dev.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <limits>
3 #include <alpaka/alpaka.hpp>
4 
6 
10 #include "TimeComputationKernels.h"
11 
12 //#define DEBUG
13 //#define ECAL_RECO_ALPAKA_DEBUG
14 
16 
17  using namespace cms::alpakatools;
18 
19  void launchKernels(Queue& queue,
20  InputProduct const& digisDevEB,
21  InputProduct const& digisDevEE,
22  OutputProduct& uncalibRecHitsDevEB,
23  OutputProduct& uncalibRecHitsDevEE,
24  EcalMultifitConditionsDevice const& conditionsDev,
25  EcalMultifitParametersDevice const& paramsDev,
26  ConfigurationParameters const& configParams) {
27  using digis_type = std::vector<uint16_t>;
28  using dids_type = std::vector<uint32_t>;
29  // according to the cpu setup //----> hardcoded
31  // according to the cpu setup //----> hardcoded
33  auto constexpr kMaxSamples = EcalDataFrame::MAXSAMPLES;
34 
35  auto const ebSize = static_cast<uint32_t>(uncalibRecHitsDevEB.const_view().metadata().size());
36  auto const totalChannels = ebSize + static_cast<uint32_t>(uncalibRecHitsDevEE.const_view().metadata().size());
37 
38  EventDataForScratchDevice scratch(configParams, totalChannels, queue);
39 
40  //
41  // 1d preparation kernel
42  //
43  uint32_t constexpr nchannels_per_block = 32;
44  auto constexpr threads_1d = kMaxSamples * nchannels_per_block;
45  auto const blocks_1d = cms::alpakatools::divide_up_by(totalChannels * kMaxSamples, threads_1d);
46  auto workDivPrep1D = cms::alpakatools::make_workdiv<Acc1D>(blocks_1d, threads_1d);
47  // Since the ::ecal::multifit::X objects are non-dynamic Eigen::Matrix types the returned pointers from the buffers
48  // and the ::ecal::multifit::X* both point to the data.
49  alpaka::exec<Acc1D>(queue,
50  workDivPrep1D,
52  digisDevEB.const_view(),
53  digisDevEE.const_view(),
54  uncalibRecHitsDevEB.view(),
55  uncalibRecHitsDevEE.view(),
56  conditionsDev.const_view(),
57  reinterpret_cast<::ecal::multifit::SampleVector*>(scratch.samplesDevBuf.data()),
58  reinterpret_cast<::ecal::multifit::SampleGainVector*>(scratch.gainsNoiseDevBuf.data()),
59  scratch.hasSwitchToGain6DevBuf.data(),
60  scratch.hasSwitchToGain1DevBuf.data(),
61  scratch.isSaturatedDevBuf.data(),
62  scratch.acStateDevBuf.data(),
63  reinterpret_cast<::ecal::multifit::BXVectorType*>(scratch.activeBXsDevBuf.data()),
66 
67  //
68  // 2d preparation kernel
69  //
70  Vec2D const blocks_2d{1u, totalChannels}; // {y, x} coordiantes
71  Vec2D const threads_2d{kMaxSamples, kMaxSamples};
72  auto workDivPrep2D = cms::alpakatools::make_workdiv<Acc2D>(blocks_2d, threads_2d);
73  alpaka::exec<Acc2D>(queue,
74  workDivPrep2D,
76  digisDevEB.const_view(),
77  digisDevEE.const_view(),
78  conditionsDev.const_view(),
79  reinterpret_cast<::ecal::multifit::SampleGainVector*>(scratch.gainsNoiseDevBuf.data()),
80  reinterpret_cast<::ecal::multifit::SampleMatrix*>(scratch.noisecovDevBuf.data()),
81  reinterpret_cast<::ecal::multifit::PulseMatrixType*>(scratch.pulse_matrixDevBuf.data()),
82  scratch.hasSwitchToGain6DevBuf.data(),
83  scratch.hasSwitchToGain1DevBuf.data(),
84  scratch.isSaturatedDevBuf.data());
85 
86  // run minimization kernels
88  digisDevEB,
89  digisDevEE,
90  uncalibRecHitsDevEB,
91  uncalibRecHitsDevEE,
92  scratch,
93  conditionsDev,
94  configParams,
95  totalChannels);
96 
97  if (configParams.shouldRunTimingComputation) {
98  //
99  // TODO: this guy can run concurrently with other kernels,
100  // there is no dependence on the order of execution
101  //
102  auto const blocks_time_init = blocks_1d;
103  auto const threads_time_init = threads_1d;
104  auto workDivTimeCompInit1D = cms::alpakatools::make_workdiv<Acc1D>(blocks_time_init, threads_time_init);
105  alpaka::exec<Acc1D>(queue,
106  workDivTimeCompInit1D,
108  digisDevEB.const_view(),
109  digisDevEE.const_view(),
110  conditionsDev.const_view(),
111  scratch.sample_valuesDevBuf.value().data(),
112  scratch.sample_value_errorsDevBuf.value().data(),
113  scratch.ampMaxErrorDevBuf.value().data(),
114  scratch.useless_sample_valuesDevBuf.value().data(),
115  scratch.pedestal_numsDevBuf.value().data());
116 
117  //
118  // TODO: small kernel only for EB. It needs to be checked if
120  //
121  // we are running only over EB digis
122  // therefore we need to create threads/blocks only for that
123  auto const threadsFixMGPA = threads_1d;
124  auto const blocksFixMGPA = cms::alpakatools::divide_up_by(kMaxSamples * ebSize, threadsFixMGPA);
125  auto workDivTimeFixMGPAslew1D = cms::alpakatools::make_workdiv<Acc1D>(blocksFixMGPA, threadsFixMGPA);
126  alpaka::exec<Acc1D>(queue,
127  workDivTimeFixMGPAslew1D,
129  digisDevEB.const_view(),
130  digisDevEE.const_view(),
131  conditionsDev.const_view(),
132  scratch.sample_valuesDevBuf.value().data(),
133  scratch.sample_value_errorsDevBuf.value().data(),
134  scratch.useless_sample_valuesDevBuf.value().data());
135 
136  auto const threads_nullhypot = threads_1d;
137  auto const blocks_nullhypot = blocks_1d;
138  auto workDivTimeNullhypot1D = cms::alpakatools::make_workdiv<Acc1D>(blocks_nullhypot, threads_nullhypot);
139  alpaka::exec<Acc1D>(queue,
140  workDivTimeNullhypot1D,
142  scratch.sample_valuesDevBuf.value().data(),
143  scratch.sample_value_errorsDevBuf.value().data(),
144  scratch.useless_sample_valuesDevBuf.value().data(),
145  scratch.chi2sNullHypotDevBuf.value().data(),
146  scratch.sum0sNullHypotDevBuf.value().data(),
147  scratch.sumAAsNullHypotDevBuf.value().data(),
148  totalChannels);
149 
150  constexpr uint32_t nchannels_per_block_makeratio = kMaxSamples;
151  constexpr auto nthreads_per_channel =
152  nchannels_per_block_makeratio * (nchannels_per_block_makeratio - 1) / 2; // n(n-1)/2
153  constexpr auto threads_makeratio = nthreads_per_channel * nchannels_per_block_makeratio;
154  auto const blocks_makeratio =
155  cms::alpakatools::divide_up_by(nthreads_per_channel * totalChannels, threads_makeratio);
156  auto workDivTimeMakeRatio1D = cms::alpakatools::make_workdiv<Acc1D>(blocks_makeratio, threads_makeratio);
157  alpaka::exec<Acc1D>(queue,
158  workDivTimeMakeRatio1D,
160  digisDevEB.const_view(),
161  digisDevEE.const_view(),
162  scratch.sample_valuesDevBuf.value().data(),
163  scratch.sample_value_errorsDevBuf.value().data(),
164  scratch.useless_sample_valuesDevBuf.value().data(),
165  scratch.pedestal_numsDevBuf.value().data(),
166  scratch.sumAAsNullHypotDevBuf.value().data(),
167  scratch.sum0sNullHypotDevBuf.value().data(),
168  scratch.tMaxAlphaBetasDevBuf.value().data(),
169  scratch.tMaxErrorAlphaBetasDevBuf.value().data(),
170  scratch.accTimeMaxDevBuf.value().data(),
171  scratch.accTimeWgtDevBuf.value().data(),
172  scratch.tcStateDevBuf.value().data(),
173  paramsDev.const_view(),
174  configParams.timeFitLimitsFirstEB,
175  configParams.timeFitLimitsFirstEE,
176  configParams.timeFitLimitsSecondEB,
177  configParams.timeFitLimitsSecondEE);
178 
179  auto const threads_findamplchi2 = threads_1d;
180  auto const blocks_findamplchi2 = blocks_1d;
181  auto workDivTimeFindAmplChi21D = cms::alpakatools::make_workdiv<Acc1D>(blocks_findamplchi2, threads_findamplchi2);
182  alpaka::exec<Acc1D>(queue,
183  workDivTimeFindAmplChi21D,
185  digisDevEB.const_view(),
186  digisDevEE.const_view(),
187  scratch.sample_valuesDevBuf.value().data(),
188  scratch.sample_value_errorsDevBuf.value().data(),
189  scratch.useless_sample_valuesDevBuf.value().data(),
190  scratch.tMaxAlphaBetasDevBuf.value().data(),
191  scratch.tMaxErrorAlphaBetasDevBuf.value().data(),
192  scratch.accTimeMaxDevBuf.value().data(),
193  scratch.accTimeWgtDevBuf.value().data(),
194  scratch.sumAAsNullHypotDevBuf.value().data(),
195  scratch.sum0sNullHypotDevBuf.value().data(),
196  scratch.chi2sNullHypotDevBuf.value().data(),
197  scratch.tcStateDevBuf.value().data(),
198  scratch.ampMaxAlphaBetaDevBuf.value().data(),
199  scratch.ampMaxErrorDevBuf.value().data(),
200  scratch.timeMaxDevBuf.value().data(),
201  scratch.timeErrorDevBuf.value().data(),
202  paramsDev.const_view());
203 
204  auto const threads_timecorr = 32;
205  auto const blocks_timecorr = cms::alpakatools::divide_up_by(totalChannels, threads_timecorr);
206  auto workDivCorrFinal1D = cms::alpakatools::make_workdiv<Acc1D>(blocks_timecorr, threads_timecorr);
207  alpaka::exec<Acc1D>(queue,
208  workDivCorrFinal1D,
210  digisDevEB.const_view(),
211  digisDevEE.const_view(),
212  uncalibRecHitsDevEB.view(),
213  uncalibRecHitsDevEE.view(),
214  conditionsDev.const_view(),
215  scratch.timeMaxDevBuf.value().data(),
216  scratch.timeErrorDevBuf.value().data(),
217  configParams.timeConstantTermEB,
218  configParams.timeConstantTermEE,
219  configParams.timeNconstEB,
220  configParams.timeNconstEE,
221  configParams.amplitudeThreshEB,
222  configParams.amplitudeThreshEE,
223  configParams.outOfTimeThreshG12pEB,
224  configParams.outOfTimeThreshG12pEE,
225  configParams.outOfTimeThreshG12mEB,
226  configParams.outOfTimeThreshG12mEE,
227  configParams.outOfTimeThreshG61pEB,
228  configParams.outOfTimeThreshG61pEE,
229  configParams.outOfTimeThreshG61mEB,
230  configParams.outOfTimeThreshG61mEE);
231  }
232  }
233 
234 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::ecal::multifit
void minimization_procedure(Queue &queue, InputProduct const &digisDevEB, InputProduct const &digisDevEE, OutputProduct &uncalibRecHitsDevEB, OutputProduct &uncalibRecHitsDevEE, EventDataForScratchDevice &scratch, EcalMultifitConditionsDevice const &conditionsDev, ConfigurationParameters const &configParams, uint32_t const totalChannels)
constexpr Idx divide_up_by(Idx value, Idx divisor)
Definition: workdivision.h:19
PortableCollection< EcalMultifitParametersSoA > EcalMultifitParametersDevice
Vec< Dim2D > Vec2D
Definition: config.h:26
PortableCollection< EcalMultifitConditionsSoA > EcalMultifitConditionsDevice
Eigen::Matrix< data_type, SampleVectorSize, 1 > SampleVector
Eigen::Matrix< char, SampleVectorSize, 1 > SampleGainVector
void launchKernels(Queue &queue, InputProduct const &digisDevEB, InputProduct const &digisDevEE, OutputProduct &uncalibRecHitsDevEB, OutputProduct &uncalibRecHitsDevEE, EcalMultifitConditionsDevice const &conditionsDev, EcalMultifitParametersDevice const &paramsDev, ConfigurationParameters const &configParams)
EcalUncalibratedRecHitDeviceCollection OutputProduct
Eigen::Matrix< char, SampleVectorSize, 1 > BXVectorType
static constexpr int MAXSAMPLES
Definition: EcalDataFrame.h:48