CMS 3D CMS Logo

fitVertices.h
Go to the documentation of this file.
1 #ifndef RecoVertex_PixelVertexFinding_plugins_alpaka_fitVertices_h
2 #define RecoVertex_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  ALPAKA_FN_ACC ALPAKA_FN_INLINE __attribute__((always_inline)) void fitVertices(Acc1D const& 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__ trkdata = ptrkdata;
28  auto& __restrict__ ws = pws;
29  auto nt = ws.ntrks();
30  float const* __restrict__ zt = ws.zt();
31  float const* __restrict__ ezt2 = ws.ezt2();
32  float* __restrict__ zv = data.zv();
33  float* __restrict__ wv = data.wv();
34  float* __restrict__ chi2 = data.chi2();
35  uint32_t& nvFinal = data.nvFinal();
36  uint32_t& nvIntermediate = ws.nvIntermediate();
37 
38  int32_t* __restrict__ nn = trkdata.ndof();
39  int32_t* __restrict__ iv = ws.iv();
40 
44 
45  // zero
47  zv[i] = 0;
48  wv[i] = 0;
49  chi2[i] = 0;
50  }
51 
52  // only for test
53  auto& noise = alpaka::declareSharedVar<int, __COUNTER__>(acc);
54 
57  noise = 0;
58  }
59  alpaka::syncBlockThreads(acc);
60 
61  // compute cluster location
63  if (iv[i] > 9990) {
64  if constexpr (verbose)
65  alpaka::atomicAdd(acc, &noise, 1, alpaka::hierarchy::Threads{});
66  continue;
67  }
68  ALPAKA_ASSERT_ACC(iv[i] >= 0);
70  auto w = 1.f / ezt2[i];
71  alpaka::atomicAdd(acc, &zv[iv[i]], zt[i] * w, alpaka::hierarchy::Threads{});
72  alpaka::atomicAdd(acc, &wv[iv[i]], w, alpaka::hierarchy::Threads{});
73  }
74 
75  alpaka::syncBlockThreads(acc);
76  // reuse nn
78  bool const wv_cond = (wv[i] > 0.f);
79  if (not wv_cond) {
80  printf("ERROR: wv[%d] (%f) > 0.f failed\n", i, wv[i]);
81  // printing info on tracks associated to this vertex
82  for (auto trk_i = 0u; trk_i < nt; ++trk_i) {
83  if (iv[trk_i] != int(i)) {
84  continue;
85  }
86  printf(" iv[%d]=%d zt[%d]=%f ezt2[%d]=%f\n", trk_i, iv[trk_i], trk_i, zt[trk_i], trk_i, ezt2[trk_i]);
87  }
88  ALPAKA_ASSERT_ACC(false);
89  }
90 
91  zv[i] /= wv[i];
92  nn[i] = -1; // ndof
93  }
94  alpaka::syncBlockThreads(acc);
95 
96  // compute chi2
97  for (auto i : cms::alpakatools::uniform_elements(acc, nt)) {
98  if (iv[i] > 9990)
99  continue;
100 
101  auto c2 = zv[iv[i]] - zt[i];
102  c2 *= c2 / ezt2[i];
103  if (c2 > chi2Max) {
104  iv[i] = 9999;
105  continue;
106  }
109  }
110  alpaka::syncBlockThreads(acc);
111 
113  if (nn[i] > 0) {
114  wv[i] *= float(nn[i]) / chi2[i];
115  }
116  }
117  if constexpr (verbose) {
119  printf("found %d proto clusters ", foundClusters);
120  printf("and %d noise\n", noise);
121  }
122  }
123  }
124 
126  public:
127  ALPAKA_FN_ACC void operator()(Acc1D const& acc,
130  WsSoAView pws,
131  float chi2Max // for outlier rejection
132  ) const {
134  }
135  };
136 
137 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder
138 
139 #endif // RecoVertex_PixelVertexFinding_plugins_alpaka_fitVertices_h
ALPAKA_FN_ACC ALPAKA_FN_INLINE VtxSoAView & pdata
Definition: fitVertices.h:19
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:311
ALPAKA_ASSERT_ACC(nvFinal<=nvIntermediate)
T w() const
bool verbose
fitVertices(pdata, pws, maxChi2ForFirstFit)
ALPAKA_FN_ACC ALPAKA_FN_INLINE __attribute__((always_inline)) void fitVertices(Acc1D const &acc
ALPAKA_FN_ACC ALPAKA_FN_INLINE VtxSoAView TrkSoAView WsSoAView & pws
Definition: fitVertices.h:19
std::vector< Block > Blocks
Definition: Block.h:99
ALPAKA_FN_ACC void operator()(Acc1D const &acc, VtxSoAView pdata, TrkSoAView ptrkdata, WsSoAView pws, float chi2Max) const
Definition: fitVertices.h:127
::vertexFinder::PixelVertexWorkSpaceSoAView WsSoAView
Definition: sortByPt2.h:24
ALPAKA_FN_ACC ALPAKA_FN_INLINE VtxSoAView TrkSoAView & ptrkdata
Definition: fitVertices.h:19
ALPAKA_FN_ACC ALPAKA_FN_INLINE VtxSoAView TrkSoAView WsSoAView float chi2Max
Definition: fitVertices.h:23
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:61