Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-03-29 07:44:12

0001 #ifndef CUDADataFormats_Track_TrackHeterogeneousT_H
0002 #define CUDADataFormats_Track_TrackHeterogeneousT_H
0003 
0004 #include <string>
0005 #include <algorithm>
0006 
0007 #include "CUDADataFormats/Track/interface/TrajectoryStateSoAT.h"
0008 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0009 #include "HeterogeneousCore/CUDAUtilities/interface/HistoContainer.h"
0010 
0011 #include "CUDADataFormats/Common/interface/HeterogeneousSoA.h"
0012 
0013 namespace pixelTrack {
0014   enum class Quality : uint8_t { bad = 0, edup, dup, loose, strict, tight, highPurity, notQuality };
0015   constexpr uint32_t qualitySize{uint8_t(Quality::notQuality)};
0016   const std::string qualityName[qualitySize]{"bad", "edup", "dup", "loose", "strict", "tight", "highPurity"};
0017   inline Quality qualityByName(std::string const &name) {
0018     auto qp = std::find(qualityName, qualityName + qualitySize, name) - qualityName;
0019     return static_cast<Quality>(qp);
0020   }
0021 }  // namespace pixelTrack
0022 
0023 template <int32_t S>
0024 class TrackSoAHeterogeneousT {
0025 public:
0026   static constexpr int32_t stride() { return S; }
0027 
0028   using Quality = pixelTrack::Quality;
0029   using hindex_type = uint32_t;
0030   using HitContainer = cms::cuda::OneToManyAssoc<hindex_type, S + 1, 5 * S>;
0031 
0032   // Always check quality is at least loose!
0033   // CUDA does not support enums  in __lgc ...
0034 private:
0035   eigenSoA::ScalarSoA<uint8_t, S> quality_;
0036 
0037 public:
0038   constexpr Quality quality(int32_t i) const { return (Quality)(quality_(i)); }
0039   constexpr Quality &quality(int32_t i) { return (Quality &)(quality_(i)); }
0040   constexpr Quality const *qualityData() const { return (Quality const *)(quality_.data()); }
0041   constexpr Quality *qualityData() { return (Quality *)(quality_.data()); }
0042 
0043   // this is chi2/ndof as not necessarely all hits are used in the fit
0044   eigenSoA::ScalarSoA<float, S> chi2;
0045 
0046   eigenSoA::ScalarSoA<int8_t, S> nLayers;
0047 
0048   constexpr int nTracks() const { return nTracks_; }
0049   constexpr void setNTracks(int n) { nTracks_ = n; }
0050 
0051   constexpr int nHits(int i) const { return detIndices.size(i); }
0052 
0053   constexpr bool isTriplet(int i) const { return nLayers(i) == 3; }
0054 
0055   constexpr int computeNumberOfLayers(int32_t i) const {
0056     // layers are in order and we assume tracks are either forward or backward
0057     auto pdet = detIndices.begin(i);
0058     int nl = 1;
0059     auto ol = phase1PixelTopology::getLayer(*pdet);
0060     for (; pdet < detIndices.end(i); ++pdet) {
0061       auto il = phase1PixelTopology::getLayer(*pdet);
0062       if (il != ol)
0063         ++nl;
0064       ol = il;
0065     }
0066     return nl;
0067   }
0068 
0069   // State at the Beam spot
0070   // phi,tip,1/pt,cotan(theta),zip
0071   TrajectoryStateSoAT<S> stateAtBS;
0072   eigenSoA::ScalarSoA<float, S> eta;
0073   eigenSoA::ScalarSoA<float, S> pt;
0074   constexpr float charge(int32_t i) const { return std::copysign(1.f, stateAtBS.state(i)(2)); }
0075   constexpr float phi(int32_t i) const { return stateAtBS.state(i)(0); }
0076   constexpr float tip(int32_t i) const { return stateAtBS.state(i)(1); }
0077   constexpr float zip(int32_t i) const { return stateAtBS.state(i)(4); }
0078 
0079   // state at the detector of the outermost hit
0080   // representation to be decided...
0081   // not yet filled on GPU
0082   // TrajectoryStateSoA<S> stateAtOuterDet;
0083 
0084   HitContainer hitIndices;
0085   HitContainer detIndices;
0086 
0087 private:
0088   int nTracks_;
0089 };
0090 
0091 namespace pixelTrack {
0092 
0093 #ifdef GPU_SMALL_EVENTS
0094   // kept for testing and debugging
0095   constexpr uint32_t maxNumber() { return 2 * 1024; }
0096 #else
0097   // tested on MC events with 55-75 pileup events
0098   constexpr uint32_t maxNumber() { return 32 * 1024; }
0099 #endif
0100 
0101   using TrackSoA = TrackSoAHeterogeneousT<maxNumber()>;
0102   using TrajectoryState = TrajectoryStateSoAT<maxNumber()>;
0103   using HitContainer = TrackSoA::HitContainer;
0104 
0105 }  // namespace pixelTrack
0106 
0107 #endif  // CUDADataFormats_Track_TrackHeterogeneousT_H