CMS 3D CMS Logo

HTXSRivetProducer.cc
Go to the documentation of this file.
1 /* \class HTXSRivetProducer
2  *
3  * \author David Sperka, University of Florida
4  *
5  * $Id: HTXSRivetProducer.cc,v 1.1 2016/09/27 13:07:29 dsperka Exp $
6  *
7  */
8 
14 
17 
18 #include "Rivet/AnalysisHandler.hh"
21 
22 #include <vector>
23 #include <cstdio>
24 #include <cstring>
25 
26 using namespace Rivet;
27 using namespace edm;
28 using namespace std;
29 
30 class HTXSRivetProducer : public edm::one::EDProducer<edm::one::WatchRuns, edm::one::SharedResources> {
31 public:
33  : _hepmcCollection(consumes<HepMCProduct>(cfg.getParameter<edm::InputTag>("HepMCCollection"))),
34  _lheRunInfo(consumes<LHERunInfoProduct, edm::InRun>(cfg.getParameter<edm::InputTag>("LHERunInfo"))) {
35  usesResource("Rivet");
36 
38 
39  _isFirstEvent = true;
40  _prodMode = cfg.getParameter<string>("ProductionMode");
41  m_HiggsProdMode = HTXS::UNKNOWN;
42 
43  produces<HTXS::HiggsClassification>("HiggsClassification").setBranchAlias("HiggsClassification");
44  }
45 
46 private:
47  void produce(edm::Event&, const edm::EventSetup&) override;
48 
49  void beginRun(edm::Run const& iRun, edm::EventSetup const& es) override;
50  void endRun(edm::Run const& iRun, edm::EventSetup const& es) override;
51 
54 
55  Rivet::AnalysisHandler _analysisHandler;
57 
61 
63 };
64 
66  //get the hepmc product from the event
68 
69  bool product_exists = iEvent.getByToken(_hepmcCollection, evt);
70  if (product_exists) {
71  // get HepMC GenEvent
72  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
73 
74  if (_prodMode == "AUTO") {
75  // for these prod modes, don't change what is set in BeginRun
76  if (m_HiggsProdMode != HTXS::GGF && m_HiggsProdMode != HTXS::VBF && m_HiggsProdMode != HTXS::GG2ZH) {
77  unsigned nWs = 0;
78  unsigned nZs = 0;
79  unsigned nTs = 0;
80  unsigned nBs = 0;
81  unsigned nHs = 0;
82 
83  HepMC::GenVertex* HSvtx = myGenEvent->signal_process_vertex();
84 
85  if (HSvtx) {
86  for (auto ptcl : particles(HSvtx, HepMC::children)) {
87  if (std::abs(ptcl->pdg_id()) == 24)
88  ++nWs;
89  if (ptcl->pdg_id() == 23)
90  ++nZs;
91  if (abs(ptcl->pdg_id()) == 6)
92  ++nTs;
93  if (abs(ptcl->pdg_id()) == 5)
94  ++nBs;
95  if (ptcl->pdg_id() == 25)
96  ++nHs;
97  }
98  }
99 
100  if (nZs == 1 && nHs == 1 && (nWs + nTs) == 0) {
101  m_HiggsProdMode = HTXS::QQ2ZH;
102  } else if (nWs == 1 && nHs == 1 && (nZs + nTs) == 0) {
103  m_HiggsProdMode = HTXS::WH;
104  } else if (nTs == 2 && nHs == 1 && nZs == 0) {
105  m_HiggsProdMode = HTXS::TTH;
106  } else if (nTs == 1 && nHs == 1 && nZs == 0) {
107  m_HiggsProdMode = HTXS::TH;
108  } else if (nBs == 2 && nHs == 1 && nZs == 0) {
109  m_HiggsProdMode = HTXS::BBH;
110  }
111 
112  _HTXS->setHiggsProdMode(m_HiggsProdMode);
113  }
114  }
115 
116  if (_isFirstEvent) {
117  _analysisHandler.addAnalysis(_HTXS);
118 
119  // set the production mode if not done already
120  if (_prodMode == "GGF")
121  m_HiggsProdMode = HTXS::GGF;
122  else if (_prodMode == "VBF")
123  m_HiggsProdMode = HTXS::VBF;
124  else if (_prodMode == "WH")
125  m_HiggsProdMode = HTXS::WH;
126  else if (_prodMode == "ZH")
127  m_HiggsProdMode = HTXS::QQ2ZH;
128  else if (_prodMode == "QQ2ZH")
129  m_HiggsProdMode = HTXS::QQ2ZH;
130  else if (_prodMode == "GG2ZH")
131  m_HiggsProdMode = HTXS::GG2ZH;
132  else if (_prodMode == "TTH")
133  m_HiggsProdMode = HTXS::TTH;
134  else if (_prodMode == "BBH")
135  m_HiggsProdMode = HTXS::BBH;
136  else if (_prodMode == "TH")
137  m_HiggsProdMode = HTXS::TH;
138  else if (_prodMode == "AUTO") {
139  edm::LogInfo("HTXSRivetProducer")
140  << "Using AUTO for HiggsProdMode, found it to be: " << m_HiggsProdMode << "\n";
141  edm::LogInfo("HTXSRivetProducer")
142  << "(UNKNOWN=0, GGF=1, VBF=2, WH=3, QQ2ZH=4, GG2ZH=5, TTH=6, BBH=7, TH=8)" << endl;
143  } else {
144  throw cms::Exception("HTXSRivetProducer")
145  << "ProductionMode must be one of: GGF,VBF,WH,ZH,QQ2ZH,GG2ZH,TTH,BBH,TH,AUTO ";
146  }
147  _HTXS->setHiggsProdMode(m_HiggsProdMode);
148 
149  // at this point the production mode must be known
150  if (m_HiggsProdMode == HTXS::UNKNOWN) {
151  edm::LogInfo("HTXSRivetProducer") << "HTXSRivetProducer WARNING: HiggsProduction mode is UNKNOWN" << endl;
152  }
153 
154  // initialize rivet analysis
155  _analysisHandler.init(*myGenEvent);
156  _isFirstEvent = false;
157  }
158 
159  // classify the event
160  Rivet::HiggsClassification rivet_cat = _HTXS->classifyEvent(*myGenEvent, m_HiggsProdMode);
161  cat_ = HTXS::Rivet2Root(rivet_cat);
162 
163  unique_ptr<HTXS::HiggsClassification> cat(new HTXS::HiggsClassification(cat_));
164 
165  iEvent.put(std::move(cat), "HiggsClassification");
166  }
167 }
168 
169 void HTXSRivetProducer::endRun(edm::Run const& iRun, edm::EventSetup const& es) { _HTXS->printClassificationSummary(); }
170 
172  if (_prodMode == "AUTO") {
174  bool product_exists = iRun.getByLabel(edm::InputTag("externalLHEProducer"), run);
175  if (product_exists) {
176  typedef std::vector<LHERunInfoProduct::Header>::const_iterator headers_const_iterator;
177  LHERunInfoProduct myLHERunInfoProduct = *(run.product());
178  for (headers_const_iterator iter = myLHERunInfoProduct.headers_begin(); iter != myLHERunInfoProduct.headers_end();
179  iter++) {
180  std::vector<std::string> lines = iter->lines();
181  for (unsigned int iLine = 0; iLine < lines.size(); iLine++) {
182  std::string line = lines.at(iLine);
183  // POWHEG
184  if (line.find("gg_H_quark-mass-effects") != std::string::npos) {
185  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
186  m_HiggsProdMode = HTXS::GGF;
187  break;
188  }
189  if (line.find("Process: HJ") != std::string::npos) {
190  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
191  m_HiggsProdMode = HTXS::GGF;
192  break;
193  }
194  if (line.find("Process: HJJ") != std::string::npos) {
195  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
196  m_HiggsProdMode = HTXS::GGF;
197  break;
198  }
199  if (line.find("VBF_H") != std::string::npos) {
200  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
201  m_HiggsProdMode = HTXS::VBF;
202  break;
203  }
204  if (line.find("HZJ") != std::string::npos) {
205  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
206  m_HiggsProdMode = HTXS::QQ2ZH;
207  break;
208  }
209  if (line.find("ggHZ") != std::string::npos) {
210  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
211  m_HiggsProdMode = HTXS::GG2ZH;
212  break;
213  }
214  // MC@NLO
215  if (line.find("ggh012j") != std::string::npos) {
216  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
217  m_HiggsProdMode = HTXS::GGF;
218  break;
219  }
220  if (line.find("vbfh") != std::string::npos) {
221  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
222  m_HiggsProdMode = HTXS::VBF;
223  break;
224  }
225  if (line.find("zh012j") != std::string::npos) {
226  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
227  m_HiggsProdMode = HTXS::QQ2ZH;
228  break;
229  }
230  if (line.find("ggzh01j") != std::string::npos) {
231  edm::LogInfo("HTXSRivetProducer") << iLine << " " << line << std::endl;
232  m_HiggsProdMode = HTXS::GG2ZH;
233  break;
234  }
235  }
236 
237  if (m_HiggsProdMode != HTXS::UNKNOWN)
238  break;
239  }
240  }
241  }
242 }
243 
245 
T getParameter(std::string const &) const
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:280
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
headers_const_iterator headers_end() const
HiggsProdMode
Higgs production modes, corresponding to input sample.
HTXS::HiggsProdMode m_HiggsProdMode
Rivet::AnalysisHandler _analysisHandler
Rivet::HiggsTemplateCrossSections * _HTXS
int iEvent
Definition: GenABIO.cc:224
def cat(path)
Definition: eostools.py:401
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
headers_const_iterator headers_begin() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void produce(edm::Event &, const edm::EventSetup &) override
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
void beginRun(edm::Run const &iRun, edm::EventSetup const &es) override
T const * product() const
Definition: Handle.h:74
Rivet routine for classifying MC events according to the Higgs template cross section categories...
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
HLT enums.
HTXS::HiggsClassification Rivet2Root(category const &htxs_cat_rivet)
HTXS::HiggsClassification cat_
void endRun(edm::Run const &iRun, edm::EventSetup const &es) override
HTXSRivetProducer(const edm::ParameterSet &cfg)
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< LHERunInfoProduct > _lheRunInfo
Definition: Run.h:45