Skip to content

Draft: Parse LHCb ID multi event containers

This MR provides a flexible and scalable manner to parse LHCb IDs from filtered objects in selection algorithms.

The problem this solves

Selection algorithms provide as outputs dev_decisions_t and dev_decisions_offsets_t. With these two outputs it is possible to traverse the selection output and determine which objects triggered the decision. The binary decision is stored as booleans with the following nested structure:

  • There are n events.
  • For each event, there are m objects (be it tracks, secondary vertices, etc.).

Even though this list of decisions is traversable, the LHCb IDs of the objects were not traversable in a unified manner. That's why the selreports added two algorithms to have a plain-ified list of LHCb IDs (count_long_track_hits, make_hits_container) and an additional output per selection algorithm that indicated which LHCb ID container in the plain list to refer to (host_lhcbid_container_t).

That solution has several problems. Mainly, the list of LHCb IDs is statically set to a group of containers and thus is not flexible without manually changing the algorithm, and the index had to be set manually as well to a growing number of containers. For "complex sequences" involving several LHCb ID containers of the same subdetectors this would become especially confusing. In other words, it is not scalable.

The proposed solution

With the creation of the new Event Model views in Allen, this MR proposes creating a LHCb ID container interface that should be implemented by each container with accessible LHCb IDs. As an example, this MR implements that for the VELO Event Model (which was recently merged into Allen), and for a sample velo_example_line. Thus, this is a realization of that proposal only for the VELO event model, and the UT, Forward, SVs, etc. containers would have to be developed as well if this proof-of-concept is agreed upon.

The ConsolidatedTypes.cuh now defines three interfaces:

  struct ILHCbIDSequence {
    virtual __host__ __device__ unsigned number_of_ids() const = 0;
    virtual __host__ __device__ unsigned id(const unsigned) const = 0;
    virtual __host__ __device__ ~ILHCbIDSequence() {}
  };

  struct ILHCbIDContainer {
    virtual __host__ __device__ unsigned number_of_id_sequences() const = 0;
    virtual __host__ __device__ const ILHCbIDSequence& id_sequence(const unsigned) const = 0;
    virtual __host__ __device__ ~ILHCbIDContainer() {}
  };

  struct IMultiEventLHCbIDContainer {
    virtual __host__ __device__ unsigned number_of_id_containers() const = 0;
    virtual __host__ __device__ const ILHCbIDContainer& id_container(const unsigned) const = 0;
    virtual __host__ __device__ ~IMultiEventLHCbIDContainer() {}
  };

An example of a specialization of ILHCbIDSequence is Allen::Views::Velo::Consolidated::Track:

      struct Track : Allen::ILHCbIDSequence {
        private:
          const Hits* m_hits = nullptr;
          unsigned m_track_index = 0;
          unsigned m_offset = 0;
          unsigned m_number_of_hits = 0;

        public:
          __host__ __device__ Track(
            const Hits* hits,
            const unsigned* offset_tracks,
            const unsigned* offset_track_hit_number,
            const unsigned track_index,
            const unsigned event_number) :
            m_hits(hits + event_number),
            m_track_index(track_index)
          {
            const auto offset_event = offset_track_hit_number + offset_tracks[event_number];
            m_offset = offset_event[track_index] - offset_event[0];
            m_number_of_hits = offset_event[track_index + 1] - offset_event[track_index];
          }

          __host__ __device__ unsigned track_index() const { return m_track_index; }

          __host__ __device__ unsigned number_of_hits() const { return m_number_of_hits; }

          __host__ __device__ Hit hit(const unsigned index) const {
            assert(m_hits != nullptr);
            assert(index < m_number_of_hits);
            return m_hits->hit(m_offset + index);
          }

          __host__ __device__ State state(const States& states_view) const { return states_view.state(m_track_index); }

          __host__ __device__ unsigned number_of_ids() const override { return number_of_hits(); }

          __host__ __device__ unsigned id(const unsigned index) const override
          {
            return hit(index).id();
          }
        };

The relevant methods above are number_of_ids and id.

