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