Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-07-03 00:42:42

0001 #ifndef RecoTracker_PixelSeeding_plugins_alpaka_CAHitNtupletGeneratorKernelsImpl_h
0002 #define RecoTracker_PixelSeeding_plugins_alpaka_CAHitNtupletGeneratorKernelsImpl_h
0003 
0004 // #define GPU_DEBUG
0005 // #define NTUPLE_DEBUG
0006 // #define CA_DEBUG
0007 // #define CA_WARNINGS
0008 
0009 // C++ includes
0010 #include <cmath>
0011 #include <cstdint>
0012 #include <cstdio>
0013 #include <limits>
0014 #include <type_traits>
0015 
0016 // Alpaka includes
0017 #include <alpaka/alpaka.hpp>
0018 
0019 // CMSSW includes
0020 #include "DataFormats/TrackSoA/interface/TrackDefinitions.h"
0021 #include "DataFormats/TrackSoA/interface/TracksSoA.h"
0022 #include "DataFormats/TrackSoA/interface/alpaka/TrackUtilities.h"
0023 #include "HeterogeneousCore/AlpakaInterface/interface/AtomicPairCounter.h"
0024 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0025 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0026 #include "FWCore/Utilities/interface/isFinite.h"
0027 #include "RecoTracker/PixelSeeding/interface/CAPairSoA.h"
0028 
0029 // local includes
0030 #include "CACell.h"
0031 #include "CAHitNtupletGeneratorKernels.h"
0032 #include "CAStructures.h"
0033 
0034 namespace ALPAKA_ACCELERATOR_NAMESPACE::caHitNtupletGeneratorKernels {
0035 
0036   using namespace ::caStructures;
0037 
0038   constexpr uint32_t tkNotFound = std::numeric_limits<uint32_t>::max();
0039   constexpr float maxScore = std::numeric_limits<float>::max();
0040   constexpr float nSigma2 = 25.f;
0041 
0042   // all of these below are mostly to avoid carrying around the relative namespace
0043 
0044   using Quality = ::pixelTrack::Quality;
0045   using TkSoAView = ::reco::TrackSoAView;
0046   using TkHitSoAView = ::reco::TrackHitSoAView;
0047 
0048   template <typename TrackerTraits>
0049   using QualityCuts = ::pixelTrack::QualityCutsT<TrackerTraits>;
0050 
0051   using Counters = caHitNtupletGenerator::Counters;
0052   using HitToTuple = caStructures::GenericContainer;
0053   using HitContainer = caStructures::SequentialContainer;
0054   using TupleMultiplicity = caStructures::GenericContainer;
0055   using HitToCell = caStructures::GenericContainer;
0056   using CellToCell = caStructures::GenericContainer;
0057   using CellToTrack = caStructures::GenericContainer;
0058 
0059   using namespace cms::alpakatools;
0060 
0061   class SetHitsLayerStart {
0062   public:
0063     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0064                                   const reco::HitModuleSoAConstView &mm,
0065                                   const reco::CALayersSoAConstView &ll,
0066                                   uint32_t *__restrict__ hitsLayerStart) const {
0067       ALPAKA_ASSERT_ACC(0 == mm.moduleStart()[0]);
0068 
0069       for (int32_t i : cms::alpakatools::uniform_elements(acc, ll.metadata().size())) {
0070         hitsLayerStart[i] = mm.moduleStart()[ll.layerStarts()[i]];
0071 #ifdef GPU_DEBUG
0072         int old = i == 0 ? 0 : mm.moduleStart()[ll.layerStarts()[i - 1]];
0073         printf("LayerStart %d/%d at module %d: %d - %d\n",
0074                i,
0075                ll.metadata().size() - 1,
0076                ll.layerStarts()[i],
0077                hitsLayerStart[i],
0078                hitsLayerStart[i] - old);
0079 #endif
0080       }
0081     }
0082   };
0083 
0084   class Kernel_printSizes {
0085   public:
0086     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0087                                   HitsConstView hh,
0088                                   TkSoAView tt,
0089                                   uint32_t const *__restrict__ nCells,
0090                                   uint32_t const *__restrict__ nTrips,
0091                                   uint32_t const *__restrict__ nCellTracks) const {
0092       if (cms::alpakatools::once_per_grid(acc))
0093         printf("nSizes:%d;%d;%d;%d;%d;%d;%d\n",
0094                hh.metadata().size(),
0095                hh.metadata().size() - hh.offsetBPIX2(),
0096                *nCells,
0097                *nTrips,
0098                *nCellTracks,
0099                tt.nTracks(),
0100                tt.metadata().size());
0101     }
0102   };
0103 
0104   template <typename TrackerTraits>
0105   class Kernel_checkOverflows {
0106   public:
0107     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0108                                   TkSoAView tracks_view,
0109                                   HitContainer const *__restrict__ foundNtuplets,
0110                                   TupleMultiplicity const *tupleMultiplicity,
0111                                   HitToTuple const *hitToTuple,
0112                                   cms::alpakatools::AtomicPairCounter *apc,
0113                                   CACell<TrackerTraits> const *__restrict__ cells,
0114                                   uint32_t const *__restrict__ nCells,
0115                                   uint32_t const *__restrict__ nTrips,
0116                                   uint32_t const *__restrict__ nCellTracks,
0117                                   caStructures::CAPairSoAConstView cellCell,
0118                                   caStructures::CAPairSoAConstView cellTrack,
0119                                   int32_t nHits,
0120                                   uint32_t maxNumberOfDoublets,
0121                                   AlgoParams const &params,
0122                                   Counters *counters) const {
0123       auto &c = *counters;
0124       // counters once per event
0125       if (cms::alpakatools::once_per_grid(acc)) {
0126         alpaka::atomicAdd(acc, &c.nEvents, 1ull, alpaka::hierarchy::Blocks{});
0127         alpaka::atomicAdd(acc, &c.nHits, static_cast<unsigned long long>(nHits), alpaka::hierarchy::Blocks{});
0128         alpaka::atomicAdd(acc, &c.nCells, static_cast<unsigned long long>(*nCells), alpaka::hierarchy::Blocks{});
0129         alpaka::atomicAdd(
0130             acc, &c.nTuples, static_cast<unsigned long long>(apc->get().first), alpaka::hierarchy::Blocks{});
0131         alpaka::atomicAdd(acc,
0132                           &c.nFitTracks,
0133                           static_cast<unsigned long long>(tupleMultiplicity->size()),
0134                           alpaka::hierarchy::Blocks{});
0135       }
0136 
0137 #ifdef NTUPLE_DEBUGS
0138       if (cms::alpakatools::once_per_grid(acc)) {
0139         printf("number of found cells %d \n found tuples %d with total hits %d out of %d\n",
0140                *nCells,
0141                apc->get().first,
0142                apc->get().second,
0143                nHits);
0144         if (apc->get().first < tracks_view.metadata().size()) {
0145           ALPAKA_ASSERT_ACC(foundNtuplets->size(apc->get().first) == 0);
0146           ALPAKA_ASSERT_ACC(foundNtuplets->size() == apc->get().second);
0147         }
0148       }
0149 
0150       for (auto idx : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes())) {
0151         if (foundNtuplets->size(idx) > TrackerTraits::maxHitsOnTrack)  // current real limit
0152           printf("ERROR %d, %d\n", idx, foundNtuplets->size(idx));
0153         ALPAKA_ASSERT_ACC(foundNtuplets->size(idx) <= TrackerTraits::maxHitsOnTrack);
0154         for (auto ih = foundNtuplets->begin(idx); ih != foundNtuplets->end(idx); ++ih)
0155           ALPAKA_ASSERT_ACC(int(*ih) < nHits);
0156       }
0157 #endif
0158 
0159       if (cms::alpakatools::once_per_grid(acc)) {
0160         if (apc->get().first >= uint32_t(tracks_view.metadata().size()))
0161           printf("Tuples overflow\n");
0162         if (*nCells >= maxNumberOfDoublets)
0163           printf("Cells overflow\n");
0164         if (*nTrips >= uint32_t(cellCell.metadata().size()))
0165           printf("Triplets overflow\n");
0166         if (*nCellTracks >= uint32_t(cellTrack.metadata().size()))
0167           printf("TracksToCell overflow\n");
0168       }
0169 
0170       for (auto idx : cms::alpakatools::uniform_elements(acc, *nCells)) {
0171         auto const &thisCell = cells[idx];
0172         if (thisCell.hasFishbone() && !thisCell.isKilled())
0173           alpaka::atomicAdd(acc, &c.nFishCells, 1ull, alpaka::hierarchy::Blocks{});
0174         if (thisCell.isKilled())
0175           alpaka::atomicAdd(acc, &c.nKilledCells, 1ull, alpaka::hierarchy::Blocks{});
0176         if (!thisCell.unused())
0177           alpaka::atomicAdd(acc, &c.nEmptyCells, 1ull, alpaka::hierarchy::Blocks{});
0178         if ((0 == hitToTuple->size(thisCell.inner_hit_id())) && (0 == hitToTuple->size(thisCell.outer_hit_id())))
0179           alpaka::atomicAdd(acc, &c.nZeroTrackCells, 1ull, alpaka::hierarchy::Blocks{});
0180       }
0181     }
0182   };
0183 
0184   template <typename TrackerTraits>
0185   class Kernel_fishboneCleaner {
0186   public:
0187     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0188                                   CACell<TrackerTraits> const *cells,
0189                                   uint32_t const *__restrict__ nCells,
0190                                   CellToTrack const *__restrict__ cellTracksHisto,
0191                                   TkSoAView tracks_view) const {
0192       constexpr auto reject = Quality::dup;
0193 
0194       for (auto idx : cms::alpakatools::uniform_elements(acc, *nCells)) {
0195         auto const &thisCell = cells[idx];
0196         if (!thisCell.isKilled())
0197           continue;
0198 
0199         auto const *__restrict__ tracksOfCell = cellTracksHisto->begin(idx);
0200         for (auto i = 0u; i < cellTracksHisto->size(idx); i++)
0201           tracks_view[tracksOfCell[i]].quality() = reject;
0202       }
0203     }
0204   };
0205 
0206   // remove shorter tracks if sharing a cell
0207   // It does not seem to affect efficiency in any way!
0208   template <typename TrackerTraits>
0209   class Kernel_earlyDuplicateRemover {
0210   public:
0211     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0212                                   CACell<TrackerTraits> const *cells,
0213                                   uint32_t const *__restrict__ nCells,
0214                                   CellToTrack const *__restrict__ cellTracksHisto,
0215                                   TkSoAView tracks_view,
0216                                   bool dupPassThrough) const {
0217       // quality to mark rejected
0218       constexpr auto reject = Quality::edup;  /// cannot be loose
0219       ALPAKA_ASSERT_ACC(nCells);
0220       for (auto idx : cms::alpakatools::uniform_elements(acc, *nCells)) {
0221 #ifdef CA_SIZES
0222         printf("cellTracksSizes;%d;%d;%d\n", idx, cT.size(), cT.capacity());
0223 #endif
0224         if (cellTracksHisto->size(idx) < 2)
0225           continue;
0226 
0227         int8_t maxNl = 0;
0228         auto const *__restrict__ tracksOfCell = cellTracksHisto->begin(idx);
0229 
0230         // find maxNl
0231         for (auto i = 0u; i < cellTracksHisto->size(idx); i++) {
0232           if (int(tracksOfCell[i]) > tracks_view.metadata().size())
0233             printf(">WARNING: %d %d %d %d\n", idx, i, int(tracksOfCell[i]), tracks_view.metadata().size());
0234           auto nl = tracks_view[tracksOfCell[i]].nLayers();
0235           maxNl = std::max(nl, maxNl);
0236         }
0237 
0238         // if (maxNl<4) continue;
0239         // quad pass through (leave it here for tests)
0240         //  maxNl = std::min(4, maxNl);
0241 
0242         for (auto i = 0u; i < cellTracksHisto->size(idx); i++) {
0243           auto it = tracksOfCell[i];
0244 
0245           if (int(it) > tracks_view.metadata().size())
0246             printf(">WARNING: %d %d %d\n", i, it, tracks_view.metadata().size());
0247           if (tracks_view[it].nLayers() < maxNl)
0248             tracks_view[it].quality() = reject;  // no race: simple assignment of the same constant
0249         }
0250       }
0251     }
0252   };
0253 
0254   // assume the above (so, short tracks already removed)
0255   template <typename TrackerTraits>
0256   class Kernel_fastDuplicateRemover {
0257   public:
0258     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0259                                   CACell<TrackerTraits> const *__restrict__ cells,
0260                                   uint32_t const *__restrict__ nCells,
0261                                   CellToTrack const *__restrict__ cellTracksHisto,
0262                                   TkSoAView tracks_view,
0263                                   bool dupPassThrough) const {
0264       // quality to mark rejected
0265       auto const reject = dupPassThrough ? Quality::loose : Quality::dup;
0266       constexpr auto loose = Quality::loose;
0267 
0268       ALPAKA_ASSERT_ACC(nCells);
0269       const auto ntNCells = (*nCells);
0270 
0271       for (auto idx : cms::alpakatools::uniform_elements(acc, ntNCells)) {
0272         if (cellTracksHisto->size(idx) < 2)
0273           continue;
0274 
0275         float mc = maxScore;
0276         uint32_t im = tkNotFound;
0277 
0278         auto score = [&](auto it) { return std::abs(reco::tip(tracks_view, it)); };
0279 
0280         // full crazy combinatorics
0281         auto const *__restrict__ thisCellTracks = cellTracksHisto->begin(idx);
0282         int ntr = cellTracksHisto->size(idx);
0283         for (int i = 0; i < ntr - 1; i++) {
0284           auto it = thisCellTracks[i];
0285           auto qi = tracks_view[it].quality();
0286           if (qi <= reject)
0287             continue;
0288           auto opi = tracks_view[it].state()(2);
0289           auto e2opi = tracks_view[it].covariance()(9);
0290           auto cti = tracks_view[it].state()(3);
0291           auto e2cti = tracks_view[it].covariance()(12);
0292           for (int j = i + 1; j < ntr; ++j) {
0293             auto jt = thisCellTracks[j];
0294             auto qj = tracks_view[jt].quality();
0295             if (qj <= reject)
0296               continue;
0297             auto opj = tracks_view[jt].state()(2);
0298             auto ctj = tracks_view[jt].state()(3);
0299             auto dct = nSigma2 * (tracks_view[jt].covariance()(12) + e2cti);
0300             if ((cti - ctj) * (cti - ctj) > dct)
0301               continue;
0302             auto dop = nSigma2 * (tracks_view[jt].covariance()(9) + e2opi);
0303             if ((opi - opj) * (opi - opj) > dop)
0304               continue;
0305             if ((qj < qi) || (qj == qi && score(it) < score(jt)))
0306               tracks_view[jt].quality() = reject;
0307             else {
0308               tracks_view[it].quality() = reject;
0309               break;
0310             }
0311           }
0312         }
0313 
0314         // find maxQual
0315         auto maxQual = reject;  // no duplicate!
0316         for (int i = 0; i < ntr; i++) {
0317           auto it = thisCellTracks[i];
0318           if (tracks_view[it].quality() > maxQual)
0319             maxQual = tracks_view[it].quality();
0320         }
0321 
0322         if (maxQual <= loose)
0323           continue;
0324 
0325         // find min score
0326         for (int i = 0; i < ntr; i++) {
0327           auto it = thisCellTracks[i];
0328           if (tracks_view[it].quality() == maxQual && score(it) < mc) {
0329             mc = score(it);
0330             im = it;
0331           }
0332         }
0333 
0334         if (tkNotFound == im)
0335           continue;
0336 
0337         // mark all other duplicates  (not yet, keep it loose)
0338         for (int i = 0; i < ntr; i++) {
0339           auto it = thisCellTracks[i];
0340           if (tracks_view[it].quality() > loose && it != im)
0341             tracks_view[it].quality() = loose;  //no race:  simple assignment of the same constant
0342         }
0343       }
0344     }
0345   };
0346 
0347   template <typename TrackerTraits>
0348   class Kernel_connect {
0349   public:
0350     ALPAKA_FN_ACC void operator()(Acc2D const &acc,
0351                                   cms::alpakatools::AtomicPairCounter *apc,  // just to zero them
0352                                   HitsConstView hh,
0353                                   reco::CALayersSoAConstView ll,
0354                                   caStructures::CAPairSoAView cn,
0355                                   CACell<TrackerTraits> *cells,
0356                                   uint32_t const *nCells,
0357                                   uint32_t *nTrips,
0358                                   HitToCell const *__restrict__ outerHitHisto,
0359                                   CellToCell *cellNeighborsHisto,
0360                                   AlgoParams const &params) const {
0361       using Cell = CACell<TrackerTraits>;
0362       uint32_t maxTriplets = cn.metadata().size();
0363 
0364       if (cms::alpakatools::once_per_grid(acc)) {
0365         *apc = 0;
0366       }  // ready for next kernel
0367 
0368       // loop on outer cells
0369       for (uint32_t cellIndex : cms::alpakatools::uniform_elements_y(acc, *nCells)) {
0370         auto &thisCell = cells[cellIndex];
0371         auto innerHitId = thisCell.inner_hit_id() - hh.offsetBPIX2();
0372 
0373         if (int(innerHitId) < 0)
0374           continue;
0375 
0376         auto const *__restrict__ outerHitCells = outerHitHisto->begin(innerHitId);
0377         auto const numberOfPossibleNeighbors = outerHitHisto->size(innerHitId);
0378 
0379 #ifdef CA_DEBUG
0380         printf("numberOfPossibleFromHisto;%d;%d;%d;%d;%d\n",
0381                *nCells,
0382                innerHitId,
0383                cellIndex,
0384                thisCell.innerLayer(),
0385                numberOfPossibleNeighbors);
0386 #endif
0387         auto ri = thisCell.inner_r(hh);
0388         auto zi = thisCell.inner_z(hh);
0389         auto ro = thisCell.outer_r(hh);
0390         auto zo = thisCell.outer_z(hh);
0391         auto thetaCut = ll[thisCell.innerLayer()].caThetaCut();
0392 
0393         // loop on inner cells
0394         for (uint32_t j : cms::alpakatools::independent_group_elements_x(acc, numberOfPossibleNeighbors)) {
0395           auto otherCell = outerHitCells[j];
0396           auto &oc = cells[otherCell];
0397           auto r1 = oc.inner_r(hh);
0398           auto z1 = oc.inner_z(hh);
0399           auto dcaCut = ll[oc.innerLayer()].caDCACut();
0400           bool aligned = Cell::areAlignedRZ(r1, z1, ri, zi, ro, zo, params.ptmin_, thetaCut);
0401           if (aligned && thisCell.dcaCut(hh, oc, dcaCut, params.hardCurvCut_)) {
0402             auto t_ind = alpaka::atomicAdd(acc, nTrips, 1u, alpaka::hierarchy::Blocks{});
0403 #ifdef CA_DEBUG
0404             printf("Triplet no. %d %.5f %.5f (%d %d) - %d %d -> (%d, %d, %d, %d) \n",
0405                    t_ind,
0406                    thetaCut,
0407                    dcaCut,
0408                    thisCell.layerPairId(),
0409                    oc.layerPairId(),
0410                    otherCell,
0411                    cellIndex,
0412                    thisCell.inner_hit_id(),
0413                    thisCell.outer_hit_id(),
0414                    oc.inner_hit_id(),
0415                    oc.outer_hit_id());
0416 #endif
0417 
0418 #ifdef CA_DEBUG
0419             printf("filling cell no. %d %d: %d -> %d\n", t_ind, cellNeighborsHisto->size(), otherCell, cellIndex);
0420 #endif
0421 
0422             if (t_ind >= maxTriplets) {
0423 #ifdef CA_WARNINGS
0424               printf("Warning!!!! Too many cell->cell (triplets) associations (limit = %d)!\n", cn.metadata().size());
0425 #endif
0426               alpaka::atomicSub(acc, nTrips, 1u, alpaka::hierarchy::Blocks{});
0427               break;
0428             }
0429 
0430             cellNeighborsHisto->count(acc, otherCell);
0431 
0432             cn[t_ind].inner() = otherCell;
0433             cn[t_ind].outer() = cellIndex;
0434             thisCell.setStatusBits(Cell::StatusBit::kUsed);
0435             oc.setStatusBits(Cell::StatusBit::kUsed);
0436           }
0437 
0438         }  // loop on inner cells
0439       }  // loop on outer cells
0440     }
0441   };
0442 
0443   template <typename TrackerTraits>
0444   class FillDoubletsHisto {
0445   public:
0446     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0447                                   CACell<TrackerTraits> const *__restrict__ cells,
0448                                   uint32_t *nCells,
0449                                   uint32_t offsetBPIX2,
0450                                   HitToCell *outerHitHisto) const {
0451       for (auto cellIndex : cms::alpakatools::uniform_elements(acc, *nCells)) {
0452 #ifdef DOUBLETS_DEBUG
0453         printf("outerHitHisto;%d;%d\n", cellIndex, cells[cellIndex].outer_hit_id());
0454 #endif
0455         outerHitHisto->fill(acc, cells[cellIndex].outer_hit_id() - offsetBPIX2, cellIndex);
0456       }
0457     }
0458   };
0459 
0460   class Kernel_fillGenericPair {
0461   public:
0462     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0463                                   caStructures::CAPairSoAConstView cn,
0464                                   uint32_t const *nElements,
0465                                   GenericContainer *genericHisto) const {
0466       for (uint32_t index : cms::alpakatools::uniform_elements(acc, *nElements)) {
0467         genericHisto->fill(acc, cn[index].inner(), cn[index].outer());
0468       }
0469     }
0470   };
0471 
0472   template <typename TrackerTraits>
0473   class Kernel_find_ntuplets {
0474   public:
0475     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0476                                   const ::reco::CAGraphSoAConstView &cc,
0477                                   TkSoAView tracks_view,
0478                                   HitContainer *foundNtuplets,
0479                                   CellToCell const *__restrict__ cellNeighborsHisto,
0480                                   CellToTrack *cellTracksHisto,
0481                                   caStructures::CAPairSoAView ct,
0482                                   CACell<TrackerTraits> *__restrict__ cells,
0483                                   uint32_t *nCellTracks,
0484                                   uint32_t const *nTriplets,
0485                                   uint32_t const *nCells,
0486                                   cms::alpakatools::AtomicPairCounter *apc,
0487                                   AlgoParams const &params) const {
0488       using Cell = CACell<TrackerTraits>;
0489 
0490 #ifdef GPU_DEBUG
0491       if (cms::alpakatools::once_per_grid(acc))
0492         printf("starting producing ntuplets from %d cells and %d triplets \n", *nCells, *nTriplets);
0493 #endif
0494 
0495       for (auto idx : cms::alpakatools::uniform_elements(acc, (*nCells))) {
0496         auto const &thisCell = cells[idx];
0497 
0498         // cut by earlyFishbone
0499         if (thisCell.isKilled())
0500           continue;
0501 
0502         // we require at least three hits
0503         if (cellNeighborsHisto->size(idx) == 0)
0504           continue;
0505 
0506         auto pid = thisCell.layerPairId();
0507         bool doit = cc[pid].startingPair();
0508 
0509         constexpr uint32_t maxDepth = TrackerTraits::maxDepth;
0510 #ifdef CA_DEBUG
0511         printf("LayerPairId %d doit ? %d From cell %d with nNeighbors = %d \n",
0512                pid,
0513                doit,
0514                idx,
0515                cellNeighborsHisto->size(idx));
0516 #endif
0517 
0518         if (doit) {
0519           typename Cell::TmpTuple stack;
0520 
0521           stack.reset();
0522           thisCell.template find_ntuplets<maxDepth>(acc,
0523                                                     cc,
0524                                                     cells,
0525                                                     *foundNtuplets,
0526                                                     cellNeighborsHisto,
0527                                                     cellTracksHisto,
0528                                                     nCellTracks,
0529                                                     ct,
0530                                                     *apc,
0531                                                     tracks_view.quality(),
0532                                                     stack,
0533                                                     params.minHitsPerNtuplet_);
0534           ALPAKA_ASSERT_ACC(stack.empty());
0535         }
0536       }
0537     }
0538   };
0539 
0540   template <typename TrackerTraits>
0541   class Kernel_mark_used {
0542   public:
0543     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0544                                   CACell<TrackerTraits> *__restrict__ cells,
0545                                   CellToTrack const *__restrict__ cellTracksHisto,
0546                                   uint32_t const *nCells) const {
0547       using Cell = CACell<TrackerTraits>;
0548       for (auto idx : cms::alpakatools::uniform_elements(acc, (*nCells))) {
0549         auto &thisCell = cells[idx];
0550         if (cellTracksHisto->size(idx) > 0)
0551           thisCell.setStatusBits(Cell::StatusBit::kInTrack);
0552       }
0553     }
0554   };
0555 
0556   template <typename TrackerTraits>
0557   class Kernel_countMultiplicity {
0558   public:
0559     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0560                                   TkSoAView tracks_view,
0561                                   HitContainer const *__restrict__ foundNtuplets,
0562                                   TupleMultiplicity *tupleMultiplicity) const {
0563       for (auto it : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes())) {
0564         auto nhits = foundNtuplets->size(it);
0565         // printf("it: %d nhits: %d \n",it,nhits);
0566         if (nhits < 3)
0567           continue;
0568         if (tracks_view[it].quality() == Quality::edup)
0569           continue;
0570         ALPAKA_ASSERT_ACC(tracks_view[it].quality() == Quality::bad);
0571         if (nhits > TrackerTraits::maxHitsOnTrack)  // current limit
0572           printf("wrong mult %d %d\n", it, nhits);
0573         ALPAKA_ASSERT_ACC(nhits <= TrackerTraits::maxHitsOnTrack);
0574         tupleMultiplicity->count(acc, nhits);
0575       }
0576     }
0577   };
0578 
0579   template <typename TrackerTraits>
0580   class Kernel_fillMultiplicity {
0581   public:
0582     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0583                                   TkSoAView tracks_view,
0584                                   HitContainer const *__restrict__ foundNtuplets,
0585                                   TupleMultiplicity *tupleMultiplicity) const {
0586       for (auto it : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes())) {
0587         auto nhits = foundNtuplets->size(it);
0588 
0589         if (nhits < 3)
0590           continue;
0591         if (tracks_view[it].quality() == Quality::edup)
0592           continue;
0593         ALPAKA_ASSERT_ACC(tracks_view[it].quality() == Quality::bad);
0594         if (nhits > TrackerTraits::maxHitsOnTrack)
0595           printf("wrong mult %d %d\n", it, nhits);
0596         ALPAKA_ASSERT_ACC(nhits <= TrackerTraits::maxHitsOnTrack);
0597         tupleMultiplicity->fill(acc, nhits, it);
0598       }
0599     }
0600   };
0601 
0602   template <typename TrackerTraits>
0603   class Kernel_classifyTracks {
0604   public:
0605     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0606                                   TkSoAView tracks_view,
0607                                   HitContainer const *__restrict__ foundNtuplets,
0608                                   QualityCuts<TrackerTraits> cuts) const {
0609       for (auto it : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes())) {
0610         auto nhits = foundNtuplets->size(it);
0611         if (nhits == 0)
0612           break;  // guard
0613 
0614         // if duplicate: not even fit
0615         if (tracks_view[it].quality() == Quality::edup)
0616           continue;
0617 
0618         ALPAKA_ASSERT_ACC(tracks_view[it].quality() == Quality::bad);
0619 
0620         // mark doublets as bad
0621         if (nhits < 3)
0622           continue;
0623 
0624         // if the fit has any invalid parameters, mark it as bad
0625         bool isNaN = false;
0626         for (int i = 0; i < 5; ++i) {
0627           isNaN |= edm::isNotFinite(tracks_view[it].state()(i));
0628         }
0629         if (isNaN) {
0630 #ifdef NTUPLE_DEBUG
0631           printf("NaN in fit %d size %d chi2 %f\n", it, foundNtuplets->size(it), tracks_view[it].chi2());
0632 #endif
0633           continue;
0634         }
0635 
0636         tracks_view[it].quality() = Quality::strict;
0637 
0638         if (cuts.strictCut(tracks_view, it))
0639           continue;
0640 
0641         tracks_view[it].quality() = Quality::tight;
0642 
0643         if (cuts.isHP(tracks_view, nhits, it))
0644           tracks_view[it].quality() = Quality::highPurity;
0645       }
0646     }
0647   };
0648 
0649   template <typename TrackerTraits>
0650   class Kernel_doStatsForTracks {
0651   public:
0652     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0653                                   TkSoAView tracks_view,
0654                                   HitContainer const *__restrict__ foundNtuplets,
0655                                   Counters *counters) const {
0656       for (auto idx : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes())) {
0657         if (foundNtuplets->size(idx) == 0)
0658           break;  //guard
0659         if (tracks_view[idx].quality() < Quality::loose)
0660           continue;
0661         alpaka::atomicAdd(acc, &(counters->nLooseTracks), 1ull, alpaka::hierarchy::Blocks{});
0662         if (tracks_view[idx].quality() < Quality::strict)
0663           continue;
0664         alpaka::atomicAdd(acc, &(counters->nGoodTracks), 1ull, alpaka::hierarchy::Blocks{});
0665       }
0666     }
0667   };
0668 
0669   template <typename TrackerTraits>
0670   class Kernel_countHitInTracks {
0671   public:
0672     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0673                                   TkSoAView tracks_view,
0674                                   HitContainer const *__restrict__ foundNtuplets,
0675                                   HitToTuple *hitToTuple) const {
0676       for (auto idx : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes())) {
0677         if (foundNtuplets->size(idx) == 0)
0678           break;  // guard
0679         for (auto h = foundNtuplets->begin(idx); h != foundNtuplets->end(idx); ++h)
0680           hitToTuple->count(acc, *h);
0681       }
0682     }
0683   };
0684 
0685   template <typename TrackerTraits>
0686   class Kernel_fillHitInTracks {
0687   public:
0688     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0689                                   TkSoAView tracks_view,
0690                                   HitContainer const *__restrict__ foundNtuplets,
0691                                   HitToTuple *hitToTuple) const {
0692       for (auto idx : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes())) {
0693         if (foundNtuplets->size(idx) == 0)
0694           break;  // guard
0695         for (auto h = foundNtuplets->begin(idx); h != foundNtuplets->end(idx); ++h)
0696           hitToTuple->fill(acc, *h, idx);
0697       }
0698     }
0699   };
0700 
0701   template <typename TrackerTraits>
0702   class Kernel_fillHitDetIndices {
0703   public:
0704     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0705                                   TkSoAView tracks_view,
0706                                   TkHitSoAView track_hits_view,
0707                                   HitContainer const *__restrict__ foundNtuplets,
0708                                   HitsConstView hh) const {
0709       // copy offsets
0710       for (auto idx : cms::alpakatools::uniform_elements(acc, foundNtuplets->nOnes() - 1)) {
0711         tracks_view[idx].hitOffsets() = foundNtuplets->off[idx + 1];  // offset for track 0 is always 0
0712       }
0713       // fill hit indices
0714       for (auto idx : cms::alpakatools::uniform_elements(acc, foundNtuplets->size())) {
0715         ALPAKA_ASSERT_ACC(foundNtuplets->content[idx] < (uint32_t)hh.metadata().size());
0716         track_hits_view[idx].id() = foundNtuplets->content[idx];
0717         track_hits_view[idx].detId() = hh[foundNtuplets->content[idx]].detectorIndex();
0718 #ifdef CA_DEBUG
0719         printf("Kernel_fillHitDetIndices %d %d %d \n",
0720                idx,
0721                foundNtuplets->content[idx],
0722                track_hits_view.metadata().size());
0723 #endif
0724       }
0725     }
0726   };
0727 
0728   template <typename TrackerTraits>
0729   class Kernel_fillNLayers {
0730   public:
0731     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0732                                   TkSoAView tracks_view,
0733                                   TkHitSoAView track_hits_view,
0734                                   uint32_t const *__restrict__ layerStarts,
0735                                   uint16_t maxLayers,
0736                                   cms::alpakatools::AtomicPairCounter *apc) const {
0737       // clamp the number of tracks to the capacity of the SoA
0738       auto ntracks = std::min<int>(apc->get().first, tracks_view.metadata().size() - 1);
0739 
0740       if (cms::alpakatools::once_per_grid(acc))
0741         tracks_view.nTracks() = ntracks;
0742       for (auto idx : cms::alpakatools::uniform_elements(acc, ntracks)) {
0743         ALPAKA_ASSERT_ACC(reco::nHits(tracks_view, idx) >= 3);
0744         tracks_view[idx].nLayers() = reco::nLayers(tracks_view, track_hits_view, maxLayers, layerStarts, idx);
0745 #ifdef CA_DEBUG
0746         printf("Kernel_fillNLayers %d %d %d - %d %d\n",
0747                idx,
0748                ntracks,
0749                tracks_view[idx].nLayers(),
0750                apc->get().first,
0751                tracks_view.metadata().size() - 1);
0752 #endif
0753       }
0754     }
0755   };
0756 
0757   template <typename TrackerTraits>
0758   class Kernel_doStatsForHitInTracks {
0759   public:
0760     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0761                                   HitToTuple const *__restrict__ hitToTuple,
0762                                   Counters *counters) const {
0763       auto &c = *counters;
0764       for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple->nOnes())) {
0765         if (hitToTuple->size(idx) == 0)
0766           continue;  // SHALL NOT BE break
0767         alpaka::atomicAdd(acc, &c.nUsedHits, 1ull, alpaka::hierarchy::Blocks{});
0768         if (hitToTuple->size(idx) > 1)
0769           alpaka::atomicAdd(acc, &c.nDupHits, 1ull, alpaka::hierarchy::Blocks{});
0770       }
0771     }
0772   };
0773 
0774   template <typename TrackerTraits>
0775   class Kernel_countSharedHit {
0776   public:
0777     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0778                                   int *__restrict__ nshared,
0779                                   HitContainer const *__restrict__ ptuples,
0780                                   Quality const *__restrict__ quality,
0781                                   HitToTuple const *__restrict__ phitToTuple) const {
0782       constexpr auto loose = Quality::loose;
0783 
0784       auto &hitToTuple = *phitToTuple;
0785       auto const &foundNtuplets = *ptuples;
0786       for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
0787         if (hitToTuple.size(idx) < 2)
0788           continue;
0789 
0790         int nt = 0;
0791 
0792         // count "good" tracks
0793         for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
0794           if (quality[*it] < loose)
0795             continue;
0796           ++nt;
0797         }
0798 
0799         if (nt < 2)
0800           continue;
0801 
0802         // now mark  each track triplet as sharing a hit
0803         for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
0804           if (foundNtuplets.size(*it) > 3)
0805             continue;
0806           alpaka::atomicAdd(acc, &nshared[*it], 1, alpaka::hierarchy::Blocks{});
0807         }
0808 
0809       }  //  hit loop
0810     }
0811   };
0812 
0813   template <typename TrackerTraits>
0814   class Kernel_markSharedHit {
0815     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0816                                   int const *__restrict__ nshared,
0817                                   HitContainer const *__restrict__ tuples,
0818                                   Quality *__restrict__ quality,
0819                                   bool dupPassThrough) const {
0820       // constexpr auto bad = Quality::bad;
0821       constexpr auto dup = Quality::dup;
0822       constexpr auto loose = Quality::loose;
0823       // constexpr auto strict = Quality::strict;
0824 
0825       // quality to mark rejected
0826       auto const reject = dupPassThrough ? loose : dup;
0827       for (auto idx : cms::alpakatools::uniform_elements(acc, tuples->nOnes())) {
0828         if (tuples->size(idx) == 0)
0829           break;  //guard
0830         if (quality[idx] <= reject)
0831           continue;
0832         if (nshared[idx] > 2)
0833           quality[idx] = reject;
0834       }
0835     }
0836   };
0837 
0838   // mostly for very forward triplets.....
0839   template <typename TrackerTraits>
0840   class Kernel_rejectDuplicate {
0841   public:
0842     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0843                                   TkSoAView tracks_view,
0844                                   bool dupPassThrough,
0845                                   HitToTuple const *__restrict__ phitToTuple) const {
0846       // quality to mark rejected
0847       auto const reject = dupPassThrough ? Quality::loose : Quality::dup;
0848 
0849       auto &hitToTuple = *phitToTuple;
0850 
0851       for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
0852         if (hitToTuple.size(idx) < 2)
0853           continue;
0854 
0855         auto score = [&](auto it, auto nl) { return std::abs(reco::tip(tracks_view, it)); };
0856 
0857         // full combinatorics
0858         for (auto ip = hitToTuple.begin(idx); ip < hitToTuple.end(idx) - 1; ++ip) {
0859           auto const it = *ip;
0860           auto qi = tracks_view[it].quality();
0861           if (qi <= reject)
0862             continue;
0863           auto opi = tracks_view[it].state()(2);
0864           auto e2opi = tracks_view[it].covariance()(9);
0865           auto cti = tracks_view[it].state()(3);
0866           auto e2cti = tracks_view[it].covariance()(12);
0867           auto nli = tracks_view[it].nLayers();
0868           for (auto jp = ip + 1; jp < hitToTuple.end(idx); ++jp) {
0869             auto const jt = *jp;
0870             auto qj = tracks_view[jt].quality();
0871             if (qj <= reject)
0872               continue;
0873             auto opj = tracks_view[jt].state()(2);
0874             auto ctj = tracks_view[jt].state()(3);
0875             auto dct = nSigma2 * (tracks_view[jt].covariance()(12) + e2cti);
0876             if ((cti - ctj) * (cti - ctj) > dct)
0877               continue;
0878             auto dop = nSigma2 * (tracks_view[jt].covariance()(9) + e2opi);
0879             if ((opi - opj) * (opi - opj) > dop)
0880               continue;
0881             auto nlj = tracks_view[jt].nLayers();
0882             if (nlj < nli || (nlj == nli && (qj < qi || (qj == qi && score(it, nli) < score(jt, nlj)))))
0883               tracks_view[jt].quality() = reject;
0884             else {
0885               tracks_view[it].quality() = reject;
0886               break;
0887             }
0888           }
0889         }
0890       }
0891     }
0892   };
0893 
0894   template <typename TrackerTraits>
0895   class Kernel_sharedHitCleaner {
0896   public:
0897     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0898                                   HitsConstView hh,
0899                                   uint32_t const *__restrict__ layerStarts,
0900                                   TkSoAView tracks_view,
0901                                   int nmin,
0902                                   bool dupPassThrough,
0903                                   HitToTuple const *__restrict__ phitToTuple) const {
0904       // quality to mark rejected
0905       auto const reject = dupPassThrough ? Quality::loose : Quality::dup;
0906       // quality of longest track
0907       auto const longTqual = Quality::highPurity;
0908 
0909       auto &hitToTuple = *phitToTuple;
0910 
0911       uint32_t l1end = layerStarts[1];
0912 
0913       for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
0914         if (hitToTuple.size(idx) < 2)
0915           continue;
0916 
0917         int8_t maxNl = 0;
0918 
0919         // find maxNl
0920         for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
0921           if (tracks_view[*it].quality() < longTqual)
0922             continue;
0923           // if (tracks_view[*it].nHits()==3) continue;
0924           auto nl = tracks_view[*it].nLayers();
0925           maxNl = std::max(nl, maxNl);
0926         }
0927 
0928         if (maxNl < 4)
0929           continue;
0930 
0931         // quad pass through (leave for tests)
0932         // maxNl = std::min(4, maxNl);
0933 
0934         // kill all tracks shorter than maxHl (only triplets???
0935         for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
0936           auto nl = tracks_view[*it].nLayers();
0937 
0938           //checking if shared hit is on bpix1 and if the tuple is short enough
0939           if (idx < l1end and nl > nmin)
0940             continue;
0941 
0942           if (nl < maxNl && tracks_view[*it].quality() > reject)
0943             tracks_view[*it].quality() = reject;
0944         }
0945       }
0946     }
0947   };
0948   template <typename TrackerTraits>
0949   class Kernel_tripletCleaner {
0950   public:
0951     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
0952                                   TkSoAView tracks_view,
0953                                   bool dupPassThrough,
0954                                   HitToTuple const *__restrict__ phitToTuple) const {
0955       // quality to mark rejected
0956       auto const reject = Quality::loose;
0957       /// min quality of good
0958       auto const good = Quality::strict;
0959 
0960       auto &hitToTuple = *phitToTuple;
0961 
0962       for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
0963         if (hitToTuple.size(idx) < 2)
0964           continue;
0965 
0966         float mc = maxScore;
0967         uint32_t im = tkNotFound;
0968         bool onlyTriplets = true;
0969 
0970         // check if only triplets
0971         for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
0972           if (tracks_view[*it].quality() <= good)
0973             continue;
0974           onlyTriplets &= reco::isTriplet(tracks_view, *it);
0975           if (!onlyTriplets)
0976             break;
0977         }
0978 
0979         // only triplets
0980         if (!onlyTriplets)
0981           continue;
0982 
0983         // for triplets choose best tip!  (should we first find best quality???)
0984         for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
0985           auto const it = *ip;
0986           if (tracks_view[it].quality() >= good && std::abs(reco::tip(tracks_view, it)) < mc) {
0987             mc = std::abs(reco::tip(tracks_view, it));
0988             im = it;
0989           }
0990         }
0991 
0992         if (tkNotFound == im)
0993           continue;
0994 
0995         // mark worse ambiguities
0996         for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
0997           auto const it = *ip;
0998           if (tracks_view[it].quality() > reject && it != im)
0999             tracks_view[it].quality() = reject;  //no race:  simple assignment of the same constant
1000         }
1001 
1002       }  // loop over hits
1003     }
1004   };
1005 
1006   template <typename TrackerTraits>
1007   class Kernel_simpleTripletCleaner {
1008   public:
1009     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
1010                                   TkSoAView tracks_view,
1011                                   bool dupPassThrough,
1012                                   HitToTuple const *__restrict__ phitToTuple) const {
1013       // quality to mark rejected
1014       auto const reject = Quality::loose;
1015       /// min quality of good
1016       auto const good = Quality::loose;
1017 
1018       auto &hitToTuple = *phitToTuple;
1019 
1020       for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
1021         if (hitToTuple.size(idx) < 2)
1022           continue;
1023 
1024         float mc = maxScore;
1025         uint32_t im = tkNotFound;
1026 
1027         // choose best tip!  (should we first find best quality???)
1028         for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
1029           auto const it = *ip;
1030           if (tracks_view[it].quality() >= good && std::abs(reco::tip(tracks_view, it)) < mc) {
1031             mc = std::abs(reco::tip(tracks_view, it));
1032             im = it;
1033           }
1034         }
1035 
1036         if (tkNotFound == im)
1037           continue;
1038 
1039         // mark worse ambiguities
1040         for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
1041           auto const it = *ip;
1042           if (tracks_view[it].quality() > reject && reco::isTriplet(tracks_view, it) && it != im)
1043             tracks_view[it].quality() = reject;  //no race:  simple assignment of the same constant
1044         }
1045 
1046       }  // loop over hits
1047     }
1048   };
1049 
1050   template <typename TrackerTraits>
1051   class Kernel_print_found_ntuplets {
1052   public:
1053     ALPAKA_FN_ACC void operator()(Acc1D const &acc,
1054                                   HitsConstView hh,
1055                                   TkSoAView tracks_view,
1056                                   HitContainer const *__restrict__ foundNtuplets,
1057                                   HitToTuple const *__restrict__ phitToTuple,
1058                                   int32_t firstPrint,
1059                                   int32_t lastPrint,
1060                                   int iev) const {
1061       constexpr auto loose = Quality::loose;
1062 
1063       for (auto i : cms::alpakatools::uniform_elements(acc, firstPrint, std::min(lastPrint, foundNtuplets->nOnes()))) {
1064         auto nh = foundNtuplets->size(i);
1065         if (nh < 3)
1066           continue;
1067         if (tracks_view[i].quality() < loose)
1068           continue;
1069         printf("TK: %d %d %d %d %f %f %f %f %f %f %f %.3f %.3f %.3f %.3f %.3f %.3f %.3f\n",
1070                10000 * iev + i,
1071                int(tracks_view[i].quality()),
1072                nh,
1073                tracks_view[i].nLayers(),
1074                reco::charge(tracks_view, i),
1075                tracks_view[i].pt(),
1076                tracks_view[i].eta(),
1077                reco::phi(tracks_view, i),
1078                reco::tip(tracks_view, i),
1079                reco::zip(tracks_view, i),
1080                tracks_view[i].chi2(),
1081                hh[*foundNtuplets->begin(i)].zGlobal(),
1082                hh[*(foundNtuplets->begin(i) + 1)].zGlobal(),
1083                hh[*(foundNtuplets->begin(i) + 2)].zGlobal(),
1084                nh > 3 ? hh[int(*(foundNtuplets->begin(i) + 3))].zGlobal() : 0,
1085                nh > 4 ? hh[int(*(foundNtuplets->begin(i) + 4))].zGlobal() : 0,
1086                nh > 5 ? hh[int(*(foundNtuplets->begin(i) + 5))].zGlobal() : 0,
1087                nh > 6 ? hh[int(*(foundNtuplets->begin(i) + nh - 1))].zGlobal() : 0);
1088       }
1089     }
1090   };
1091 
1092   class Kernel_printCounters {
1093   public:
1094     ALPAKA_FN_ACC void operator()(Acc1D const &acc, Counters const *counters) const {
1095       auto const &c = *counters;
1096       printf(
1097           "||Counters | nEvents | nHits | nCells | nTuples | nFitTacks  |  nLooseTracks  |  nGoodTracks | nUsedHits | "
1098           "nDupHits | nFishCells | nKilledCells | nUsedCells | nZeroTrackCells ||\n");
1099       printf("Counters Raw %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld\n",
1100              c.nEvents,
1101              c.nHits,
1102              c.nCells,
1103              c.nTuples,
1104              c.nFitTracks,
1105              c.nLooseTracks,
1106              c.nGoodTracks,
1107              c.nUsedHits,
1108              c.nDupHits,
1109              c.nFishCells,
1110              c.nKilledCells,
1111              c.nEmptyCells,
1112              c.nZeroTrackCells);
1113       printf(
1114           "Counters Norm %lld ||  %.1f|  %.1f|  %.1f|  %.1f|  %.1f|  %.1f|  %.1f|  %.1f|  %.3f|  %.3f|  %.3f|  "
1115           "%.3f||\n",
1116           c.nEvents,
1117           c.nHits / double(c.nEvents),
1118           c.nCells / double(c.nEvents),
1119           c.nTuples / double(c.nEvents),
1120           c.nFitTracks / double(c.nEvents),
1121           c.nLooseTracks / double(c.nEvents),
1122           c.nGoodTracks / double(c.nEvents),
1123           c.nUsedHits / double(c.nEvents),
1124           c.nDupHits / double(c.nEvents),
1125           c.nFishCells / double(c.nCells),
1126           c.nKilledCells / double(c.nCells),
1127           c.nEmptyCells / double(c.nCells),
1128           c.nZeroTrackCells / double(c.nCells));
1129     }
1130   };
1131 
1132 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE::caHitNtupletGeneratorKernels
1133 
1134 #endif  // RecoTracker_PixelSeeding_plugins_alpaka_CAHitNtupletGeneratorKernelsImpl_h