CMS 3D CMS Logo

fitVertices.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelVertexFinding_plugins_alpaka_fitVertices_h
2 #define RecoTracker_PixelVertexFinding_plugins_alpaka_fitVertices_h
3 
4 #include <algorithm>
5 #include <cmath>
6 #include <cstdint>
7 
8 #include <alpaka/alpaka.hpp>
9 
13 
14 #include "vertexFinder.h"
15 
17 
18  template <typename TAcc>
19  ALPAKA_FN_ACC ALPAKA_FN_INLINE __attribute__((always_inline)) void fitVertices(const TAcc& acc,
21  WsSoAView& pws,
22  float chi2Max // for outlier rejection
23  ) {
24  constexpr bool verbose = false; // in principle the compiler should optmize out if false
25 
26  auto& __restrict__ data = pdata;
27  auto& __restrict__ ws = pws;
28  auto nt = ws.ntrks();
29  float const* __restrict__ zt = ws.zt();
30  float const* __restrict__ ezt2 = ws.ezt2();
31  float* __restrict__ zv = data.zv();
32  float* __restrict__ wv = data.wv();
33  float* __restrict__ chi2 = data.chi2();
34  uint32_t& nvFinal = data.nvFinal();
35  uint32_t& nvIntermediate = ws.nvIntermediate();
36 
37  int32_t* __restrict__ nn = data.ndof();
38  int32_t* __restrict__ iv = ws.iv();
39 
42  auto foundClusters = nvFinal;
43 
44  // zero
46  zv[i] = 0;
47  wv[i] = 0;
48  chi2[i] = 0;
49  }
50 
51  // only for test
52  auto& noise = alpaka::declareSharedVar<int, __COUNTER__>(acc);
53 
54  if constexpr (verbose) {
56  noise = 0;
57  }
58  alpaka::syncBlockThreads(acc);
59 
60  // compute cluster location
61  for (auto i : cms::alpakatools::uniform_elements(acc, nt)) {
62  if (iv[i] > 9990) {
63  if constexpr (verbose)
64  alpaka::atomicAdd(acc, &noise, 1, alpaka::hierarchy::Threads{});
65  continue;
66  }
67  ALPAKA_ASSERT_ACC(iv[i] >= 0);
69  auto w = 1.f / ezt2[i];
70  alpaka::atomicAdd(acc, &zv[iv[i]], zt[i] * w, alpaka::hierarchy::Threads{});
71  alpaka::atomicAdd(acc, &wv[iv[i]], w, alpaka::hierarchy::Threads{});
72  }
73 
74  alpaka::syncBlockThreads(acc);
75  // reuse nn
77  ALPAKA_ASSERT_ACC(wv[i] > 0.f);
78  zv[i] /= wv[i];
79  nn[i] = -1; // ndof
80  }
81  alpaka::syncBlockThreads(acc);
82 
83  // compute chi2
84  for (auto i : cms::alpakatools::uniform_elements(acc, nt)) {
85  if (iv[i] > 9990)
86  continue;
87 
88  auto c2 = zv[iv[i]] - zt[i];
89  c2 *= c2 / ezt2[i];
90  if (c2 > chi2Max) {
91  iv[i] = 9999;
92  continue;
93  }
96  }
97  alpaka::syncBlockThreads(acc);
98 
100  if (nn[i] > 0) {
101  wv[i] *= float(nn[i]) / chi2[i];
102  }
103  }
104  if constexpr (verbose) {
106  printf("found %d proto clusters ", foundClusters);
107  printf("and %d noise\n", noise);
108  }
109  }
110  }
111 
113  public:
114  template <typename TAcc>
115  ALPAKA_FN_ACC void operator()(const TAcc& acc,
117  WsSoAView pws,
118  float chi2Max // for outlier rejection
119  ) const {
120  fitVertices(acc, pdata, pws, chi2Max);
121  }
122  };
123 
124 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder
125 
126 #endif // RecoTracker_PixelVertexFinding_plugins_alpaka_fitVertices_h
ALPAKA_FN_ACC constexpr bool once_per_block(TAcc const &acc)
ALPAKA_FN_ACC auto uniform_elements(TAcc const &acc, TArgs... args)
Definition: workdivision.h:303
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView & pdata
T w() const
fitVertices(pdata, pws, maxChi2ForFirstFit)
ALPAKA_FN_ACC ALPAKA_FN_INLINE void VtxSoAView WsSoAView & pws
std::vector< Block > Blocks
Definition: Block.h:99
double f[11][100]
::vertexFinder::PixelVertexWorkSpaceSoAView WsSoAView
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:61
ALPAKA_FN_ACC ALPAKA_FN_INLINE void __attribute__((always_inline)) clusterTracksByDensity(const TAcc &acc
ALPAKA_FN_ACC ALPAKA_FN_INLINE VtxSoAView WsSoAView float chi2Max
Definition: fitVertices.h:23
ALPAKA_FN_ACC void operator()(const TAcc &acc, VtxSoAView pdata, WsSoAView pws, float chi2Max) const
Definition: fitVertices.h:115