CMS 3D CMS Logo

TrackUtilities.h
Go to the documentation of this file.
1 #ifndef DataFormats_Track_interface_alpaka_TrackUtilities_h
2 #define DataFormats_Track_interface_alpaka_TrackUtilities_h
3 
7 
8 // Methods that operate on View and ConstView of the TrackSoA, and cannot be class methods.
9 template <typename TrackerTraits>
10 struct TracksUtilities {
12  using TrackSoAConstView = typename reco::TrackSoA<TrackerTraits>::template Layout<>::ConstView;
14 
15  // state at the beam spot: { phi, tip, 1/pt, cotan(theta), zip }
16 
17  template <typename V3, typename M3, typename V2, typename M2>
18  ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE static constexpr void copyFromCircle(
19  TrackSoAView &tracks, V3 const &cp, M3 const &ccov, V2 const &lp, M2 const &lcov, float b, int32_t i) {
20  tracks[i].state() << cp.template cast<float>(), lp.template cast<float>();
21 
22  tracks[i].state()(2) = tracks[i].state()(2) * b;
23  auto cov = tracks[i].covariance();
24  cov(0) = ccov(0, 0);
25  cov(1) = ccov(0, 1);
26  cov(2) = b * float(ccov(0, 2));
27  cov(4) = cov(3) = 0;
28  cov(5) = ccov(1, 1);
29  cov(6) = b * float(ccov(1, 2));
30  cov(8) = cov(7) = 0;
31  cov(9) = b * b * float(ccov(2, 2));
32  cov(11) = cov(10) = 0;
33  cov(12) = lcov(0, 0);
34  cov(13) = lcov(0, 1);
35  cov(14) = lcov(1, 1);
36  }
37 
38  template <typename V5, typename M5>
39  ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE static constexpr void copyFromDense(TrackSoAView &tracks,
40  V5 const &v,
41  M5 const &cov,
42  int32_t i) {
43  tracks[i].state() = v.template cast<float>();
44  for (int j = 0, ind = 0; j < 5; ++j)
45  for (auto k = j; k < 5; ++k)
46  tracks[i].covariance()(ind++) = cov(j, k);
47  }
48 
49  template <typename V5, typename M5>
50  ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE static constexpr void copyToDense(const TrackSoAConstView &tracks,
51  V5 &v,
52  M5 &cov,
53  int32_t i) {
54  v = tracks[i].state().template cast<typename V5::Scalar>();
55  for (int j = 0, ind = 0; j < 5; ++j) {
56  cov(j, j) = tracks[i].covariance()(ind++);
57  for (auto k = j + 1; k < 5; ++k)
58  cov(k, j) = cov(j, k) = tracks[i].covariance()(ind++);
59  }
60  }
61 
62  ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE static constexpr int computeNumberOfLayers(const TrackSoAConstView &tracks,
63  int32_t i) {
64  auto pdet = tracks.detIndices().begin(i);
65  int nl = 1;
66  auto ol = pixelTopology::getLayer<TrackerTraits>(*pdet);
67  for (; pdet < tracks.detIndices().end(i); ++pdet) {
68  auto il = pixelTopology::getLayer<TrackerTraits>(*pdet);
69  if (il != ol)
70  ++nl;
71  ol = il;
72  }
73  return nl;
74  }
75 
76  ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE static constexpr int nHits(const TrackSoAConstView &tracks, int i) {
77  return tracks.detIndices().size(i);
78  }
79 };
80 
81 namespace pixelTrack {
82 
83  template <typename TrackerTraits, typename Enable = void>
84  struct QualityCutsT {};
85 
86  template <typename TrackerTraits>
87  struct QualityCutsT<TrackerTraits, pixelTopology::isPhase1Topology<TrackerTraits>> {
89  using TrackSoAConstView = typename reco::TrackSoA<TrackerTraits>::template Layout<>::ConstView;
90  float chi2Coeff[4];
91  float chi2MaxPt; // GeV
92  float chi2Scale;
93 
94  struct Region {
95  float maxTip; // cm
96  float minPt; // GeV
97  float maxZip; // cm
98  };
99 
100  Region triplet;
101  Region quadruplet;
102 
103  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isHP(const TrackSoAConstView &tracks, int nHits, int it) const {
104  // impose "region cuts" based on the fit results (phi, Tip, pt, cotan(theta)), Zip)
105  // default cuts:
106  // - for triplets: |Tip| < 0.3 cm, pT > 0.5 GeV, |Zip| < 12.0 cm
107  // - for quadruplets: |Tip| < 0.5 cm, pT > 0.3 GeV, |Zip| < 12.0 cm
108  // (see CAHitNtupletGeneratorGPU.cc)
109  auto const &region = (nHits > 3) ? quadruplet : triplet;
110  return (std::abs(reco::tip(tracks, it)) < region.maxTip) and (tracks.pt(it) > region.minPt) and
111  (std::abs(reco::zip(tracks, it)) < region.maxZip);
112  }
113 
114  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool strictCut(const TrackSoAConstView &tracks, int it) const {
115  auto roughLog = [](float x) {
116  // max diff [0.5,12] at 1.25 0.16143
117  // average diff 0.0662998
118  union IF {
119  uint32_t i;
120  float f;
121  };
122  IF z;
123  z.f = x;
124  uint32_t lsb = 1 < 21;
125  z.i += lsb;
126  z.i >>= 21;
127  auto f = z.i & 3;
128  int ex = int(z.i >> 2) - 127;
129 
130  // log2(1+0.25*f)
131  // averaged over bins
132  const float frac[4] = {0.160497f, 0.452172f, 0.694562f, 0.901964f};
133  return float(ex) + frac[f];
134  };
135 
136  float pt = std::min<float>(tracks.pt(it), chi2MaxPt);
137  float chi2Cut = chi2Scale * (chi2Coeff[0] + roughLog(pt) * chi2Coeff[1]);
138  if (tracks.chi2(it) >= chi2Cut) {
139 #ifdef NTUPLE_FIT_DEBUG
140  printf("Bad chi2 %d pt %f eta %f chi2 %f\n", it, tracks.pt(it), tracks.eta(it), tracks.chi2(it));
141 #endif
142  return true;
143  }
144  return false;
145  }
146  };
147 
148  template <typename TrackerTraits>
149  struct QualityCutsT<TrackerTraits, pixelTopology::isPhase2Topology<TrackerTraits>> {
151  using TrackSoAConstView = typename reco::TrackSoA<TrackerTraits>::template Layout<>::ConstView;
152 
153  float maxChi2;
154  float minPt;
155  float maxTip;
156  float maxZip;
157 
158  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isHP(const TrackSoAConstView &tracks, int nHits, int it) const {
159  return (std::abs(reco::tip(tracks, it)) < maxTip) and (tracks.pt(it) > minPt) and
160  (std::abs(reco::zip(tracks, it)) < maxZip);
161  }
162  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool strictCut(const TrackSoAConstView &tracks, int it) const {
163  return tracks.chi2(it) >= maxChi2;
164  }
165  };
166 
167 } // namespace pixelTrack
168 
169 // TODO: Should those be placed in the ALPAKA_ACCELERATOR_NAMESPACE
172 
173 #endif
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isHP(const TrackSoAConstView &tracks, int nHits, int it) const
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::ConstView TrackSoAConstView
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool strictCut(const TrackSoAConstView &tracks, int it) const
typename std::enable_if< std::is_base_of< Phase2, T >::value >::type isPhase2Topology
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::View TrackSoAView
minPt
Definition: PV_cfg.py:223
ALPAKA_FN_HOST_ACC static ALPAKA_FN_INLINE constexpr void copyToDense(const TrackSoAConstView &tracks, V5 &v, M5 &cov, int32_t i)
ALPAKA_FN_HOST_ACC static ALPAKA_FN_INLINE constexpr void copyFromDense(TrackSoAView &tracks, V5 const &v, M5 const &cov, int32_t i)
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::View TrackSoAView
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::ConstView TrackSoAConstView
ALPAKA_FN_HOST_ACC static ALPAKA_FN_INLINE constexpr int computeNumberOfLayers(const TrackSoAConstView &tracks, int32_t i)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float tip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:85
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::View TrackSoAView
typename std::enable_if< std::is_base_of< Phase1, T >::value >::type isPhase1Topology
ALPAKA_FN_HOST_ACC static ALPAKA_FN_INLINE constexpr void copyFromCircle(TrackSoAView &tracks, V3 const &cp, M3 const &ccov, V2 const &lp, M2 const &lcov, float b, int32_t i)
uint32_t hindex_type
Definition: TracksSoA.h:25
typename TrackSoA< TrackerTraits >::hindex_type hindex_type
double b
Definition: hdecay.h:120
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isHP(const TrackSoAConstView &tracks, int nHits, int it) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool strictCut(const TrackSoAConstView &tracks, int it) const
ALPAKA_FN_HOST_ACC static ALPAKA_FN_INLINE constexpr int nHits(const TrackSoAConstView &tracks, int i)
float x
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
typename TrackSoA< TrackerTraits >::template TrackSoALayout<>::ConstView TrackSoAConstView