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  // columns
55  assert(view.metadata().addressOf_x() == view.x());
56  assert(view.metadata().addressOf_x() == &view.x(0));
57  assert(view.metadata().addressOf_x() == &view[0].x());
58  assert(view.metadata().addressOf_y() == view.y());
59  assert(view.metadata().addressOf_y() == &view.y(0));
60  assert(view.metadata().addressOf_y() == &view[0].y());
61  assert(view.metadata().addressOf_z() == view.z());
62  assert(view.metadata().addressOf_z() == &view.z(0));
63  assert(view.metadata().addressOf_z() == &view[0].z());
64  assert(view.metadata().addressOf_id() == view.id());
65  assert(view.metadata().addressOf_id() == &view.id(0));
66  assert(view.metadata().addressOf_id() == &view[0].id());
67  // scalars
68  assert(view.metadata().addressOf_r() == &view.r());
69  //assert(view.metadata().addressOf_r() == &view.r(0)); // cannot access a scalar with an index
70  //assert(view.metadata().addressOf_r() == &view[0].r()); // cannot access a scalar via a SoA row-like accessor
71  // columns of arrays
72  assert(view.metadata().addressOf_flags() == view.flags());
73  assert(view.metadata().addressOf_flags() == &view.flags(0));
74  assert(view.metadata().addressOf_flags() == &view[0].flags());
75  // columns of Eigen matrices
76  assert(view.metadata().addressOf_m() == view.m());
77  assert(view.metadata().addressOf_m() == &view.m(0).coeffRef(0, 0));
78  assert(view.metadata().addressOf_m() == &view[0].m().coeffRef(0, 0));
79  }
80 
81 } // namespace
82 
84 public:
86  : source_{config.getParameter<edm::InputTag>("source")},
88  expectSize_{config.getParameter<int>("expectSize")} {
89  if (std::string const& eb = config.getParameter<std::string>("expectBackend"); not eb.empty()) {
92  }
93  }
94 
95  void analyze(edm::StreamID sid, edm::Event const& event, edm::EventSetup const&) const override {
96  portabletest::TestHostCollection const& product = event.get(token_);
97  auto const& view = product.const_view();
98  auto& mview = product.view();
99  auto const& cmview = product.view();
100 
101  if (expectSize_ >= 0 and expectSize_ != view.metadata().size()) {
102  throw cms::Exception("Assert") << "Expected input collection size " << expectSize_ << ", got "
103  << view.metadata().size();
104  }
105 
106  {
107  edm::LogInfo msg("TestAlpakaAnalyzer");
108  msg << source_.encode() << ".size() = " << view.metadata().size() << '\n';
109  msg << " data @ " << product.buffer().data() << ",\n"
110  << " x @ " << view.metadata().addressOf_x() << " = " << Column(view.x(), view.metadata().size()) << ",\n"
111  << " y @ " << view.metadata().addressOf_y() << " = " << Column(view.y(), view.metadata().size()) << ",\n"
112  << " z @ " << view.metadata().addressOf_z() << " = " << Column(view.z(), view.metadata().size()) << ",\n"
113  << " id @ " << view.metadata().addressOf_id() << " = " << Column(view.id(), view.metadata().size())
114  << ",\n"
115  << " r @ " << view.metadata().addressOf_r() << " = " << view.r() << '\n'
116  << " flags @ " << view.metadata().addressOf_flags() << " = " << Column(view.flags(), view.metadata().size())
117  << ",\n"
118  << " m @ " << view.metadata().addressOf_m() << " = { ... {" << view[1].m()(1, Eigen::indexing::all)
119  << " } ... } \n";
120  msg << std::hex << " [y - x] = 0x"
121  << reinterpret_cast<intptr_t>(view.metadata().addressOf_y()) -
122  reinterpret_cast<intptr_t>(view.metadata().addressOf_x())
123  << " [z - y] = 0x"
124  << reinterpret_cast<intptr_t>(view.metadata().addressOf_z()) -
125  reinterpret_cast<intptr_t>(view.metadata().addressOf_y())
126  << " [id - z] = 0x"
127  << reinterpret_cast<intptr_t>(view.metadata().addressOf_id()) -
128  reinterpret_cast<intptr_t>(view.metadata().addressOf_z())
129  << " [r - id] = 0x"
130  << reinterpret_cast<intptr_t>(view.metadata().addressOf_r()) -
131  reinterpret_cast<intptr_t>(view.metadata().addressOf_id())
132  << " [flags - r] = 0x"
133  << reinterpret_cast<intptr_t>(view.metadata().addressOf_flags()) -
134  reinterpret_cast<intptr_t>(view.metadata().addressOf_r())
135  << " [m - flags] = 0x"
136  << reinterpret_cast<intptr_t>(view.metadata().addressOf_m()) -
137  reinterpret_cast<intptr_t>(view.metadata().addressOf_flags());
138  }
139 
140  checkViewAddresses(view);
141  checkViewAddresses(mview);
142  checkViewAddresses(cmview);
143 
144  const portabletest::Matrix matrix{{1, 2, 3, 4, 5, 6}, {2, 4, 6, 8, 10, 12}, {3, 6, 9, 12, 15, 18}};
145  const portabletest::Array flags{{6, 4, 2, 0}};
146  assert(view.r() == 1.);
147  for (int32_t i = 0; i < view.metadata().size(); ++i) {
148  auto vi = view[i];
149  assert(vi.x() == 0.);
150  assert(vi.y() == 0.);
151  assert(vi.z() == 0.);
152  assert(vi.id() == i);
153  assert(vi.flags() == flags);
154  assert(vi.m() == matrix * i);
155  }
156 
157  if (expectBackend_) {
158  auto backend = static_cast<cms::alpakatools::Backend>(event.get(backendToken_));
159  if (expectBackend_ != backend) {
160  throw cms::Exception("Assert") << "Expected input backend " << cms::alpakatools::toString(*expectBackend_)
161  << ", got " << cms::alpakatools::toString(backend);
162  }
163  }
164  }
165 
168  desc.add<edm::InputTag>("source");
169  desc.add<int>("expectSize", -1)
170  ->setComment("Expected size of the input collection. Values < 0 mean the check is not performed. Default: -1");
171  desc.add<std::string>("expectBackend", "")
172  ->setComment(
173  "Expected backend of the input collection. Empty value means to not perform the check. Default: empty "
174  "string");
175  descriptions.addWithDefaultLabel(desc);
176  }
177 
178 private:
182  std::optional<cms::alpakatools::Backend> expectBackend_;
183  const int expectSize_;
184 };
185 
size
Write out results.
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< portabletest::TestHostCollection > token_
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:19
void analyze(edm::StreamID sid, edm::Event const &event, edm::EventSetup const &) const override
std::optional< cms::alpakatools::Backend > expectBackend_
Definition: config.py:1
std::string const & label() const
Definition: InputTag.h:36
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)
std::string_view toString(Backend backend)
Definition: Backend.cc:24
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
edm::EDGetTokenT< unsigned short > backendToken_
Backend toBackend(std::string_view name)
Definition: Backend.cc:13
#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:80
TestAlpakaAnalyzer(edm::ParameterSet const &config)
std::string const & process() const
Definition: InputTag.h:40
long double T
Definition: event.py:1