CMS 3D CMS Logo

LHETablesProducer.cc
Go to the documentation of this file.
9 #include "TLorentzVector.h"
10 
11 #include <vector>
12 #include <iostream>
13 
15 public:
17  : lheTag_(edm::vector_transform(params.getParameter<std::vector<edm::InputTag>>("lheInfo"),
18  [this](const edm::InputTag& tag) { return mayConsume<LHEEventProduct>(tag); })),
19  precision_(params.getParameter<int>("precision")),
20  storeLHEParticles_(params.getParameter<bool>("storeLHEParticles")) {
21  produces<nanoaod::FlatTable>("LHE");
23  produces<nanoaod::FlatTable>("LHEPart");
24  }
25 
26  ~LHETablesProducer() override {}
27 
28  void produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const override {
29  auto lheTab = std::make_unique<nanoaod::FlatTable>(1, "LHE", true);
30 
32  for (const auto& lheTag : lheTag_) {
33  iEvent.getByToken(lheTag, lheInfo);
34  if (lheInfo.isValid()) {
35  break;
36  }
37  }
38  if (lheInfo.isValid()) {
39  auto lhePartTab = fillLHEObjectTable(*lheInfo, *lheTab);
41  iEvent.put(std::move(lhePartTab), "LHEPart");
42  } else {
43  if (storeLHEParticles_) { // need to store a dummy table anyway to make the framework happy
44  auto lhePartTab = std::make_unique<nanoaod::FlatTable>(1, "LHEPart", true);
45  iEvent.put(std::move(lhePartTab), "LHEPart");
46  }
47  }
48  iEvent.put(std::move(lheTab), "LHE");
49  }
50 
51  std::unique_ptr<nanoaod::FlatTable> fillLHEObjectTable(const LHEEventProduct& lheProd,
52  nanoaod::FlatTable& out) const {
53  double lheHT = 0, lheHTIncoming = 0;
54  unsigned int lheNj = 0, lheNb = 0, lheNc = 0, lheNuds = 0, lheNglu = 0;
55  double lheVpt = 0;
56 
57  const auto& hepeup = lheProd.hepeup();
58  const auto& pup = hepeup.PUP;
59  int lep = -1, lepBar = -1, nu = -1, nuBar = -1;
60  std::vector<float> vals_pt;
61  std::vector<float> vals_eta;
62  std::vector<float> vals_phi;
63  std::vector<float> vals_mass;
64  std::vector<int> vals_pid;
65  for (unsigned int i = 0, n = pup.size(); i < n; ++i) {
66  int status = hepeup.ISTUP[i];
67  int idabs = std::abs(hepeup.IDUP[i]);
68  if (status == 1) {
69  TLorentzVector p4(pup[i][0], pup[i][1], pup[i][2], pup[i][3]); // x,y,z,t
70  vals_pt.push_back(p4.Pt());
71  vals_eta.push_back(p4.Eta());
72  vals_phi.push_back(p4.Phi());
73  vals_mass.push_back(p4.M());
74  vals_pid.push_back(hepeup.IDUP[i]);
75  }
76  if ((status == 1) && ((idabs == 21) || (idabs > 0 && idabs < 7))) { //# gluons and quarks
77  // object counters
78  lheNj++;
79  if (idabs == 5)
80  lheNb++;
81  else if (idabs == 4)
82  lheNc++;
83  else if (idabs <= 3)
84  lheNuds++;
85  else if (idabs == 21)
86  lheNglu++;
87  // HT
88  double pt = std::hypot(pup[i][0], pup[i][1]); // first entry is px, second py
89  lheHT += pt;
90  int mothIdx = std::max(
91  hepeup.MOTHUP[i].first - 1,
92  0); //first and last mother as pair; first entry has index 1 in LHE; incoming particles return motherindex 0
93  int mothIdxTwo = std::max(hepeup.MOTHUP[i].second - 1, 0);
94  int mothStatus = hepeup.ISTUP[mothIdx];
95  int mothStatusTwo = hepeup.ISTUP[mothIdxTwo];
96  bool hasIncomingAsMother = mothStatus < 0 || mothStatusTwo < 0;
97  if (hasIncomingAsMother)
98  lheHTIncoming += pt;
99  } else if (idabs == 12 || idabs == 14 || idabs == 16) {
100  (hepeup.IDUP[i] > 0 ? nu : nuBar) = i;
101  } else if (idabs == 11 || idabs == 13 || idabs == 15) {
102  (hepeup.IDUP[i] > 0 ? lep : lepBar) = i;
103  }
104  }
105  std::pair<int, int> v(0, 0);
106  if (lep != -1 && lepBar != -1)
107  v = std::make_pair(lep, lepBar);
108  else if (lep != -1 && nuBar != -1)
109  v = std::make_pair(lep, nuBar);
110  else if (nu != -1 && lepBar != -1)
111  v = std::make_pair(nu, lepBar);
112  else if (nu != -1 && nuBar != -1)
113  v = std::make_pair(nu, nuBar);
114  if (v.first != -1 && v.second != -1) {
115  lheVpt = std::hypot(pup[v.first][0] + pup[v.second][0], pup[v.first][1] + pup[v.second][1]);
116  }
117 
118  out.addColumnValue<uint8_t>(
119  "Njets", lheNj, "Number of jets (partons) at LHE step", nanoaod::FlatTable::UInt8Column);
120  out.addColumnValue<uint8_t>("Nb", lheNb, "Number of b partons at LHE step", nanoaod::FlatTable::UInt8Column);
121  out.addColumnValue<uint8_t>("Nc", lheNc, "Number of c partons at LHE step", nanoaod::FlatTable::UInt8Column);
122  out.addColumnValue<uint8_t>(
123  "Nuds", lheNuds, "Number of u,d,s partons at LHE step", nanoaod::FlatTable::UInt8Column);
124  out.addColumnValue<uint8_t>(
125  "Nglu", lheNglu, "Number of gluon partons at LHE step", nanoaod::FlatTable::UInt8Column);
126  out.addColumnValue<float>("HT", lheHT, "HT, scalar sum of parton pTs at LHE step", nanoaod::FlatTable::FloatColumn);
127  out.addColumnValue<float>("HTIncoming",
128  lheHTIncoming,
129  "HT, scalar sum of parton pTs at LHE step, restricted to partons",
131  out.addColumnValue<float>("Vpt", lheVpt, "pT of the W or Z boson at LHE step", nanoaod::FlatTable::FloatColumn);
132  out.addColumnValue<uint8_t>("NpNLO", lheProd.npNLO(), "number of partons at NLO", nanoaod::FlatTable::UInt8Column);
133  out.addColumnValue<uint8_t>("NpLO", lheProd.npLO(), "number of partons at LO", nanoaod::FlatTable::UInt8Column);
134 
135  auto outPart = std::make_unique<nanoaod::FlatTable>(vals_pt.size(), "LHEPart", false);
136  outPart->addColumn<float>("pt", vals_pt, "Pt of LHE particles", nanoaod::FlatTable::FloatColumn, this->precision_);
137  outPart->addColumn<float>(
138  "eta", vals_eta, "Pseodorapidity of LHE particles", nanoaod::FlatTable::FloatColumn, this->precision_);
139  outPart->addColumn<float>(
140  "phi", vals_phi, "Phi of LHE particles", nanoaod::FlatTable::FloatColumn, this->precision_);
141  outPart->addColumn<float>(
142  "mass", vals_mass, "Mass of LHE particles", nanoaod::FlatTable::FloatColumn, this->precision_);
143  outPart->addColumn<int>("pdgId", vals_pid, "PDG ID of LHE particles", nanoaod::FlatTable::IntColumn);
144 
145  return outPart;
146  }
147 
150  desc.add<std::vector<edm::InputTag>>("lheInfo", std::vector<edm::InputTag>{{"externalLHEProducer"}, {"source"}})
151  ->setComment("tag(s) for the LHE information (LHEEventProduct)");
152  desc.add<int>("precision", -1)->setComment("precision on the 4-momenta of the LHE particles");
153  desc.add<bool>("storeLHEParticles", false)
154  ->setComment("Whether we want to store the 4-momenta of the status 1 particles at LHE level");
155  descriptions.add("lheInfoTable", desc);
156  }
157 
158 protected:
159  const std::vector<edm::EDGetTokenT<LHEEventProduct>> lheTag_;
160  const unsigned int precision_;
161  const bool storeLHEParticles_;
162 };
163 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
const lhef::HEPEUP & hepeup() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
LHETablesProducer(edm::ParameterSet const &params)
void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
const unsigned int precision_
std::unique_ptr< nanoaod::FlatTable > fillLHEObjectTable(const LHEEventProduct &lheProd, nanoaod::FlatTable &out) const
const bool storeLHEParticles_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~LHETablesProducer() override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double p4[4]
Definition: TauolaWrapper.h:92
std::vector< FiveVector > PUP
Definition: LesHouches.h:246
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const std::vector< edm::EDGetTokenT< LHEEventProduct > > lheTag_
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:124
int npLO() const
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
int npNLO() const
def move(src, dest)
Definition: eostools.py:511