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  // Temporary fix: add variables to perform STXS 1.3 classification with nanoAOD on-the-fly
306  // Vector-boson pt for VH production modes
307  if (isVH(prodMode)) {
308  cat.V_pt = cat.V.pt();
309  } else {
310  cat.V_pt = -999;
311  }
312  // Dijet variables using jets30 collection
313  if (cat.jets30.size() >= 2) {
314  cat.Mjj = (cat.jets30[0].mom() + cat.jets30[1].mom()).mass();
315  cat.ptHjj = (cat.jets30[0].mom() + cat.jets30[1].mom() + cat.higgs.momentum()).pt();
316  cat.dPhijj = cat.jets30[0].mom().phi() - cat.jets30[1].mom().phi();
317  // Return phi angle in the interval [-PI,PI)
318  if (cat.dPhijj >= Rivet::pi)
319  cat.dPhijj -= 2 * Rivet::pi;
320  else if (cat.dPhijj < -1 * Rivet::pi)
321  cat.dPhijj += 2 * Rivet::pi;
322  } else {
323  cat.Mjj = -999;
324  cat.ptHjj = -999;
325  cat.dPhijj = -999;
326  }
327 
328  // check that four mometum sum of all stable particles satisfies momentum consevation
329  /*
330  if ( sum.pt()>0.1 )
331  return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
332  std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
333 */
334  // check if V-boson was not included in the event record but decay particles were
335  // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
336  if (is_uncatdV)
337  return error(cat, HTXS::VH_IDENTIFICATION, "Failed to identify associated V-boson!");
338 
339  /*****
340  * Step 4.
341  * Classify and save output
342  */
343 
344  // Apply the categorization categorization
345  cat.isZ2vvDecay = false;
346  if ((prodMode == HTXS::GG2ZH || prodMode == HTXS::QQ2ZH) && !quarkDecay(cat.V) && !ChLeptonDecay(cat.V))
347  cat.isZ2vvDecay = true;
348  cat.stage0_cat = getStage0Category(prodMode, cat.higgs, cat.V);
349  cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode, cat.higgs, cat.jets25, cat.V);
350  cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode, cat.higgs, cat.jets30, cat.V);
351  cat.stage1_1_cat_pTjet25GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets25, cat.V);
352  cat.stage1_1_cat_pTjet30GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets30, cat.V);
353  cat.stage1_1_fine_cat_pTjet25GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
354  cat.stage1_1_fine_cat_pTjet30GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
355  cat.stage1_2_cat_pTjet25GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets25, cat.V);
356  cat.stage1_2_cat_pTjet30GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets30, cat.V);
357  cat.stage1_2_fine_cat_pTjet25GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
358  cat.stage1_2_fine_cat_pTjet30GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
359 
360  cat.errorCode = HTXS::SUCCESS;
362  ++m_sumevents;
363 
364  return cat;
365  }
366 
372 
374  int getBin(double x, std::vector<double> bins) {
375  for (size_t i = 1; i < bins.size(); ++i)
376  if (x < bins[i])
377  return i - 1;
378  return bins.size() - 1;
379  }
380 
384  int vbfTopology(const Jets &jets, const Particle &higgs) {
385  if (jets.size() < 2)
386  return 0;
387  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
388  bool VBFtopo = (j1 + j2).mass() > 400.0 && std::abs(j1.rapidity() - j2.rapidity()) > 2.8;
389  return VBFtopo ? (j1 + j2 + higgs.momentum()).pt() < 25 ? 2 : 1 : 0;
390  }
391 
396  int vbfTopology_Stage1_X(const Jets &jets, const Particle &higgs) {
397  if (jets.size() < 2)
398  return 0;
399  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
400  double mjj = (j1 + j2).mass();
401  if (mjj > 350 && mjj <= 700)
402  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
403  else if (mjj > 700)
404  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
405  else
406  return 0;
407  }
408 
415  int vbfTopology_Stage1_X_Fine(const Jets &jets, const Particle &higgs) {
416  if (jets.size() < 2)
417  return 0;
418  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
419  double mjj = (j1 + j2).mass();
420  if (mjj > 350 && mjj <= 700)
421  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
422  else if (mjj > 700 && mjj <= 1000)
423  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
424  else if (mjj > 1000 && mjj <= 1500)
425  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 5 : 6;
426  else if (mjj > 1500)
427  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 7 : 8;
428  else
429  return 0;
430  }
431 
433  bool isVH(HTXS::HiggsProdMode p) { return p == HTXS::WH || p == HTXS::QQ2ZH || p == HTXS::GG2ZH; }
434 
437  const Particle &higgs,
438  const Particle &V) {
439  using namespace HTXS::Stage0;
440  int ctrlHiggs = std::abs(higgs.rapidity()) < 2.5;
441  // Special cases first, qq→Hqq
442  if ((prodMode == HTXS::WH || prodMode == HTXS::QQ2ZH) && quarkDecay(V)) {
443  return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
444  } else if (prodMode == HTXS::GG2ZH && quarkDecay(V)) {
445  return Category(HTXS::GGF * 10 + ctrlHiggs);
446  }
447  // General case after
448  return Category(prodMode * 10 + ctrlHiggs);
449  }
450 
453  const Particle &higgs,
454  const Jets &jets,
455  const Particle &V) {
456  using namespace HTXS::Stage1;
457  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
458  double pTj1 = !jets.empty() ? jets[0].momentum().pt() : 0;
459  int vbfTopo = vbfTopology(jets, higgs);
460 
461  // 1. GGF Stage 1 categories
462  // Following YR4 write-up: XXXXX
463  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
464  if (fwdHiggs)
465  return GG2H_FWDH;
466  if (Njets == 0)
467  return GG2H_0J;
468  else if (Njets == 1)
469  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
470  else if (Njets >= 2) {
471  // events with pT_H>200 get priority over VBF cuts
472  if (higgs.pt() <= 200) {
473  if (vbfTopo == 2)
474  return GG2H_VBFTOPO_JET3VETO;
475  else if (vbfTopo == 1)
476  return GG2H_VBFTOPO_JET3;
477  }
478  // Njets >= 2jets without VBF topology
479  return Category(GG2H_GE2J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
480  }
481  }
482  // 2. Electroweak qq->Hqq Stage 1 categories
483  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
484  if (std::abs(higgs.rapidity()) > 2.5)
485  return QQ2HQQ_FWDH;
486  if (pTj1 > 200)
487  return QQ2HQQ_PTJET1_GT200;
488  if (vbfTopo == 2)
490  if (vbfTopo == 1)
491  return QQ2HQQ_VBFTOPO_JET3;
492  double mjj = jets.size() > 1 ? (jets[0].mom() + jets[1].mom()).mass() : 0;
493  if (60 < mjj && mjj < 120)
494  return QQ2HQQ_VH2JET;
495  return QQ2HQQ_REST;
496  }
497  // 3. WH->Hlv categories
498  else if (prodMode == HTXS::WH) {
499  if (fwdHiggs)
500  return QQ2HLNU_FWDH;
501  else if (V.pt() < 150)
502  return QQ2HLNU_PTV_0_150;
503  else if (V.pt() > 250)
504  return QQ2HLNU_PTV_GT250;
505  // 150 < pTV/GeV < 250
507  }
508  // 4. qq->ZH->llH categories
509  else if (prodMode == HTXS::QQ2ZH) {
510  if (fwdHiggs)
511  return QQ2HLL_FWDH;
512  else if (V.pt() < 150)
513  return QQ2HLL_PTV_0_150;
514  else if (V.pt() > 250)
515  return QQ2HLL_PTV_GT250;
516  // 150 < pTV/GeV < 250
518  }
519  // 5. gg->ZH->llH categories
520  else if (prodMode == HTXS::GG2ZH) {
521  if (fwdHiggs)
522  return GG2HLL_FWDH;
523  if (V.pt() < 150)
524  return GG2HLL_PTV_0_150;
525  else if (jets.empty())
526  return GG2HLL_PTV_GT150_0J;
527  return GG2HLL_PTV_GT150_GE1J;
528  }
529  // 6.ttH,bbH,tH categories
530  else if (prodMode == HTXS::TTH)
531  return Category(TTH_FWDH + ctrlHiggs);
532  else if (prodMode == HTXS::BBH)
533  return Category(BBH_FWDH + ctrlHiggs);
534  else if (prodMode == HTXS::TH)
535  return Category(TH_FWDH + ctrlHiggs);
536  return UNKNOWN;
537  }
538 
541  const Particle &higgs,
542  const Jets &jets,
543  const Particle &V) {
544  using namespace HTXS::Stage1_1;
545  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
546  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
547 
548  // 1. GGF Stage 1 categories
549  // Following YR4 write-up: XXXXX
550  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
551  if (fwdHiggs)
552  return GG2H_FWDH;
553  if (higgs.pt() > 200)
554  return GG2H_PTH_GT200;
555  if (Njets == 0)
556  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
557  if (Njets == 1)
558  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
559  if (Njets > 1) {
560  //VBF topology
561  if (vbfTopo)
562  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
563  //Njets >= 2jets without VBF topology (mjj<350)
564  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
565  }
566  }
567 
568  // 2. Electroweak qq->Hqq Stage 1.1 categories
569  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
570  if (std::abs(higgs.rapidity()) > 2.5)
571  return QQ2HQQ_FWDH;
572  int Njets = jets.size();
573  if (Njets == 0)
574  return QQ2HQQ_0J;
575  else if (Njets == 1)
576  return QQ2HQQ_1J;
577  else if (Njets >= 2) {
578  double mjj = (jets[0].mom() + jets[1].mom()).mass();
579  if (mjj < 60)
580  return QQ2HQQ_MJJ_0_60;
581  else if (60 < mjj && mjj < 120)
582  return QQ2HQQ_MJJ_60_120;
583  else if (120 < mjj && mjj < 350)
584  return QQ2HQQ_MJJ_120_350;
585  else if (mjj > 350) {
586  if (higgs.pt() > 200)
588  if (vbfTopo)
589  return Category(QQ2HQQ_MJJ_GT350_PTH_GT200 + vbfTopo);
590  }
591  }
592  }
593  // 3. WH->Hlv categories
594  else if (prodMode == HTXS::WH) {
595  if (fwdHiggs)
596  return QQ2HLNU_FWDH;
597  else if (V.pt() < 75)
598  return QQ2HLNU_PTV_0_75;
599  else if (V.pt() < 150)
600  return QQ2HLNU_PTV_75_150;
601  else if (V.pt() > 250)
602  return QQ2HLNU_PTV_GT250;
603  // 150 < pTV/GeV < 250
605  }
606  // 4. qq->ZH->llH categories
607  else if (prodMode == HTXS::QQ2ZH) {
608  if (fwdHiggs)
609  return QQ2HLL_FWDH;
610  else if (V.pt() < 75)
611  return QQ2HLL_PTV_0_75;
612  else if (V.pt() < 150)
613  return QQ2HLL_PTV_75_150;
614  else if (V.pt() > 250)
615  return QQ2HLL_PTV_GT250;
616  // 150 < pTV/GeV < 250
618  }
619  // 5. gg->ZH->llH categories
620  else if (prodMode == HTXS::GG2ZH) {
621  if (fwdHiggs)
622  return GG2HLL_FWDH;
623  else if (V.pt() < 75)
624  return GG2HLL_PTV_0_75;
625  else if (V.pt() < 150)
626  return GG2HLL_PTV_75_150;
627  else if (V.pt() > 250)
628  return GG2HLL_PTV_GT250;
630  }
631  // 6.ttH,bbH,tH categories
632  else if (prodMode == HTXS::TTH)
633  return Category(TTH_FWDH + ctrlHiggs);
634  else if (prodMode == HTXS::BBH)
635  return Category(BBH_FWDH + ctrlHiggs);
636  else if (prodMode == HTXS::TH)
637  return Category(TH_FWDH + ctrlHiggs);
638  return UNKNOWN;
639  }
640 
643  const Particle &higgs,
644  const Jets &jets,
645  const Particle &V) {
646  using namespace HTXS::Stage1_1_Fine;
647  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
648  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
649 
650  // 1. GGF Stage 1.1 categories
651  // Following YR4 write-up: XXXXX
652  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
653  if (fwdHiggs)
654  return GG2H_FWDH;
655  if (higgs.pt() > 200)
656  return GG2H_PTH_GT200;
657  if (Njets == 0)
658  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
659  if (Njets == 1)
660  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
661  if (Njets > 1) {
662  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
663  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
664  //VBF topology
665  if (vbfTopo)
666  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
667  //Njets >= 2jets without VBF topology (mjj<350)
668  if (pTHjj < 25)
669  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
670  else
671  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
672  }
673  }
674 
675  // 2. Electroweak qq->Hqq Stage 1.1 categories
676  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
677  if (std::abs(higgs.rapidity()) > 2.5)
678  return QQ2HQQ_FWDH;
679  int Njets = jets.size();
680  if (Njets == 0)
681  return QQ2HQQ_0J;
682  else if (Njets == 1)
683  return QQ2HQQ_1J;
684  else if (Njets >= 2) {
685  double mjj = (jets[0].mom() + jets[1].mom()).mass();
686  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
687  if (mjj < 350) {
688  if (pTHjj < 25)
689  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
690  else
691  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
692  } else { //mjj>350 GeV
693  if (higgs.pt() < 200) {
694  return Category(QQ2HQQ_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
695  } else {
697  }
698  }
699  }
700  }
701  // 3. WH->Hlv categories
702  else if (prodMode == HTXS::WH) {
703  if (fwdHiggs)
704  return QQ2HLNU_FWDH;
705  int Njets = jets.size();
706  if (Njets == 0)
707  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
708  if (Njets == 1)
709  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
710  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
711  }
712  // 4. qq->ZH->llH categories
713  else if (prodMode == HTXS::QQ2ZH) {
714  if (fwdHiggs)
715  return QQ2HLL_FWDH;
716  int Njets = jets.size();
717  if (Njets == 0)
718  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
719  if (Njets == 1)
720  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
721  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
722  }
723  // 5. gg->ZH->llH categories
724  else if (prodMode == HTXS::GG2ZH) {
725  if (fwdHiggs)
726  return GG2HLL_FWDH;
727  int Njets = jets.size();
728  if (Njets == 0)
729  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
730  if (Njets == 1)
731  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
732  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
733  }
734  // 6.ttH,bbH,tH categories
735  else if (prodMode == HTXS::TTH)
736  return Category(TTH_FWDH + ctrlHiggs);
737  else if (prodMode == HTXS::BBH)
738  return Category(BBH_FWDH + ctrlHiggs);
739  else if (prodMode == HTXS::TH)
740  return Category(TH_FWDH + ctrlHiggs);
741  return UNKNOWN;
742  }
743 
746  const Particle &higgs,
747  const Jets &jets,
748  const Particle &V) {
749  using namespace HTXS::Stage1_2;
750  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
751  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
752 
753  // 1. GGF Stage 1 categories
754  // Following YR4 write-up: XXXXX
755  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
756  if (fwdHiggs)
757  return GG2H_FWDH;
758  if (higgs.pt() > 200)
759  return Category(GG2H_PTH_200_300 + getBin(higgs.pt(), {200, 300, 450, 650}));
760  if (Njets == 0)
761  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
762  if (Njets == 1)
763  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
764  if (Njets > 1) {
765  //VBF topology
766  if (vbfTopo)
768  //Njets >= 2jets without VBF topology (mjj<350)
769  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
770  }
771  }
772 
773  // 2. Electroweak qq->Hqq Stage 1.2 categories
774  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
775  if (std::abs(higgs.rapidity()) > 2.5)
776  return QQ2HQQ_FWDH;
777  int Njets = jets.size();
778  if (Njets == 0)
779  return QQ2HQQ_0J;
780  else if (Njets == 1)
781  return QQ2HQQ_1J;
782  else if (Njets >= 2) {
783  double mjj = (jets[0].mom() + jets[1].mom()).mass();
784  if (mjj < 60)
785  return QQ2HQQ_GE2J_MJJ_0_60;
786  else if (60 < mjj && mjj < 120)
787  return QQ2HQQ_GE2J_MJJ_60_120;
788  else if (120 < mjj && mjj < 350)
790  else if (mjj > 350) {
791  if (higgs.pt() > 200)
793  if (vbfTopo)
794  return Category(QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200 + vbfTopo);
795  }
796  }
797  }
798  // 3. WH->Hlv categories
799  else if (prodMode == HTXS::WH) {
800  if (fwdHiggs)
801  return QQ2HLNU_FWDH;
802  else if (V.pt() < 75)
803  return QQ2HLNU_PTV_0_75;
804  else if (V.pt() < 150)
805  return QQ2HLNU_PTV_75_150;
806  else if (V.pt() > 250)
807  return QQ2HLNU_PTV_GT250;
808  // 150 < pTV/GeV < 250
810  }
811  // 4. qq->ZH->llH categories
812  else if (prodMode == HTXS::QQ2ZH) {
813  if (fwdHiggs)
814  return QQ2HLL_FWDH;
815  else if (V.pt() < 75)
816  return QQ2HLL_PTV_0_75;
817  else if (V.pt() < 150)
818  return QQ2HLL_PTV_75_150;
819  else if (V.pt() > 250)
820  return QQ2HLL_PTV_GT250;
821  // 150 < pTV/GeV < 250
823  }
824  // 5. gg->ZH->llH categories
825  else if (prodMode == HTXS::GG2ZH) {
826  if (fwdHiggs)
827  return GG2HLL_FWDH;
828  else if (V.pt() < 75)
829  return GG2HLL_PTV_0_75;
830  else if (V.pt() < 150)
831  return GG2HLL_PTV_75_150;
832  else if (V.pt() > 250)
833  return GG2HLL_PTV_GT250;
835  }
836  // 6.ttH,bbH,tH categories
837  else if (prodMode == HTXS::TTH) {
838  if (fwdHiggs)
839  return TTH_FWDH;
840  else
841  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300}));
842  } else if (prodMode == HTXS::BBH)
843  return Category(BBH_FWDH + ctrlHiggs);
844  else if (prodMode == HTXS::TH)
845  return Category(TH_FWDH + ctrlHiggs);
846  return UNKNOWN;
847  }
848 
851  const Particle &higgs,
852  const Jets &jets,
853  const Particle &V) {
854  using namespace HTXS::Stage1_2_Fine;
855  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
856  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
857 
858  // 1. GGF Stage 1.2 categories
859  // Following YR4 write-up: XXXXX
860  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
861  if (fwdHiggs)
862  return GG2H_FWDH;
863  if (higgs.pt() > 200) {
864  if (Njets > 0) {
865  double pTHj = (jets[0].momentum() + higgs.momentum()).pt();
866  if (pTHj / higgs.pt() > 0.15)
867  return Category(GG2H_PTH_200_300_PTHJoverPTH_GT15 + getBin(higgs.pt(), {200, 300, 450, 650}));
868  else
869  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
870  } else
871  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
872  }
873  if (Njets == 0)
874  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
875  if (Njets == 1)
876  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
877  if (Njets > 1) {
878  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
879  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
880  //VBF topology
881  if (vbfTopo)
883  //Njets >= 2jets without VBF topology (mjj<350)
884  if (pTHjj < 25)
885  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
886  else
887  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
888  }
889  }
890 
891  // 2. Electroweak qq->Hqq Stage 1.2 categories
892  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
893  if (std::abs(higgs.rapidity()) > 2.5)
894  return QQ2HQQ_FWDH;
895  int Njets = jets.size();
896  if (Njets == 0)
897  return QQ2HQQ_0J;
898  else if (Njets == 1)
899  return QQ2HQQ_1J;
900  else if (Njets >= 2) {
901  double mjj = (jets[0].mom() + jets[1].mom()).mass();
902  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
903  if (mjj < 350) {
904  if (pTHjj < 25)
905  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
906  else
907  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
908  } else { //mjj>350 GeV
909  if (higgs.pt() < 200) {
911  } else {
913  }
914  }
915  }
916  }
917  // 3. WH->Hlv categories
918  else if (prodMode == HTXS::WH) {
919  if (fwdHiggs)
920  return QQ2HLNU_FWDH;
921  int Njets = jets.size();
922  if (Njets == 0)
923  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
924  if (Njets == 1)
925  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
926  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
927  }
928  // 4. qq->ZH->llH categories
929  else if (prodMode == HTXS::QQ2ZH) {
930  if (fwdHiggs)
931  return QQ2HLL_FWDH;
932  int Njets = jets.size();
933  if (Njets == 0)
934  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
935  if (Njets == 1)
936  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
937  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
938  }
939  // 5. gg->ZH->llH categories
940  else if (prodMode == HTXS::GG2ZH) {
941  if (fwdHiggs)
942  return GG2HLL_FWDH;
943  int Njets = jets.size();
944  if (Njets == 0)
945  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
946  if (Njets == 1)
947  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
948  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
949  }
950  // 6.ttH,bbH,tH categories
951  else if (prodMode == HTXS::TTH) {
952  if (fwdHiggs)
953  return TTH_FWDH;
954  else
955  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300, 450}));
956  } else if (prodMode == HTXS::BBH)
957  return Category(BBH_FWDH + ctrlHiggs);
958  else if (prodMode == HTXS::TH)
959  return Category(TH_FWDH + ctrlHiggs);
960  return UNKNOWN;
961  }
962 
964 
967 
969  void setHiggsProdMode(HTXS::HiggsProdMode prodMode) { m_HiggsProdMode = prodMode; }
970 
974  void init() override {
975  printf("==============================================================\n");
976  printf("======== HiggsTemplateCrossSections Initialization =========\n");
977  printf("==============================================================\n");
978  // check that the production mode has been set
979  // if running in standalone Rivet the production mode is set through an env variable
981  char *pm_env = std::getenv("HIGGSPRODMODE");
982  string pm(pm_env == nullptr ? "" : pm_env);
983  if (pm == "GGF")
985  else if (pm == "VBF")
987  else if (pm == "WH")
989  else if (pm == "ZH")
991  else if (pm == "QQ2ZH")
993  else if (pm == "GG2ZH")
995  else if (pm == "TTH")
997  else if (pm == "BBH")
999  else if (pm == "TH")
1001  else {
1002  MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
1003  }
1004  }
1005 
1006  // Projections for final state particles
1007  const FinalState FS;
1008  declare(FS, "FS");
1009 
1010  // initialize the histograms with for each of the stages
1011  initializeHistos();
1012  m_sumw = 0.0;
1013  m_sumevents = 0;
1014  printf("==============================================================\n");
1015  printf("======== Higgs prod mode %d =========\n", m_HiggsProdMode);
1016  printf("======== Sucessful Initialization =========\n");
1017  printf("==============================================================\n");
1018  }
1019 
1020  // Perform the per-event analysis
1021  void analyze(const Event &event) override {
1022  // get the classification
1023  HiggsClassification cat = classifyEvent(event, m_HiggsProdMode);
1024 
1025  // Fill histograms: categorization --> linerize the categories
1026  const double weight = 1.0;
1027  m_sumw += weight;
1028 
1029  int F = cat.stage0_cat % 10, P = cat.stage1_cat_pTjet30GeV / 100;
1030  hist_stage0->fill(cat.stage0_cat / 10 * 2 + F, weight);
1031 
1032  // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
1033  vector<int> offset({0, 1, 13, 19, 24, 29, 33, 35, 37, 39});
1034  int off = offset[P];
1035  hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV % 100 + off, weight);
1036  hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV % 100 + off, weight);
1037 
1038  // 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
1039  static vector<int> offset1_2({0, 1, 18, 29, 35, 41, 47, 53, 55, 57});
1040  int off1_2 = offset1_2[P];
1041  // 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
1042  static vector<int> offset1_2f({0, 1, 29, 54, 70, 86, 102, 109, 111, 113});
1043  int off1_2f = offset1_2f[P];
1044  hist_stage1_2_pTjet25->fill(cat.stage1_2_cat_pTjet25GeV % 100 + off1_2, weight);
1045  hist_stage1_2_pTjet30->fill(cat.stage1_2_cat_pTjet30GeV % 100 + off1_2, weight);
1046  hist_stage1_2_fine_pTjet25->fill(cat.stage1_2_fine_cat_pTjet25GeV % 100 + off1_2f, weight);
1047  hist_stage1_2_fine_pTjet30->fill(cat.stage1_2_fine_cat_pTjet30GeV % 100 + off1_2f, weight);
1048 
1049  // Fill histograms: variables used in the categorization
1050  hist_pT_Higgs->fill(cat.higgs.pT(), weight);
1051  hist_y_Higgs->fill(cat.higgs.rapidity(), weight);
1052  hist_pT_V->fill(cat.V.pT(), weight);
1053 
1054  hist_Njets25->fill(cat.jets25.size(), weight);
1055  hist_Njets30->fill(cat.jets30.size(), weight);
1056 
1057  hist_isZ2vv->fill(cat.isZ2vvDecay, weight);
1058 
1059  // Jet variables. Use jet collection with pT threshold at 30 GeV
1060  if (!cat.jets30.empty())
1061  hist_pT_jet1->fill(cat.jets30[0].pt(), weight);
1062  if (cat.jets30.size() >= 2) {
1063  const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
1064  hist_deltay_jj->fill(std::abs(j1.rapidity() - j2.rapidity()), weight);
1065  hist_dijet_mass->fill((j1 + j2).mass(), weight);
1066  hist_pT_Hjj->fill((j1 + j2 + cat.higgs.momentum()).pt(), weight);
1067  }
1068  }
1069 
1071  MSG_INFO(" ====================================================== ");
1072  MSG_INFO(" Higgs Template X-Sec Categorization Tool ");
1073  MSG_INFO(" Status Code Summary ");
1074  MSG_INFO(" ====================================================== ");
1075  bool allSuccess = (m_sumevents == m_errorCount[HTXS::SUCCESS]);
1076  if (allSuccess)
1077  MSG_INFO(" >>>> All " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized!");
1078  else {
1079  MSG_INFO(" >>>> " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized");
1080  MSG_INFO(" >>>> --> the following errors occured:");
1081  MSG_INFO(" >>>> " << m_errorCount[HTXS::PRODMODE_DEFINED] << " had an undefined Higgs production mode.");
1082  MSG_INFO(" >>>> " << m_errorCount[HTXS::MOMENTUM_CONSERVATION] << " failed momentum conservation.");
1083  MSG_INFO(" >>>> " << m_errorCount[HTXS::HIGGS_IDENTIFICATION]
1084  << " failed to identify a valid Higgs boson.");
1085  MSG_INFO(" >>>> " << m_errorCount[HTXS::HS_VTX_IDENTIFICATION]
1086  << " failed to identify the hard scatter vertex.");
1087  MSG_INFO(" >>>> " << m_errorCount[HTXS::VH_IDENTIFICATION] << " VH: to identify a valid V-boson.");
1088  MSG_INFO(" >>>> " << m_errorCount[HTXS::TOP_W_IDENTIFICATION]
1089  << " failed to identify valid Ws from top decay.");
1090  }
1091  MSG_INFO(" ====================================================== ");
1092  MSG_INFO(" ====================================================== ");
1093  }
1094 
1095  void finalize() override {
1097  double sf = m_sumw > 0 ? 1.0 / m_sumw : 1.0;
1098  scale({hist_stage0,
1105  hist_Njets25,
1106  hist_Njets30,
1107  hist_pT_Higgs,
1108  hist_y_Higgs,
1109  hist_pT_V,
1110  hist_pT_jet1,
1113  hist_pT_Hjj},
1114  sf);
1115  }
1116 
1117  /*
1118  * initialize histograms
1119  */
1120 
1122  book(hist_stage0, "HTXS_stage0", 20, 0, 20);
1123  book(hist_stage1_pTjet25, "HTXS_stage1_pTjet25", 40, 0, 40);
1124  book(hist_stage1_pTjet30, "HTXS_stage1_pTjet30", 40, 0, 40);
1125  book(hist_stage1_2_pTjet25, "HTXS_stage1_2_pTjet25", 57, 0, 57);
1126  book(hist_stage1_2_pTjet30, "HTXS_stage1_2_pTjet30", 57, 0, 57);
1127  book(hist_stage1_2_fine_pTjet25, "HTXS_stage1_2_fine_pTjet25", 113, 0, 113);
1128  book(hist_stage1_2_fine_pTjet30, "HTXS_stage1_2_fine_pTjet30", 113, 0, 113);
1129  book(hist_pT_Higgs, "pT_Higgs", 80, 0, 400);
1130  book(hist_y_Higgs, "y_Higgs", 80, -4, 4);
1131  book(hist_pT_V, "pT_V", 80, 0, 400);
1132  book(hist_pT_jet1, "pT_jet1", 80, 0, 400);
1133  book(hist_deltay_jj, "deltay_jj", 50, 0, 10);
1134  book(hist_dijet_mass, "m_jj", 50, 0, 2000);
1135  book(hist_pT_Hjj, "pT_Hjj", 50, 0, 250);
1136  book(hist_Njets25, "Njets25", 10, 0, 10);
1137  book(hist_Njets30, "Njets30", 10, 0, 10);
1138  book(hist_isZ2vv, "isZ2vv", 2, 0, 2);
1139  }
1141 
1142  /*
1143  * initialize private members used in the classification procedure
1144  */
1145 
1146  private:
1147  double m_sumw;
1148  size_t m_sumevents;
1150  std::map<HTXS::ErrorCode, size_t> m_errorCount;
1151  Histo1DPtr hist_stage0;
1159  Histo1DPtr hist_isZ2vv;
1160  };
1161 
1162  // the PLUGIN only needs to be decleared when running standalone Rivet
1163  // and causes compilation / linking issues if included in Athena / RootCore
1164  //check for Rivet environment variable RIVET_ANALYSIS_PATH
1165 #ifdef RIVET_ANALYSIS_PATH
1166  // The hook for the plugin system
1167  DECLARE_RIVET_PLUGIN(HiggsTemplateCrossSections);
1168 #endif
1169 
1170 } // namespace Rivet
1171 
1172 #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)
const Double_t pi
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