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 
27 using namespace Rivet;
28 using namespace edm;
29 using namespace std;
30 
31 class HTXSRivetProducer : public edm::one::EDProducer<edm::one::WatchRuns,edm::one::SharedResources> {
32 public:
33 
35  _hepmcCollection(consumes<HepMCProduct>(cfg.getParameter<edm::InputTag>("HepMCCollection"))),
36  _lheRunInfo(consumes<LHERunInfoProduct,edm::InRun>(cfg.getParameter<edm::InputTag>("LHERunInfo")))
37  {
38  usesResource("Rivet");
39 
41 
42  _isFirstEvent = true;
43  _prodMode = cfg.getParameter<string>("ProductionMode");
44  m_HiggsProdMode = HTXS::UNKNOWN;
45 
46  produces<HTXS::HiggsClassification>("HiggsClassification").setBranchAlias("HiggsClassification");
47 
48  }
50 
51 private:
52 
53  void produce( edm::Event&, const edm::EventSetup&) ;
54 
55  void beginRun(edm::Run const& iRun, edm::EventSetup const& es) override ;
56  void endRun(edm::Run const& iRun, edm::EventSetup const& es) override ;
57 
60 
61  Rivet::AnalysisHandler _analysisHandler;
63 
67 
69 
70 };
71 
73 }
74 
76 
77  //get the hepmc product from the event
79 
80  bool product_exists = iEvent.getByToken(_hepmcCollection, evt);
81  if( product_exists ){
82 
83  // get HepMC GenEvent
84  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
85 
86  if (_prodMode == "AUTO") {
87 
88  // for these prod modes, don't change what is set in BeginRun
89  if (m_HiggsProdMode != HTXS::GGF && m_HiggsProdMode != HTXS::VBF && m_HiggsProdMode != HTXS::GG2ZH) {
90 
91  unsigned nWs = 0;
92  unsigned nZs = 0;
93  unsigned nTs = 0;
94  unsigned nBs = 0;
95  unsigned nHs = 0;
96 
97  HepMC::GenVertex *HSvtx = myGenEvent->signal_process_vertex();
98 
99  if(HSvtx){
100  for (auto ptcl:particles(HSvtx,HepMC::children)) {
101  if (std::abs(ptcl->pdg_id()) == 24) ++nWs;
102  if (ptcl->pdg_id() == 23) ++nZs;
103  if (abs(ptcl->pdg_id()) == 6) ++nTs;
104  if (abs(ptcl->pdg_id()) == 5) ++nBs;
105  if (ptcl->pdg_id() == 25) ++nHs;
106  }
107  }
108 
109  if (nZs==1 && nHs==1 && (nWs+nTs)==0) {
110  m_HiggsProdMode = HTXS::QQ2ZH;
111  } else if (nWs==1 && nHs==1 && (nZs+nTs)==0) {
112  m_HiggsProdMode = HTXS::WH;
113  } else if (nTs==2 && nHs==1 && nZs==0) {
114  m_HiggsProdMode = HTXS::TTH;
115  } else if (nTs==1 && nHs==1 && nZs==0) {
116  m_HiggsProdMode = HTXS::TH;
117  } else if (nBs==2 && nHs==1 && nZs==0) {
118  m_HiggsProdMode = HTXS::BBH;
119  }
120 
121  _HTXS->setHiggsProdMode(m_HiggsProdMode);
122 
123  }
124  }
125 
126 
127  if (_isFirstEvent){
128 
129  _analysisHandler.addAnalysis(_HTXS);
130 
131  // set the production mode if not done already
132  if ( _prodMode == "GGF" ) m_HiggsProdMode = HTXS::GGF;
133  else if ( _prodMode == "VBF" ) m_HiggsProdMode = HTXS::VBF;
134  else if ( _prodMode == "WH" ) m_HiggsProdMode = HTXS::WH;
135  else if ( _prodMode == "ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
136  else if ( _prodMode == "QQ2ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
137  else if ( _prodMode == "GG2ZH" ) m_HiggsProdMode = HTXS::GG2ZH;
138  else if ( _prodMode == "TTH" ) m_HiggsProdMode = HTXS::TTH;
139  else if ( _prodMode == "BBH" ) m_HiggsProdMode = HTXS::BBH;
140  else if ( _prodMode == "TH" ) m_HiggsProdMode = HTXS::TH;
141  else if ( _prodMode == "AUTO" ) {
142  edm::LogInfo ("HTXSRivetProducer") <<"Using AUTO for HiggsProdMode, found it to be: "<<m_HiggsProdMode<< "\n";
143  edm::LogInfo ("HTXSRivetProducer") <<"(UNKNOWN=0, GGF=1, VBF=2, WH=3, QQ2ZH=4, GG2ZH=5, TTH=6, BBH=7, TH=8)"<<endl;
144  } else {
145  throw cms::Exception("HTXSRivetProducer") << "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 
160  // classify the event
161  Rivet::HiggsClassification rivet_cat = _HTXS->classifyEvent(*myGenEvent,m_HiggsProdMode);
162  cat_ = HTXS::Rivet2Root(rivet_cat);
163 
164  unique_ptr<HTXS::HiggsClassification> cat( new HTXS::HiggsClassification( cat_ ) );
165 
166  iEvent.put(std::move(cat),"HiggsClassification");
167  }
168 }
169 
171 {
172  _HTXS->printClassificationSummary();
173 }
174 
176 
177  if (_prodMode == "AUTO") {
178 
180  bool product_exists = iRun.getByLabel( edm::InputTag("externalLHEProducer"), run );
181  if( product_exists ){
182 
183  typedef std::vector<LHERunInfoProduct::Header>::const_iterator headers_const_iterator;
184  LHERunInfoProduct myLHERunInfoProduct = *(run.product());
185  for (headers_const_iterator iter=myLHERunInfoProduct.headers_begin(); iter!=myLHERunInfoProduct.headers_end(); iter++){
186  std::vector<std::string> lines = iter->lines();
187  for (unsigned int iLine = 0; iLine<lines.size(); iLine++) {
188  std::string line=lines.at(iLine);
189  // POWHEG
190  if (line.find("gg_H_quark-mass-effects") != std::string::npos) {
191  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
192  m_HiggsProdMode = HTXS::GGF;
193  break;
194  }
195  if (line.find("Process: HJ") != std::string::npos) {
196  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
197  m_HiggsProdMode = HTXS::GGF;
198  break;
199  }
200  if (line.find("Process: HJJ") != std::string::npos) {
201  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
202  m_HiggsProdMode = HTXS::GGF;
203  break;
204  }
205  if (line.find("VBF_H") != std::string::npos) {
206  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
207  m_HiggsProdMode = HTXS::VBF;
208  break;
209  }
210  if (line.find("HZJ") != std::string::npos) {
211  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
212  m_HiggsProdMode = HTXS::QQ2ZH;
213  break;
214  }
215  if (line.find("ggHZ") != std::string::npos) {
216  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
217  m_HiggsProdMode = HTXS::GG2ZH;
218  break;
219  }
220  // MC@NLO
221  if (line.find("ggh012j") != std::string::npos) {
222  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
223  m_HiggsProdMode = HTXS::GGF;
224  break;
225  }
226  if (line.find("vbfh") != std::string::npos) {
227  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
228  m_HiggsProdMode = HTXS::VBF;
229  break;
230  }
231  if (line.find("zh012j") != std::string::npos) {
232  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
233  m_HiggsProdMode = HTXS::QQ2ZH;
234  break;
235  }
236  if (line.find("ggzh01j") != std::string::npos) {
237  edm::LogInfo ("HTXSRivetProducer") <<iLine<< " "<<line<<std::endl;
238  m_HiggsProdMode = HTXS::GG2ZH;
239  break;
240  }
241  }
242 
243  if ( m_HiggsProdMode != HTXS::UNKNOWN) break;
244  }
245  }
246  }
247 }
248 
249 
251 
T getParameter(std::string const &) const
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:303
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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:230
def cat(path)
Definition: eostools.py:400
headers_const_iterator headers_begin() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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:81
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:510
edm::EDGetTokenT< LHERunInfoProduct > _lheRunInfo
Definition: Run.h:44
void produce(edm::Event &, const edm::EventSetup &)