CMS 3D CMS Logo

RiemannFit.dev.cc
Go to the documentation of this file.
1 #include <cstdint>
2 
3 #include <alpaka/alpaka.hpp>
4 
12 
13 #include "HelixFit.h"
14 #include "CAStructures.h"
15 
16 template <typename TrackerTraits>
18 template <typename TrackerTraits>
20 template <typename TrackerTraits>
22 
24  using namespace alpaka;
25  using namespace cms::alpakatools;
26 
27  template <int N, typename TrackerTraits>
29  public:
30  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
31  ALPAKA_FN_ACC void operator()(TAcc const &acc,
32  Tuples<TrackerTraits> const *__restrict__ foundNtuplets,
34  uint32_t nHits,
36  pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const *__restrict__ cpeParams,
37  double *__restrict__ phits,
38  float *__restrict__ phits_ge,
39  double *__restrict__ pfast_fit,
40  uint32_t offset) const {
41  constexpr uint32_t hitsInFit = N;
42 
43  ALPAKA_ASSERT_ACC(hitsInFit <= nHits);
44 
48 
49  // look in bin for this hit multiplicity
50 
51 #ifdef RIEMANN_DEBUG
52  const uint32_t threadIdx(alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc)[0u]);
54  printf("%d Ntuple of size %d for %d hits to fit\n", tupleMultiplicity->size(nHits), nHits, hitsInFit);
55 #endif
56 
58  for (auto local_idx : cms::alpakatools::uniform_elements(acc, nt)) {
59  auto tuple_idx = local_idx + offset;
60  if (tuple_idx >= tupleMultiplicity->size(nHits))
61  break;
62 
63  // get it from the ntuple container (one to one to helix)
64  auto tkid = *(tupleMultiplicity->begin(nHits) + tuple_idx);
65  ALPAKA_ASSERT_ACC(static_cast<int>(tkid) < foundNtuplets->nOnes());
66 
67  ALPAKA_ASSERT_ACC(foundNtuplets->size(tkid) == nHits);
68 
69  riemannFit::Map3xNd<N> hits(phits + local_idx);
71  riemannFit::Map6xNf<N> hits_ge(phits_ge + local_idx);
72 
73  // Prepare data structure
74  auto const *hitId = foundNtuplets->begin(tkid);
75  for (unsigned int i = 0; i < hitsInFit; ++i) {
76  auto hit = hitId[i];
77  float ge[6];
78  cpeParams->detParams(hh[hit].detectorIndex()).frame.toGlobal(hh[hit].xerrLocal(), 0, hh[hit].yerrLocal(), ge);
79 
80  hits.col(i) << hh[hit].xGlobal(), hh[hit].yGlobal(), hh[hit].zGlobal();
81  hits_ge.col(i) << ge[0], ge[1], ge[2], ge[3], ge[4], ge[5];
82  }
84 
85 #ifdef RIEMANN_DEBUG
86  // any NaN value should cause the track to be rejected at a later stage
91 #endif
92  }
93  }
94  };
95 
96  template <int N, typename TrackerTraits>
98  public:
99  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
100  ALPAKA_FN_ACC void operator()(TAcc const &acc,
102  uint32_t nHits,
103  double bField,
104  double *__restrict__ phits,
105  float *__restrict__ phits_ge,
106  double *__restrict__ pfast_fit_input,
108  uint32_t offset) const {
111 
112  // same as above...
113 
114  // look in bin for this hit multiplicity
116  for (auto local_idx : cms::alpakatools::uniform_elements(acc, nt)) {
117  auto tuple_idx = local_idx + offset;
118  if (tuple_idx >= tupleMultiplicity->size(nHits))
119  break;
120 
121  riemannFit::Map3xNd<N> hits(phits + local_idx);
123  riemannFit::Map6xNf<N> hits_ge(phits_ge + local_idx);
124 
125  riemannFit::VectorNd<N> rad = (hits.block(0, 0, 2, N).colwise().norm());
126 
128  riemannFit::loadCovariance2D(acc, hits_ge, hits_cov);
129 
130  circle_fit[local_idx] =
131  riemannFit::circleFit(acc, hits.block(0, 0, 2, N), hits_cov, fast_fit, rad, bField, true);
132 
133 #ifdef RIEMANN_DEBUG
134 // auto tkid = *(tupleMultiplicity->begin(nHits) + tuple_idx);
135 // printf("kernelCircleFit circle.par(0,1,2): %d %f,%f,%f\n", tkid,
136 // circle_fit[local_idx].par(0), circle_fit[local_idx].par(1), circle_fit[local_idx].par(2));
137 #endif
138  }
139  }
140  };
141 
142  template <int N, typename TrackerTraits>
144  public:
145  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
146  ALPAKA_FN_ACC void operator()(TAcc const &acc,
148  uint32_t nHits,
149  double bField,
151  double *__restrict__ phits,
152  float *__restrict__ phits_ge,
153  double *__restrict__ pfast_fit_input,
154  riemannFit::CircleFit *__restrict__ circle_fit,
155  uint32_t offset) const {
158 
159  // same as above...
160 
161  // look in bin for this hit multiplicity
163  for (auto local_idx : cms::alpakatools::uniform_elements(acc, nt)) {
164  auto tuple_idx = local_idx + offset;
165  if (tuple_idx >= tupleMultiplicity->size(nHits))
166  break;
167 
168  // get it for the ntuple container (one to one to helix)
169  int32_t tkid = *(tupleMultiplicity->begin(nHits) + tuple_idx);
170 
171  riemannFit::Map3xNd<N> hits(phits + local_idx);
173  riemannFit::Map6xNf<N> hits_ge(phits_ge + local_idx);
174 
175  auto const &line_fit = riemannFit::lineFit(acc, hits, hits_ge, circle_fit[local_idx], fast_fit, bField, true);
176 
178 
180  circle_fit[local_idx].par,
181  circle_fit[local_idx].cov,
182  line_fit.par,
183  line_fit.cov,
184  1.f / float(bField),
185  tkid);
186  results_view[tkid].pt() = bField / std::abs(circle_fit[local_idx].par(2));
187  results_view[tkid].eta() = asinhf(line_fit.par(0));
188  results_view[tkid].chi2() = (circle_fit[local_idx].chi2 + line_fit.chi2) / (2 * N - 5);
189 
190 #ifdef RIEMANN_DEBUG
191  printf("kernelLineFit size %d for %d hits circle.par(0,1,2): %d %f,%f,%f\n",
192  N,
193  nHits,
194  tkid,
195  circle_fit[local_idx].par(0),
196  circle_fit[local_idx].par(1),
197  circle_fit[local_idx].par(2));
198  printf("kernelLineFit line.par(0,1): %d %f,%f\n", tkid, line_fit.par(0), line_fit.par(1));
199  printf("kernelLineFit chi2 cov %f/%f %e,%e,%e,%e,%e\n",
200  circle_fit[local_idx].chi2,
201  line_fit.chi2,
202  circle_fit[local_idx].cov(0, 0),
203  circle_fit[local_idx].cov(1, 1),
204  circle_fit[local_idx].cov(2, 2),
205  line_fit.cov(0, 0),
206  line_fit.cov(1, 1));
207 #endif
208  }
209  }
210  };
211 
212  template <typename TrackerTraits>
215  uint32_t nhits,
216  uint32_t maxNumberOfTuples,
217  Queue &queue) {
218  assert(tuples_);
219 
220  auto blockSize = 64;
221  auto numberOfBlocks = (maxNumberOfConcurrentFits_ + blockSize - 1) / blockSize;
222  const auto workDivTriplets = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
223  const auto workDivQuadsPenta = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks / 4, blockSize);
224 
225  // Fit internals
226  auto hitsDevice = cms::alpakatools::make_device_buffer<double[]>(
227  queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::Matrix3xNd<4>) / sizeof(double));
228  auto hits_geDevice = cms::alpakatools::make_device_buffer<float[]>(
229  queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::Matrix6x4f) / sizeof(float));
230  auto fast_fit_resultsDevice = cms::alpakatools::make_device_buffer<double[]>(
231  queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::Vector4d) / sizeof(double));
232  auto circle_fit_resultsDevice_holder =
233  cms::alpakatools::make_device_buffer<char[]>(queue, maxNumberOfConcurrentFits_ * sizeof(riemannFit::CircleFit));
234  riemannFit::CircleFit *circle_fit_resultsDevice_ =
235  (riemannFit::CircleFit *)(circle_fit_resultsDevice_holder.data());
236 
237  for (uint32_t offset = 0; offset < maxNumberOfTuples; offset += maxNumberOfConcurrentFits_) {
238  // triplets
239  alpaka::exec<Acc1D>(queue,
240  workDivTriplets,
242  tuples_,
243  tupleMultiplicity_,
244  3,
245  hv,
246  cpeParams,
247  hitsDevice.data(),
248  hits_geDevice.data(),
249  fast_fit_resultsDevice.data(),
250  offset);
251 
252  alpaka::exec<Acc1D>(queue,
253  workDivTriplets,
255  tupleMultiplicity_,
256  3,
257  bField_,
258  hitsDevice.data(),
259  hits_geDevice.data(),
260  fast_fit_resultsDevice.data(),
261  circle_fit_resultsDevice_,
262  offset);
263 
264  alpaka::exec<Acc1D>(queue,
265  workDivTriplets,
267  tupleMultiplicity_,
268  3,
269  bField_,
270  outputSoa_,
271  hitsDevice.data(),
272  hits_geDevice.data(),
273  fast_fit_resultsDevice.data(),
274  circle_fit_resultsDevice_,
275  offset);
276 
277  // quads
278  alpaka::exec<Acc1D>(queue,
279  workDivQuadsPenta,
281  tuples_,
282  tupleMultiplicity_,
283  4,
284  hv,
285  cpeParams,
286  hitsDevice.data(),
287  hits_geDevice.data(),
288  fast_fit_resultsDevice.data(),
289  offset);
290 
291  alpaka::exec<Acc1D>(queue,
292  workDivQuadsPenta,
294  tupleMultiplicity_,
295  4,
296  bField_,
297  hitsDevice.data(),
298  hits_geDevice.data(),
299  fast_fit_resultsDevice.data(),
300  circle_fit_resultsDevice_,
301  offset);
302 
303  alpaka::exec<Acc1D>(queue,
304  workDivQuadsPenta,
306  tupleMultiplicity_,
307  4,
308  bField_,
309  outputSoa_,
310  hitsDevice.data(),
311  hits_geDevice.data(),
312  fast_fit_resultsDevice.data(),
313  circle_fit_resultsDevice_,
314  offset);
315 
316  if (fitNas4_) {
317  // penta
318  alpaka::exec<Acc1D>(queue,
319  workDivQuadsPenta,
321  tuples_,
322  tupleMultiplicity_,
323  5,
324  hv,
325  cpeParams,
326  hitsDevice.data(),
327  hits_geDevice.data(),
328  fast_fit_resultsDevice.data(),
329  offset);
330 
331  alpaka::exec<Acc1D>(queue,
332  workDivQuadsPenta,
334  tupleMultiplicity_,
335  5,
336  bField_,
337  hitsDevice.data(),
338  hits_geDevice.data(),
339  fast_fit_resultsDevice.data(),
340  circle_fit_resultsDevice_,
341  offset);
342 
343  alpaka::exec<Acc1D>(queue,
344  workDivQuadsPenta,
346  tupleMultiplicity_,
347  5,
348  bField_,
349  outputSoa_,
350  hitsDevice.data(),
351  hits_geDevice.data(),
352  fast_fit_resultsDevice.data(),
353  circle_fit_resultsDevice_,
354  offset);
355  } else {
356  // penta all 5
357  alpaka::exec<Acc1D>(queue,
358  workDivQuadsPenta,
360  tuples_,
361  tupleMultiplicity_,
362  5,
363  hv,
364  cpeParams,
365  hitsDevice.data(),
366  hits_geDevice.data(),
367  fast_fit_resultsDevice.data(),
368  offset);
369 
370  alpaka::exec<Acc1D>(queue,
371  workDivQuadsPenta,
373  tupleMultiplicity_,
374  5,
375  bField_,
376  hitsDevice.data(),
377  hits_geDevice.data(),
378  fast_fit_resultsDevice.data(),
379  circle_fit_resultsDevice_,
380  offset);
381 
382  alpaka::exec<Acc1D>(queue,
383  workDivQuadsPenta,
385  tupleMultiplicity_,
386  5,
387  bField_,
388  outputSoa_,
389  hitsDevice.data(),
390  hits_geDevice.data(),
391  fast_fit_resultsDevice.data(),
392  circle_fit_resultsDevice_,
393  offset);
394  }
395  }
396  }
397 
398  template class HelixFit<pixelTopology::Phase1>;
399  template class HelixFit<pixelTopology::Phase2>;
401 
402 } // 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