File indexing completed on 2024-04-06 12:15:39
0001 #include <array>
0002 #include <string>
0003 #include <vector>
0004
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/Framework/interface/global/EDFilter.h"
0008 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "HeterogeneousCore/AlpakaInterface/interface/Backend.h"
0012
0013 class AlpakaBackendFilter : public edm::global::EDFilter<> {
0014 public:
0015 explicit AlpakaBackendFilter(edm::ParameterSet const& config)
0016 : producer_(consumes<unsigned short>(config.getParameter<edm::InputTag>("producer"))), backends_{} {
0017 for (auto const& backend : config.getParameter<std::vector<std::string>>("backends")) {
0018 backends_[static_cast<unsigned short>(cms::alpakatools::toBackend(backend))] = true;
0019 }
0020 }
0021
0022 bool filter(edm::StreamID sid, edm::Event& event, edm::EventSetup const& setup) const final {
0023 return backends_[event.get(producer_)];
0024 }
0025
0026 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0027
0028 private:
0029 const edm::EDGetTokenT<unsigned short> producer_;
0030 std::array<bool, static_cast<short>(cms::alpakatools::Backend::size)> backends_;
0031 };
0032
0033 void AlpakaBackendFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0034 edm::ParameterSetDescription desc;
0035 desc.add<edm::InputTag>("producer", edm::InputTag("alpakaBackendProducer", "backend"))
0036 ->setComment(
0037 "Use the 'backend' instance label to read the backend indicator that is implicitly produced by every alpaka "
0038 "EDProducer.");
0039 desc.add<std::vector<std::string>>("backends", {"SerialSync"})
0040 ->setComment("Valid backends are 'SerialSync', 'CudaAsync', 'ROCmAsync', and 'TbbAsync'.");
0041 descriptions.addWithDefaultLabel(desc);
0042 descriptions.setComment(
0043 "This EDFilter accepts events if the alpaka EDProducer 'producer' was run on a backend among those listed by the "
0044 "'backends' parameter.");
0045 }
0046
0047 #include "FWCore/Framework/interface/MakerMacros.h"
0048 DEFINE_FWK_MODULE(AlpakaBackendFilter);