CMS 3D CMS Logo

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