CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackSoAHeterogeneousT.h
Go to the documentation of this file.
1 #ifndef CUDADataFormats_Track_TrackHeterogeneousT_H
2 #define CUDADataFormats_Track_TrackHeterogeneousT_H
3 
4 #include <string>
5 #include <algorithm>
6 
10 
12 
13 namespace pixelTrack {
14  enum class Quality : uint8_t { bad = 0, edup, dup, loose, strict, tight, highPurity, notQuality };
15  constexpr uint32_t qualitySize{uint8_t(Quality::notQuality)};
16  const std::string qualityName[qualitySize]{"bad", "edup", "dup", "loose", "strict", "tight", "highPurity"};
19  return static_cast<Quality>(qp);
20  }
21 } // namespace pixelTrack
22 
23 template <int32_t S>
25 public:
26  static constexpr int32_t stride() { return S; }
27 
29  using hindex_type = uint32_t;
31 
32  // Always check quality is at least loose!
33  // CUDA does not support enums in __lgc ...
34 private:
36 
37 public:
38  constexpr Quality quality(int32_t i) const { return (Quality)(quality_(i)); }
39  constexpr Quality &quality(int32_t i) { return (Quality &)(quality_(i)); }
40  constexpr Quality const *qualityData() const { return (Quality const *)(quality_.data()); }
41  constexpr Quality *qualityData() { return (Quality *)(quality_.data()); }
42 
43  // this is chi2/ndof as not necessarely all hits are used in the fit
45 
47 
48  constexpr int nHits(int i) const { return detIndices.size(i); }
49 
50  constexpr bool isTriplet(int i) const { return nLayers(i) == 3; }
51 
52  constexpr int computeNumberOfLayers(int32_t i) const {
53  // layers are in order and we assume tracks are either forward or backward
54  auto pdet = detIndices.begin(i);
55  int nl = 1;
56  auto ol = phase1PixelTopology::getLayer(*pdet);
57  for (; pdet < detIndices.end(i); ++pdet) {
58  auto il = phase1PixelTopology::getLayer(*pdet);
59  if (il != ol)
60  ++nl;
61  ol = il;
62  }
63  return nl;
64  }
65 
66  // State at the Beam spot
67  // phi,tip,1/pt,cotan(theta),zip
71  constexpr float charge(int32_t i) const { return std::copysign(1.f, stateAtBS.state(i)(2)); }
72  constexpr float phi(int32_t i) const { return stateAtBS.state(i)(0); }
73  constexpr float tip(int32_t i) const { return stateAtBS.state(i)(1); }
74  constexpr float zip(int32_t i) const { return stateAtBS.state(i)(4); }
75 
76  // state at the detector of the outermost hit
77  // representation to be decided...
78  // not yet filled on GPU
79  // TrajectoryStateSoA<S> stateAtOuterDet;
80 
83 };
84 
85 namespace pixelTrack {
86 
87 #ifdef GPU_SMALL_EVENTS
88  // kept for testing and debugging
89  constexpr uint32_t maxNumber() { return 2 * 1024; }
90 #else
91  // tested on MC events with 55-75 pileup events
92  constexpr uint32_t maxNumber() { return 32 * 1024; }
93 #endif
94 
98 
99 } // namespace pixelTrack
100 
101 #endif // CUDADataFormats_Track_TrackHeterogeneousT_H
eigenSoA::ScalarSoA< int8_t, S > nLayers
constexpr bool isTriplet(int i) const
constexpr int computeNumberOfLayers(int32_t i) const
Quality qualityByName(std::string const &name)
const std::string qualityName[qualitySize]
constexpr Quality * qualityData()
static constexpr int32_t stride()
constexpr Quality const * qualityData() const
constexpr uint32_t maxNumber()
constexpr float zip(int32_t i) const
eigenSoA::ScalarSoA< uint8_t, S > quality_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
constexpr float charge(int32_t i) const
eigenSoA::ScalarSoA< float, S > pt
eigenSoA::ScalarSoA< float, S > chi2
constexpr float tip(int32_t i) const
constexpr int nHits(int i) const
constexpr uint32_t qualitySize
TrajectoryStateSoAT< S > stateAtBS
constexpr uint8_t getLayer(uint32_t detId)
cms::cuda::OneToManyAssoc< hindex_type, S+1, 5 *S > HitContainer
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:97
constexpr float phi(int32_t i) const
constexpr Quality quality(int32_t i) const
__host__ __device__ constexpr Scalar * data()
Definition: eigenSoA.h:26
constexpr Quality & quality(int32_t i)
eigenSoA::ScalarSoA< float, S > eta