6 #include <alpaka/alpaka.hpp> 16 template <
typename TrackerTraits>
18 template <
typename TrackerTraits>
20 template <
typename TrackerTraits>
26 template <
int N,
typename TrackerTraits>
29 template <
typename TAcc,
typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
35 typename TrackerTraits::tindex_type *__restrict__
ptkids,
36 double *__restrict__
phits,
57 #ifdef BROKENLINE_DEBUG 58 const uint32_t
threadIdx(alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc)[0u]);
61 printf(
"%d Ntuple of size %d/%d for %d hits to fit\n",
totTK,
nHitsL,
nHitsH, hitsInFit);
66 auto tuple_idx = local_idx +
offset;
67 if ((
int)tuple_idx >=
totTK) {
87 auto &&
done = alpaka::declareSharedVar<int, __COUNTER__>(acc);
89 alpaka::syncBlockThreads(acc);
100 auto dx =
hh[hitId[hitsInFit - 1]].xGlobal() -
hh[hitId[0]].xGlobal();
101 auto dy =
hh[hitId[hitsInFit - 1]].yGlobal() -
hh[hitId[0]].yGlobal();
102 auto dz =
hh[hitId[hitsInFit - 1]].zGlobal() -
hh[hitId[0]].zGlobal();
108 for (uint32_t
i = 0;
i < hitsInFit; ++
i) {
110 if (hitsInFit - 1 ==
i)
124 dp.frame.rotation().multiply(
dx,
dy,
dz, ux, uy, uz);
127 int(cb * (
float(phase1PixelTopology::pixelThickess) /
float(phase1PixelTopology::pixelPitchY)) * 8.
f) - 4;
131 bin = std::clamp(
bin, low_value, high_value);
132 float yerr =
dp.sigmay[
bin] * 1.e-4
f;
133 yerr *=
dp.yfact[qbin];
136 yerr = nok ?
hh[
hit].yerrLocal() : yerr;
137 dp.frame.toGlobal(
hh[
hit].xerrLocal(), 0, yerr, ge);
145 printf(
"Track id %d %d Hit %d on %d\nGlobal: hits.col(%d) << %f,%f,%f\n",
149 hh[
hit].detectorIndex(),
154 printf(
"Error: hits_ge.col(%d) << %e,%e,%e,%e,%e,%e\n",
i, ge[0], ge[1], ge[2], ge[3], ge[4], ge[5]);
159 hits_ge.col(
i) << ge[0], ge[1], ge[2], ge[3], ge[4], ge[5];
163 #ifdef BROKENLINE_DEBUG 174 template <
int N,
typename TrackerTraits>
177 template <
typename TAcc,
typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
182 typename TrackerTraits::tindex_type
const *__restrict__
ptkids,
183 double *__restrict__
phits,
198 auto tkid =
ptkids[local_idx];
221 #ifdef BROKENLINE_DEBUG 222 if (!(circle.
chi2 >= 0) || !(
line.chi2 >= 0))
223 printf(
"kernelBLFit failed! %f/%f\n", circle.
chi2,
line.chi2);
224 printf(
"kernelBLFit size %d for %d hits circle.par(0,1,2): %d %f,%f,%f\n",
231 printf(
"kernelBLHits line.par(0,1): %d %f,%f\n", tkid,
line.par(0),
line.par(1));
232 printf(
"kernelBLHits chi2 cov %f/%f %e,%e,%e,%e,%e\n",
245 template <
typename TrackerTraits>
250 uint32_t maxNumberOfTuples,
254 uint32_t blockSize = 64;
256 const WorkDiv1D workDivTriplets = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
257 const WorkDiv1D workDivQuadsPenta = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks / 4, blockSize);
261 cms::alpakatools::make_device_buffer<typename TrackerTraits::tindex_type[]>(
queue, maxNumberOfConcurrentFits_);
262 auto hitsDevice = cms::alpakatools::make_device_buffer<double[]>(
264 auto hits_geDevice = cms::alpakatools::make_device_buffer<float[]>(
266 auto fast_fit_resultsDevice = cms::alpakatools::make_device_buffer<double[]>(
269 for (uint32_t
offset = 0;
offset < maxNumberOfTuples;
offset += maxNumberOfConcurrentFits_) {
272 alpaka::exec<Acc1D>(
queue,
281 hits_geDevice.data(),
282 fast_fit_resultsDevice.data(),
287 alpaka::exec<Acc1D>(
queue,
295 hits_geDevice.data(),
296 fast_fit_resultsDevice.data());
300 riemannFit::rolling_fits<4, TrackerTraits::maxHitsOnTrack, 1>([
this,
306 &fast_fit_resultsDevice,
309 &workDivQuadsPenta](
auto i) {
310 alpaka::exec<Acc1D>(
queue,
319 hits_geDevice.data(),
320 fast_fit_resultsDevice.data(),
325 alpaka::exec<Acc1D>(
queue,
333 hits_geDevice.data(),
334 fast_fit_resultsDevice.data());
338 riemannFit::rolling_fits<4, TrackerTraits::maxHitsOnTrackForFullFit, 1>([
this,
344 &fast_fit_resultsDevice,
347 &workDivQuadsPenta](
auto i) {
348 alpaka::exec<Acc1D>(
queue,
357 hits_geDevice.data(),
358 fast_fit_resultsDevice.data(),
363 alpaka::exec<Acc1D>(
queue,
371 hits_geDevice.data(),
372 fast_fit_resultsDevice.data());
375 static_assert(TrackerTraits::maxHitsOnTrackForFullFit < TrackerTraits::maxHitsOnTrack);
378 alpaka::exec<Acc1D>(
queue,
387 hits_geDevice.data(),
388 fast_fit_resultsDevice.data(),
389 TrackerTraits::maxHitsOnTrackForFullFit,
390 TrackerTraits::maxHitsOnTrack - 1,
393 alpaka::exec<Acc1D>(
queue,
401 hits_geDevice.data(),
402 fast_fit_resultsDevice.data());
ALPAKA_FN_ACC ALPAKA_FN_INLINE void fastFit(const TAcc &acc, const M3xN &hits, V4 &result)
A very fast helix fit.
constexpr DetParams const &__restrict__ detParams(int i) const
void launchBrokenLineKernels(const HitConstView &hv, ParamsOnDevice const *cpeParams, uint32_t nhits, uint32_t maxNumberOfTuples, Queue &queue)
constexpr int kGenErrorQBins
Vector3d par
parameter: (X0,Y0,R)
typename reco::TrackSoA< TrackerTraits >::HitContainer Tuples
TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity
ALPAKA_FN_ACC void operator()(TAcc const &acc, TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity, double bField, OutputSoAView< TrackerTraits > results_view, typename TrackerTraits::tindex_type const *__restrict__ ptkids, double *__restrict__ phits, float *__restrict__ phits_ge, double *__restrict__ pfast_fit) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE void lineFit(const TAcc &acc, const M6xN &hits_ge, const V4 &fast_fit, const double bField, const PreparedBrokenLineData< n > &data, riemannFit::LineFit &line_results)
Performs the Broken Line fit in the straight track case (that is, the fit parameters are only the int...
__host__ __device__ void prepareBrokenLineData(const M3xN &hits, const V4 &fast_fit, const double bField, PreparedBrokenLineData< n > &results)
Computes the data needed for the Broken Line fit procedure that are mainly common for the circle and ...
Eigen::Matrix< float, 6, N > Matrix6xNf
WorkDiv< Dim1D > WorkDiv1D
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ float *__restrict__ double *__restrict__ uint32_t nHitsL
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ phits
Eigen::Map< Vector4d, 0, Eigen::InnerStride< stride > > Map4d
Eigen::Map< Matrix3xNd< N >, 0, Eigen::Stride< 3 *stride, stride > > Map3xNd
ALPAKA_FN_ACC void operator()(TAcc const &acc, Tuples< TrackerTraits > const *__restrict__ foundNtuplets, TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity, TrackingRecHitSoAConstView< TrackerTraits > hh, pixelCPEforDevice::ParamsOnDeviceT< TrackerTraits > const *__restrict__ cpeParams, typename TrackerTraits::tindex_type *__restrict__ ptkids, double *__restrict__ phits, float *__restrict__ phits_ge, double *__restrict__ pfast_fit, uint32_t nHitsL, uint32_t nHitsH, int32_t offset) const
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ ptkids
std::vector< Block > Blocks
Abs< T >::type abs(const T &t)
Eigen::Matrix< double, 3, N > Matrix3xNd
double OutputSoAView< TrackerTraits > results_view
Eigen::Map< Matrix6xNf< N >, 0, Eigen::Stride< 6 *stride, stride > > Map6xNf
ALPAKA_FN_ACC ALPAKA_FN_INLINE void uint32_t const uint32_t CACellT< TrackerTraits > uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > hh
ALPAKA_FN_ACC ALPAKA_FN_INLINE void const M3xN const V4 & fast_fit
constexpr auto invalidTkId
typename TrackingRecHitSoA< TrackerTraits >::template TrackingRecHitSoALayout<>::ConstView TrackingRecHitSoAConstView
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ float *__restrict__ double *__restrict__ uint32_t uint32_t nHitsH
auto const & foundNtuplets
char data[epos_bytes_allocation]
data needed for the Broken Line fit procedure.
static constexpr __host__ __device__ void copyFromCircle(TrackSoAView &tracks, V3 const &cp, M3 const &ccov, V2 const &lp, M2 const &lcov, float b, int32_t i)
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
ALPAKA_FN_ACC ALPAKA_FN_INLINE void circleFit(const TAcc &acc, const M3xN &hits, const M6xN &hits_ge, const V4 &fast_fit, const double bField, PreparedBrokenLineData< n > &data, karimaki_circle_fit &circle_results)
Performs the Broken Line fit in the curved track case (that is, the fit parameters are the intercepti...
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ float *__restrict__ phits_ge
ALPAKA_ASSERT_ACC(offsets)
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ float *__restrict__ double *__restrict__ pfast_fit
T1 atomicAdd(T1 *a, T2 b)
constexpr int kNumErrorBins
reco::TrackSoAView< TrackerTraits > OutputSoAView
constexpr uint32_t maxNumberOfConcurrentFits
typename reco::TrackSoA< TrackerTraits >::template Layout<>::View TrackSoAView