CMS 3D CMS Logo

gpuFishbone.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelSeeding_plugins_gpuFishbone_h
2 #define RecoTracker_PixelSeeding_plugins_gpuFishbone_h
3 
4 #include <algorithm>
5 #include <cmath>
6 #include <cstdint>
7 #include <cstdio>
8 #include <limits>
9 
13 
14 #include "GPUCACell.h"
15 #include "CAStructures.h"
16 
17 namespace gpuPixelDoublets {
18 
19  template <typename TrackerTraits>
21  template <typename TrackerTraits>
23  template <typename TrackerTraits>
25  template <typename TrackerTraits>
27  template <typename TrackerTraits>
29  template <typename TrackerTraits>
31 
32  template <typename TrackerTraits>
35  uint32_t const* __restrict__ nCells,
37  int32_t nHits,
38  bool checkTrack) {
40 
41  auto const isOuterHitOfCell = isOuterHitOfCellWrap.container;
42  int32_t offset = isOuterHitOfCellWrap.offset;
43 
44  // x runs faster...
45  auto firstY = threadIdx.y + blockIdx.y * blockDim.y;
46  auto firstX = threadIdx.x;
47 
48  float x[maxCellsPerHit], y[maxCellsPerHit], z[maxCellsPerHit], n[maxCellsPerHit];
49  uint32_t cc[maxCellsPerHit];
50  uint16_t d[maxCellsPerHit];
51  uint8_t l[maxCellsPerHit];
52 
53  for (int idy = firstY, nt = nHits - offset; idy < nt; idy += gridDim.y * blockDim.y) {
54  auto const& vc = isOuterHitOfCell[idy];
55  auto size = vc.size();
56  if (size < 2)
57  continue;
58  // if alligned kill one of the two.
59  // in principle one could try to relax the cut (only in r-z?) for jumping-doublets
60  auto const& c0 = cells[vc[0]];
61  auto xo = c0.outer_x(hh);
62  auto yo = c0.outer_y(hh);
63  auto zo = c0.outer_z(hh);
64  auto sg = 0;
65  for (int32_t ic = 0; ic < size; ++ic) {
66  auto& ci = cells[vc[ic]];
67  if (ci.unused())
68  continue; // for triplets equivalent to next
69  if (checkTrack && ci.tracks().empty())
70  continue;
71  cc[sg] = vc[ic];
72  l[sg] = ci.layerPairId();
73  d[sg] = ci.inner_detIndex(hh);
74  x[sg] = ci.inner_x(hh) - xo;
75  y[sg] = ci.inner_y(hh) - yo;
76  z[sg] = ci.inner_z(hh) - zo;
77  n[sg] = x[sg] * x[sg] + y[sg] * y[sg] + z[sg] * z[sg];
78  ++sg;
79  }
80  if (sg < 2)
81  continue;
82  // here we parallelize
83  for (int32_t ic = firstX; ic < sg - 1; ic += blockDim.x) {
84  auto& ci = cells[cc[ic]];
85  for (auto jc = ic + 1; jc < sg; ++jc) {
86  auto& cj = cells[cc[jc]];
87  // must be different detectors
88  // if (d[ic]==d[jc]) continue;
89  auto cos12 = x[ic] * x[jc] + y[ic] * y[jc] + z[ic] * z[jc];
90  if (d[ic] != d[jc] && cos12 * cos12 >= 0.99999f * (n[ic] * n[jc])) {
91  // alligned: kill farthest (prefer consecutive layers)
92  // if same layer prefer farthest (longer level arm) and make space for intermediate hit
93  bool sameLayer = l[ic] == l[jc];
94  if (n[ic] > n[jc]) {
95  if (sameLayer) {
96  cj.kill(); // closest
97  ci.setFishbone(cj.inner_hit_id(), cj.inner_z(hh), hh);
98  } else {
99  ci.kill(); // farthest
100  // break; // removed to improve reproducibility. keep it for reference and tests
101  }
102  } else {
103  if (!sameLayer) {
104  cj.kill(); // farthest
105  } else {
106  ci.kill(); // closest
107  cj.setFishbone(ci.inner_hit_id(), ci.inner_z(hh), hh);
108  // break; // removed to improve reproducibility. keep it for reference and tests
109  }
110  }
111  }
112  } // cj
113  } // ci
114  } // hits
115  }
116 
117 } // namespace gpuPixelDoublets
118 
119 #endif // RecoTracker_PixelSeeding_plugins_gpuFishbone_h
const dim3 threadIdx
Definition: cudaCompat.h:29
float z[maxCellsPerHit]
Definition: gpuFishbone.h:48
float n[maxCellsPerHit]
Definition: gpuFishbone.h:48
TrackingRecHitSoAConstView< TrackerTraits > HitsConstView
Definition: GPUCACell.h:34
const dim3 gridDim
Definition: cudaCompat.h:33
uint16_t d[maxCellsPerHit]
Definition: gpuFishbone.h:50
float x[maxCellsPerHit]
Definition: gpuFishbone.h:48
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
#define __global__
Definition: cudaCompat.h:19
GPUCACellT< TrackerTraits > uint32_t const *__restrict__ OuterHitOfCell< TrackerTraits > const isOuterHitOfCellWrap
Definition: gpuFishbone.h:34
const dim3 blockDim
Definition: cudaCompat.h:30
float y[maxCellsPerHit]
Definition: gpuFishbone.h:48
GPUCACellT< TrackerTraits > * cells
Definition: gpuFishbone.h:34
uint8_t l[maxCellsPerHit]
Definition: gpuFishbone.h:51
auto const isOuterHitOfCell
Definition: gpuFishbone.h:41
double f[11][100]
GPUCACellT< TrackerTraits > uint32_t const *__restrict__ nCells
Definition: gpuFishbone.h:34
const dim3 blockIdx
Definition: cudaCompat.h:32
int nt
Definition: AMPTWrapper.h:42
uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > hh
typename GPUCACellT< TrackerTraits >::HitsConstView HitsConstView
Definition: gpuFishbone.h:30
GPUCACellT< TrackerTraits > uint32_t const *__restrict__ OuterHitOfCell< TrackerTraits > const int32_t bool checkTrack
Definition: gpuFishbone.h:38
GPUCACellT< TrackerTraits > uint32_t const *__restrict__ OuterHitOfCell< TrackerTraits > const int32_t nHits
Definition: gpuFishbone.h:34