Energy

Eta

ID

JetType1

JetType2

JetType3

JetType4

JetType5

Label

Phi

Px

Py

Pz

testTableFilling

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188
#include "FWCore/SOA/interface/Table.h"
#include "FWCore/SOA/interface/TableView.h"
#include "FWCore/SOA/interface/Column.h"

#include <cppunit/extensions/HelperMacros.h>

class testTableFilling : public CppUnit::TestFixture {
  CPPUNIT_TEST_SUITE(testTableFilling);

  CPPUNIT_TEST(soaDeclareDefaultTest1);
  CPPUNIT_TEST(soaDeclareDefaultTest2);
  CPPUNIT_TEST(soaDeclareDefaultTest3);
  CPPUNIT_TEST(soaDeclareDefaultTest4);
  CPPUNIT_TEST(soaDeclareDefaultTest5);
  CPPUNIT_TEST_SUITE_END();

public:
  void setUp() {}
  void tearDown() {}

  void soaDeclareDefaultTest1();
  void soaDeclareDefaultTest2();
  void soaDeclareDefaultTest3();
  void soaDeclareDefaultTest4();
  void soaDeclareDefaultTest5();
};

namespace ts {

  namespace col {
    SOA_DECLARE_COLUMN(Eta, float, "eta");
    SOA_DECLARE_COLUMN(Phi, float, "phi");
    SOA_DECLARE_COLUMN(Energy, float, "energy");
    SOA_DECLARE_COLUMN(ID, int, "id");
    SOA_DECLARE_COLUMN(Label, std::string, "label");

    SOA_DECLARE_COLUMN(Px, double, "p_x");
    SOA_DECLARE_COLUMN(Py, double, "p_y");
    SOA_DECLARE_COLUMN(Pz, double, "p_z");
  }  // namespace col

  using EtaPhiTable = edm::soa::Table<ts::col::Eta, ts::col::Phi>;
  using EtaPhiTableView = edm::soa::ViewFromTable_t<EtaPhiTable>;
}  // namespace ts

namespace edm::soa {

  // The SOA_DECLARE_DEFAULT macro should be used in the edm::soa namespace
  SOA_DECLARE_DEFAULT(ts::col::Eta, eta());
  SOA_DECLARE_DEFAULT(ts::col::Phi, phi());
  SOA_DECLARE_DEFAULT(ts::col::Energy, energy());
  SOA_DECLARE_DEFAULT(ts::col::ID, id());
  SOA_DECLARE_DEFAULT(ts::col::Label, label());
  SOA_DECLARE_DEFAULT(ts::col::Px, px());
  SOA_DECLARE_DEFAULT(ts::col::Py, py());
  SOA_DECLARE_DEFAULT(ts::col::Pz, pz());

}  // namespace edm::soa

///registration of the test so that the runner can find it
CPPUNIT_TEST_SUITE_REGISTRATION(testTableFilling);

namespace {
  template <class Object>
  void validateEtaPhiTable(ts::EtaPhiTableView table, std::vector<Object> const& objects) {
    int iRow = 0;
    for (auto const& row : table) {
      CPPUNIT_ASSERT(row.get<ts::col::Eta>() == objects[iRow].eta_);
      CPPUNIT_ASSERT(row.get<ts::col::Phi>() == objects[iRow].phi_);
      ++iRow;
    }
  }
}  // namespace

struct JetType1 {
  // jet type that is fully compatible with the
  // SOA_DECLARE_DEFAULT-generated value_for_column functions
  auto eta() const { return eta_; }
  auto phi() const { return phi_; }

  const float eta_;
  const float phi_;
};

void testTableFilling::soaDeclareDefaultTest1() {
  std::vector<JetType1> jets = {{1., 3.14}, {2., 0.}, {4., 1.3}};
  std::vector<std::string> labels = {{"jet0", "jet1", "jet2"}};

  ts::EtaPhiTable table(jets);

  validateEtaPhiTable(table, jets);
}

struct JetType2 {
  // jet type that is NOT compatible with the SOA_DECLARE_DEFAULT-generated
  // value_for_column functions...
  auto eta() const { return eta_; }

  const float eta_;
  const float phi_;
};

// ...so we need to write our own value_for_column function for JetType2
double value_for_column(JetType2 const& iJ, ts::col::Phi*) { return iJ.phi_; }

void testTableFilling::soaDeclareDefaultTest2() {
  std::vector<JetType2> jets = {{1., 3.14}, {2., 0.}, {4., 1.3}};
  std::vector<std::string> labels = {{"jet0", "jet1", "jet2"}};

  ts::EtaPhiTable table(jets);

  validateEtaPhiTable(table, jets);
}

namespace ts::reco {

  struct JetType3 {
    auto eta() const { return eta_; }

    const float eta_;
    const float phi_;
  };

  double value_for_column(ts::reco::JetType3 const& iJ, ts::col::Phi*) { return iJ.phi_; }
}  // namespace ts::reco

// same tests as for JetType2 but with the type and the value_for_column
// defined in their own test namespace
void testTableFilling::soaDeclareDefaultTest3() {
  std::vector<ts::reco::JetType3> jets = {{1., 3.14}, {2., 0.}, {4., 1.3}};
  std::vector<std::string> labels = {{"jet0", "jet1", "jet2"}};

  ts::EtaPhiTable table(jets);

  validateEtaPhiTable(table, jets);
}

struct JetType4 {
  // jet type that is compatible with the SOA_DECLARE_DEFAULT-generated
  // value_for_column functions, but it would not give the right result.
  // This is to check that the custom value_for_column defined below is
  // prioritized.
  auto eta() const { return eta_; }
  auto phi() const { return 0.f; }

  const float eta_;
  const float phi_;
};

double value_for_column(JetType4 const& iJ, ts::col::Phi*) { return iJ.phi_; }

void testTableFilling::soaDeclareDefaultTest4() {
  std::vector<JetType4> jets = {{1., 3.14}, {2., 0.}, {4., 1.3}};
  std::vector<std::string> labels = {{"jet0", "jet1", "jet2"}};

  ts::EtaPhiTable table(jets);

  validateEtaPhiTable(table, jets);
}

namespace ts::reco {

  struct JetType5 {
    // jet type that is compatible with the SOA_DECLARE_DEFAULT-generated
    // value_for_column functions, but it would not give the right result.
    // This is to check that the custom value_for_column defined below is
    // prioritized.
    auto eta() const { return eta_; }
    auto phi() const { return 0.f; }

    const float eta_;
    const float phi_;
  };

  double value_for_column(ts::reco::JetType5 const& iJ, ts::col::Phi*) { return iJ.phi_; }

}  // namespace ts::reco

void testTableFilling::soaDeclareDefaultTest5() {
  std::vector<ts::reco::JetType5> jets = {{1., 3.14}, {2., 0.}, {4., 1.3}};
  std::vector<std::string> labels = {{"jet0", "jet1", "jet2"}};

  ts::EtaPhiTable table(jets);

  validateEtaPhiTable(table, jets);
}

#include <Utilities/Testing/interface/CppUnit_testdriver.icpp>