CMS 3D CMS Logo

gpuPixelDoublets.h
Go to the documentation of this file.
1 #ifndef RecoPixelVertexing_PixelTriplets_plugins_gpuPixelDoublets_h
2 #define RecoPixelVertexing_PixelTriplets_plugins_gpuPixelDoublets_h
3 
5 
6 #define CONSTANT_VAR __constant__
7 
8 namespace gpuPixelDoublets {
9 
10  constexpr int nPairsForQuadruplets = 13; // quadruplets require hits in all layers
11  constexpr int nPairsForTriplets = nPairsForQuadruplets + 2; // include barrel "jumping" layer pairs
12  constexpr int nPairs = nPairsForTriplets + 4; // include forward "jumping" layer pairs
14 
15  // start constants
16  // clang-format off
17 
18  CONSTANT_VAR const uint8_t layerPairs[2 * nPairs] = {
19  0, 1, 0, 4, 0, 7, // BPIX1 (3)
20  1, 2, 1, 4, 1, 7, // BPIX2 (6)
21  4, 5, 7, 8, // FPIX1 (8)
22  2, 3, 2, 4, 2, 7, 5, 6, 8, 9, // BPIX3 & FPIX2 (13)
23  0, 2, 1, 3, // Jumping Barrel (15)
24  0, 5, 0, 8, // Jumping Forward (BPIX1,FPIX2)
25  4, 6, 7, 9 // Jumping Forward (19)
26  };
27 
28  constexpr int16_t phi0p05 = 522; // round(521.52189...) = phi2short(0.05);
29  constexpr int16_t phi0p06 = 626; // round(625.82270...) = phi2short(0.06);
30  constexpr int16_t phi0p07 = 730; // round(730.12648...) = phi2short(0.07);
31 
33  phi0p07,
34  phi0p07,
35  phi0p05,
36  phi0p06,
37  phi0p06,
38  phi0p05,
39  phi0p05,
40  phi0p06,
41  phi0p06,
42  phi0p06,
43  phi0p05,
44  phi0p05,
45  phi0p05,
46  phi0p05,
47  phi0p05,
48  phi0p05,
49  phi0p05,
50  phi0p05};
51  // phi0p07, phi0p07, phi0p06,phi0p06, phi0p06,phi0p06}; // relaxed cuts
52 
53  CONSTANT_VAR float const minz[nPairs] = {
54  -20., 0., -30., -22., 10., -30., -70., -70., -22., 15., -30, -70., -70., -20., -22., 0, -30., -70., -70.};
55  CONSTANT_VAR float const maxz[nPairs] = {
56  20., 30., 0., 22., 30., -10., 70., 70., 22., 30., -15., 70., 70., 20., 22., 30., 0., 70., 70.};
57  CONSTANT_VAR float const maxr[nPairs] = {
58  20., 9., 9., 20., 7., 7., 5., 5., 20., 6., 6., 5., 5., 20., 20., 9., 9., 9., 9.};
59 
60  // end constants
61  // clang-format on
62 
67 
69  int nHits,
71  CellNeighbors* cellNeighborsContainer,
73  CellTracks* cellTracksContainer) {
74  assert(isOuterHitOfCell.container);
75  int first = blockIdx.x * blockDim.x + threadIdx.x;
76  for (int i = first; i < nHits - isOuterHitOfCell.offset; i += gridDim.x * blockDim.x)
77  isOuterHitOfCell.container[i].reset();
78 
79  if (0 == first) {
82  auto i = cellNeighbors->extend();
83  assert(0 == i);
84  (*cellNeighbors)[0].reset();
85  i = cellTracks->extend();
86  assert(0 == i);
87  (*cellTracks)[0].reset();
88  }
89  }
90 
91  constexpr auto getDoubletsFromHistoMaxBlockSize = 64; // for both x and y
93 
95 #ifdef __CUDACC__
97 #endif
98  void getDoubletsFromHisto(GPUCACell* cells,
99  uint32_t* nCells,
102  TrackingRecHit2DSOAView const* __restrict__ hhp,
104  int nActualPairs,
105  bool ideal_cond,
106  bool doClusterCut,
107  bool doZ0Cut,
108  bool doPtCut,
109  uint32_t maxNumOfDoublets) {
110  auto const& __restrict__ hh = *hhp;
112  nActualPairs,
113  cells,
114  nCells,
116  cellTracks,
117  hh,
119  phicuts,
120  minz,
121  maxz,
122  maxr,
123  ideal_cond,
124  doClusterCut,
125  doZ0Cut,
126  doPtCut,
128  }
129 
130 } // namespace gpuPixelDoublets
131 
132 #endif // RecoPixelVertexing_PixelTriplets_plugins_gpuPixelDoublets_h
const dim3 threadIdx
Definition: cudaCompat.h:29
uint32_t CellNeighborsVector CellTracksVector * cellTracks
cms::cuda::VecArray< tindex_type, maxCellTracks > CellTracks
Definition: CAConstants.h:72
constexpr int16_t phi0p06
__constant__ const uint8_t layerPairs[2 *nPairs]
constexpr uint32_t maxNumOfActiveDoublets
Definition: CAConstants.h:41
cms::cuda::SimpleVector< CellNeighbors > CellNeighborsVector
Definition: CAConstants.h:74
__constant__ float const maxz[nPairs]
__constant__ float const minz[nPairs]
const dim3 gridDim
Definition: cudaCompat.h:33
__device__ uint32_t GPUCACell * cells
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool bool doZ0Cut
#define __global__
Definition: cudaCompat.h:19
constexpr uint32_t maxNumberOfLayerPairs
Definition: CAConstants.h:44
const dim3 blockDim
Definition: cudaCompat.h:30
constexpr int16_t phi0p05
__constant__ float const maxr[nPairs]
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool bool bool doPtCut
constexpr int nPairs
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell isOuterHitOfCell
constexpr int16_t phi0p07
cms::cuda::SimpleVector< CellTracks > CellTracksVector
Definition: CAConstants.h:75
__device__ int extend(int size=1)
Definition: SimpleVector.h:84
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool doClusterCut
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool ideal_cond
const dim3 blockIdx
Definition: cudaCompat.h:32
constexpr int nPairsForTriplets
__constant__ const int16_t phicuts[nPairs]
constexpr auto getDoubletsFromHistoMinBlocksPerMP
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int nActualPairs
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
constexpr void reset()
Definition: SimpleVector.h:108
uint32_t CellNeighborsVector * cellNeighbors
constexpr auto getDoubletsFromHistoMaxBlockSize
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool bool bool uint32_t maxNumOfDoublets
constexpr void construct(int capacity, T *data)
Definition: SimpleVector.h:19
constexpr int nPairsForQuadruplets
#define CONSTANT_VAR
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ hhp
__device__ uint32_t GPUCACell uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const &__restrict__ hh
doubletsFromHisto(layerPairs, nActualPairs, cells, nCells, cellNeighbors, cellTracks, hh, isOuterHitOfCell, phicuts, minz, maxz, maxr, ideal_cond, doClusterCut, doZ0Cut, doPtCut, maxNumOfDoublets)
cms::cuda::VecArray< uint32_t, maxCellNeighbors > CellNeighbors
Definition: CAConstants.h:71