|
|
Go to the documentation of this file. 1 #ifndef RecoPixelVertexing_PixelTriplets_plugins_gpuPixelDoubletsAlgos_h
2 #define RecoPixelVertexing_PixelTriplets_plugins_gpuPixelDoubletsAlgos_h
33 int16_t
const* __restrict__
phicuts,
34 float const* __restrict__
minz,
35 float const* __restrict__
maxz,
36 float const* __restrict__
maxr,
44 constexpr
int minYsizeB1 = 36;
56 uint32_t
const* __restrict__
offsets =
hh.hitsLayerStart();
107 auto mi =
hh.detectorIndex(
i);
117 auto mez =
hh.zGlobal(
i);
133 if (mes > 0 && mes < minYsizeB1)
139 auto mep =
hh.iphi(
i);
140 auto mer =
hh.rGlobal(
i);
143 constexpr
float z0cut = 12.f;
144 constexpr
float hardPtCut = 0.5f;
146 constexpr
float minRadius = hardPtCut * 87.78f;
147 constexpr
float minRadius2T4 = 4.f * minRadius * minRadius;
148 auto ptcut = [&](
int j, int16_t idphi) {
149 auto r2t4 = minRadius2T4;
151 auto ro =
hh.rGlobal(
j);
153 return dphi * dphi * (r2t4 - ri * ro) > (ro - ri) * (ro - ri);
155 auto z0cutoff = [&](
int j) {
156 auto zo =
hh.zGlobal(
j);
157 auto ro =
hh.rGlobal(
j);
162 auto zsizeCut = [&](
int j) {
163 auto onlyBarrel =
outer < 4;
164 auto so =
hh.clusterSizeY(
j);
169 auto zo =
hh.zGlobal(
j);
170 auto ro =
hh.rGlobal(
j);
171 return onlyBarrel ? mes > 0 && so > 0 &&
std::abs(so - mes) >
dy
172 : (inner < 4) && mes > 0 &&
190 for (
auto kk = kl;
kk != khh; incr(
kk)) {
192 if (
kk != kl &&
kk != kh)
195 auto const* __restrict__
p =
phiBinner.begin(
kk + hoff);
202 auto mo =
hh.detectorIndex(oi);
209 auto mop =
hh.iphi(oi);
236 printf(
"OuterHitOfCell full for %d in layer %d/%d, %d,%d %d\n",
i,
inner,
outer,
nmin, tot, tooMany);
243 #endif // RecoPixelVertexing_PixelTriplets_plugins_gpuPixelDoubletsAlgos_h
T1 atomicSub(T1 *a, T2 b)
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool bool bool doPtCut
uint32_t CellNeighborsVector * cellNeighbors
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool ideal_cond
cms::cuda::VecArray< tindex_type, maxCellTracks > CellTracks
cms::cuda::VecArray< uint32_t, maxCellNeighbors > CellNeighbors
__device__ uint32_t GPUCACell * cells
__constant__ const float maxz[nPairs]
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool bool bool uint32_t maxNumOfDoublets
__device__ uint32_t GPUCACell uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const &__restrict__ hh
T1 atomicAdd(T1 *a, T2 b)
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool bool doZ0Cut
constexpr uint16_t maxNumModules
caConstants::CellNeighborsVector CellNeighborsVector
const HitContainer *__restrict__ const TkSoA *__restrict__ Quality *__restrict__ uint16_t nmin
uint32_t CellNeighborsVector CellTracksVector * cellTracks
const __constant__ uint8_t layerPairs[2 *nPairs]
cms::cuda::SimpleVector< CellNeighbors > CellNeighborsVector
doubletsFromHisto(layerPairs, nActualPairs, cells, nCells, cellNeighbors, cellTracks, hh, isOuterHitOfCell, phicuts, minz, maxz, maxr, ideal_cond, doClusterCut, doZ0Cut, doPtCut, maxNumOfDoublets)
__constant__ const float maxr[nPairs]
auto const &__restrict__ phiBinner
const __constant__ int16_t phicuts[nPairs]
uint32_t const *__restrict__ offsets
static constexpr uint32_t nbins()
caConstants::CellTracks CellTracks
__device__ int push_back(const T &element)
static constexpr UT bin(T t)
constexpr int maxDYsize12
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell * isOuterHitOfCell
caConstants::CellNeighbors CellNeighbors
constexpr uint32_t maxNumberOfLayerPairs
Abs< T >::type abs(const T &t)
uint32_t CellNeighborsVector CellTracksVector TrackingRecHit2DSOAView const *__restrict__ GPUCACell::OuterHitOfCell int bool bool doClusterCut
cms::cuda::HistoContainer< int16_t, 128, gpuClustering::maxNumClusters, 8 *sizeof(int16_t), hindex_type, 10 > PhiBinner
constexpr float short2phi(short x)
caConstants::CellTracksVector CellTracksVector
__constant__ const float minz[nPairs]
static constexpr auto histOff(uint32_t nh)
__shared__ uint32_t innerLayerCumulativeSize[nPairsMax]
cms::cuda::SimpleVector< CellTracks > CellTracksVector