CMS 3D CMS Logo

CAHitNtupletGeneratorKernels.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelSeeding_plugins_CAHitNtupletGeneratorKernels_h
2 #define RecoTracker_PixelSeeding_plugins_CAHitNtupletGeneratorKernels_h
3 
4 // #define GPU_DEBUG
5 
6 #include "GPUCACell.h"
7 #include "gpuPixelDoublets.h"
8 
14 
15 // #define DUMP_GPU_TK_TUPLES
16 
18 
19  //Configuration params common to all topologies, for the algorithms
20  struct AlgoParams {
21  const bool onGPU_;
22  const uint32_t minHitsForSharingCut_;
23  const bool useRiemannFit_;
24  const bool fitNas4_;
26  const bool earlyFishbone_;
27  const bool lateFishbone_;
28  const bool doStats_;
29  const bool doSharedHitCut_;
30  const bool dupPassThrough_;
32  };
33 
34  //CAParams
35  struct CACommon {
36  const uint32_t maxNumberOfDoublets_;
37  const uint32_t minHitsPerNtuplet_;
38  const float ptmin_;
39  const float CAThetaCutBarrel_;
40  const float CAThetaCutForward_;
41  const float hardCurvCut_;
42  const float dcaCutInnerTriplet_;
43  const float dcaCutOuterTriplet_;
44  };
45 
46  template <typename TrackerTraits, typename Enable = void>
47  struct CAParamsT : public CACommon {
48  __device__ __forceinline__ bool startingLayerPair(int16_t pid) const { return false; };
49  __device__ __forceinline__ bool startAt0(int16_t pid) const { return false; };
50  };
51 
52  template <typename TrackerTraits>
53  struct CAParamsT<TrackerTraits, pixelTopology::isPhase1Topology<TrackerTraits>> : public CACommon {
55  __device__ __forceinline__ bool startingLayerPair(int16_t pid) const {
56  return minHitsPerNtuplet_ > 3 ? pid < 3 : pid < 8 || pid > 12;
57  }
58 
60  __device__ __forceinline__ bool startAt0(int16_t pid) const {
61  assert((pixelTopology::Phase1::layerPairs[pid * 2] == 0) ==
62  (pid < 3 || pid == 13 || pid == 15 || pid == 16)); // to be 100% sure it's working, may be removed
63  return pixelTopology::Phase1::layerPairs[pid * 2] == 0;
64  }
65  };
66 
67  template <typename TrackerTraits>
68  struct CAParamsT<TrackerTraits, pixelTopology::isPhase2Topology<TrackerTraits>> : public CACommon {
69  const bool includeFarForwards_;
71  __device__ __forceinline__ bool startingLayerPair(int16_t pid) const {
72  return pid < 33; // in principle one could remove 5,6,7 23, 28 and 29
73  }
74 
76  __device__ __forceinline__ bool startAt0(int16_t pid) const {
77  assert((pixelTopology::Phase2::layerPairs[pid * 2] == 0) == ((pid < 3) | (pid >= 23 && pid < 28)));
78  return pixelTopology::Phase2::layerPairs[pid * 2] == 0;
79  }
80  };
81 
82  //Full list of params = algo params + ca params + cell params + quality cuts
83  //Generic template
84  template <typename TrackerTraits, typename Enable = void>
85  struct ParamsT : public AlgoParams {
86  // one should define the params for its own pixelTopology
87  // not defining anything here
88  inline uint32_t nPairs() const { return 0; }
89  };
90 
91  template <typename TrackerTraits>
92  struct ParamsT<TrackerTraits, pixelTopology::isPhase1Topology<TrackerTraits>> : public AlgoParams {
93  using TT = TrackerTraits;
94  using QualityCuts = pixelTrack::QualityCutsT<TT>; //track quality cuts
95  using CellCuts = gpuPixelDoublets::CellCutsT<TT>; //cell building cuts
96  using CAParams = CAParamsT<TT>; //params to be used on device
97 
98  ParamsT(AlgoParams const& commonCuts,
99  CellCuts const& cellCuts,
100  QualityCuts const& cutsCuts,
101  CAParams const& caParams)
102  : AlgoParams(commonCuts), cellCuts_(cellCuts), qualityCuts_(cutsCuts), caParams_(caParams) {}
103 
105  const QualityCuts qualityCuts_{// polynomial coefficients for the pT-dependent chi2 cut
106  {0.68177776, 0.74609577, -0.08035491, 0.00315399},
107  // max pT used to determine the chi2 cut
108  10.,
109  // chi2 scale factor: 30 for broken line fit, 45 for Riemann fit
110  30.,
111  // regional cuts for triplets
112  {
113  0.3, // |Tip| < 0.3 cm
114  0.5, // pT > 0.5 GeV
115  12.0 // |Zip| < 12.0 cm
116  },
117  // regional cuts for quadruplets
118  {
119  0.5, // |Tip| < 0.5 cm
120  0.3, // pT > 0.3 GeV
121  12.0 // |Zip| < 12.0 cm
122  }};
125  inline uint32_t nPairs() const {
126  // take all layer pairs into account
127  uint32_t nActualPairs = TT::nPairs;
129  // exclude forward "jumping" layer pairs
130  nActualPairs = TT::nPairsForTriplets;
131  }
132  if (caParams_.minHitsPerNtuplet_ > 3) {
133  // for quadruplets, exclude all "jumping" layer pairs
134  nActualPairs = TT::nPairsForQuadruplets;
135  }
136 
137  return nActualPairs;
138  }
139 
140  }; // Params Phase1
141 
142  template <typename TrackerTraits>
143  struct ParamsT<TrackerTraits, pixelTopology::isPhase2Topology<TrackerTraits>> : public AlgoParams {
144  using TT = TrackerTraits;
148 
149  ParamsT(AlgoParams const& commonCuts,
150  CellCuts const& cellCuts,
151  QualityCuts const& qualityCuts,
152  CAParams const& caParams)
153  : AlgoParams(commonCuts), cellCuts_(cellCuts), qualityCuts_(qualityCuts), caParams_(caParams) {}
154 
155  // quality cuts
157  const QualityCuts qualityCuts_{5.0f, /*chi2*/ 0.9f, /* pT in Gev*/ 0.4f, /*zip in cm*/ 12.0f /*tip in cm*/};
159 
160  inline uint32_t nPairs() const {
161  // take all layer pairs into account
162  uint32_t nActualPairs = TT::nPairsMinimal;
163  if (caParams_.includeFarForwards_) {
164  // considera far forwards (> 11 & > 23)
165  nActualPairs = TT::nPairsFarForwards;
166  }
168  // include jumping forwards
170  }
171 
172  return nActualPairs;
173  }
174 
175  }; // Params Phase1
176 
177  // counters
178  struct Counters {
179  unsigned long long nEvents;
180  unsigned long long nHits;
181  unsigned long long nCells;
182  unsigned long long nTuples;
183  unsigned long long nFitTracks;
184  unsigned long long nLooseTracks;
185  unsigned long long nGoodTracks;
186  unsigned long long nUsedHits;
187  unsigned long long nDupHits;
188  unsigned long long nFishCells;
189  unsigned long long nKilledCells;
190  unsigned long long nEmptyCells;
191  unsigned long long nZeroTrackCells;
192  };
193 
195 
196 } // namespace caHitNtupletGenerator
197 
198 template <typename TTraits, typename TTTraits>
200 public:
201  using Traits = TTraits;
202  using TrackerTraits = TTTraits;
208 
209  template <typename T>
211 
215 
224 
226 
229 
231  : params_(params), paramsMaxDoubletes3Quarters_(3 * params.caParams_.maxNumberOfDoublets_ / 4) {}
232 
233  ~CAHitNtupletGeneratorKernels() = default;
234 
236 
237  void launchKernels(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
238 
239  void classifyTuples(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
240 
241  void buildDoublets(const HitsConstView& hh, cudaStream_t stream);
242  void allocateOnGPU(int32_t nHits, cudaStream_t stream);
243  void cleanup(cudaStream_t cudaStream);
244 
245  static void printCounters(Counters const* counters);
247 
248 protected:
249  Counters* counters_ = nullptr;
250 
251  // workspace
257 
261  uint32_t* device_nCells_ = nullptr;
262 
266 
268 
270 
272 
274 
276 
277  // params
282  inline uint32_t nDoubletBlocks(uint32_t blockSize) {
283  // We want (3 * params_.maxNumberOfDoublets_ / 4 + blockSize - 1) / blockSize, but first part is pre-computed.
284  return (paramsMaxDoubletes3Quarters_ + blockSize - 1) / blockSize;
285  }
286 
288  inline uint32_t nQuadrupletBlocks(uint32_t blockSize) {
289  // pixelTopology::maxNumberOfQuadruplets is a constexpr, so the compiler will pre compute the 3*max/4
290  return (3 * TrackerTraits::maxNumberOfQuadruplets / 4 + blockSize - 1) / blockSize;
291  }
292 };
293 
294 template <typename TrackerTraits>
295 class CAHitNtupletGeneratorKernelsGPU : public CAHitNtupletGeneratorKernels<cms::cudacompat::GPUTraits, TrackerTraits> {
297 
300 
302 
307 
310 
312 
313 public:
314  void launchKernels(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
315  void classifyTuples(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
316  void buildDoublets(const HitsConstView& hh, int32_t offsetBPIX2, cudaStream_t stream);
317  void allocateOnGPU(int32_t nHits, cudaStream_t stream);
318  static void printCounters(Counters const* counters);
319 };
320 
321 template <typename TrackerTraits>
322 class CAHitNtupletGeneratorKernelsCPU : public CAHitNtupletGeneratorKernels<cms::cudacompat::CPUTraits, TrackerTraits> {
324 
327 
329 
334 
337 
339 
340 public:
341  void launchKernels(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
342  void classifyTuples(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
343  void buildDoublets(const HitsConstView& hh, int32_t offsetBPIX2, cudaStream_t stream);
344  void allocateOnGPU(int32_t nHits, cudaStream_t stream);
345  static void printCounters(Counters const* counters);
346 };
347 
348 #endif // RecoTracker_PixelSeeding_plugins_CAHitNtupletGeneratorKernels_h
void classifyTuples(const HitsConstView &hh, TkSoAView &track_view, cudaStream_t cudaStream)
CAHitNtupletGeneratorKernels(Params const &params)
TrackingRecHitSoAConstView< TrackerTraits > HitsConstView
static void printCounters(Counters const *counters)
#define __forceinline__
Definition: cudaCompat.h:22
unique_ptr< OuterHitOfCellContainer[]> device_isOuterHitOfCell_
void launchKernels(const HitsConstView &hh, TkSoAView &track_view, cudaStream_t cudaStream)
void classifyTuples(const HitsConstView &hh, TkSoAView &track_view, cudaStream_t cudaStream)
unique_ptr< CellTracksVector > device_theCellTracks_
typename std::enable_if< std::is_base_of< Phase2, T >::value >::type isPhase2Topology
caHitNtupletGenerator::Counters Counters
cms::cuda::AtomicPairCounter * device_hitToTuple_apc_
TrackSoAView< TrackerTraits > TkSoAView
const uint32_t paramsMaxDoubletes3Quarters_
Intermediate result avoiding repeated computations.
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
unique_ptr< HitToTuple > device_hitToTuple_
assert(be >=bs)
unique_ptr< uint32_t[]> device_hitToTupleStorage_
void classifyTuples(const HitsConstView &hh, TkSoAView &track_view, cudaStream_t cudaStream)
void allocateOnGPU(int32_t nHits, cudaStream_t stream)
ParamsT(AlgoParams const &commonCuts, CellCuts const &cellCuts, QualityCuts const &cutsCuts, CAParams const &caParams)
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ uint32_t const *__restrict__ CellNeighborsVector< TrackerTraits > const CellTracksVector< TrackerTraits > const OuterHitOfCell< TrackerTraits > const int32_t uint32_t Counters * counters
static void printCounters(Counters const *counters)
unique_ptr< TupleMultiplicity > device_tupleMultiplicity_
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
void buildDoublets(const HitsConstView &hh, cudaStream_t stream)
typename std::enable_if< std::is_base_of< Phase1, T >::value >::type isPhase1Topology
TupleMultiplicity const * tupleMultiplicity() const
ParamsT(AlgoParams const &commonCuts, CellCuts const &cellCuts, QualityCuts const &qualityCuts, CAParams const &caParams)
uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > OuterHitOfCell< TrackerTraits > int nActualPairs
typename TrackSoA< TrackerTraits >::HitContainer HitContainer
uint32_t nQuadrupletBlocks(uint32_t blockSize)
Compute the number of quadruplet blocks for block size.
void launchKernels(const HitsConstView &hh, TkSoAView &track_view, cudaStream_t cudaStream)
caHitNtupletGenerator::Counters Counters
unique_ptr< unsigned char[]> cellStorage_
uint32_t nDoubletBlocks(uint32_t blockSize)
Compute the number of doublet blocks for block size.
unique_ptr< CellNeighborsVector > device_theCellNeighbors_
void allocateOnGPU(int32_t nHits, cudaStream_t stream)
typename TrackSoA< TrackerTraits >::HitContainer HitContainer
TrackingRecHitSoAConstView< TrackerTraits > HitsConstView
typename TrackingRecHitSoA< TrackerTraits >::template TrackingRecHitSoALayout<>::ConstView TrackingRecHitSoAConstView
void allocateOnGPU(int32_t nHits, cudaStream_t stream)
void launchKernels(const HitsConstView &hh, TkSoAView &track_view, cudaStream_t cudaStream)
void buildDoublets(const HitsConstView &hh, int32_t offsetBPIX2, cudaStream_t stream)
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::View TrackSoAView
void buildDoublets(const HitsConstView &hh, int32_t offsetBPIX2, cudaStream_t stream)
typename TrackingRecHitSoA< TrackerTraits >::template TrackingRecHitSoALayout<>::View TrackingRecHitSoAView
unique_ptr< cms::cuda::AtomicPairCounter::c_type[]> device_storage_
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > hh
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
#define __device__
static constexpr uint8_t const * layerPairs
TrackSoAView< TrackerTraits > TkSoAView
static void printCounters(Counters const *counters)
static constexpr uint8_t const * layerPairs
cms::cuda::AtomicPairCounter * device_hitTuple_apc_
void cleanup(cudaStream_t cudaStream)