CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
gpuFishbone.h
Go to the documentation of this file.
1 #ifndef RecoPixelVertexing_PixelTriplets_plugins_gpuFishbone_h
2 #define RecoPixelVertexing_PixelTriplets_plugins_gpuFishbone_h
3 
4 #include <algorithm>
5 #include <cmath>
6 #include <cstdint>
7 #include <cstdio>
8 #include <limits>
9 
14 
15 #include "GPUCACell.h"
16 
17 namespace gpuPixelDoublets {
18 
19  __global__ void fishbone(GPUCACell::Hits const* __restrict__ hhp,
21  uint32_t const* __restrict__ nCells,
22  GPUCACell::OuterHitOfCell const isOuterHitOfCellWrap,
23  int32_t nHits,
24  bool checkTrack) {
26 
27  auto const& hh = *hhp;
28 
29  auto const isOuterHitOfCell = isOuterHitOfCellWrap.container;
30  int32_t offset = isOuterHitOfCellWrap.offset;
31 
32  // x runs faster...
33  auto firstY = threadIdx.y + blockIdx.y * blockDim.y;
34  auto firstX = threadIdx.x;
35 
37  uint32_t cc[maxCellsPerHit];
38  uint16_t d[maxCellsPerHit];
39  uint8_t l[maxCellsPerHit];
40 
41  for (int idy = firstY, nt = nHits - offset; idy < nt; idy += gridDim.y * blockDim.y) {
42  auto const& vc = isOuterHitOfCell[idy];
43  auto size = vc.size();
44  if (size < 2)
45  continue;
46  // if alligned kill one of the two.
47  // in principle one could try to relax the cut (only in r-z?) for jumping-doublets
48  auto const& c0 = cells[vc[0]];
49  auto xo = c0.outer_x(hh);
50  auto yo = c0.outer_y(hh);
51  auto zo = c0.outer_z(hh);
52  auto sg = 0;
53  for (int32_t ic = 0; ic < size; ++ic) {
54  auto& ci = cells[vc[ic]];
55  if (ci.unused())
56  continue; // for triplets equivalent to next
57  if (checkTrack && ci.tracks().empty())
58  continue;
59  cc[sg] = vc[ic];
60  l[sg] = ci.layerPairId();
61  d[sg] = ci.inner_detIndex(hh);
62  x[sg] = ci.inner_x(hh) - xo;
63  y[sg] = ci.inner_y(hh) - yo;
64  z[sg] = ci.inner_z(hh) - zo;
65  n[sg] = x[sg] * x[sg] + y[sg] * y[sg] + z[sg] * z[sg];
66  ++sg;
67  }
68  if (sg < 2)
69  continue;
70  // here we parallelize
71  for (int32_t ic = firstX; ic < sg - 1; ic += blockDim.x) {
72  auto& ci = cells[cc[ic]];
73  for (auto jc = ic + 1; jc < sg; ++jc) {
74  auto& cj = cells[cc[jc]];
75  // must be different detectors
76  // if (d[ic]==d[jc]) continue;
77  auto cos12 = x[ic] * x[jc] + y[ic] * y[jc] + z[ic] * z[jc];
78  if (d[ic] != d[jc] && cos12 * cos12 >= 0.99999f * n[ic] * n[jc]) {
79  // alligned: kill farthest (prefer consecutive layers)
80  // if same layer prefer farthest (longer level arm) and make space for intermediate hit
81  bool sameLayer = l[ic] == l[jc];
82  if (n[ic] > n[jc]) {
83  if (sameLayer) {
84  cj.kill(); // closest
85  ci.setFishbone(cj.inner_hit_id());
86  } else {
87  ci.kill(); // farthest
88  break;
89  }
90  } else {
91  if (!sameLayer) {
92  cj.kill(); // farthest
93  } else {
94  ci.kill(); // closest
95  cj.setFishbone(ci.inner_hit_id());
96  break;
97  }
98  }
99  }
100  } // cj
101  } // ci
102  } // hits
103  }
104 
105 } // namespace gpuPixelDoublets
106 
107 #endif // RecoPixelVertexing_PixelTriplets_plugins_gpuFishbone_h
const dim3 threadIdx
Definition: cudaCompat.h:29
const dim3 gridDim
Definition: cudaCompat.h:33
__device__ uint32_t GPUCACell * cells
#define __global__
Definition: cudaCompat.h:19
const dim3 blockDim
Definition: cudaCompat.h:30
tuple d
Definition: ztail.py:151
const dim3 blockIdx
Definition: cudaCompat.h:32
int nt
Definition: AMPTWrapper.h:42
static constexpr auto maxCellsPerHit
Definition: GPUCACell.h:24
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell isOuterHitOfCell
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 uint32_t maxCellsPerHit
Definition: CAConstants.h:38
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ hhp
__device__ CellTracksVector Hits const int layerPairId
Definition: GPUCACell.h:46
tuple size
Write out results.
__device__ uint32_t GPUCACell uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const &__restrict__ hh
OuterHitOfCellContainer * container
Definition: CAConstants.h:83