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
24  HiggsTemplateCrossSections() : Analysis("HiggsTemplateCrossSections"), m_HiggsProdMode(HTXS::UNKNOWN) {}
25 
26  public:
31 
34  if (ptcl.genParticle()->end_vertex()) {
35  if (!hasChild(ptcl.genParticle(), ptcl.pid()))
36  return ptcl;
37  else
38  return getLastInstance(ptcl.children()[0]);
39  }
40  return ptcl;
41  }
42 
44  bool originateFrom(const Particle &p, const Particles &ptcls) {
45  const ConstGenVertexPtr prodVtx = p.genParticle()->production_vertex();
46  if (prodVtx == nullptr)
47  return false;
48  // for each ancestor, check if it matches any of the input particles
49  for (const auto &ancestor : HepMCUtils::particles(prodVtx, Relatives::ANCESTORS)) {
50  for (const auto &part : ptcls)
51  if (ancestor == part.genParticle())
52  return true;
53  }
54  // if we get here, no ancetor matched any input particle
55  return false;
56  }
57 
59  bool originateFrom(const Particle &p, const Particle &p2) {
60  Particles ptcls = {p2};
61  return originateFrom(p, ptcls);
62  }
63 
65  bool hasChild(const ConstGenParticlePtr &ptcl, int pdgID) {
66  for (const auto &child : Particle(*ptcl).children()) {
67  if (child.pid() == pdgID) {
68  return true;
69  }
70  }
71  return false;
72  }
73 
75  bool hasParent(const ConstGenParticlePtr &ptcl, int pdgID) {
76  for (auto parent : HepMCUtils::particles(ptcl->production_vertex(), Relatives::PARENTS))
77  if (parent->pdg_id() == pdgID)
78  return true;
79  return false;
80  }
81 
83  bool quarkDecay(const Particle &p) {
84  for (const auto &child : p.children()) {
85  if (PID::isQuark(child.pid())) {
86  return true;
87  }
88  }
89  return false;
90  }
91 
93  bool ChLeptonDecay(const Particle &p) {
94  for (const auto &child : p.children()) {
95  if (PID::isChargedLepton(child.pid())) {
96  return true;
97  }
98  }
99  return false;
100  }
101 
104  HiggsClassification error(HiggsClassification &cat,
106  std::string msg = "",
107  int NmaxWarnings = 20) {
108  // Set the error, and keep statistics
109  cat.errorCode = err;
110  ++m_errorCount[err];
111 
112  // Print warning message to the screen/log
113  static std::atomic<int> Nwarnings{0};
114  if (!msg.empty() && ++Nwarnings < NmaxWarnings)
115  MSG_WARNING(msg);
116 
117  return cat;
118  }
120 
122  HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode) {
124  m_HiggsProdMode = prodMode;
125 
126  // the classification object
127  HiggsClassification cat;
128  cat.prodMode = prodMode;
129  cat.errorCode = HTXS::UNDEFINED;
130  cat.stage0_cat = HTXS::Stage0::UNKNOWN;
131  cat.stage1_cat_pTjet25GeV = HTXS::Stage1::UNKNOWN;
132  cat.stage1_cat_pTjet30GeV = HTXS::Stage1::UNKNOWN;
133  cat.stage1_1_cat_pTjet25GeV = HTXS::Stage1_1::UNKNOWN;
134  cat.stage1_1_cat_pTjet30GeV = HTXS::Stage1_1::UNKNOWN;
135  cat.stage1_1_fine_cat_pTjet25GeV = HTXS::Stage1_1_Fine::UNKNOWN;
136  cat.stage1_1_fine_cat_pTjet30GeV = HTXS::Stage1_1_Fine::UNKNOWN;
137  cat.stage1_2_cat_pTjet25GeV = HTXS::Stage1_2::UNKNOWN;
138  cat.stage1_2_cat_pTjet30GeV = HTXS::Stage1_2::UNKNOWN;
139  cat.stage1_2_fine_cat_pTjet25GeV = HTXS::Stage1_2_Fine::UNKNOWN;
140  cat.stage1_2_fine_cat_pTjet30GeV = HTXS::Stage1_2_Fine::UNKNOWN;
141 
142  if (prodMode == HTXS::UNKNOWN)
143  return error(cat,
145  "Unkown Higgs production mechanism. Cannot classify event."
146  " Classification for all events will most likely fail.");
147 
148  /*****
149  * Step 1.
150  * Idenfify the Higgs boson and the hard scatter vertex
151  * There should be only one of each.
152  */
153 
154  ConstGenVertexPtr HSvtx = nullptr;
155  int Nhiggs = 0;
156  for (const ConstGenParticlePtr &ptcl : HepMCUtils::particles(event.genEvent())) {
157  // a) Reject all non-Higgs particles
158  if (!PID::isHiggs(ptcl->pdg_id()))
159  continue;
160  // b) select only the final Higgs boson copy, prior to decay
161  if (ptcl->end_vertex() && !hasChild(ptcl, PID::HIGGS)) {
162  cat.higgs = Particle(ptcl);
163  ++Nhiggs;
164  }
165  // c) HepMC3 does not provide signal_proces_vertex anymore
166  // set hard-scatter vertex based on first Higgs boson
167  if (HSvtx == nullptr && ptcl->production_vertex() && !hasParent(ptcl, PID::HIGGS))
168  HSvtx = ptcl->production_vertex();
169  }
170 
171  // Make sure things are in order so far
172  if (Nhiggs != 1)
173  return error(cat,
175  "Current event has " + std::to_string(Nhiggs) + " Higgs bosons. There must be only one.");
176  if (cat.higgs.children().size() < 2)
177  return error(cat, HTXS::HIGGS_DECAY_IDENTIFICATION, "Could not identify Higgs boson decay products.");
178 
179  if (HSvtx == nullptr)
180  return error(cat, HTXS::HS_VTX_IDENTIFICATION, "Cannot find hard-scatter vertex of current event.");
181 
182  /*****
183  * Step 2.
184  * Identify associated vector bosons
185  */
186 
187  // Find associated vector bosons
188  bool is_uncatdV = false;
189  Particles uncatV_decays;
190  FourMomentum uncatV_p4(0, 0, 0, 0);
191  FourVector uncatV_v4(0, 0, 0, 0);
192  int nWs = 0, nZs = 0;
193  if (isVH(prodMode)) {
194  for (auto ptcl : HepMCUtils::particles(HSvtx, Relatives::CHILDREN)) {
195  if (PID::isW(ptcl->pdg_id())) {
196  ++nWs;
197  cat.V = Particle(ptcl);
198  }
199  if (PID::isZ(ptcl->pdg_id())) {
200  ++nZs;
201  cat.V = Particle(ptcl);
202  }
203  }
204  if (nWs + nZs > 0)
205  cat.V = getLastInstance(cat.V);
206  else {
207  for (auto ptcl : HepMCUtils::particles(HSvtx, Relatives::CHILDREN)) {
208  if (!PID::isHiggs(ptcl->pdg_id())) {
209  uncatV_decays += Particle(ptcl);
210  uncatV_p4 += Particle(ptcl).momentum();
211  // uncatV_v4 += Particle(ptcl).origin();
212  }
213  }
214  // is_uncatdV = true; cat.V = Particle(24,uncatV_p4,uncatV_v4);
215  is_uncatdV = true;
216  cat.V = Particle(24, uncatV_p4);
217  }
218  }
219 
220  if (!is_uncatdV) {
221  if (isVH(prodMode) && !cat.V.genParticle()->end_vertex())
222  return error(cat, HTXS::VH_DECAY_IDENTIFICATION, "Vector boson does not decay!");
223 
224  if (isVH(prodMode) && cat.V.children().size() < 2)
225  return error(cat, HTXS::VH_DECAY_IDENTIFICATION, "Vector boson does not decay!");
226 
227  if ((prodMode == HTXS::WH && (nZs > 0 || nWs != 1)) ||
228  ((prodMode == HTXS::QQ2ZH || prodMode == HTXS::GG2ZH) && (nZs != 1 || nWs > 0)))
229  return error(cat,
231  "Found " + std::to_string(nWs) + " W-bosons and " + std::to_string(nZs) +
232  " Z-bosons. Inconsitent with VH expectation.");
233  }
234 
235  // Find and store the W-bosons from ttH->WbWbH
236  Particles Ws;
237  if (prodMode == HTXS::TTH || prodMode == HTXS::TH) {
238  // loop over particles produced in hard-scatter vertex
239  for (auto ptcl : HepMCUtils::particles(HSvtx, Relatives::CHILDREN)) {
240  if (!PID::isTop(ptcl->pdg_id()))
241  continue;
242  Particle top = getLastInstance(Particle(ptcl));
243  if (top.genParticle()->end_vertex())
244  for (const auto &child : top.children())
245  if (PID::isW(child.pid()))
246  Ws += getLastInstance(child);
247  }
248  }
249 
250  // Make sure result make sense
251  if ((prodMode == HTXS::TTH && Ws.size() < 2) || (prodMode == HTXS::TH && Ws.empty()))
252  return error(cat, HTXS::TOP_W_IDENTIFICATION, "Failed to identify W-boson(s) from t-decay!");
253 
254  /*****
255  * Step 3.
256  * Build jets
257  * Make sure all stable particles are present
258  */
259 
260  // Create a list of the vector bosons that decay leptonically
261  // Either the vector boson produced in association with the Higgs boson,
262  // or the ones produced from decays of top quarks produced with the Higgs
263  Particles leptonicVs;
264  if (!is_uncatdV) {
265  if (isVH(prodMode) && !quarkDecay(cat.V))
266  leptonicVs += cat.V;
267  } else
268  leptonicVs = uncatV_decays;
269  for (const auto &W : Ws)
270  if (W.genParticle()->end_vertex() && !quarkDecay(W))
271  leptonicVs += W;
272 
273  // Obtain all stable, final-state particles
274  const Particles FS = apply<FinalState>(event, "FS").particles();
275  Particles hadrons;
276  FourMomentum sum(0, 0, 0, 0), vSum(0, 0, 0, 0), hSum(0, 0, 0, 0);
277  for (const Particle &p : FS) {
278  // Add up the four momenta of all stable particles as a cross check
279  sum += p.momentum();
280  // ignore particles from the Higgs boson
281  if (originateFrom(p, cat.higgs)) {
282  hSum += p.momentum();
283  continue;
284  }
285  // Cross-check the V decay products for VH
286  if (isVH(prodMode) && !is_uncatdV && originateFrom(p, Ws))
287  vSum += p.momentum();
288  // ignore final state particles from leptonic V decays
289  if (!leptonicVs.empty() && originateFrom(p, leptonicVs))
290  continue;
291  // All particles reaching here are considered hadrons and will be used to build jets
292  hadrons += p;
293  }
294 
295  cat.p4decay_higgs = hSum;
296  cat.p4decay_V = vSum;
297 
298  FinalState fps_temp;
299  FastJets jets(fps_temp, JetAlg::ANTIKT, 0.4);
300  jets.calc(hadrons);
301 
302  cat.jets25 = jets.jetsByPt(Cuts::pT > 25.0);
303  cat.jets30 = jets.jetsByPt(Cuts::pT > 30.0);
304 
305  // check that four mometum sum of all stable particles satisfies momentum consevation
306  /*
307  if ( sum.pt()>0.1 )
308  return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
309  std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
310 */
311  // check if V-boson was not included in the event record but decay particles were
312  // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
313  if (is_uncatdV)
314  return error(cat, HTXS::VH_IDENTIFICATION, "Failed to identify associated V-boson!");
315 
316  /*****
317  * Step 4.
318  * Classify and save output
319  */
320 
321  // Apply the categorization categorization
322  cat.isZ2vvDecay = false;
323  if ((prodMode == HTXS::GG2ZH || prodMode == HTXS::QQ2ZH) && !quarkDecay(cat.V) && !ChLeptonDecay(cat.V))
324  cat.isZ2vvDecay = true;
325  cat.stage0_cat = getStage0Category(prodMode, cat.higgs, cat.V);
326  cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode, cat.higgs, cat.jets25, cat.V);
327  cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode, cat.higgs, cat.jets30, cat.V);
328  cat.stage1_1_cat_pTjet25GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets25, cat.V);
329  cat.stage1_1_cat_pTjet30GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets30, cat.V);
330  cat.stage1_1_fine_cat_pTjet25GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
331  cat.stage1_1_fine_cat_pTjet30GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
332  cat.stage1_2_cat_pTjet25GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets25, cat.V);
333  cat.stage1_2_cat_pTjet30GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets30, cat.V);
334  cat.stage1_2_fine_cat_pTjet25GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
335  cat.stage1_2_fine_cat_pTjet30GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
336 
337  cat.errorCode = HTXS::SUCCESS;
339  ++m_sumevents;
340 
341  return cat;
342  }
343 
349 
351  int getBin(double x, std::vector<double> bins) {
352  for (size_t i = 1; i < bins.size(); ++i)
353  if (x < bins[i])
354  return i - 1;
355  return bins.size() - 1;
356  }
357 
361  int vbfTopology(const Jets &jets, const Particle &higgs) {
362  if (jets.size() < 2)
363  return 0;
364  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
365  bool VBFtopo = (j1 + j2).mass() > 400.0 && std::abs(j1.rapidity() - j2.rapidity()) > 2.8;
366  return VBFtopo ? (j1 + j2 + higgs.momentum()).pt() < 25 ? 2 : 1 : 0;
367  }
368 
373  int vbfTopology_Stage1_X(const Jets &jets, const Particle &higgs) {
374  if (jets.size() < 2)
375  return 0;
376  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
377  double mjj = (j1 + j2).mass();
378  if (mjj > 350 && mjj <= 700)
379  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
380  else if (mjj > 700)
381  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
382  else
383  return 0;
384  }
385 
392  int vbfTopology_Stage1_X_Fine(const Jets &jets, const Particle &higgs) {
393  if (jets.size() < 2)
394  return 0;
395  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
396  double mjj = (j1 + j2).mass();
397  if (mjj > 350 && mjj <= 700)
398  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
399  else if (mjj > 700 && mjj <= 1000)
400  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
401  else if (mjj > 1000 && mjj <= 1500)
402  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 5 : 6;
403  else if (mjj > 1500)
404  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 7 : 8;
405  else
406  return 0;
407  }
408 
410  bool isVH(HTXS::HiggsProdMode p) { return p == HTXS::WH || p == HTXS::QQ2ZH || p == HTXS::GG2ZH; }
411 
414  const Particle &higgs,
415  const Particle &V) {
416  using namespace HTXS::Stage0;
417  int ctrlHiggs = std::abs(higgs.rapidity()) < 2.5;
418  // Special cases first, qq→Hqq
419  if ((prodMode == HTXS::WH || prodMode == HTXS::QQ2ZH) && quarkDecay(V)) {
420  return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
421  } else if (prodMode == HTXS::GG2ZH && quarkDecay(V)) {
422  return Category(HTXS::GGF * 10 + ctrlHiggs);
423  }
424  // General case after
425  return Category(prodMode * 10 + ctrlHiggs);
426  }
427 
430  const Particle &higgs,
431  const Jets &jets,
432  const Particle &V) {
433  using namespace HTXS::Stage1;
434  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
435  double pTj1 = !jets.empty() ? jets[0].momentum().pt() : 0;
436  int vbfTopo = vbfTopology(jets, higgs);
437 
438  // 1. GGF Stage 1 categories
439  // Following YR4 write-up: XXXXX
440  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
441  if (fwdHiggs)
442  return GG2H_FWDH;
443  if (Njets == 0)
444  return GG2H_0J;
445  else if (Njets == 1)
446  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
447  else if (Njets >= 2) {
448  // events with pT_H>200 get priority over VBF cuts
449  if (higgs.pt() <= 200) {
450  if (vbfTopo == 2)
451  return GG2H_VBFTOPO_JET3VETO;
452  else if (vbfTopo == 1)
453  return GG2H_VBFTOPO_JET3;
454  }
455  // Njets >= 2jets without VBF topology
456  return Category(GG2H_GE2J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
457  }
458  }
459  // 2. Electroweak qq->Hqq Stage 1 categories
460  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
461  if (std::abs(higgs.rapidity()) > 2.5)
462  return QQ2HQQ_FWDH;
463  if (pTj1 > 200)
464  return QQ2HQQ_PTJET1_GT200;
465  if (vbfTopo == 2)
467  if (vbfTopo == 1)
468  return QQ2HQQ_VBFTOPO_JET3;
469  double mjj = jets.size() > 1 ? (jets[0].mom() + jets[1].mom()).mass() : 0;
470  if (60 < mjj && mjj < 120)
471  return QQ2HQQ_VH2JET;
472  return QQ2HQQ_REST;
473  }
474  // 3. WH->Hlv categories
475  else if (prodMode == HTXS::WH) {
476  if (fwdHiggs)
477  return QQ2HLNU_FWDH;
478  else if (V.pt() < 150)
479  return QQ2HLNU_PTV_0_150;
480  else if (V.pt() > 250)
481  return QQ2HLNU_PTV_GT250;
482  // 150 < pTV/GeV < 250
484  }
485  // 4. qq->ZH->llH categories
486  else if (prodMode == HTXS::QQ2ZH) {
487  if (fwdHiggs)
488  return QQ2HLL_FWDH;
489  else if (V.pt() < 150)
490  return QQ2HLL_PTV_0_150;
491  else if (V.pt() > 250)
492  return QQ2HLL_PTV_GT250;
493  // 150 < pTV/GeV < 250
495  }
496  // 5. gg->ZH->llH categories
497  else if (prodMode == HTXS::GG2ZH) {
498  if (fwdHiggs)
499  return GG2HLL_FWDH;
500  if (V.pt() < 150)
501  return GG2HLL_PTV_0_150;
502  else if (jets.empty())
503  return GG2HLL_PTV_GT150_0J;
504  return GG2HLL_PTV_GT150_GE1J;
505  }
506  // 6.ttH,bbH,tH categories
507  else if (prodMode == HTXS::TTH)
508  return Category(TTH_FWDH + ctrlHiggs);
509  else if (prodMode == HTXS::BBH)
510  return Category(BBH_FWDH + ctrlHiggs);
511  else if (prodMode == HTXS::TH)
512  return Category(TH_FWDH + ctrlHiggs);
513  return UNKNOWN;
514  }
515 
518  const Particle &higgs,
519  const Jets &jets,
520  const Particle &V) {
521  using namespace HTXS::Stage1_1;
522  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
523  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
524 
525  // 1. GGF Stage 1 categories
526  // Following YR4 write-up: XXXXX
527  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
528  if (fwdHiggs)
529  return GG2H_FWDH;
530  if (higgs.pt() > 200)
531  return GG2H_PTH_GT200;
532  if (Njets == 0)
533  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
534  if (Njets == 1)
535  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
536  if (Njets > 1) {
537  //VBF topology
538  if (vbfTopo)
539  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
540  //Njets >= 2jets without VBF topology (mjj<350)
541  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
542  }
543  }
544 
545  // 2. Electroweak qq->Hqq Stage 1.1 categories
546  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
547  if (std::abs(higgs.rapidity()) > 2.5)
548  return QQ2HQQ_FWDH;
549  int Njets = jets.size();
550  if (Njets == 0)
551  return QQ2HQQ_0J;
552  else if (Njets == 1)
553  return QQ2HQQ_1J;
554  else if (Njets >= 2) {
555  double mjj = (jets[0].mom() + jets[1].mom()).mass();
556  if (mjj < 60)
557  return QQ2HQQ_MJJ_0_60;
558  else if (60 < mjj && mjj < 120)
559  return QQ2HQQ_MJJ_60_120;
560  else if (120 < mjj && mjj < 350)
561  return QQ2HQQ_MJJ_120_350;
562  else if (mjj > 350) {
563  if (higgs.pt() > 200)
565  if (vbfTopo)
566  return Category(QQ2HQQ_MJJ_GT350_PTH_GT200 + vbfTopo);
567  }
568  }
569  }
570  // 3. WH->Hlv categories
571  else if (prodMode == HTXS::WH) {
572  if (fwdHiggs)
573  return QQ2HLNU_FWDH;
574  else if (V.pt() < 75)
575  return QQ2HLNU_PTV_0_75;
576  else if (V.pt() < 150)
577  return QQ2HLNU_PTV_75_150;
578  else if (V.pt() > 250)
579  return QQ2HLNU_PTV_GT250;
580  // 150 < pTV/GeV < 250
582  }
583  // 4. qq->ZH->llH categories
584  else if (prodMode == HTXS::QQ2ZH) {
585  if (fwdHiggs)
586  return QQ2HLL_FWDH;
587  else if (V.pt() < 75)
588  return QQ2HLL_PTV_0_75;
589  else if (V.pt() < 150)
590  return QQ2HLL_PTV_75_150;
591  else if (V.pt() > 250)
592  return QQ2HLL_PTV_GT250;
593  // 150 < pTV/GeV < 250
595  }
596  // 5. gg->ZH->llH categories
597  else if (prodMode == HTXS::GG2ZH) {
598  if (fwdHiggs)
599  return GG2HLL_FWDH;
600  else if (V.pt() < 75)
601  return GG2HLL_PTV_0_75;
602  else if (V.pt() < 150)
603  return GG2HLL_PTV_75_150;
604  else if (V.pt() > 250)
605  return GG2HLL_PTV_GT250;
607  }
608  // 6.ttH,bbH,tH categories
609  else if (prodMode == HTXS::TTH)
610  return Category(TTH_FWDH + ctrlHiggs);
611  else if (prodMode == HTXS::BBH)
612  return Category(BBH_FWDH + ctrlHiggs);
613  else if (prodMode == HTXS::TH)
614  return Category(TH_FWDH + ctrlHiggs);
615  return UNKNOWN;
616  }
617 
620  const Particle &higgs,
621  const Jets &jets,
622  const Particle &V) {
623  using namespace HTXS::Stage1_1_Fine;
624  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
625  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
626 
627  // 1. GGF Stage 1.1 categories
628  // Following YR4 write-up: XXXXX
629  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
630  if (fwdHiggs)
631  return GG2H_FWDH;
632  if (higgs.pt() > 200)
633  return GG2H_PTH_GT200;
634  if (Njets == 0)
635  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
636  if (Njets == 1)
637  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
638  if (Njets > 1) {
639  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
640  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
641  //VBF topology
642  if (vbfTopo)
643  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
644  //Njets >= 2jets without VBF topology (mjj<350)
645  if (pTHjj < 25)
646  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
647  else
648  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
649  }
650  }
651 
652  // 2. Electroweak qq->Hqq Stage 1.1 categories
653  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
654  if (std::abs(higgs.rapidity()) > 2.5)
655  return QQ2HQQ_FWDH;
656  int Njets = jets.size();
657  if (Njets == 0)
658  return QQ2HQQ_0J;
659  else if (Njets == 1)
660  return QQ2HQQ_1J;
661  else if (Njets >= 2) {
662  double mjj = (jets[0].mom() + jets[1].mom()).mass();
663  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
664  if (mjj < 350) {
665  if (pTHjj < 25)
666  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
667  else
668  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
669  } else { //mjj>350 GeV
670  if (higgs.pt() < 200) {
671  return Category(QQ2HQQ_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
672  } else {
674  }
675  }
676  }
677  }
678  // 3. WH->Hlv categories
679  else if (prodMode == HTXS::WH) {
680  if (fwdHiggs)
681  return QQ2HLNU_FWDH;
682  int Njets = jets.size();
683  if (Njets == 0)
684  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
685  if (Njets == 1)
686  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
687  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
688  }
689  // 4. qq->ZH->llH categories
690  else if (prodMode == HTXS::QQ2ZH) {
691  if (fwdHiggs)
692  return QQ2HLL_FWDH;
693  int Njets = jets.size();
694  if (Njets == 0)
695  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
696  if (Njets == 1)
697  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
698  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
699  }
700  // 5. gg->ZH->llH categories
701  else if (prodMode == HTXS::GG2ZH) {
702  if (fwdHiggs)
703  return GG2HLL_FWDH;
704  int Njets = jets.size();
705  if (Njets == 0)
706  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
707  if (Njets == 1)
708  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
709  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
710  }
711  // 6.ttH,bbH,tH categories
712  else if (prodMode == HTXS::TTH)
713  return Category(TTH_FWDH + ctrlHiggs);
714  else if (prodMode == HTXS::BBH)
715  return Category(BBH_FWDH + ctrlHiggs);
716  else if (prodMode == HTXS::TH)
717  return Category(TH_FWDH + ctrlHiggs);
718  return UNKNOWN;
719  }
720 
723  const Particle &higgs,
724  const Jets &jets,
725  const Particle &V) {
726  using namespace HTXS::Stage1_2;
727  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
728  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
729 
730  // 1. GGF Stage 1 categories
731  // Following YR4 write-up: XXXXX
732  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
733  if (fwdHiggs)
734  return GG2H_FWDH;
735  if (higgs.pt() > 200)
736  return Category(GG2H_PTH_200_300 + getBin(higgs.pt(), {200, 300, 450, 650}));
737  if (Njets == 0)
738  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
739  if (Njets == 1)
740  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
741  if (Njets > 1) {
742  //VBF topology
743  if (vbfTopo)
745  //Njets >= 2jets without VBF topology (mjj<350)
746  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
747  }
748  }
749 
750  // 2. Electroweak qq->Hqq Stage 1.2 categories
751  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
752  if (std::abs(higgs.rapidity()) > 2.5)
753  return QQ2HQQ_FWDH;
754  int Njets = jets.size();
755  if (Njets == 0)
756  return QQ2HQQ_0J;
757  else if (Njets == 1)
758  return QQ2HQQ_1J;
759  else if (Njets >= 2) {
760  double mjj = (jets[0].mom() + jets[1].mom()).mass();
761  if (mjj < 60)
762  return QQ2HQQ_GE2J_MJJ_0_60;
763  else if (60 < mjj && mjj < 120)
764  return QQ2HQQ_GE2J_MJJ_60_120;
765  else if (120 < mjj && mjj < 350)
767  else if (mjj > 350) {
768  if (higgs.pt() > 200)
770  if (vbfTopo)
771  return Category(QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200 + vbfTopo);
772  }
773  }
774  }
775  // 3. WH->Hlv categories
776  else if (prodMode == HTXS::WH) {
777  if (fwdHiggs)
778  return QQ2HLNU_FWDH;
779  else if (V.pt() < 75)
780  return QQ2HLNU_PTV_0_75;
781  else if (V.pt() < 150)
782  return QQ2HLNU_PTV_75_150;
783  else if (V.pt() > 250)
784  return QQ2HLNU_PTV_GT250;
785  // 150 < pTV/GeV < 250
787  }
788  // 4. qq->ZH->llH categories
789  else if (prodMode == HTXS::QQ2ZH) {
790  if (fwdHiggs)
791  return QQ2HLL_FWDH;
792  else if (V.pt() < 75)
793  return QQ2HLL_PTV_0_75;
794  else if (V.pt() < 150)
795  return QQ2HLL_PTV_75_150;
796  else if (V.pt() > 250)
797  return QQ2HLL_PTV_GT250;
798  // 150 < pTV/GeV < 250
800  }
801  // 5. gg->ZH->llH categories
802  else if (prodMode == HTXS::GG2ZH) {
803  if (fwdHiggs)
804  return GG2HLL_FWDH;
805  else if (V.pt() < 75)
806  return GG2HLL_PTV_0_75;
807  else if (V.pt() < 150)
808  return GG2HLL_PTV_75_150;
809  else if (V.pt() > 250)
810  return GG2HLL_PTV_GT250;
812  }
813  // 6.ttH,bbH,tH categories
814  else if (prodMode == HTXS::TTH) {
815  if (fwdHiggs)
816  return TTH_FWDH;
817  else
818  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300}));
819  } else if (prodMode == HTXS::BBH)
820  return Category(BBH_FWDH + ctrlHiggs);
821  else if (prodMode == HTXS::TH)
822  return Category(TH_FWDH + ctrlHiggs);
823  return UNKNOWN;
824  }
825 
828  const Particle &higgs,
829  const Jets &jets,
830  const Particle &V) {
831  using namespace HTXS::Stage1_2_Fine;
832  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
833  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
834 
835  // 1. GGF Stage 1.2 categories
836  // Following YR4 write-up: XXXXX
837  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
838  if (fwdHiggs)
839  return GG2H_FWDH;
840  if (higgs.pt() > 200) {
841  if (Njets > 0) {
842  double pTHj = (jets[0].momentum() + higgs.momentum()).pt();
843  if (pTHj / higgs.pt() > 0.15)
844  return Category(GG2H_PTH_200_300_PTHJoverPTH_GT15 + getBin(higgs.pt(), {200, 300, 450, 650}));
845  else
846  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
847  } else
848  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
849  }
850  if (Njets == 0)
851  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
852  if (Njets == 1)
853  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
854  if (Njets > 1) {
855  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
856  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
857  //VBF topology
858  if (vbfTopo)
860  //Njets >= 2jets without VBF topology (mjj<350)
861  if (pTHjj < 25)
862  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
863  else
864  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
865  }
866  }
867 
868  // 2. Electroweak qq->Hqq Stage 1.2 categories
869  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
870  if (std::abs(higgs.rapidity()) > 2.5)
871  return QQ2HQQ_FWDH;
872  int Njets = jets.size();
873  if (Njets == 0)
874  return QQ2HQQ_0J;
875  else if (Njets == 1)
876  return QQ2HQQ_1J;
877  else if (Njets >= 2) {
878  double mjj = (jets[0].mom() + jets[1].mom()).mass();
879  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
880  if (mjj < 350) {
881  if (pTHjj < 25)
882  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
883  else
884  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
885  } else { //mjj>350 GeV
886  if (higgs.pt() < 200) {
888  } else {
890  }
891  }
892  }
893  }
894  // 3. WH->Hlv categories
895  else if (prodMode == HTXS::WH) {
896  if (fwdHiggs)
897  return QQ2HLNU_FWDH;
898  int Njets = jets.size();
899  if (Njets == 0)
900  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
901  if (Njets == 1)
902  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
903  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
904  }
905  // 4. qq->ZH->llH categories
906  else if (prodMode == HTXS::QQ2ZH) {
907  if (fwdHiggs)
908  return QQ2HLL_FWDH;
909  int Njets = jets.size();
910  if (Njets == 0)
911  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
912  if (Njets == 1)
913  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
914  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
915  }
916  // 5. gg->ZH->llH categories
917  else if (prodMode == HTXS::GG2ZH) {
918  if (fwdHiggs)
919  return GG2HLL_FWDH;
920  int Njets = jets.size();
921  if (Njets == 0)
922  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
923  if (Njets == 1)
924  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
925  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
926  }
927  // 6.ttH,bbH,tH categories
928  else if (prodMode == HTXS::TTH) {
929  if (fwdHiggs)
930  return TTH_FWDH;
931  else
932  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300, 450}));
933  } else if (prodMode == HTXS::BBH)
934  return Category(BBH_FWDH + ctrlHiggs);
935  else if (prodMode == HTXS::TH)
936  return Category(TH_FWDH + ctrlHiggs);
937  return UNKNOWN;
938  }
939 
941 
944 
946  void setHiggsProdMode(HTXS::HiggsProdMode prodMode) { m_HiggsProdMode = prodMode; }
947 
951  void init() override {
952  printf("==============================================================\n");
953  printf("======== HiggsTemplateCrossSections Initialization =========\n");
954  printf("==============================================================\n");
955  // check that the production mode has been set
956  // if running in standalone Rivet the production mode is set through an env variable
958  char *pm_env = std::getenv("HIGGSPRODMODE");
959  string pm(pm_env == nullptr ? "" : pm_env);
960  if (pm == "GGF")
962  else if (pm == "VBF")
964  else if (pm == "WH")
966  else if (pm == "ZH")
968  else if (pm == "QQ2ZH")
970  else if (pm == "GG2ZH")
972  else if (pm == "TTH")
974  else if (pm == "BBH")
976  else if (pm == "TH")
978  else {
979  MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
980  }
981  }
982 
983  // Projections for final state particles
984  const FinalState FS;
985  declare(FS, "FS");
986 
987  // initialize the histograms with for each of the stages
989  m_sumw = 0.0;
990  m_sumevents = 0;
991  printf("==============================================================\n");
992  printf("======== Higgs prod mode %d =========\n", m_HiggsProdMode);
993  printf("======== Sucessful Initialization =========\n");
994  printf("==============================================================\n");
995  }
996 
997  // Perform the per-event analysis
998  void analyze(const Event &event) override {
999  // get the classification
1000  HiggsClassification cat = classifyEvent(event, m_HiggsProdMode);
1001 
1002  // Fill histograms: categorization --> linerize the categories
1003  const double weight = 1.0;
1004  m_sumw += weight;
1005 
1006  int F = cat.stage0_cat % 10, P = cat.stage1_cat_pTjet30GeV / 100;
1007  hist_stage0->fill(cat.stage0_cat / 10 * 2 + F, weight);
1008 
1009  // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
1010  vector<int> offset({0, 1, 13, 19, 24, 29, 33, 35, 37, 39});
1011  int off = offset[P];
1012  hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV % 100 + off, weight);
1013  hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV % 100 + off, weight);
1014 
1015  // 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
1016  static vector<int> offset1_2({0, 1, 18, 29, 35, 41, 47, 53, 55, 57});
1017  int off1_2 = offset1_2[P];
1018  // 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
1019  static vector<int> offset1_2f({0, 1, 29, 54, 70, 86, 102, 109, 111, 113});
1020  int off1_2f = offset1_2f[P];
1021  hist_stage1_2_pTjet25->fill(cat.stage1_2_cat_pTjet25GeV % 100 + off1_2, weight);
1022  hist_stage1_2_pTjet30->fill(cat.stage1_2_cat_pTjet30GeV % 100 + off1_2, weight);
1023  hist_stage1_2_fine_pTjet25->fill(cat.stage1_2_fine_cat_pTjet25GeV % 100 + off1_2f, weight);
1024  hist_stage1_2_fine_pTjet30->fill(cat.stage1_2_fine_cat_pTjet30GeV % 100 + off1_2f, weight);
1025 
1026  // Fill histograms: variables used in the categorization
1027  hist_pT_Higgs->fill(cat.higgs.pT(), weight);
1028  hist_y_Higgs->fill(cat.higgs.rapidity(), weight);
1029  hist_pT_V->fill(cat.V.pT(), weight);
1030 
1031  hist_Njets25->fill(cat.jets25.size(), weight);
1032  hist_Njets30->fill(cat.jets30.size(), weight);
1033 
1034  hist_isZ2vv->fill(cat.isZ2vvDecay, weight);
1035 
1036  // Jet variables. Use jet collection with pT threshold at 30 GeV
1037  if (!cat.jets30.empty())
1038  hist_pT_jet1->fill(cat.jets30[0].pt(), weight);
1039  if (cat.jets30.size() >= 2) {
1040  const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
1041  hist_deltay_jj->fill(std::abs(j1.rapidity() - j2.rapidity()), weight);
1042  hist_dijet_mass->fill((j1 + j2).mass(), weight);
1043  hist_pT_Hjj->fill((j1 + j2 + cat.higgs.momentum()).pt(), weight);
1044  }
1045  }
1046 
1048  MSG_INFO(" ====================================================== ");
1049  MSG_INFO(" Higgs Template X-Sec Categorization Tool ");
1050  MSG_INFO(" Status Code Summary ");
1051  MSG_INFO(" ====================================================== ");
1052  bool allSuccess = (m_sumevents == m_errorCount[HTXS::SUCCESS]);
1053  if (allSuccess)
1054  MSG_INFO(" >>>> All " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized!");
1055  else {
1056  MSG_INFO(" >>>> " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized");
1057  MSG_INFO(" >>>> --> the following errors occured:");
1058  MSG_INFO(" >>>> " << m_errorCount[HTXS::PRODMODE_DEFINED] << " had an undefined Higgs production mode.");
1059  MSG_INFO(" >>>> " << m_errorCount[HTXS::MOMENTUM_CONSERVATION] << " failed momentum conservation.");
1060  MSG_INFO(" >>>> " << m_errorCount[HTXS::HIGGS_IDENTIFICATION]
1061  << " failed to identify a valid Higgs boson.");
1062  MSG_INFO(" >>>> " << m_errorCount[HTXS::HS_VTX_IDENTIFICATION]
1063  << " failed to identify the hard scatter vertex.");
1064  MSG_INFO(" >>>> " << m_errorCount[HTXS::VH_IDENTIFICATION] << " VH: to identify a valid V-boson.");
1065  MSG_INFO(" >>>> " << m_errorCount[HTXS::TOP_W_IDENTIFICATION]
1066  << " failed to identify valid Ws from top decay.");
1067  }
1068  MSG_INFO(" ====================================================== ");
1069  MSG_INFO(" ====================================================== ");
1070  }
1071 
1072  void finalize() override {
1074  double sf = m_sumw > 0 ? 1.0 / m_sumw : 1.0;
1075  scale({hist_stage0,
1082  hist_Njets25,
1083  hist_Njets30,
1084  hist_pT_Higgs,
1085  hist_y_Higgs,
1086  hist_pT_V,
1087  hist_pT_jet1,
1090  hist_pT_Hjj},
1091  sf);
1092  }
1093 
1094  /*
1095  * initialize histograms
1096  */
1097 
1099  book(hist_stage0, "HTXS_stage0", 20, 0, 20);
1100  book(hist_stage1_pTjet25, "HTXS_stage1_pTjet25", 40, 0, 40);
1101  book(hist_stage1_pTjet30, "HTXS_stage1_pTjet30", 40, 0, 40);
1102  book(hist_stage1_2_pTjet25, "HTXS_stage1_2_pTjet25", 57, 0, 57);
1103  book(hist_stage1_2_pTjet30, "HTXS_stage1_2_pTjet30", 57, 0, 57);
1104  book(hist_stage1_2_fine_pTjet25, "HTXS_stage1_2_fine_pTjet25", 113, 0, 113);
1105  book(hist_stage1_2_fine_pTjet30, "HTXS_stage1_2_fine_pTjet30", 113, 0, 113);
1106  book(hist_pT_Higgs, "pT_Higgs", 80, 0, 400);
1107  book(hist_y_Higgs, "y_Higgs", 80, -4, 4);
1108  book(hist_pT_V, "pT_V", 80, 0, 400);
1109  book(hist_pT_jet1, "pT_jet1", 80, 0, 400);
1110  book(hist_deltay_jj, "deltay_jj", 50, 0, 10);
1111  book(hist_dijet_mass, "m_jj", 50, 0, 2000);
1112  book(hist_pT_Hjj, "pT_Hjj", 50, 0, 250);
1113  book(hist_Njets25, "Njets25", 10, 0, 10);
1114  book(hist_Njets30, "Njets30", 10, 0, 10);
1115  book(hist_isZ2vv, "isZ2vv", 2, 0, 2);
1116  }
1118 
1119  /*
1120  * initialize private members used in the classification procedure
1121  */
1122 
1123  private:
1124  double m_sumw;
1125  size_t m_sumevents;
1127  std::map<HTXS::ErrorCode, size_t> m_errorCount;
1128  Histo1DPtr hist_stage0;
1136  Histo1DPtr hist_isZ2vv;
1137  };
1138 
1139  // the PLUGIN only needs to be decleared when running standalone Rivet
1140  // and causes compilation / linking issues if included in Athena / RootCore
1141  //check for Rivet environment variable RIVET_ANALYSIS_PATH
1142 #ifdef RIVET_ANALYSIS_PATH
1143  // The hook for the plugin system
1144  DECLARE_RIVET_PLUGIN(HiggsTemplateCrossSections);
1145 #endif
1146 
1147 } // namespace Rivet
1148 
1149 #endif
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.
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
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...
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
bool isChargedLepton(const HepMC::GenParticle *part)
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.
std::map< HTXS::ErrorCode, size_t > m_errorCount
void init() override
default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Hig...
static std::string to_string(const XMLCh *ch)
void analyze(const Event &event) override
def cat(path)
Definition: eostools.py:401
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.
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t V
bool quarkDecay(const Particle &p)
Return true is particle decays to quarks.
failed momentum conservation
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.
bool hasChild(const ConstGenParticlePtr &ptcl, int pdgID)
Checks whether the input particle has a child with a given PDGID.
failed to identify associated vector boson decay products
void setHiggsProdMode(HTXS::HiggsProdMode prodMode)
Sets the Higgs production mode.
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:286
HTXS::Stage0::Category getStage0Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Particle &V)
Stage-0 HTXS categorization.
bool hasParent(const ConstGenParticlePtr &ptcl, int pdgID)
Checks whether the input particle has a parent with a given PDGID.
float x
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:163
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