CMS 3D CMS Logo

LHETablesProducer.cc
Go to the documentation of this file.
8 
9 #include <vector>
10 #include <iostream>
11 
12 
14  public:
16  lheTag_(consumes<LHEEventProduct>(params.getParameter<edm::InputTag>("lheInfo")))
17  {
18  produces<nanoaod::FlatTable>();
19  }
20 
21  ~LHETablesProducer() override {}
22 
23  void produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const override {
24  auto lheTab = std::make_unique<nanoaod::FlatTable>(1, "LHE", true);
25 
27  if (iEvent.getByToken(lheTag_, lheInfo)) {
28  fillLHEObjectTable(*lheInfo, *lheTab);
29  }
30 
31  iEvent.put(std::move(lheTab));
32  }
33 
34  void fillLHEObjectTable(const LHEEventProduct & lheProd, nanoaod::FlatTable & out) const {
35  double lheHT = 0, lheHTIncoming = 0;
36  unsigned int lheNj = 0, lheNb = 0, lheNc = 0, lheNuds = 0, lheNglu = 0;
37  double lheVpt = 0;
38 
39  const auto & hepeup = lheProd.hepeup();
40  const auto & pup = hepeup.PUP;
41  int lep = -1, lepBar = -1, nu = -1, nuBar = -1;
42  for (unsigned int i = 0, n = pup.size(); i < n; ++i) {
43  int status = hepeup.ISTUP[i];
44  int idabs = std::abs(hepeup.IDUP[i]);
45  if ( (status == 1) && ( ( idabs == 21 ) || (idabs > 0 && idabs < 7) ) ) { //# gluons and quarks
46  // object counters
47  lheNj++;
48  if (idabs==5) lheNb++;
49  else if (idabs==4) lheNc++;
50  else if (idabs <= 3 ) lheNuds++;
51  else if (idabs == 21) lheNglu++;
52  // HT
53  double pt = std::hypot( pup[i][0], pup[i][1] ); // first entry is px, second py
54  lheHT += pt;
55  int mothIdx = std::max(hepeup.MOTHUP[i].first-1,0); //first and last mother as pair; first entry has index 1 in LHE; incoming particles return motherindex 0
56  int mothIdxTwo = std::max(hepeup.MOTHUP[i].second-1,0);
57  int mothStatus = hepeup.ISTUP[mothIdx];
58  int mothStatusTwo = hepeup.ISTUP[mothIdxTwo];
59  bool hasIncomingAsMother = mothStatus<0 || mothStatusTwo<0;
60  if (hasIncomingAsMother) lheHTIncoming += pt;
61  } else if (idabs == 12 || idabs == 14 || idabs == 16) {
62  (hepeup.IDUP[i] > 0 ? nu : nuBar) = i;
63  } else if (idabs == 11 || idabs == 13 || idabs == 15) {
64  (hepeup.IDUP[i] > 0 ? lep : lepBar) = i;
65  }
66  }
67  std::pair<int,int> v(0,0);
68  if (lep != -1 && lepBar != -1) v = std::make_pair(lep,lepBar);
69  else if (lep != -1 && nuBar != -1) v = std::make_pair(lep, nuBar);
70  else if (nu != -1 && lepBar != -1) v = std::make_pair(nu ,lepBar);
71  else if (nu != -1 && nuBar != -1) v = std::make_pair(nu , nuBar);
72  if (v.first != -1 && v.second != -1) {
73  lheVpt = std::hypot( pup[v.first][0] + pup[v.second][0], pup[v.first][1] + pup[v.second][1] );
74  }
75 
76  out.addColumnValue<uint8_t>("Njets", lheNj, "Number of jets (partons) at LHE step", nanoaod::FlatTable::UInt8Column);
77  out.addColumnValue<uint8_t>("Nb", lheNb, "Number of b partons at LHE step", nanoaod::FlatTable::UInt8Column);
78  out.addColumnValue<uint8_t>("Nc", lheNc, "Number of c partons at LHE step", nanoaod::FlatTable::UInt8Column);
79  out.addColumnValue<uint8_t>("Nuds", lheNuds, "Number of u,d,s partons at LHE step", nanoaod::FlatTable::UInt8Column);
80  out.addColumnValue<uint8_t>("Nglu", lheNglu, "Number of gluon partons at LHE step", nanoaod::FlatTable::UInt8Column);
81  out.addColumnValue<float>("HT", lheHT, "HT, scalar sum of parton pTs at LHE step", nanoaod::FlatTable::FloatColumn);
82  out.addColumnValue<float>("HTIncoming", lheHTIncoming, "HT, scalar sum of parton pTs at LHE step, restricted to partons", nanoaod::FlatTable::FloatColumn);
83  out.addColumnValue<float>("Vpt", lheVpt, "pT of the W or Z boson at LHE step", nanoaod::FlatTable::FloatColumn);
84  out.addColumnValue<uint8_t>("NpNLO", lheProd.npNLO(), "number of partons at NLO", nanoaod::FlatTable::UInt8Column);
85  out.addColumnValue<uint8_t>("NpLO", lheProd.npLO(), "number of partons at LO", nanoaod::FlatTable::UInt8Column);
86  }
87 
88  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions) {
90  desc.add<edm::InputTag>("lheInfo", edm::InputTag("externalLHEProducer"))->setComment("tag for the LHE information (LHEEventProduct)");
91  descriptions.add("lheInfoTable", desc);
92  }
93 
94  protected:
96 };
97 
100 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
const lhef::HEPEUP & hepeup() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
LHETablesProducer(edm::ParameterSet const &params)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:230
~LHETablesProducer() override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< LHEEventProduct > lheTag_
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void addColumnValue(const std::string &name, const C &value, const std::string &docString, ColumnType type=defaultColumnType< T >(), int mantissaBits=-1)
Definition: FlatTable.h:101
int npLO() const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
void fillLHEObjectTable(const LHEEventProduct &lheProd, nanoaod::FlatTable &out) const
int npNLO() const
def move(src, dest)
Definition: eostools.py:510