An example algorithm that takes a IMultiEventLHCbIDContainer* and traverses it is sel_test_t. It takes a single input DEVICE_INPUT(dev_multiev_lhcb_id_container_t, Allen::IMultiEventLHCbIDContainer*) dev_multiev_lhcb_id_container; and has the following traversal implementation:

__global__ void sel_test_kernel(sel_test::Parameters parameters, const unsigned number_of_lines)
{
  if (threadIdx.x == 0) {
    for (unsigned id_container_index = 0; id_container_index < number_of_lines; ++id_container_index) {
      const Allen::IMultiEventLHCbIDContainer* multievent_lhcbidcont = parameters.dev_multiev_lhcb_id_container[id_container_index];

      printf("Number of id containers (events): %i\n", multievent_lhcbidcont->number_of_id_containers());
      for (unsigned event_number = 0; event_number < multievent_lhcbidcont->number_of_id_containers(); event_number++) {

        printf("Event %i:\n", event_number);
        const auto& id_container = multievent_lhcbidcont->id_container(event_number);
        for (unsigned cont = 0; cont < id_container.number_of_id_sequences(); ++cont) {

          printf("Sequence %u: ", cont);
          const auto& seq = id_container.id_sequence(cont);
          for (unsigned i = 0; i < seq.number_of_ids(); ++i) {
            printf("%u, ", seq.id(i));
          }
          printf("\n");
        }
      }
    }
  }
}

This, combined with the boolean decisions can be used to access the relevant LHCb IDs. The dev_multiev_lhcb_id_container_t container of sel_test_t has been generated by GatherSelections by concatenating together all the IMultiEventLHCbIDContainer* outputs from each line. For selection algorithms with a single decision for each event, the IMultiEventLHCbIDContainer* will be nullptr.

Note that this solution does not require the configuration to support assigning derived to base types. Each line is expected to have an input dev_lhcbid_container_t which implements a IMultiEventLHCbIDContainer specialization.

Where are views created

In order to avoid dynamic memory allocation (not supported in Allen) in-placement new is used. In order to avoid object slicing, the above interfaces return const &. These objects have to be guaranteed to exist to avoid unexpected behaviour, so they are created in a single __global__ function to that end:

/**
 * @brief Creates VELO views for hits, track, tracks and multieventtracks.
 */
__global__ void create_velo_views(velo_consolidate_tracks::Parameters parameters)
{
  const unsigned number_of_events = parameters.dev_number_of_events[0];
  const unsigned event_number = blockIdx.x;

  const auto event_tracks_offset = parameters.dev_offsets_all_velo_tracks[event_number];
  const auto event_number_of_tracks = parameters.dev_offsets_all_velo_tracks[event_number + 1] - event_tracks_offset;

  for (unsigned track_index = threadIdx.x; track_index < event_number_of_tracks; track_index += blockDim.x) {
    new (parameters.dev_velo_track_view + event_tracks_offset + track_index) Allen::Views::Velo::Consolidated::Track {
      parameters.dev_velo_hits_view,
      parameters.dev_offsets_all_velo_tracks,
      parameters.dev_offsets_velo_track_hit_number,
      track_index,
      event_number};
  }

  if (threadIdx.x == 0) {
    new (parameters.dev_velo_hits_view + event_number) Allen::Views::Velo::Consolidated::Hits {
      parameters.dev_velo_track_hits,
      parameters.dev_offsets_all_velo_tracks,
      parameters.dev_offsets_velo_track_hit_number,
      event_number,
      number_of_events};

    new (parameters.dev_velo_tracks_view + event_number) Allen::Views::Velo::Consolidated::Tracks {
      parameters.dev_velo_track_view, parameters.dev_offsets_all_velo_tracks, event_number};
  }

  if (blockIdx.x == 0 && threadIdx.x == 0) {
    new (parameters.dev_velo_multi_event_tracks_view)
      Allen::Views::Velo::Consolidated::MultiEventTracks {parameters.dev_velo_tracks_view, number_of_events};
  }
}
Edited by Daniel Campora Perez

Merge request reports