File indexing completed on 2024-04-06 12:05:23
0001 #ifndef DataFormats_TrackSoA_interface_TracksSoA_h
0002 #define DataFormats_TrackSoA_interface_TracksSoA_h
0003
0004 #include <alpaka/alpaka.hpp>
0005
0006 #include <Eigen/Core>
0007
0008 #include "HeterogeneousCore/AlpakaInterface/interface/OneToManyAssoc.h"
0009 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0010 #include "DataFormats/SoATemplate/interface/SoALayout.h"
0011 #include "DataFormats/TrackSoA/interface/TrackDefinitions.h"
0012
0013 namespace reco {
0014
0015 template <typename TrackerTraits>
0016 struct TrackSoA {
0017 static constexpr int32_t S = TrackerTraits::maxNumberOfTuples;
0018 static constexpr int32_t H = TrackerTraits::avgHitsPerTrack;
0019
0020
0021 using Vector5f = Eigen::Matrix<float, 5, 1>;
0022 using Vector15f = Eigen::Matrix<float, 15, 1>;
0023 using Quality = pixelTrack::Quality;
0024
0025 using hindex_type = uint32_t;
0026
0027 using HitContainer = cms::alpakatools::OneToManyAssocSequential<hindex_type, S + 1, H * S>;
0028
0029 GENERATE_SOA_LAYOUT(Layout,
0030 SOA_COLUMN(Quality, quality),
0031 SOA_COLUMN(float, chi2),
0032 SOA_COLUMN(int8_t, nLayers),
0033 SOA_COLUMN(float, eta),
0034 SOA_COLUMN(float, pt),
0035
0036 SOA_EIGEN_COLUMN(Vector5f, state),
0037 SOA_EIGEN_COLUMN(Vector15f, covariance),
0038 SOA_SCALAR(int, nTracks),
0039 SOA_SCALAR(HitContainer, hitIndices),
0040 SOA_SCALAR(HitContainer, detIndices))
0041 };
0042
0043 template <typename TrackerTraits>
0044 using TrackLayout = typename reco::TrackSoA<TrackerTraits>::template Layout<>;
0045 template <typename TrackerTraits>
0046 using TrackSoAView = typename reco::TrackSoA<TrackerTraits>::template Layout<>::View;
0047 template <typename TrackerTraits>
0048 using TrackSoAConstView = typename reco::TrackSoA<TrackerTraits>::template Layout<>::ConstView;
0049
0050
0051
0052
0053
0054
0055 template <typename T>
0056 struct IsTrackSoAConstView : std::false_type {};
0057 template <>
0058 struct IsTrackSoAConstView<TrackSoAConstView<pixelTopology::Phase1>> : std::true_type {};
0059 template <>
0060 struct IsTrackSoAConstView<TrackSoAView<pixelTopology::Phase1>> : std::true_type {};
0061 template <>
0062 struct IsTrackSoAConstView<TrackSoAConstView<pixelTopology::Phase2>> : std::true_type {};
0063 template <>
0064 struct IsTrackSoAConstView<TrackSoAView<pixelTopology::Phase2>> : std::true_type {};
0065 template <>
0066 struct IsTrackSoAConstView<TrackSoAConstView<pixelTopology::HIonPhase1>> : std::true_type {};
0067 template <>
0068 struct IsTrackSoAConstView<TrackSoAView<pixelTopology::HIonPhase1>> : std::true_type {};
0069
0070 template <typename T>
0071 constexpr bool isTrackSoAConstView = IsTrackSoAConstView<T>::value;
0072
0073 template <typename ConstView, typename = std::enable_if_t<isTrackSoAConstView<ConstView>>>
0074 ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float charge(ConstView const& tracks, int32_t i) {
0075
0076 float v = tracks[i].state()(2);
0077 return float((0.0f < v) - (v < 0.0f));
0078 }
0079
0080 template <typename ConstView, typename = std::enable_if_t<isTrackSoAConstView<ConstView>>>
0081 ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float phi(ConstView const& tracks, int32_t i) {
0082 return tracks[i].state()(0);
0083 }
0084
0085 template <typename ConstView, typename = std::enable_if_t<isTrackSoAConstView<ConstView>>>
0086 ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float tip(ConstView const& tracks, int32_t i) {
0087 return tracks[i].state()(1);
0088 }
0089
0090 template <typename ConstView, typename = std::enable_if_t<isTrackSoAConstView<ConstView>>>
0091 ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const& tracks, int32_t i) {
0092 return tracks[i].state()(4);
0093 }
0094
0095 template <typename ConstView, typename = std::enable_if_t<isTrackSoAConstView<ConstView>>>
0096 ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr bool isTriplet(ConstView const& tracks, int32_t i) {
0097 return tracks[i].nLayers() == 3;
0098 }
0099
0100 }
0101
0102 #endif