CMS 3D CMS Logo

HiggsTemplateCrossSections.cc
Go to the documentation of this file.
1 #ifndef TRUTHRIVETTOOLS_HIGGSTEMPLATECROSSSECTIONS_CC
2 #define TRUTHRIVETTOOLS_HIGGSTEMPLATECROSSSECTIONS_CC
3 
4 // -*- C++ -*-
5 #include "Rivet/Analysis.hh"
6 #include "Rivet/Particle.hh"
7 #include "Rivet/Projections/FastJets.hh"
8 
9 // Definition of the StatusCode and Category enums
10 //#include "HiggsTemplateCrossSections.h"
12 
13 #include <atomic>
14 
15 namespace Rivet {
16 
21  class HiggsTemplateCrossSections : public Analysis {
22  public:
23  // Constructor
25  : Analysis("HiggsTemplateCrossSections"),
27 
28  public:
29 
34 
37  if ( ptcl.genParticle()->end_vertex() ) {
38  if ( !hasChild(ptcl.genParticle(),ptcl.pdgId()) ) return ptcl;
39  else return getLastInstance(ptcl.children()[0]);
40  }
41  return ptcl;
42  }
43 
45  bool originateFrom(const Particle& p, const Particles& ptcls ) {
46  const GenVertex* prodVtx = p.genParticle()->production_vertex();
47  if (prodVtx == nullptr) return false;
48  // for each ancestor, check if it matches any of the input particles
49  for (const auto & ancestor:particles(prodVtx, HepMC::ancestors)){
50  for ( auto part:ptcls )
51  if ( ancestor==part.genParticle() ) return true;
52  }
53  // if we get here, no ancetor matched any input particle
54  return false;
55  }
56 
58  bool originateFrom(const Particle& p, const Particle& p2 ) {
59  Particles ptcls = {p2}; return originateFrom(p,ptcls);
60  }
61 
63  bool hasChild(const GenParticle *ptcl, int pdgID) {
64  for (auto child:Particle(*ptcl).children()) {
65  if (child.pdgId()==pdgID) return true;
66  }
67  return false;
68  }
69 
71  bool hasParent(const GenParticle *ptcl, int pdgID) {
72  for (auto parent:particles(ptcl->production_vertex(),HepMC::parents))
73  if (parent->pdg_id()==pdgID) return true;
74  return false;
75  }
76 
78  bool quarkDecay(const Particle &p) {
79  for (auto child:p.children()) {
80  if (PID::isQuark(child.pdgId())) return true;
81  }
82  return false;
83  }
84 
86  bool ChLeptonDecay(const Particle &p) {
87  for (auto child : p.children()) {
88  if (PID::isChLepton(child.pdgId())) {
89  return true;
90  }
91  }
92  return false;
93  }
94 
97  HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err,
98  std::string msg="", int NmaxWarnings=20) {
99  // Set the error, and keep statistics
100  cat.errorCode = err;
101  ++m_errorCount[err];
102 
103  // Print warning message to the screen/log
104  static std::atomic<int> Nwarnings{0};
105  if ( !msg.empty() && ++Nwarnings < NmaxWarnings )
106  MSG_WARNING(msg);
107 
108  return cat;
109  }
111 
113  HiggsClassification classifyEvent(const Event& event, const HTXS::HiggsProdMode prodMode ) {
114 
116 
117  // the classification object
118  HiggsClassification cat;
119  cat.prodMode = prodMode;
120  cat.errorCode = HTXS::UNDEFINED;
121  cat.stage0_cat = HTXS::Stage0::UNKNOWN;
122  cat.stage1_cat_pTjet25GeV = HTXS::Stage1::UNKNOWN;
123  cat.stage1_cat_pTjet30GeV = HTXS::Stage1::UNKNOWN;
124  cat.stage1_1_cat_pTjet25GeV = HTXS::Stage1_1::UNKNOWN;
125  cat.stage1_1_cat_pTjet30GeV = HTXS::Stage1_1::UNKNOWN;
126  cat.stage1_1_fine_cat_pTjet25GeV = HTXS::Stage1_1_Fine::UNKNOWN;
127  cat.stage1_1_fine_cat_pTjet30GeV = HTXS::Stage1_1_Fine::UNKNOWN;
128  cat.stage1_2_cat_pTjet25GeV = HTXS::Stage1_2::UNKNOWN;
129  cat.stage1_2_cat_pTjet30GeV = HTXS::Stage1_2::UNKNOWN;
130  cat.stage1_2_fine_cat_pTjet25GeV = HTXS::Stage1_2_Fine::UNKNOWN;
131  cat.stage1_2_fine_cat_pTjet30GeV = HTXS::Stage1_2_Fine::UNKNOWN;
132 
133  if (prodMode == HTXS::UNKNOWN)
134  return error(cat,HTXS::PRODMODE_DEFINED,
135  "Unkown Higgs production mechanism. Cannot classify event."
136  " Classification for all events will most likely fail.");
137 
138  /*****
139  * Step 1.
140  * Idenfify the Higgs boson and the hard scatter vertex
141  * There should be only one of each.
142  */
143 
144  GenVertex *HSvtx = event.genEvent()->signal_process_vertex();
145  int Nhiggs=0;
146  for ( const GenParticle *ptcl : Rivet::particles(event.genEvent()) ) {
147 
148  // a) Reject all non-Higgs particles
149  if ( !PID::isHiggs(ptcl->pdg_id()) ) continue;
150  // b) select only the final Higgs boson copy, prior to decay
151  if ( ptcl->end_vertex() && !hasChild(ptcl,PID::HIGGS) ) {
152  cat.higgs = Particle(ptcl); ++Nhiggs;
153  }
154  // c) if HepMC::signal_proces_vertex is missing
155  // set hard-scatter vertex based on first Higgs boson
156  if ( HSvtx==nullptr && ptcl->production_vertex() && !hasParent(ptcl,PID::HIGGS) )
157  HSvtx = ptcl->production_vertex();
158  }
159 
160  // Make sure things are in order so far
161  if (Nhiggs!=1)
163  "Current event has "+std::to_string(Nhiggs)+" Higgs bosons. There must be only one.");
164  if (cat.higgs.children().size()<2)
166  "Could not identify Higgs boson decay products.");
167 
168  if (HSvtx == nullptr)
169  return error(cat,HTXS::HS_VTX_IDENTIFICATION,"Cannot find hard-scatter vertex of current event.");
170 
171  /*****
172  * Step 2.
173  * Identify associated vector bosons
174  */
175 
176  // Find associated vector bosons
177  bool is_uncatdV = false;
178  Particles uncatV_decays;
179  FourMomentum uncatV_p4(0,0,0,0);
180  FourVector uncatV_v4(0,0,0,0);
181  int nWs=0, nZs=0, nTop=0;
182  if ( isVH(prodMode) ) {
183  for (auto ptcl:particles(HSvtx,HepMC::children)) {
184  if (PID::isW(ptcl->pdg_id())) { ++nWs; cat.V=Particle(ptcl); }
185  if (PID::isZ(ptcl->pdg_id())) { ++nZs; cat.V=Particle(ptcl); }
186  }
187  if(nWs+nZs>0) cat.V = getLastInstance(cat.V);
188  else {
189  for (auto ptcl:particles(HSvtx,HepMC::children)) {
190  if (!PID::isHiggs(ptcl->pdg_id())) {
191  uncatV_decays += Particle(ptcl);
192  uncatV_p4 += Particle(ptcl).momentum();
193  // uncatV_v4 += Particle(ptcl).origin();
194  }
195  }
196  // is_uncatdV = true; cat.V = Particle(24,uncatV_p4,uncatV_v4);
197  is_uncatdV = true; cat.V = Particle(24,uncatV_p4);
198  }
199  }
200 
201  if ( !is_uncatdV ){
202 
203  if ( isVH(prodMode) && !cat.V.genParticle()->end_vertex() )
204  return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
205 
206  if ( isVH(prodMode) && cat.V.children().size()<2 )
207  return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
208 
209  if ( ( prodMode==HTXS::WH && (nZs>0||nWs!=1) ) ||
210  ( (prodMode==HTXS::QQ2ZH||prodMode==HTXS::GG2ZH) && (nZs!=1||nWs>0) ) )
211  return error(cat,HTXS::VH_IDENTIFICATION,"Found "+std::to_string(nWs)+" W-bosons and "+
212  std::to_string(nZs)+" Z-bosons. Inconsitent with VH expectation.");
213  }
214 
215  // Find and store the W-bosons from ttH->WbWbH
216  Particles Ws;
217  if ( prodMode==HTXS::TTH || prodMode==HTXS::TH ){
218  // loop over particles produced in hard-scatter vertex
219  for ( auto ptcl : particles(HSvtx,HepMC::children) ) {
220  if ( !PID::isTop(ptcl->pdg_id()) ) continue;
221  ++nTop;
222  Particle top = getLastInstance(Particle(ptcl));
223  if ( top.genParticle()->end_vertex() )
224  for (auto child:top.children())
225  if ( PID::isW(child.pdgId()) ) Ws += getLastInstance(child);
226  }
227  }
228 
229  // Make sure result make sense
230  if ( (prodMode==HTXS::TTH && Ws.size()<2) || (prodMode==HTXS::TH && Ws.empty() ) )
231  return error(cat,HTXS::TOP_W_IDENTIFICATION,"Failed to identify W-boson(s) from t-decay!");
232 
233  /*****
234  * Step 3.
235  * Build jets
236  * Make sure all stable particles are present
237  */
238 
239  // Create a list of the vector bosons that decay leptonically
240  // Either the vector boson produced in association with the Higgs boson,
241  // or the ones produced from decays of top quarks produced with the Higgs
242  Particles leptonicVs;
243  if ( !is_uncatdV ){
244  if ( isVH(prodMode) && !quarkDecay(cat.V) ) leptonicVs += cat.V;
245  }else leptonicVs = uncatV_decays;
246  for ( auto W:Ws ) if ( W.genParticle()->end_vertex() && !quarkDecay(W) ) leptonicVs += W;
247 
248  // Obtain all stable, final-state particles
249  const ParticleVector FS = applyProjection<FinalState>(event, "FS").particles();
250  Particles hadrons;
251  FourMomentum sum(0,0,0,0), vSum(0,0,0,0), hSum(0,0,0,0);
252  for ( const Particle &p : FS ) {
253  // Add up the four momenta of all stable particles as a cross check
254  sum += p.momentum();
255  // ignore particles from the Higgs boson
256  if ( originateFrom(p,cat.higgs) ) { hSum += p.momentum(); continue; }
257  // Cross-check the V decay products for VH
258  if ( isVH(prodMode) && !is_uncatdV && originateFrom(p,Ws) ) vSum += p.momentum();
259  // ignore final state particles from leptonic V decays
260  if ( !leptonicVs.empty() && originateFrom(p,leptonicVs) ) continue;
261  // All particles reaching here are considered hadrons and will be used to build jets
262  hadrons += p;
263  }
264 
265  cat.p4decay_higgs = hSum;
266  cat.p4decay_V = vSum;
267 
268  FinalState fps_temp;
269  FastJets jets(fps_temp, FastJets::ANTIKT, 0.4 );
270  jets.calc(hadrons);
271 
272  cat.jets25 = jets.jetsByPt( Cuts::pT > 25.0 );
273  cat.jets30 = jets.jetsByPt( Cuts::pT > 30.0 );
274 
275  // check that four mometum sum of all stable particles satisfies momentum consevation
276 /*
277  if ( sum.pt()>0.1 )
278  return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
279  std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
280 */
281  // check if V-boson was not included in the event record but decay particles were
282  // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
283  if(is_uncatdV)
284  return error(cat,HTXS::VH_IDENTIFICATION,"Failed to identify associated V-boson!");
285 
286  /*****
287  * Step 4.
288  * Classify and save output
289  */
290 
291  // Apply the categorization categorization
292  cat.isZ2vvDecay = false;
293  if ((prodMode == HTXS::GG2ZH || prodMode == HTXS::QQ2ZH) && !quarkDecay(cat.V) && !ChLeptonDecay(cat.V))
294  cat.isZ2vvDecay = true;
295  cat.stage0_cat = getStage0Category(prodMode, cat.higgs, cat.V);
296  cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode, cat.higgs, cat.jets25, cat.V);
297  cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode, cat.higgs, cat.jets30, cat.V);
298  cat.stage1_1_cat_pTjet25GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets25, cat.V);
299  cat.stage1_1_cat_pTjet30GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets30, cat.V);
300  cat.stage1_1_fine_cat_pTjet25GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
301  cat.stage1_1_fine_cat_pTjet30GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
302  cat.stage1_2_cat_pTjet25GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets25, cat.V);
303  cat.stage1_2_cat_pTjet30GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets30, cat.V);
304  cat.stage1_2_fine_cat_pTjet25GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
305  cat.stage1_2_fine_cat_pTjet30GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
306 
307  cat.errorCode = HTXS::SUCCESS;
309 
310  return cat;
311  }
312 
318 
320  int getBin(double x, std::vector<double> bins) {
321  for (size_t i=1;i<bins.size();++i)
322  if (x<bins[i]) return i-1;
323  return bins.size()-1;
324  }
325 
329  int vbfTopology(const Jets &jets, const Particle &higgs) {
330  if (jets.size()<2) return 0;
331  const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
332  bool VBFtopo = (j1+j2).mass() > 400.0 && std::abs(j1.rapidity()-j2.rapidity()) > 2.8;
333  return VBFtopo ? (j1+j2+higgs.momentum()).pt()<25 ? 2 : 1 : 0;
334  }
335 
340  int vbfTopology_Stage1_X(const Jets &jets, const Particle &higgs) {
341  if (jets.size() < 2)
342  return 0;
343  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
344  double mjj = (j1 + j2).mass();
345  if (mjj > 350 && mjj <= 700)
346  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
347  else if (mjj > 700)
348  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
349  else
350  return 0;
351  }
352 
359  int vbfTopology_Stage1_X_Fine(const Jets &jets, const Particle &higgs) {
360  if (jets.size() < 2)
361  return 0;
362  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
363  double mjj = (j1 + j2).mass();
364  if (mjj > 350 && mjj <= 700)
365  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
366  else if (mjj > 700 && mjj <= 1000)
367  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
368  else if (mjj > 1000 && mjj <= 1500)
369  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 5 : 6;
370  else if (mjj > 1500)
371  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 7 : 8;
372  else
373  return 0;
374  }
375 
377  bool isVH(HTXS::HiggsProdMode p) { return p==HTXS::WH || p==HTXS::QQ2ZH || p==HTXS::GG2ZH; }
378 
381  const Particle &higgs,
382  const Particle &V) {
383  using namespace HTXS::Stage0;
384  int ctrlHiggs = std::abs(higgs.rapidity())<2.5;
385  // Special cases first, qq→Hqq
386  if ( (prodMode==HTXS::WH||prodMode==HTXS::QQ2ZH) && quarkDecay(V) ) {
387  return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
388  } else if ( prodMode==HTXS::GG2ZH && quarkDecay(V) ) {
389  return Category(HTXS::GGF*10 + ctrlHiggs);
390  }
391  // General case after
392  return Category(prodMode*10 + ctrlHiggs);
393  }
394 
397  const Particle &higgs,
398  const Jets &jets,
399  const Particle &V) {
400  using namespace HTXS::Stage1;
401  int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
402  double pTj1 = !jets.empty() ? jets[0].momentum().pt() : 0;
403  int vbfTopo = vbfTopology(jets,higgs);
404 
405  // 1. GGF Stage 1 categories
406  // Following YR4 write-up: XXXXX
407  if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
408  if (fwdHiggs) return GG2H_FWDH;
409  if (Njets==0) return GG2H_0J;
410  else if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
411  else if (Njets>=2) {
412  // events with pT_H>200 get priority over VBF cuts
413  if(higgs.pt()<=200){
414  if (vbfTopo==2) return GG2H_VBFTOPO_JET3VETO;
415  else if (vbfTopo==1) return GG2H_VBFTOPO_JET3;
416  }
417  // Njets >= 2jets without VBF topology
418  return Category(GG2H_GE2J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
419  }
420  }
421  // 2. Electroweak qq->Hqq Stage 1 categories
422  else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
423  if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
424  if (pTj1>200) return QQ2HQQ_PTJET1_GT200;
425  if (vbfTopo==2) return QQ2HQQ_VBFTOPO_JET3VETO;
426  if (vbfTopo==1) return QQ2HQQ_VBFTOPO_JET3;
427  double mjj = jets.size()>1 ? (jets[0].mom()+jets[1].mom()).mass():0;
428  if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_VH2JET;
429  return QQ2HQQ_REST;
430  }
431  // 3. WH->Hlv categories
432  else if (prodMode==HTXS::WH) {
433  if (fwdHiggs) return QQ2HLNU_FWDH;
434  else if (V.pt()<150) return QQ2HLNU_PTV_0_150;
435  else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
436  // 150 < pTV/GeV < 250
437  return jets.empty() ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
438  }
439  // 4. qq->ZH->llH categories
440  else if (prodMode==HTXS::QQ2ZH) {
441  if (fwdHiggs) return QQ2HLL_FWDH;
442  else if (V.pt()<150) return QQ2HLL_PTV_0_150;
443  else if (V.pt()>250) return QQ2HLL_PTV_GT250;
444  // 150 < pTV/GeV < 250
445  return jets.empty() ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
446  }
447  // 5. gg->ZH->llH categories
448  else if (prodMode==HTXS::GG2ZH ) {
449  if (fwdHiggs) return GG2HLL_FWDH;
450  if (V.pt()<150) return GG2HLL_PTV_0_150;
451  else if (jets.empty()) return GG2HLL_PTV_GT150_0J;
452  return GG2HLL_PTV_GT150_GE1J;
453  }
454  // 6.ttH,bbH,tH categories
455  else if (prodMode==HTXS::TTH) return Category(TTH_FWDH+ctrlHiggs);
456  else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
457  else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
458  return UNKNOWN;
459  }
460 
463  const Particle &higgs,
464  const Jets &jets,
465  const Particle &V) {
466  using namespace HTXS::Stage1_1;
467  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
468  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
469 
470  // 1. GGF Stage 1 categories
471  // Following YR4 write-up: XXXXX
472  if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
473  if (fwdHiggs) return GG2H_FWDH;
474  if ( higgs.pt()>200 ) return GG2H_PTH_GT200;
475  if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
476  if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
477  if (Njets>1){
478  //VBF topology
479  if(vbfTopo) return Category(GG2H_MJJ_350_700_PTHJJ_0_25+vbfTopo-1);
480  //Njets >= 2jets without VBF topology (mjj<350)
481  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
482  }
483  }
484 
485  // 2. Electroweak qq->Hqq Stage 1.1 categories
486  else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
487  if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
488  int Njets=jets.size();
489  if (Njets==0) return QQ2HQQ_0J;
490  else if (Njets==1) return QQ2HQQ_1J;
491  else if (Njets>=2) {
492  double mjj = (jets[0].mom()+jets[1].mom()).mass();
493  if ( mjj < 60 ) return QQ2HQQ_MJJ_0_60;
494  else if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_MJJ_60_120;
495  else if ( 120 < mjj && mjj < 350 ) return QQ2HQQ_MJJ_120_350;
496  else if ( mjj > 350 ) {
497  if (higgs.pt()>200) return QQ2HQQ_MJJ_GT350_PTH_GT200;
498  if(vbfTopo) return Category(QQ2HQQ_MJJ_GT350_PTH_GT200+vbfTopo);
499  }
500  }
501  }
502  // 3. WH->Hlv categories
503  else if (prodMode==HTXS::WH) {
504  if (fwdHiggs) return QQ2HLNU_FWDH;
505  else if (V.pt()<75) return QQ2HLNU_PTV_0_75;
506  else if (V.pt()<150) return QQ2HLNU_PTV_75_150;
507  else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
508  // 150 < pTV/GeV < 250
509  return jets.empty() ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
510  }
511  // 4. qq->ZH->llH categories
512  else if (prodMode==HTXS::QQ2ZH) {
513  if (fwdHiggs) return QQ2HLL_FWDH;
514  else if (V.pt()<75) return QQ2HLL_PTV_0_75;
515  else if (V.pt()<150) return QQ2HLL_PTV_75_150;
516  else if (V.pt()>250) return QQ2HLL_PTV_GT250;
517  // 150 < pTV/GeV < 250
518  return jets.empty() ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
519  }
520  // 5. gg->ZH->llH categories
521  else if (prodMode==HTXS::GG2ZH ) {
522  if (fwdHiggs) return GG2HLL_FWDH;
523  else if (V.pt()<75) return GG2HLL_PTV_0_75;
524  else if (V.pt()<150) return GG2HLL_PTV_75_150;
525  else if (V.pt()>250) return GG2HLL_PTV_GT250;
526  return jets.empty() ? GG2HLL_PTV_150_250_0J : GG2HLL_PTV_150_250_GE1J;
527  }
528  // 6.ttH,bbH,tH categories
529  else if (prodMode==HTXS::TTH) return Category(TTH_FWDH+ctrlHiggs);
530  else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
531  else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
532  return UNKNOWN;
533  }
534 
537  const Particle &higgs,
538  const Jets &jets,
539  const Particle &V) {
540  using namespace HTXS::Stage1_1_Fine;
541  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
542  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
543 
544  // 1. GGF Stage 1.1 categories
545  // Following YR4 write-up: XXXXX
546  if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
547  if (fwdHiggs) return GG2H_FWDH;
548  if ( higgs.pt()>200 ) return GG2H_PTH_GT200;
549  if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
550  if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
551  if (Njets>1){
552  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
553  double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
554  //VBF topology
555  if(vbfTopo) return Category(GG2H_MJJ_350_700_PTHJJ_0_25+vbfTopo-1);
556  //Njets >= 2jets without VBF topology (mjj<350)
557  if (pTHjj<25) return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25+getBin(higgs.pt(),{0,60,120,200}));
558  else return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25+getBin(higgs.pt(),{0,60,120,200}));
559  }
560  }
561 
562  // 2. Electroweak qq->Hqq Stage 1.1 categories
563  else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
564  if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
565  int Njets=jets.size();
566  if (Njets==0) return QQ2HQQ_0J;
567  else if (Njets==1) return QQ2HQQ_1J;
568  else if (Njets>=2) {
569  double mjj = (jets[0].mom()+jets[1].mom()).mass();
570  double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
571  if (mjj<350){
572  if (pTHjj<25) return Category(QQ2HQQ_MJJ_0_60_PTHJJ_0_25+getBin(mjj,{0,60,120,350}));
573  else return Category(QQ2HQQ_MJJ_0_60_PTHJJ_GT25+getBin(mjj,{0,60,120,350}));
574  } else { //mjj>350 GeV
575  if (higgs.pt()<200){
576  return Category(QQ2HQQ_MJJ_350_700_PTHJJ_0_25+vbfTopo-1);
577  } else {
579  }
580  }
581  }
582  }
583  // 3. WH->Hlv categories
584  else if (prodMode==HTXS::WH) {
585  if (fwdHiggs) return QQ2HLNU_FWDH;
586  int Njets=jets.size();
587  if (Njets==0) return Category(QQ2HLNU_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
588  if (Njets==1) return Category(QQ2HLNU_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
589  return Category(QQ2HLNU_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
590  }
591  // 4. qq->ZH->llH categories
592  else if (prodMode==HTXS::QQ2ZH) {
593  if (fwdHiggs) return QQ2HLL_FWDH;
594  int Njets=jets.size();
595  if (Njets==0) return Category(QQ2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
596  if (Njets==1) return Category(QQ2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
597  return Category(QQ2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
598  }
599  // 5. gg->ZH->llH categories
600  else if (prodMode==HTXS::GG2ZH ) {
601  if (fwdHiggs) return GG2HLL_FWDH;
602  int Njets=jets.size();
603  if (Njets==0) return Category(GG2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
604  if (Njets==1) return Category(GG2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
605  return Category(GG2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
606  }
607  // 6.ttH,bbH,tH categories
608  else if (prodMode==HTXS::TTH) return Category(TTH_FWDH+ctrlHiggs);
609  else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
610  else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
611  return UNKNOWN;
612  }
613 
616  const Particle &higgs,
617  const Jets &jets,
618  const Particle &V) {
619  using namespace HTXS::Stage1_2;
620  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
621  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
622 
623  // 1. GGF Stage 1 categories
624  // Following YR4 write-up: XXXXX
625  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
626  if (fwdHiggs)
627  return GG2H_FWDH;
628  if (higgs.pt() > 200)
629  return Category(GG2H_PTH_200_300 + getBin(higgs.pt(), {200, 300, 450, 650}));
630  if (Njets == 0)
631  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
632  if (Njets == 1)
633  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
634  if (Njets > 1) {
635  //VBF topology
636  if (vbfTopo)
638  //Njets >= 2jets without VBF topology (mjj<350)
639  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
640  }
641  }
642 
643  // 2. Electroweak qq->Hqq Stage 1.2 categories
644  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
645  if (std::abs(higgs.rapidity()) > 2.5)
646  return QQ2HQQ_FWDH;
647  int Njets = jets.size();
648  if (Njets == 0)
649  return QQ2HQQ_0J;
650  else if (Njets == 1)
651  return QQ2HQQ_1J;
652  else if (Njets >= 2) {
653  double mjj = (jets[0].mom() + jets[1].mom()).mass();
654  if (mjj < 60)
655  return QQ2HQQ_GE2J_MJJ_0_60;
656  else if (60 < mjj && mjj < 120)
657  return QQ2HQQ_GE2J_MJJ_60_120;
658  else if (120 < mjj && mjj < 350)
660  else if (mjj > 350) {
661  if (higgs.pt() > 200)
663  if (vbfTopo)
664  return Category(QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200 + vbfTopo);
665  }
666  }
667  }
668  // 3. WH->Hlv categories
669  else if (prodMode == HTXS::WH) {
670  if (fwdHiggs)
671  return QQ2HLNU_FWDH;
672  else if (V.pt() < 75)
673  return QQ2HLNU_PTV_0_75;
674  else if (V.pt() < 150)
675  return QQ2HLNU_PTV_75_150;
676  else if (V.pt() > 250)
677  return QQ2HLNU_PTV_GT250;
678  // 150 < pTV/GeV < 250
679  return jets.empty() ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
680  }
681  // 4. qq->ZH->llH categories
682  else if (prodMode == HTXS::QQ2ZH) {
683  if (fwdHiggs)
684  return QQ2HLL_FWDH;
685  else if (V.pt() < 75)
686  return QQ2HLL_PTV_0_75;
687  else if (V.pt() < 150)
688  return QQ2HLL_PTV_75_150;
689  else if (V.pt() > 250)
690  return QQ2HLL_PTV_GT250;
691  // 150 < pTV/GeV < 250
692  return jets.empty() ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
693  }
694  // 5. gg->ZH->llH categories
695  else if (prodMode == HTXS::GG2ZH) {
696  if (fwdHiggs)
697  return GG2HLL_FWDH;
698  else if (V.pt() < 75)
699  return GG2HLL_PTV_0_75;
700  else if (V.pt() < 150)
701  return GG2HLL_PTV_75_150;
702  else if (V.pt() > 250)
703  return GG2HLL_PTV_GT250;
704  return jets.empty() ? GG2HLL_PTV_150_250_0J : GG2HLL_PTV_150_250_GE1J;
705  }
706  // 6.ttH,bbH,tH categories
707  else if (prodMode == HTXS::TTH) {
708  if (fwdHiggs)
709  return TTH_FWDH;
710  else
711  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300}));
712  } else if (prodMode == HTXS::BBH)
713  return Category(BBH_FWDH + ctrlHiggs);
714  else if (prodMode == HTXS::TH)
715  return Category(TH_FWDH + ctrlHiggs);
716  return UNKNOWN;
717  }
718 
721  const Particle &higgs,
722  const Jets &jets,
723  const Particle &V) {
724  using namespace HTXS::Stage1_2_Fine;
725  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
726  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
727 
728  // 1. GGF Stage 1.2 categories
729  // Following YR4 write-up: XXXXX
730  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
731  if (fwdHiggs)
732  return GG2H_FWDH;
733  if (higgs.pt() > 200) {
734  if (Njets > 0) {
735  double pTHj = (jets[0].momentum() + higgs.momentum()).pt();
736  if (pTHj / higgs.pt() > 0.15)
737  return Category(GG2H_PTH_200_300_PTHJoverPTH_GT15 + getBin(higgs.pt(), {200, 300, 450, 650}));
738  else
739  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
740  } else
741  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
742  }
743  if (Njets == 0)
744  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
745  if (Njets == 1)
746  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
747  if (Njets > 1) {
748  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
749  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
750  //VBF topology
751  if (vbfTopo)
753  //Njets >= 2jets without VBF topology (mjj<350)
754  if (pTHjj < 25)
755  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
756  else
757  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
758  }
759  }
760 
761  // 2. Electroweak qq->Hqq Stage 1.2 categories
762  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
763  if (std::abs(higgs.rapidity()) > 2.5)
764  return QQ2HQQ_FWDH;
765  int Njets = jets.size();
766  if (Njets == 0)
767  return QQ2HQQ_0J;
768  else if (Njets == 1)
769  return QQ2HQQ_1J;
770  else if (Njets >= 2) {
771  double mjj = (jets[0].mom() + jets[1].mom()).mass();
772  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
773  if (mjj < 350) {
774  if (pTHjj < 25)
775  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
776  else
777  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
778  } else { //mjj>350 GeV
779  if (higgs.pt() < 200) {
781  } else {
783  }
784  }
785  }
786  }
787  // 3. WH->Hlv categories
788  else if (prodMode == HTXS::WH) {
789  if (fwdHiggs)
790  return QQ2HLNU_FWDH;
791  int Njets = jets.size();
792  if (Njets == 0)
793  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
794  if (Njets == 1)
795  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
796  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
797  }
798  // 4. qq->ZH->llH categories
799  else if (prodMode == HTXS::QQ2ZH) {
800  if (fwdHiggs)
801  return QQ2HLL_FWDH;
802  int Njets = jets.size();
803  if (Njets == 0)
804  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
805  if (Njets == 1)
806  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
807  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
808  }
809  // 5. gg->ZH->llH categories
810  else if (prodMode == HTXS::GG2ZH) {
811  if (fwdHiggs)
812  return GG2HLL_FWDH;
813  int Njets = jets.size();
814  if (Njets == 0)
815  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
816  if (Njets == 1)
817  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
818  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
819  }
820  // 6.ttH,bbH,tH categories
821  else if (prodMode == HTXS::TTH) {
822  if (fwdHiggs)
823  return TTH_FWDH;
824  else
825  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300, 450}));
826  } else if (prodMode == HTXS::BBH)
827  return Category(BBH_FWDH + ctrlHiggs);
828  else if (prodMode == HTXS::TH)
829  return Category(TH_FWDH + ctrlHiggs);
830  return UNKNOWN;
831  }
832 
834 
835 
838 
840  void setHiggsProdMode( HTXS::HiggsProdMode prodMode ){ m_HiggsProdMode = prodMode; }
841 
845  void init() override {
846  printf("==============================================================\n");
847  printf("======== HiggsTemplateCrossSections Initialization =========\n");
848  printf("==============================================================\n");
849  // check that the production mode has been set
850  // if running in standalone Rivet the production mode is set through an env variable
852  char *pm_env = getenv("HIGGSPRODMODE");
853  string pm(pm_env==nullptr?"":pm_env);
854  if ( pm == "GGF" ) m_HiggsProdMode = HTXS::GGF;
855  else if ( pm == "VBF" ) m_HiggsProdMode = HTXS::VBF;
856  else if ( pm == "WH" ) m_HiggsProdMode = HTXS::WH;
857  else if ( pm == "ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
858  else if ( pm == "QQ2ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
859  else if ( pm == "GG2ZH" ) m_HiggsProdMode = HTXS::GG2ZH;
860  else if ( pm == "TTH" ) m_HiggsProdMode = HTXS::TTH;
861  else if ( pm == "BBH" ) m_HiggsProdMode = HTXS::BBH;
862  else if ( pm == "TH" ) m_HiggsProdMode = HTXS::TH;
863  else {
864  MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
865  }
866  }
867 
868  // Projections for final state particles
869  const FinalState FS;
870  addProjection(FS,"FS");
871 
872  // initialize the histograms with for each of the stages
874  m_sumw = 0.0;
875  printf("==============================================================\n");
876  printf("======== Higgs prod mode %d =========\n",m_HiggsProdMode);
877  printf("======== Sucessful Initialization =========\n");
878  printf("==============================================================\n");
879  }
880 
881  // Perform the per-event analysis
882  void analyze(const Event& event) override {
883 
884  // get the classification
885  HiggsClassification cat = classifyEvent(event,m_HiggsProdMode);
886 
887  // Fill histograms: categorization --> linerize the categories
888  const double weight = event.weight();
889  m_sumw += weight;
890 
891  int F=cat.stage0_cat%10, P=cat.stage1_cat_pTjet30GeV/100;
892  hist_stage0->fill( cat.stage0_cat/10*2+F, weight );
893 
894  // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
895  vector<int> offset({0,1,13,19,24,29,33,35,37,39});
896  int off = offset[P];
897  hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV%100 + off, weight);
898  hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV%100 + off, weight);
899 
900  // Stage 1_2 enum offsets for each production mode: GGF=17, VBF=11, WH= 6, QQ2ZH=6, GG2ZH=6, TTH=6, BBH=2, TH=2
901  static vector<int> offset1_2({0, 1, 18, 29, 35, 41, 47, 53, 55, 57});
902  int off1_2 = offset1_2[P];
903  // Stage 1_2 Fine enum offsets for each production mode: GGF=28, VBF=25, WH= 16, QQ2ZH=16, GG2ZH=16, TTH=7, BBH=2, TH=2
904  static vector<int> offset1_2f({0, 1, 29, 54, 70, 86, 102, 109, 111, 113});
905  int off1_2f = offset1_2f[P];
906  hist_stage1_2_pTjet25->fill(cat.stage1_2_cat_pTjet25GeV % 100 + off1_2, weight);
907  hist_stage1_2_pTjet30->fill(cat.stage1_2_cat_pTjet30GeV % 100 + off1_2, weight);
908  hist_stage1_2_fine_pTjet25->fill(cat.stage1_2_fine_cat_pTjet25GeV % 100 + off1_2f, weight);
909  hist_stage1_2_fine_pTjet30->fill(cat.stage1_2_fine_cat_pTjet30GeV % 100 + off1_2f, weight);
910 
911  // Fill histograms: variables used in the categorization
912  hist_pT_Higgs->fill(cat.higgs.pT(),weight);
913  hist_y_Higgs->fill(cat.higgs.rapidity(),weight);
914  hist_pT_V->fill(cat.V.pT(),weight);
915 
916  hist_Njets25->fill(cat.jets25.size(),weight);
917  hist_Njets30->fill(cat.jets30.size(),weight);
918 
919  hist_isZ2vv->fill(cat.isZ2vvDecay, weight);
920 
921  // Jet variables. Use jet collection with pT threshold at 30 GeV
922  if (!cat.jets30.empty()) hist_pT_jet1->fill(cat.jets30[0].pt(),weight);
923  if (cat.jets30.size()>=2) {
924  const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
925  hist_deltay_jj->fill(std::abs(j1.rapidity()-j2.rapidity()),weight);
926  hist_dijet_mass->fill((j1+j2).mass(),weight);
927  hist_pT_Hjj->fill((j1+j2+cat.higgs.momentum()).pt(),weight);
928  }
929  }
930 
932  MSG_INFO (" ====================================================== ");
933  MSG_INFO (" Higgs Template X-Sec Categorization Tool ");
934  MSG_INFO (" Status Code Summary ");
935  MSG_INFO (" ====================================================== ");
936  bool allSuccess = (numEvents()==m_errorCount[HTXS::SUCCESS]);
937  if ( allSuccess ) MSG_INFO (" >>>> All "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized!");
938  else{
939  MSG_INFO (" >>>> "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized");
940  MSG_INFO (" >>>> --> the following errors occured:");
941  MSG_INFO (" >>>> "<< m_errorCount[HTXS::PRODMODE_DEFINED] <<" had an undefined Higgs production mode.");
942  MSG_INFO (" >>>> "<< m_errorCount[HTXS::MOMENTUM_CONSERVATION] <<" failed momentum conservation.");
943  MSG_INFO (" >>>> "<< m_errorCount[HTXS::HIGGS_IDENTIFICATION] <<" failed to identify a valid Higgs boson.");
944  MSG_INFO (" >>>> "<< m_errorCount[HTXS::HS_VTX_IDENTIFICATION] <<" failed to identify the hard scatter vertex.");
945  MSG_INFO (" >>>> "<< m_errorCount[HTXS::VH_IDENTIFICATION] <<" VH: to identify a valid V-boson.");
946  MSG_INFO (" >>>> "<< m_errorCount[HTXS::TOP_W_IDENTIFICATION] <<" failed to identify valid Ws from top decay.");
947  }
948  MSG_INFO (" ====================================================== ");
949  MSG_INFO (" ====================================================== ");
950  }
951 
952 
953  void finalize() override {
955  double sf = m_sumw > 0 ? 1.0 / m_sumw : 1.0;
956  for (auto hist : {hist_stage0,
963  hist_Njets25,
964  hist_Njets30,
966  hist_y_Higgs,
967  hist_pT_V,
968  hist_pT_jet1,
971  hist_pT_Hjj})
972  scale(hist, sf);
973  }
974 
975  /*
976  * initialize histograms
977  */
978 
980  hist_stage0 = bookHisto1D("HTXS_stage0",20,0,20);
981  hist_stage1_pTjet25 = bookHisto1D("HTXS_stage1_pTjet25",40,0,40);
982  hist_stage1_pTjet30 = bookHisto1D("HTXS_stage1_pTjet30",40,0,40);
983  hist_stage1_2_pTjet25 = bookHisto1D("HTXS_stage1_2_pTjet25",57,0,57);
984  hist_stage1_2_pTjet30 = bookHisto1D("HTXS_stage1_2_pTjet30",57,0,57);
985  hist_stage1_2_fine_pTjet25 = bookHisto1D("HTXS_stage1_2_fine_pTjet25",113,0,113);
986  hist_stage1_2_fine_pTjet30 = bookHisto1D("HTXS_stage1_2_fine_pTjet30",113,0,113);
987  hist_pT_Higgs = bookHisto1D("pT_Higgs",80,0,400);
988  hist_y_Higgs = bookHisto1D("y_Higgs",80,-4,4);
989  hist_pT_V = bookHisto1D("pT_V",80,0,400);
990  hist_pT_jet1 = bookHisto1D("pT_jet1",80,0,400);
991  hist_deltay_jj = bookHisto1D("deltay_jj",50,0,10);
992  hist_dijet_mass = bookHisto1D("m_jj",50,0,2000);
993  hist_pT_Hjj = bookHisto1D("pT_Hjj",50,0,250);
994  hist_Njets25 = bookHisto1D("Njets25",10,0,10);
995  hist_Njets30 = bookHisto1D("Njets30",10,0,10);
996  hist_isZ2vv = bookHisto1D("isZ2vv",2,0,2);
997  }
999 
1000  /*
1001  * initialize private members used in the classification procedure
1002  */
1003 
1004  private:
1005  double m_sumw;
1007  std::map<HTXS::ErrorCode,size_t> m_errorCount;
1008  Histo1DPtr hist_stage0;
1016  Histo1DPtr hist_isZ2vv;
1017  };
1018 
1019  // the PLUGIN only needs to be decleared when running standalone Rivet
1020  // and causes compilation / linking issues if included in Athena / RootCore
1021  //check for Rivet environment variable RIVET_ANALYSIS_PATH
1022 #ifdef RIVET_ANALYSIS_PATH
1023  // The hook for the plugin system
1025 #endif
1026 
1027 }
1028 
1029 #endif
1030 
std::map< HTXS::ErrorCode, size_t > m_errorCount
int vbfTopology(const Jets &jets, const Particle &higgs)
VBF topolog selection 0 = fail loose selction: m_jj > 400 GeV and Dy_jj > 2.8 1 pass loose...
bool isVH(HTXS::HiggsProdMode p)
Whether the Higgs is produced in association with a vector boson (VH)
failed to identify associated vector boson
int getBin(double x, std::vector< double > bins)
Return bin index of x given the provided bin edges. 0=first bin, -1=underflow bin.
ErrorCode
Error code: whether the classification was successful or failed.
TPRegexp parents
Definition: eve_filter.cc:21
HTXS::Stage1_2_Fine::Category getStage1_2_Fine_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1.2 Fine categorization.
bool originateFrom(const Particle &p, const Particle &p2)
Whether particle p originates from p2.
successful classification
bool hasParent(const GenParticle *ptcl, int pdgID)
Checks whether the input particle has a parent with a given PDGID.
int vbfTopology_Stage1_X_Fine(const Jets &jets, const Particle &higgs)
VBF topology selection for Stage1.1 and Stage 1.2 Fine 0 = fail loose selection: m_jj > 350 GeV 1 pas...
0: Unidentified isolated particle
Definition: ParticleCode.h:19
HTXS::Stage1::Category getStage1Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1 categorization.
failed to identify Higgs boson decay products
Definition: weight.py:1
HiggsProdMode
Higgs production modes, corresponding to input sample.
HTXS::Stage1_1::Category getStage1_1_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1.1 categorization.
bool ChLeptonDecay(const Particle &p)
Return true if particle decays to charged leptons.
void init() override
default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Hig...
void analyze(const Event &event) override
def cat(path)
Definition: eostools.py:400
HTXS::Stage1_1_Fine::Category getStage1_1_Fine_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1_1 categorization.
bool hasChild(const GenParticle *ptcl, int pdgID)
Checks whether the input particle has a child with a given PDGID.
bool quarkDecay(const Particle &p)
Return true is particle decays to quarks.
failed momentum conservation
vector< PseudoJet > jets
DECLARE_RIVET_PLUGIN(CMS_2013_I1224539_DIJET)
Particle getLastInstance(Particle ptcl)
follow a "propagating" particle and return its last instance
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Namespace for Stage0 categorization.
Higgs Template Cross Section namespace.
failed to identify associated vector boson decay products
void setHiggsProdMode(HTXS::HiggsProdMode prodMode)
Sets the Higgs production mode.
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
double p2[4]
Definition: TauolaWrapper.h:90
HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode)
Main classificaion method.
failed to identify Higgs boson
part
Definition: HCALResponse.h:20
Rivet routine for classifying MC events according to the Higgs template cross section categories...
failed to identify hard scatter vertex
std::pair< OmniClusterRef, TrackingParticleRef > P
tuple msg
Definition: mps_check.py:277
HTXS::Stage0::Category getStage0Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Particle &V)
Stage-0 HTXS categorization.
HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err, std::string msg="", int NmaxWarnings=20)
Returns the classification object with the error code set. Prints an warning message, and keeps track of number of errors.
bool originateFrom(const Particle &p, const Particles &ptcls)
Whether particle p originate from any of the ptcls.
HTXS::Stage1_2::Category getStage1_2_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1.2 categorization.
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
production mode not defined
int vbfTopology_Stage1_X(const Jets &jets, const Particle &higgs)
VBF topology selection Stage1.1 and Stage1.2 0 = fail loose selection: m_jj > 350 GeV 1 pass loose...
Definition: event.py:1
failed to identify top decay