Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-09-07 02:51:59

0001 #include <cassert>
0002 
0003 #include "DataFormats/PortableTestObjects/interface/TestHostCollection.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/stream/EDAnalyzer.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0012 #include "FWCore/Utilities/interface/EDGetToken.h"
0013 #include "FWCore/Utilities/interface/InputTag.h"
0014 
0015 class TestAlpakaAnalyzer : public edm::stream::EDAnalyzer<> {
0016 public:
0017   TestAlpakaAnalyzer(edm::ParameterSet const& config)
0018       : source_{config.getParameter<edm::InputTag>("source")}, token_{consumes(source_)} {}
0019 
0020   void analyze(edm::Event const& event, edm::EventSetup const&) override {
0021     portabletest::TestHostCollection const& product = event.get(token_);
0022 
0023     auto const& view = product.const_view();
0024     for (int32_t i = 0; i < view.metadata().size(); ++i) {
0025       assert(view[i].id() == i);
0026     }
0027 
0028     edm::LogInfo msg("TestAlpakaAnalyzer");
0029     msg << source_.encode() << ".size() = " << view.metadata().size() << '\n';
0030     msg << "  data = " << product.buffer().data() << ",\n"
0031         << "  x    = " << view.metadata().addressOf_x() << ",\n"
0032         << "  y    = " << view.metadata().addressOf_y() << ",\n"
0033         << "  z    = " << view.metadata().addressOf_z() << ",\n"
0034         << "  id   = " << view.metadata().addressOf_id() << ",\n"
0035         << "  r    = " << view.metadata().addressOf_r() << '\n';
0036     msg << std::hex << "  [y - x] = 0x"
0037         << reinterpret_cast<intptr_t>(view.metadata().addressOf_y()) -
0038                reinterpret_cast<intptr_t>(view.metadata().addressOf_x())
0039         << "  [z - y] = 0x"
0040         << reinterpret_cast<intptr_t>(view.metadata().addressOf_z()) -
0041                reinterpret_cast<intptr_t>(view.metadata().addressOf_y())
0042         << "  [id - z] = 0x"
0043         << reinterpret_cast<intptr_t>(view.metadata().addressOf_id()) -
0044                reinterpret_cast<intptr_t>(view.metadata().addressOf_z())
0045         << "  [r - id] = 0x"
0046         << reinterpret_cast<intptr_t>(view.metadata().addressOf_r()) -
0047                reinterpret_cast<intptr_t>(view.metadata().addressOf_id());
0048   }
0049 
0050   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0051     edm::ParameterSetDescription desc;
0052     desc.add<edm::InputTag>("source");
0053     descriptions.addWithDefaultLabel(desc);
0054   }
0055 
0056 private:
0057   const edm::InputTag source_;
0058   const edm::EDGetTokenT<portabletest::TestHostCollection> token_;
0059 };
0060 
0061 #include "FWCore/Framework/interface/MakerMacros.h"
0062 DEFINE_FWK_MODULE(TestAlpakaAnalyzer);