CMS 3D CMS Logo

RiemannFit.dev.cc
Go to the documentation of this file.
1 #include <cstdint>
2 
3 #include <alpaka/alpaka.hpp>
4 
11 
12 #include "HelixFit.h"
13 #include "CAStructures.h"
14 
15 template <typename TrackerTraits>
17 template <typename TrackerTraits>
19 template <typename TrackerTraits>
21 
23  using namespace alpaka;
24  using namespace cms::alpakatools;
25 
26  template <int N, typename TrackerTraits>
28  public:
29  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
30  ALPAKA_FN_ACC void operator()(TAcc const &acc,
31  Tuples<TrackerTraits> const *__restrict__ foundNtuplets,
33  uint32_t nHits,
35  pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const *__restrict__ cpeParams,
36  double *__restrict__ phits,
37  float *__restrict__ phits_ge,
38  double *__restrict__ pfast_fit,
39  uint32_t offset) const {
40  constexpr uint32_t hitsInFit = N;
41 
42  ALPAKA_ASSERT_ACC(hitsInFit <= nHits);
43 
47 
48  // look in bin for this hit multiplicity
49 
50 #ifdef RIEMANN_DEBUG
51  const uint32_t threadIdx(alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc)[0u]);
53  printf("%d Ntuple of size %d for %d hits to fit\n", tupleMultiplicity->size(nHits), nHits, hitsInFit);
54 #endif
55 
57  for (auto local_idx : cms::alpakatools::uniform_elements(acc, nt)) {
58  auto tuple_idx = local_idx + offset;
59  if (tuple_idx >= tupleMultiplicity->size(nHits))
60  break;
61 
62  // get it from the ntuple container (one to one to helix)
63  auto tkid = *(tupleMultiplicity->begin(nHits) + tuple_idx);
64  ALPAKA_ASSERT_ACC(static_cast<int>(tkid) < foundNtuplets->nOnes());
65 
66  ALPAKA_ASSERT_ACC(foundNtuplets->size(tkid) == nHits);
67 
68  riemannFit::Map3xNd<N> hits(phits + local_idx);
70  riemannFit::Map6xNf<N> hits_ge(phits_ge + local_idx);
71 
72  // Prepare data structure
73  auto const *hitId = foundNtuplets->begin(tkid);
74  for (unsigned int i = 0; i < hitsInFit; ++i) {
75  auto hit = hitId[i];
76  float ge[6];
77  cpeParams->detParams(hh[hit].detectorIndex()).frame.toGlobal(hh[hit].xerrLocal(), 0, hh[hit].yerrLocal(), ge);
78 
79  hits.col(i) << hh[hit].xGlobal(), hh[hit].yGlobal(), hh[hit].zGlobal();
80  hits_ge.col(i) << ge[0], ge[1], ge[2], ge[3], ge[4], ge[5];
81  }
83 
84 #ifdef RIEMANN_DEBUG
85  // any NaN value should cause the track to be rejected at a later stage
90 #endif
91  }
92  }
93  };
94 
95  template <int N, typename TrackerTraits>
97  public:
98  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
99  ALPAKA_FN_ACC void operator()(TAcc const &acc,
101  uint32_t nHits,
102  double bField,
103  double *__restrict__ phits,
104  float *__restrict__ phits_ge,
105  double *__restrict__ pfast_fit_input,
107  uint32_t offset) const {
110 
111  // same as above...
112 
113  // look in bin for this hit multiplicity
115  for (auto local_idx : cms::alpakatools::uniform_elements(acc, nt)) {
116  auto tuple_idx = local_idx + offset;
117  if (tuple_idx >= tupleMultiplicity->size(nHits))
118  break;
119 
120  riemannFit::Map3xNd<N> hits(phits + local_idx);
122  riemannFit::Map6xNf<N> hits_ge(phits_ge + local_idx);
123 
124  riemannFit::VectorNd<N> rad = (hits.block(0, 0, 2, N).colwise().norm());
125 
127  riemannFit::loadCovariance2D(acc, hits_ge, hits_cov);
128 
129  circle_fit[local_idx] =
130  riemannFit::circleFit(acc, hits.block(0, 0, 2, N), hits_cov, fast_fit, rad, bField, true);
131 
132 #ifdef RIEMANN_DEBUG
133 // auto tkid = *(tupleMultiplicity->begin(nHits) + tuple_idx);
134 // printf("kernelCircleFit circle.par(0,1,2): %d %f,%f,%f\n", tkid,
135 // circle_fit[local_idx].par(0), circle_fit[local_idx].par(1), circle_fit[local_idx].par(2));
136 #endif
137  }
138  }
139  };
140 
141  template <int N, typename TrackerTraits>
143  public:
144  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
145  ALPAKA_FN_ACC void operator()(TAcc const &acc,
147  uint32_t nHits,
148  double bField,
150  double *__restrict__ phits,
151  float *__restrict__ phits_ge,
152  double *__restrict__ pfast_fit_input,
153  riemannFit::CircleFit *__restrict__ circle_fit,
154  uint32_t offset) const {
157 
158  // same as above...
159 
160  // look in bin for this hit multiplicity
162  for (auto local_idx : cms::alpakatools::uniform_elements(acc, nt)) {
163  auto tuple_idx = local_idx + offset;
164  if (tuple_idx >= tupleMultiplicity->size(nHits))
165  break;
166 
167  // get it for the ntuple container (one to one to helix)
168  int32_t tkid = *(tupleMultiplicity->begin(nHits) + tuple_idx);
169 
170  riemannFit::Map3xNd<N> hits(phits + local_idx);
172  riemannFit::Map6xNf<N> hits_ge(phits_ge + local_idx);
173 
174  auto const &line_fit = riemannFit::lineFit(acc, hits, hits_ge, circle_fit[local_idx], fast_fit, bField, true);
175 
177 
179  circle_fit[local_idx].par,
180  circle_fit[local_idx].cov,
181  line_fit.par,
182  line_fit.cov,
183  1.f / float(bField),
184  tkid);
185  results_view[tkid].pt() = bField / std::abs(circle_fit[local_idx].par(2));
186  results_view[tkid].eta() = asinhf(line_fit.par(0));
187  results_view[tkid].chi2() = (circle_fit[local_idx].chi2 + line_fit.chi2) / (2 * N - 5);
188 
189 #ifdef RIEMANN_DEBUG
190  printf("kernelLineFit size %d for %d hits circle.par(0,1,2): %d %f,%f,%f\n",
191  N,
192  nHits,
193  tkid,
194  circle_fit[local_idx].par(0),
195  circle_fit[local_idx].par(1),
196  circle_fit[local_idx].par(2));
197  printf("kernelLineFit line.par(0,1): %d %f,%f\n", tkid, line_fit.par(0), line_fit.par(1));
198  printf("kernelLineFit chi2 cov %f/%f %e,%e,%e,%e,%e\n",
199  circle_fit[local_idx].chi2,
200  line_fit.chi2,
201  circle_fit[local_idx].cov(0, 0),
202  circle_fit[local_idx].cov(1, 1),
203  circle_fit[local_idx].cov(2, 2),
204  line_fit.cov(0, 0),
205  line_fit.cov(1, 1));
206 #endif
207  }
208  }
209  };
210 
211  template <typename TrackerTraits>
214  uint32_t nhits,
215  uint32_t maxNumberOfTuples,
216  Queue &queue) {
217  assert(tuples_);
218 
219  auto blockSize = 64;
220  auto numberOfBlocks = (maxNumberOfConcurrentFits_ + blockSize - 1) / blockSize;
221  const auto workDivTriplets = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
222  const auto workDivQuadsPenta = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks / 4, blockSize);
223 
224  // Fit internals
225  auto hitsDevice = cms::alpakatools::make_device_buffer<double[]>(
226  queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::Matrix3xNd<4>) / sizeof(double));
227  auto hits_geDevice = cms::alpakatools::make_device_buffer<float[]>(
228  queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::Matrix6x4f) / sizeof(float));
229  auto fast_fit_resultsDevice = cms::alpakatools::make_device_buffer<double[]>(
230  queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::Vector4d) / sizeof(double));
231  auto circle_fit_resultsDevice_holder =
232  cms::alpakatools::make_device_buffer<char[]>(queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::CircleFit));
233  riemannFit::CircleFit *circle_fit_resultsDevice_ =
234  (riemannFit::CircleFit *)(circle_fit_resultsDevice_holder.data());
235 
236  for (uint32_t offset = 0; offset < maxNumberOfTuples; offset += maxNumberOfConcurrentFits_) {
237  // triplets
238  alpaka::exec<Acc1D>(queue,
239  workDivTriplets,
241  tuples_,
242  tupleMultiplicity_,
243  3,
244  hv,
245  cpeParams,
246  hitsDevice.data(),
247  hits_geDevice.data(),
248  fast_fit_resultsDevice.data(),
249  offset);
250 
251  alpaka::exec<Acc1D>(queue,
252  workDivTriplets,
254  tupleMultiplicity_,
255  3,
256  bField_,
257  hitsDevice.data(),
258  hits_geDevice.data(),
259  fast_fit_resultsDevice.data(),
260  circle_fit_resultsDevice_,
261  offset);
262 
263  alpaka::exec<Acc1D>(queue,
264  workDivTriplets,
266  tupleMultiplicity_,
267  3,
268  bField_,
269  outputSoa_,
270  hitsDevice.data(),
271  hits_geDevice.data(),
272  fast_fit_resultsDevice.data(),
273  circle_fit_resultsDevice_,
274  offset);
275 
276  // quads
277  alpaka::exec<Acc1D>(queue,
278  workDivQuadsPenta,
280  tuples_,
281  tupleMultiplicity_,
282  4,
283  hv,
284  cpeParams,
285  hitsDevice.data(),
286  hits_geDevice.data(),
287  fast_fit_resultsDevice.data(),
288  offset);
289 
290  alpaka::exec<Acc1D>(queue,
291  workDivQuadsPenta,
293  tupleMultiplicity_,
294  4,
295  bField_,
296  hitsDevice.data(),
297  hits_geDevice.data(),
298  fast_fit_resultsDevice.data(),
299  circle_fit_resultsDevice_,
300  offset);
301 
302  alpaka::exec<Acc1D>(queue,
303  workDivQuadsPenta,
305  tupleMultiplicity_,
306  4,
307  bField_,
308  outputSoa_,
309  hitsDevice.data(),
310  hits_geDevice.data(),
311  fast_fit_resultsDevice.data(),
312  circle_fit_resultsDevice_,
313  offset);
314 
315  if (fitNas4_) {
316  // penta
317  alpaka::exec<Acc1D>(queue,
318  workDivQuadsPenta,
320  tuples_,
321  tupleMultiplicity_,
322  5,
323  hv,
324  cpeParams,
325  hitsDevice.data(),
326  hits_geDevice.data(),
327  fast_fit_resultsDevice.data(),
328  offset);
329 
330  alpaka::exec<Acc1D>(queue,
331  workDivQuadsPenta,
333  tupleMultiplicity_,
334  5,
335  bField_,
336  hitsDevice.data(),
337  hits_geDevice.data(),
338  fast_fit_resultsDevice.data(),
339  circle_fit_resultsDevice_,
340  offset);
341 
342  alpaka::exec<Acc1D>(queue,
343  workDivQuadsPenta,
345  tupleMultiplicity_,
346  5,
347  bField_,
348  outputSoa_,
349  hitsDevice.data(),
350  hits_geDevice.data(),
351  fast_fit_resultsDevice.data(),
352  circle_fit_resultsDevice_,
353  offset);
354  } else {
355  // penta all 5
356  alpaka::exec<Acc1D>(queue,
357  workDivQuadsPenta,
359  tuples_,
360  tupleMultiplicity_,
361  5,
362  hv,
363  cpeParams,
364  hitsDevice.data(),
365  hits_geDevice.data(),
366  fast_fit_resultsDevice.data(),
367  offset);
368 
369  alpaka::exec<Acc1D>(queue,
370  workDivQuadsPenta,
372  tupleMultiplicity_,
373  5,
374  bField_,
375  hitsDevice.data(),
376  hits_geDevice.data(),
377  fast_fit_resultsDevice.data(),
378  circle_fit_resultsDevice_,
379  offset);
380 
381  alpaka::exec<Acc1D>(queue,
382  workDivQuadsPenta,
384  tupleMultiplicity_,
385  5,
386  bField_,
387  outputSoa_,
388  hitsDevice.data(),
389  hits_geDevice.data(),
390  fast_fit_resultsDevice.data(),
391  circle_fit_resultsDevice_,
392  offset);
393  }
394  }
395  }
396 
397  template class HelixFit<pixelTopology::Phase1>;
398  template class HelixFit<pixelTopology::Phase2>;
400 
401 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
const dim3 threadIdx
Definition: cudaCompat.h:29
Eigen::Vector4d Vector4d
Definition: FitResult.h:12
def isnan(num)
constexpr DetParams const &__restrict__ detParams(int i) const
ALPAKA_FN_ACC auto uniform_elements(TAcc const &acc, TArgs... args)
Definition: workdivision.h:311
ALPAKA_FN_ACC ALPAKA_FN_INLINE LineFit lineFit(const TAcc &acc, const M3xN &hits, const M6xN &hits_ge, const CircleFit &circle, const V4 &fast_fit, const double bField, const bool error)
Perform an ordinary least square fit in the s-z plane to compute the parameters cotTheta and Zip...
Definition: RiemannFit.h:796
ALPAKA_FN_ACC void loadCovariance2D(const TAcc &acc, M6xNf const &ge, M2Nd &hits_cov)
Definition: FitUtils.h:126
ALPAKA_FN_ACC void operator()(TAcc const &acc, Tuples< TrackerTraits > const *__restrict__ foundNtuplets, TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity, uint32_t nHits, TrackingRecHitSoAConstView< TrackerTraits > hh, pixelCPEforDevice::ParamsOnDeviceT< TrackerTraits > const *__restrict__ cpeParams, double *__restrict__ phits, float *__restrict__ phits_ge, double *__restrict__ pfast_fit, uint32_t offset) const
Eigen::Matrix< float, 6, 4 > Matrix6x4f
Definition: HelixFit.h:25
Vector3d par
parameter: (X0,Y0,R)
Definition: FitResult.h:24
assert(be >=bs)
typename reco::TrackSoA< TrackerTraits >::HitContainer Tuples
TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity
ALPAKA_FN_ACC ALPAKA_FN_INLINE void fastFit(const TAcc &acc, const M3xN &hits, V4 &result)
A very fast helix fit: it fits a circle by three points (first, middle and last point) and a line by ...
Definition: RiemannFit.h:385
Eigen::Matrix< double, N, 1 > VectorNd
Definition: FitUtils.h:37
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ phits
Eigen::Map< Vector4d, 0, Eigen::InnerStride< stride > > Map4d
Definition: HelixFit.h:39
Eigen::Map< Matrix3xNd< N >, 0, Eigen::Stride< 3 *stride, stride > > Map3xNd
Definition: HelixFit.h:32
uint32_t double double *__restrict__ float *__restrict__ double *__restrict__ pfast_fit_input
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Eigen::Matrix< double, 3, N > Matrix3xNd
Definition: HelixFit.h:30
ALPAKA_FN_ACC void operator()(TAcc const &acc, TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity, uint32_t nHits, double bField, OutputSoAView< TrackerTraits > results_view, double *__restrict__ phits, float *__restrict__ phits_ge, double *__restrict__ pfast_fit_input, riemannFit::CircleFit *__restrict__ circle_fit, uint32_t offset) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE void fromCircleToPerigee(const TAcc &acc, CircleFit &circle)
Transform circle parameter from (X0,Y0,R) to (phi,Tip,q/R) and consequently covariance matrix...
Definition: FitUtils.h:239
double OutputSoAView< TrackerTraits > results_view
Eigen::Map< Matrix6xNf< N >, 0, Eigen::Stride< 6 *stride, stride > > Map6xNf
Definition: HelixFit.h:37
ALPAKA_FN_ACC ALPAKA_FN_INLINE void uint32_t const uint32_t CACellT< TrackerTraits > uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > hh
int nt
Definition: AMPTWrapper.h:42
ALPAKA_FN_ACC void operator()(TAcc const &acc, TupleMultiplicity< TrackerTraits > const *__restrict__ tupleMultiplicity, uint32_t nHits, double bField, double *__restrict__ phits, float *__restrict__ phits_ge, double *__restrict__ pfast_fit_input, riemannFit::CircleFit *circle_fit, uint32_t offset) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE void const M3xN const V4 & fast_fit
Definition: BrokenLine.h:158
Eigen::Matrix< double, 2 *N, 2 *N > Matrix2Nd
Definition: FitUtils.h:25
#define N
Definition: blowfish.cc:9
typename TrackingRecHitSoA< TrackerTraits >::template TrackingRecHitSoALayout<>::ConstView TrackingRecHitSoAConstView
ALPAKA_FN_ACC constexpr bool once_per_grid(TAcc const &acc)
uint32_t double double *__restrict__ float *__restrict__ double *__restrict__ riemannFit::CircleFit * circle_fit
void launchRiemannKernels(const HitConstView &hv, ParamsOnDevice const *cpeParams, uint32_t nhits, uint32_t maxNumberOfTuples, Queue &queue)
ALPAKA_FN_ACC ALPAKA_FN_INLINE CircleFit circleFit(const TAcc &acc, const M2xN &hits2D, const Matrix2Nd< N > &hits_cov2D, const V4 &fast_fit, const VectorNd< N > &rad, const double bField, const bool error)
Fit a generic number of 2D points with a circle using Riemann-Chernov algorithm. Covariance matrix of...
Definition: RiemannFit.h:465
static constexpr __host__ __device__ void copyFromCircle(TrackSoAView &tracks, V3 const &cp, M3 const &ccov, V2 const &lp, M2 const &lcov, float b, int32_t i)
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ float *__restrict__ phits_ge
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ float *__restrict__ double *__restrict__ pfast_fit
reco::TrackSoAView< TrackerTraits > OutputSoAView
constexpr uint32_t maxNumberOfConcurrentFits
Definition: HelixFit.h:21
typename reco::TrackSoA< TrackerTraits >::template Layout<>::View TrackSoAView
Definition: TracksSoA.h:45