CMS 3D CMS Logo

TestAlpakaAnalyzer.cc
Go to the documentation of this file.
1 #include <cassert>
2 
3 #include "DataFormats/PortableTestObjects/interface/TestHostCollection.h"
14 
15 namespace {
16 
17  template <typename T>
18  class Column {
19  public:
20  Column(T const* data, size_t size) : data_(data), size_(size) {}
21 
22  void print(std::ostream& out) const {
23  std::stringstream buffer;
24  buffer << "{ ";
25  if (size_ > 0) {
26  buffer << data_[0];
27  }
28  if (size_ > 1) {
29  buffer << ", " << data_[1];
30  }
31  if (size_ > 2) {
32  buffer << ", " << data_[2];
33  }
34  if (size_ > 3) {
35  buffer << ", ...";
36  }
37  buffer << '}';
38  out << buffer.str();
39  }
40 
41  private:
42  T const* const data_;
43  size_t const size_;
44  };
45 
46  template <typename T>
47  std::ostream& operator<<(std::ostream& out, Column<T> const& column) {
48  column.print(out);
49  return out;
50  }
51 
52  template <typename T>
53  void checkViewAddresses(T const& view) {
54  assert(view.metadata().addressOf_x() == view.x());
55  assert(view.metadata().addressOf_x() == &view.x(0));
56  assert(view.metadata().addressOf_x() == &view[0].x());
57  assert(view.metadata().addressOf_y() == view.y());
58  assert(view.metadata().addressOf_y() == &view.y(0));
59  assert(view.metadata().addressOf_y() == &view[0].y());
60  assert(view.metadata().addressOf_z() == view.z());
61  assert(view.metadata().addressOf_z() == &view.z(0));
62  assert(view.metadata().addressOf_z() == &view[0].z());
63  assert(view.metadata().addressOf_id() == view.id());
64  assert(view.metadata().addressOf_id() == &view.id(0));
65  assert(view.metadata().addressOf_id() == &view[0].id());
66  assert(view.metadata().addressOf_m() == view.m());
67  assert(view.metadata().addressOf_m() == &view.m(0).coeffRef(0, 0));
68  assert(view.metadata().addressOf_m() == &view[0].m().coeffRef(0, 0));
69  assert(view.metadata().addressOf_r() == &view.r());
70  //assert(view.metadata().addressOf_r() == &view.r(0)); // cannot access a scalar with an index
71  //assert(view.metadata().addressOf_r() == &view[0].r()); // cannot access a scalar via a SoA row-like accessor
72  }
73 
74 } // namespace
75 
77 public:
79  : source_{config.getParameter<edm::InputTag>("source")}, token_{consumes(source_)} {}
80 
81  void analyze(edm::Event const& event, edm::EventSetup const&) override {
82  portabletest::TestHostCollection const& product = event.get(token_);
83  auto const& view = product.const_view();
84  auto& mview = product.view();
85  auto const& cmview = product.view();
86 
87  {
88  edm::LogInfo msg("TestAlpakaAnalyzer");
89  msg << source_.encode() << ".size() = " << view.metadata().size() << '\n';
90  msg << " data @ " << product.buffer().data() << ",\n"
91  << " x @ " << view.metadata().addressOf_x() << " = " << Column(view.x(), view.metadata().size()) << ",\n"
92  << " y @ " << view.metadata().addressOf_y() << " = " << Column(view.y(), view.metadata().size()) << ",\n"
93  << " z @ " << view.metadata().addressOf_z() << " = " << Column(view.z(), view.metadata().size()) << ",\n"
94  << " id @ " << view.metadata().addressOf_id() << " = " << Column(view.id(), view.metadata().size())
95  << ",\n"
96  << " r @ " << view.metadata().addressOf_r() << " = " << view.r() << '\n'
97  << " m @ " << view.metadata().addressOf_m() << " = { ... {" << view[1].m()(1, Eigen::all)
98  << " } ... } \n";
99  msg << std::hex << " [y - x] = 0x"
100  << reinterpret_cast<intptr_t>(view.metadata().addressOf_y()) -
101  reinterpret_cast<intptr_t>(view.metadata().addressOf_x())
102  << " [z - y] = 0x"
103  << reinterpret_cast<intptr_t>(view.metadata().addressOf_z()) -
104  reinterpret_cast<intptr_t>(view.metadata().addressOf_y())
105  << " [id - z] = 0x"
106  << reinterpret_cast<intptr_t>(view.metadata().addressOf_id()) -
107  reinterpret_cast<intptr_t>(view.metadata().addressOf_z())
108  << " [r - id] = 0x"
109  << reinterpret_cast<intptr_t>(view.metadata().addressOf_r()) -
110  reinterpret_cast<intptr_t>(view.metadata().addressOf_id())
111  << " [m - r] = 0x"
112  << reinterpret_cast<intptr_t>(view.metadata().addressOf_m()) -
113  reinterpret_cast<intptr_t>(view.metadata().addressOf_r());
114  }
115 
116  checkViewAddresses(view);
117  checkViewAddresses(mview);
118  checkViewAddresses(cmview);
119 
120  const portabletest::Matrix matrix{{1, 2, 3, 4, 5, 6}, {2, 4, 6, 8, 10, 12}, {3, 6, 9, 12, 15, 18}};
121  assert(view.r() == 1.);
122  for (int32_t i = 0; i < view.metadata().size(); ++i) {
123  auto vi = view[i];
124  assert(vi.x() == 0.);
125  assert(vi.y() == 0.);
126  assert(vi.z() == 0.);
127  assert(vi.id() == i);
128  assert(vi.m() == matrix * i);
129  }
130  }
131 
134  desc.add<edm::InputTag>("source");
135  descriptions.addWithDefaultLabel(desc);
136  }
137 
138 private:
141 };
142 
size
Write out results.
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< portabletest::TestHostCollection > token_
void analyze(edm::Event const &event, edm::EventSetup const &) override
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
std::string encode() const
Definition: InputTag.cc:159
Eigen::Matrix< double, 3, 6 > Matrix
Definition: TestSoA.h:13
Definition: config.py:1
assert(be >=bs)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
ConstView const & const_view() const
const edm::InputTag source_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
tuple msg
Definition: mps_check.py:286
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
TestAlpakaAnalyzer(edm::ParameterSet const &config)
long double T
Definition: event.py:1