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, HepMC::ancestors)) {
50  for (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 (auto child : Particle(*ptcl).children())
67  if (child.pid() == pdgID)
68  return true;
69  return false;
70  }
71 
73  bool hasParent(const ConstGenParticlePtr &ptcl, int pdgID) {
74  for (auto parent : HepMCUtils::particles(ptcl->production_vertex(), HepMC::parents))
75  if (parent->pdg_id() == pdgID)
76  return true;
77  return false;
78  }
79 
81  bool quarkDecay(const Particle &p) {
82  for (auto child : p.children())
83  if (PID::isQuark(child.pid()))
84  return true;
85  return false;
86  }
87 
90  HiggsClassification error(HiggsClassification &cat,
92  std::string msg = "",
93  int NmaxWarnings = 20) {
94  // Set the error, and keep statistics
95  cat.errorCode = err;
96  ++m_errorCount[err];
97 
98  // Print warning message to the screen/log
99  static std::atomic<int> Nwarnings{0};
100  if (!msg.empty() && ++Nwarnings < NmaxWarnings)
101  MSG_WARNING(msg);
102 
103  return cat;
104  }
106 
108  HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode) {
110  m_HiggsProdMode = prodMode;
111 
112  // the classification object
113  HiggsClassification cat;
114  cat.prodMode = prodMode;
115  cat.errorCode = HTXS::UNDEFINED;
116  cat.stage0_cat = HTXS::Stage0::UNKNOWN;
117  cat.stage1_cat_pTjet25GeV = HTXS::Stage1::UNKNOWN;
118  cat.stage1_cat_pTjet30GeV = HTXS::Stage1::UNKNOWN;
119  cat.stage1_1_cat_pTjet25GeV = HTXS::Stage1_1::UNKNOWN;
120  cat.stage1_1_cat_pTjet30GeV = HTXS::Stage1_1::UNKNOWN;
121  cat.stage1_1_fine_cat_pTjet25GeV = HTXS::Stage1_1_Fine::UNKNOWN;
122  cat.stage1_1_fine_cat_pTjet30GeV = HTXS::Stage1_1_Fine::UNKNOWN;
123 
124  if (prodMode == HTXS::UNKNOWN)
125  return error(cat,
127  "Unkown Higgs production mechanism. Cannot classify event."
128  " Classification for all events will most likely fail.");
129 
130  /*****
131  * Step 1.
132  * Idenfify the Higgs boson and the hard scatter vertex
133  * There should be only one of each.
134  */
135 
136  ConstGenVertexPtr HSvtx = event.genEvent()->signal_process_vertex();
137  int Nhiggs = 0;
138  for (const ConstGenParticlePtr ptcl : HepMCUtils::particles(event.genEvent())) {
139  // a) Reject all non-Higgs particles
140  if (!PID::isHiggs(ptcl->pdg_id()))
141  continue;
142  // b) select only the final Higgs boson copy, prior to decay
143  if (ptcl->end_vertex() && !hasChild(ptcl, PID::HIGGS)) {
144  cat.higgs = Particle(ptcl);
145  ++Nhiggs;
146  }
147  // c) if HepMC::signal_proces_vertex is missing
148  // set hard-scatter vertex based on first Higgs boson
149  if (HSvtx == nullptr && ptcl->production_vertex() && !hasParent(ptcl, PID::HIGGS))
150  HSvtx = ptcl->production_vertex();
151  }
152 
153  // Make sure things are in order so far
154  if (Nhiggs != 1)
155  return error(cat,
157  "Current event has " + std::to_string(Nhiggs) + " Higgs bosons. There must be only one.");
158  if (cat.higgs.children().size() < 2)
159  return error(cat, HTXS::HIGGS_DECAY_IDENTIFICATION, "Could not identify Higgs boson decay products.");
160 
161  if (HSvtx == nullptr)
162  return error(cat, HTXS::HS_VTX_IDENTIFICATION, "Cannot find hard-scatter vertex of current event.");
163 
164  /*****
165  * Step 2.
166  * Identify associated vector bosons
167  */
168 
169  // Find associated vector bosons
170  bool is_uncatdV = false;
171  Particles uncatV_decays;
172  FourMomentum uncatV_p4(0, 0, 0, 0);
173  FourVector uncatV_v4(0, 0, 0, 0);
174  int nWs = 0, nZs = 0, nTop = 0;
175  if (isVH(prodMode)) {
176  for (auto ptcl : HepMCUtils::particles(HSvtx, HepMC::children)) {
177  if (PID::isW(ptcl->pdg_id())) {
178  ++nWs;
179  cat.V = Particle(ptcl);
180  }
181  if (PID::isZ(ptcl->pdg_id())) {
182  ++nZs;
183  cat.V = Particle(ptcl);
184  }
185  }
186  if (nWs + nZs > 0)
187  cat.V = getLastInstance(cat.V);
188  else {
189  for (auto ptcl : HepMCUtils::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;
198  cat.V = Particle(24, uncatV_p4);
199  }
200  }
201 
202  if (!is_uncatdV) {
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,
213  "Found " + std::to_string(nWs) + " W-bosons and " + std::to_string(nZs) +
214  " Z-bosons. Inconsitent with VH expectation.");
215  }
216 
217  // Find and store the W-bosons from ttH->WbWbH
218  Particles Ws;
219  if (prodMode == HTXS::TTH || prodMode == HTXS::TH) {
220  // loop over particles produced in hard-scatter vertex
221  for (auto ptcl : HepMCUtils::particles(HSvtx, HepMC::children)) {
222  if (!PID::isTop(ptcl->pdg_id()))
223  continue;
224  ++nTop;
225  Particle top = getLastInstance(Particle(ptcl));
226  if (top.genParticle()->end_vertex())
227  for (auto child : top.children())
228  if (PID::isW(child.pid()))
229  Ws += child;
230  }
231  }
232 
233  // Make sure result make sense
234  if ((prodMode == HTXS::TTH && Ws.size() < 2) || (prodMode == HTXS::TH && Ws.empty()))
235  return error(cat, HTXS::TOP_W_IDENTIFICATION, "Failed to identify W-boson(s) from t-decay!");
236 
237  /*****
238  * Step 3.
239  * Build jets
240  * Make sure all stable particles are present
241  */
242 
243  // Create a list of the vector bosons that decay leptonically
244  // Either the vector boson produced in association with the Higgs boson,
245  // or the ones produced from decays of top quarks produced with the Higgs
246  Particles leptonicVs;
247  if (!is_uncatdV) {
248  if (isVH(prodMode) && !quarkDecay(cat.V))
249  leptonicVs += cat.V;
250  } else
251  leptonicVs = uncatV_decays;
252  for (auto W : Ws)
253  if (W.genParticle()->end_vertex() && !quarkDecay(W))
254  leptonicVs += W;
255 
256  // Obtain all stable, final-state particles
257  const Particles FS = apply<FinalState>(event, "FS").particles();
258  Particles hadrons;
259  FourMomentum sum(0, 0, 0, 0), vSum(0, 0, 0, 0), hSum(0, 0, 0, 0);
260  for (const Particle &p : FS) {
261  // Add up the four momenta of all stable particles as a cross check
262  sum += p.momentum();
263  // ignore particles from the Higgs boson
264  if (originateFrom(p, cat.higgs)) {
265  hSum += p.momentum();
266  continue;
267  }
268  // Cross-check the V decay products for VH
269  if (isVH(prodMode) && !is_uncatdV && originateFrom(p, Ws))
270  vSum += p.momentum();
271  // ignore final state particles from leptonic V decays
272  if (!leptonicVs.empty() && originateFrom(p, leptonicVs))
273  continue;
274  // All particles reaching here are considered hadrons and will be used to build jets
275  hadrons += p;
276  }
277 
278  cat.p4decay_higgs = hSum;
279  cat.p4decay_V = vSum;
280 
281  FinalState fps_temp;
282  FastJets jets(fps_temp, FastJets::ANTIKT, 0.4);
283  jets.calc(hadrons);
284 
285  cat.jets25 = jets.jetsByPt(Cuts::pT > 25.0);
286  cat.jets30 = jets.jetsByPt(Cuts::pT > 30.0);
287 
288  // check that four mometum sum of all stable particles satisfies momentum consevation
289  /*
290  if ( sum.pt()>0.1 )
291  return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
292  std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
293 */
294  // check if V-boson was not included in the event record but decay particles were
295  // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
296  if (is_uncatdV)
297  return error(cat, HTXS::VH_IDENTIFICATION, "Failed to identify associated V-boson!");
298 
299  /*****
300  * Step 4.
301  * Classify and save output
302  */
303 
304  // Apply the categorization categorization
305  cat.stage0_cat = getStage0Category(prodMode, cat.higgs, cat.V);
306  cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode, cat.higgs, cat.jets25, cat.V);
307  cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode, cat.higgs, cat.jets30, cat.V);
308  cat.stage1_1_cat_pTjet25GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets25, cat.V);
309  cat.stage1_1_cat_pTjet30GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets30, cat.V);
310  cat.stage1_1_fine_cat_pTjet25GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
311  cat.stage1_1_fine_cat_pTjet30GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
312 
313  cat.errorCode = HTXS::SUCCESS;
315  ++m_sumevents;
316 
317  return cat;
318  }
319 
325 
327  int getBin(double x, std::vector<double> bins) {
328  for (size_t i = 1; i < bins.size(); ++i)
329  if (x < bins[i])
330  return i - 1;
331  return bins.size() - 1;
332  }
333 
337  int vbfTopology(const Jets &jets, const Particle &higgs) {
338  if (jets.size() < 2)
339  return 0;
340  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
341  bool VBFtopo = (j1 + j2).mass() > 400.0 && std::abs(j1.rapidity() - j2.rapidity()) > 2.8;
342  return VBFtopo ? (j1 + j2 + higgs.momentum()).pt() < 25 ? 2 : 1 : 0;
343  }
344 
349  int vbfTopology_Stage1_1(const Jets &jets, const Particle &higgs) {
350  if (jets.size() < 2)
351  return 0;
352  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
353  double mjj = (j1 + j2).mass();
354  if (mjj > 350 && mjj <= 700)
355  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
356  else if (mjj > 700)
357  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
358  else
359  return 0;
360  }
361 
368  int vbfTopology_Stage1_1_Fine(const Jets &jets, const Particle &higgs) {
369  if (jets.size() < 2)
370  return 0;
371  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
372  double mjj = (j1 + j2).mass();
373  if (mjj > 350 && mjj <= 700)
374  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
375  else if (mjj > 700 && mjj <= 1000)
376  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
377  else if (mjj > 1000 && mjj <= 1500)
378  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 5 : 6;
379  else if (mjj > 1500)
380  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 7 : 8;
381  else
382  return 0;
383  }
384 
386  bool isVH(HTXS::HiggsProdMode p) { return p == HTXS::WH || p == HTXS::QQ2ZH || p == HTXS::GG2ZH; }
387 
390  const Particle &higgs,
391  const Particle &V) {
392  using namespace HTXS::Stage0;
393  int ctrlHiggs = std::abs(higgs.rapidity()) < 2.5;
394  // Special cases first, qq→Hqq
395  if ((prodMode == HTXS::WH || prodMode == HTXS::QQ2ZH) && quarkDecay(V)) {
396  return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
397  } else if (prodMode == HTXS::GG2ZH && quarkDecay(V)) {
398  return Category(HTXS::GGF * 10 + ctrlHiggs);
399  }
400  // General case after
401  return Category(prodMode * 10 + ctrlHiggs);
402  }
403 
406  const Particle &higgs,
407  const Jets &jets,
408  const Particle &V) {
409  using namespace HTXS::Stage1;
410  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
411  double pTj1 = !jets.empty() ? jets[0].momentum().pt() : 0;
412  int vbfTopo = vbfTopology(jets, higgs);
413 
414  // 1. GGF Stage 1 categories
415  // Following YR4 write-up: XXXXX
416  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
417  if (fwdHiggs)
418  return GG2H_FWDH;
419  if (Njets == 0)
420  return GG2H_0J;
421  else if (Njets == 1)
422  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
423  else if (Njets >= 2) {
424  // events with pT_H>200 get priority over VBF cuts
425  if (higgs.pt() <= 200) {
426  if (vbfTopo == 2)
427  return GG2H_VBFTOPO_JET3VETO;
428  else if (vbfTopo == 1)
429  return GG2H_VBFTOPO_JET3;
430  }
431  // Njets >= 2jets without VBF topology
432  return Category(GG2H_GE2J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
433  }
434  }
435  // 2. Electroweak qq->Hqq Stage 1 categories
436  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
437  if (std::abs(higgs.rapidity()) > 2.5)
438  return QQ2HQQ_FWDH;
439  if (pTj1 > 200)
440  return QQ2HQQ_PTJET1_GT200;
441  if (vbfTopo == 2)
443  if (vbfTopo == 1)
444  return QQ2HQQ_VBFTOPO_JET3;
445  double mjj = jets.size() > 1 ? (jets[0].mom() + jets[1].mom()).mass() : 0;
446  if (60 < mjj && mjj < 120)
447  return QQ2HQQ_VH2JET;
448  return QQ2HQQ_REST;
449  }
450  // 3. WH->Hlv categories
451  else if (prodMode == HTXS::WH) {
452  if (fwdHiggs)
453  return QQ2HLNU_FWDH;
454  else if (V.pt() < 150)
455  return QQ2HLNU_PTV_0_150;
456  else if (V.pt() > 250)
457  return QQ2HLNU_PTV_GT250;
458  // 150 < pTV/GeV < 250
459  return jets.empty() ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
460  }
461  // 4. qq->ZH->llH categories
462  else if (prodMode == HTXS::QQ2ZH) {
463  if (fwdHiggs)
464  return QQ2HLL_FWDH;
465  else if (V.pt() < 150)
466  return QQ2HLL_PTV_0_150;
467  else if (V.pt() > 250)
468  return QQ2HLL_PTV_GT250;
469  // 150 < pTV/GeV < 250
470  return jets.empty() ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
471  }
472  // 5. gg->ZH->llH categories
473  else if (prodMode == HTXS::GG2ZH) {
474  if (fwdHiggs)
475  return GG2HLL_FWDH;
476  if (V.pt() < 150)
477  return GG2HLL_PTV_0_150;
478  else if (jets.empty())
479  return GG2HLL_PTV_GT150_0J;
480  return GG2HLL_PTV_GT150_GE1J;
481  }
482  // 6.ttH,bbH,tH categories
483  else if (prodMode == HTXS::TTH)
484  return Category(TTH_FWDH + ctrlHiggs);
485  else if (prodMode == HTXS::BBH)
486  return Category(BBH_FWDH + ctrlHiggs);
487  else if (prodMode == HTXS::TH)
488  return Category(TH_FWDH + ctrlHiggs);
489  return UNKNOWN;
490  }
491 
494  const Particle &higgs,
495  const Jets &jets,
496  const Particle &V) {
497  using namespace HTXS::Stage1_1;
498  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
499  int vbfTopo = vbfTopology_Stage1_1(jets, higgs);
500 
501  // 1. GGF Stage 1 categories
502  // Following YR4 write-up: XXXXX
503  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
504  if (fwdHiggs)
505  return GG2H_FWDH;
506  if (higgs.pt() > 200)
507  return GG2H_PTH_GT200;
508  if (Njets == 0)
509  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
510  if (Njets == 1)
511  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
512  if (Njets > 1) {
513  //VBF topology
514  if (vbfTopo)
515  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
516  //Njets >= 2jets without VBF topology (mjj<350)
517  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
518  }
519  }
520 
521  // 2. Electroweak qq->Hqq Stage 1.1 categories
522  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
523  if (std::abs(higgs.rapidity()) > 2.5)
524  return QQ2HQQ_FWDH;
525  int Njets = jets.size();
526  if (Njets == 0)
527  return QQ2HQQ_0J;
528  else if (Njets == 1)
529  return QQ2HQQ_1J;
530  else if (Njets >= 2) {
531  double mjj = (jets[0].mom() + jets[1].mom()).mass();
532  if (mjj < 60)
533  return QQ2HQQ_MJJ_0_60;
534  else if (60 < mjj && mjj < 120)
535  return QQ2HQQ_MJJ_60_120;
536  else if (120 < mjj && mjj < 350)
537  return QQ2HQQ_MJJ_120_350;
538  else if (mjj > 350) {
539  if (higgs.pt() > 200)
541  if (vbfTopo)
542  return Category(QQ2HQQ_MJJ_GT350_PTH_GT200 + vbfTopo);
543  }
544  }
545  }
546  // 3. WH->Hlv categories
547  else if (prodMode == HTXS::WH) {
548  if (fwdHiggs)
549  return QQ2HLNU_FWDH;
550  else if (V.pt() < 75)
551  return QQ2HLNU_PTV_0_75;
552  else if (V.pt() < 150)
553  return QQ2HLNU_PTV_75_150;
554  else if (V.pt() > 250)
555  return QQ2HLNU_PTV_GT250;
556  // 150 < pTV/GeV < 250
557  return jets.empty() ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
558  }
559  // 4. qq->ZH->llH categories
560  else if (prodMode == HTXS::QQ2ZH) {
561  if (fwdHiggs)
562  return QQ2HLL_FWDH;
563  else if (V.pt() < 75)
564  return QQ2HLL_PTV_0_75;
565  else if (V.pt() < 150)
566  return QQ2HLL_PTV_75_150;
567  else if (V.pt() > 250)
568  return QQ2HLL_PTV_GT250;
569  // 150 < pTV/GeV < 250
570  return jets.empty() ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
571  }
572  // 5. gg->ZH->llH categories
573  else if (prodMode == HTXS::GG2ZH) {
574  if (fwdHiggs)
575  return GG2HLL_FWDH;
576  else if (V.pt() < 75)
577  return GG2HLL_PTV_0_75;
578  else if (V.pt() < 150)
579  return GG2HLL_PTV_75_150;
580  else if (V.pt() > 250)
581  return GG2HLL_PTV_GT250;
582  return jets.empty() ? GG2HLL_PTV_150_250_0J : GG2HLL_PTV_150_250_GE1J;
583  }
584  // 6.ttH,bbH,tH categories
585  else if (prodMode == HTXS::TTH)
586  return Category(TTH_FWDH + ctrlHiggs);
587  else if (prodMode == HTXS::BBH)
588  return Category(BBH_FWDH + ctrlHiggs);
589  else if (prodMode == HTXS::TH)
590  return Category(TH_FWDH + ctrlHiggs);
591  return UNKNOWN;
592  }
593 
596  const Particle &higgs,
597  const Jets &jets,
598  const Particle &V) {
599  using namespace HTXS::Stage1_1_Fine;
600  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
601  int vbfTopo = vbfTopology_Stage1_1_Fine(jets, higgs);
602 
603  // 1. GGF Stage 1.1 categories
604  // Following YR4 write-up: XXXXX
605  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
606  if (fwdHiggs)
607  return GG2H_FWDH;
608  if (higgs.pt() > 200)
609  return GG2H_PTH_GT200;
610  if (Njets == 0)
611  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
612  if (Njets == 1)
613  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
614  if (Njets > 1) {
615  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
616  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
617  //VBF topology
618  if (vbfTopo)
619  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
620  //Njets >= 2jets without VBF topology (mjj<350)
621  if (pTHjj < 25)
622  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
623  else
624  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
625  }
626  }
627 
628  // 2. Electroweak qq->Hqq Stage 1.1 categories
629  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
630  if (std::abs(higgs.rapidity()) > 2.5)
631  return QQ2HQQ_FWDH;
632  int Njets = jets.size();
633  if (Njets == 0)
634  return QQ2HQQ_0J;
635  else if (Njets == 1)
636  return QQ2HQQ_1J;
637  else if (Njets >= 2) {
638  double mjj = (jets[0].mom() + jets[1].mom()).mass();
639  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
640  if (mjj < 350) {
641  if (pTHjj < 25)
642  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
643  else
644  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
645  } else { //mjj>350 GeV
646  if (higgs.pt() < 200) {
647  return Category(QQ2HQQ_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
648  } else {
650  }
651  }
652  }
653  }
654  // 3. WH->Hlv categories
655  else if (prodMode == HTXS::WH) {
656  if (fwdHiggs)
657  return QQ2HLNU_FWDH;
658  int Njets = jets.size();
659  if (Njets == 0)
660  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
661  if (Njets == 1)
662  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
663  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
664  }
665  // 4. qq->ZH->llH categories
666  else if (prodMode == HTXS::QQ2ZH) {
667  if (fwdHiggs)
668  return QQ2HLL_FWDH;
669  int Njets = jets.size();
670  if (Njets == 0)
671  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
672  if (Njets == 1)
673  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
674  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
675  }
676  // 5. gg->ZH->llH categories
677  else if (prodMode == HTXS::GG2ZH) {
678  if (fwdHiggs)
679  return GG2HLL_FWDH;
680  int Njets = jets.size();
681  if (Njets == 0)
682  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
683  if (Njets == 1)
684  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
685  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
686  }
687  // 6.ttH,bbH,tH categories
688  else if (prodMode == HTXS::TTH)
689  return Category(TTH_FWDH + ctrlHiggs);
690  else if (prodMode == HTXS::BBH)
691  return Category(BBH_FWDH + ctrlHiggs);
692  else if (prodMode == HTXS::TH)
693  return Category(TH_FWDH + ctrlHiggs);
694  return UNKNOWN;
695  }
696 
698 
701 
703  void setHiggsProdMode(HTXS::HiggsProdMode prodMode) { m_HiggsProdMode = prodMode; }
704 
708  void init() override {
709  printf("==============================================================\n");
710  printf("======== HiggsTemplateCrossSections Initialization =========\n");
711  printf("==============================================================\n");
712  // check that the production mode has been set
713  // if running in standalone Rivet the production mode is set through an env variable
715  char *pm_env = std::getenv("HIGGSPRODMODE");
716  string pm(pm_env == nullptr ? "" : pm_env);
717  if (pm == "GGF")
719  else if (pm == "VBF")
721  else if (pm == "WH")
723  else if (pm == "ZH")
725  else if (pm == "QQ2ZH")
727  else if (pm == "GG2ZH")
729  else if (pm == "TTH")
731  else if (pm == "BBH")
733  else if (pm == "TH")
735  else {
736  MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
737  }
738  }
739 
740  // Projections for final state particles
741  const FinalState FS;
742  declare(FS, "FS");
743 
744  // initialize the histograms with for each of the stages
746  m_sumw = 0.0;
747  m_sumevents = 0;
748  printf("==============================================================\n");
749  printf("======== Higgs prod mode %d =========\n", m_HiggsProdMode);
750  printf("======== Sucessful Initialization =========\n");
751  printf("==============================================================\n");
752  }
753 
754  // Perform the per-event analysis
755  void analyze(const Event &event) override {
756  // get the classification
757  HiggsClassification cat = classifyEvent(event, m_HiggsProdMode);
758 
759  // Fill histograms: categorization --> linerize the categories
760  const double weight = 1.0;
761  m_sumw += weight;
762 
763  int F = cat.stage0_cat % 10, P = cat.stage1_cat_pTjet30GeV / 100;
764  hist_stage0->fill(cat.stage0_cat / 10 * 2 + F, weight);
765 
766  // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
767  vector<int> offset({0, 1, 13, 19, 24, 29, 33, 35, 37, 39});
768  int off = offset[P];
769  hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV % 100 + off, weight);
770  hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV % 100 + off, weight);
771 
772  // Fill histograms: variables used in the categorization
773  hist_pT_Higgs->fill(cat.higgs.pT(), weight);
774  hist_y_Higgs->fill(cat.higgs.rapidity(), weight);
775  hist_pT_V->fill(cat.V.pT(), weight);
776 
777  hist_Njets25->fill(cat.jets25.size(), weight);
778  hist_Njets30->fill(cat.jets30.size(), weight);
779 
780  // Jet variables. Use jet collection with pT threshold at 30 GeV
781  if (!cat.jets30.empty())
782  hist_pT_jet1->fill(cat.jets30[0].pt(), weight);
783  if (cat.jets30.size() >= 2) {
784  const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
785  hist_deltay_jj->fill(std::abs(j1.rapidity() - j2.rapidity()), weight);
786  hist_dijet_mass->fill((j1 + j2).mass(), weight);
787  hist_pT_Hjj->fill((j1 + j2 + cat.higgs.momentum()).pt(), weight);
788  }
789  }
790 
792  MSG_INFO(" ====================================================== ");
793  MSG_INFO(" Higgs Template X-Sec Categorization Tool ");
794  MSG_INFO(" Status Code Summary ");
795  MSG_INFO(" ====================================================== ");
796  bool allSuccess = (m_sumevents == m_errorCount[HTXS::SUCCESS]);
797  if (allSuccess)
798  MSG_INFO(" >>>> All " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized!");
799  else {
800  MSG_INFO(" >>>> " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized");
801  MSG_INFO(" >>>> --> the following errors occured:");
802  MSG_INFO(" >>>> " << m_errorCount[HTXS::PRODMODE_DEFINED] << " had an undefined Higgs production mode.");
803  MSG_INFO(" >>>> " << m_errorCount[HTXS::MOMENTUM_CONSERVATION] << " failed momentum conservation.");
804  MSG_INFO(" >>>> " << m_errorCount[HTXS::HIGGS_IDENTIFICATION]
805  << " failed to identify a valid Higgs boson.");
806  MSG_INFO(" >>>> " << m_errorCount[HTXS::HS_VTX_IDENTIFICATION]
807  << " failed to identify the hard scatter vertex.");
808  MSG_INFO(" >>>> " << m_errorCount[HTXS::VH_IDENTIFICATION] << " VH: to identify a valid V-boson.");
809  MSG_INFO(" >>>> " << m_errorCount[HTXS::TOP_W_IDENTIFICATION]
810  << " failed to identify valid Ws from top decay.");
811  }
812  MSG_INFO(" ====================================================== ");
813  MSG_INFO(" ====================================================== ");
814  }
815 
816  void finalize() override {
818  double sf = m_sumw > 0 ? 1.0 / m_sumw : 1.0;
819  for (auto hist : {hist_stage0,
822  hist_Njets25,
823  hist_Njets30,
825  hist_y_Higgs,
826  hist_pT_V,
827  hist_pT_jet1,
830  hist_pT_Hjj})
831  scale(hist, sf);
832  }
833 
834  /*
835  * initialize histograms
836  */
837 
839  book(hist_stage0, "HTXS_stage0", 20, 0, 20);
840  book(hist_stage1_pTjet25, "HTXS_stage1_pTjet25", 40, 0, 40);
841  book(hist_stage1_pTjet30, "HTXS_stage1_pTjet30", 40, 0, 40);
842  book(hist_pT_Higgs, "pT_Higgs", 80, 0, 400);
843  book(hist_y_Higgs, "y_Higgs", 80, -4, 4);
844  book(hist_pT_V, "pT_V", 80, 0, 400);
845  book(hist_pT_jet1, "pT_jet1", 80, 0, 400);
846  book(hist_deltay_jj, "deltay_jj", 50, 0, 10);
847  book(hist_dijet_mass, "m_jj", 50, 0, 2000);
848  book(hist_pT_Hjj, "pT_Hjj", 50, 0, 250);
849  book(hist_Njets25, "Njets25", 10, 0, 10);
850  book(hist_Njets30, "Njets30", 10, 0, 10);
851  }
853 
854  /*
855  * initialize private members used in the classification procedure
856  */
857 
858  private:
859  double m_sumw;
860  size_t m_sumevents;
862  std::map<HTXS::ErrorCode, size_t> m_errorCount;
863  Histo1DPtr hist_stage0;
866  Histo1DPtr hist_pT_V, hist_pT_jet1;
869  };
870 
871  // the PLUGIN only needs to be decleared when running standalone Rivet
872  // and causes compilation / linking issues if included in Athena / RootCore
873  //check for Rivet environment variable RIVET_ANALYSIS_PATH
874 #ifdef RIVET_ANALYSIS_PATH
875  // The hook for the plugin system
877 #endif
878 
879 } // namespace Rivet
880 
881 #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.
TPRegexp parents
Definition: eve_filter.cc:21
int vbfTopology_Stage1_1(const Jets &jets, const Particle &higgs)
VBF topology selection Stage1.1 0 = fail loose selection: m_jj > 350 GeV 1 pass loose, but fail additional cut pT(Hjj)<25. 2 pass pT(Hjj)>25 selection 3 pass tight (m_jj>700 GeV), but fail additional cut pT(Hjj)<25. 4 pass pT(Hjj)>25 selection.
bool originateFrom(const Particle &p, const Particle &p2)
Whether particle p originates from p2.
successful classification
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.
int vbfTopology_Stage1_1_Fine(const Jets &jets, const Particle &higgs)
VBF topology selection for Stage1.1 Fine 0 = fail loose selection: m_jj > 350 GeV 1 pass loose...
std::map< HTXS::ErrorCode, size_t > m_errorCount
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: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.
bool quarkDecay(const Particle &p)
Return true is particle decays to quarks.
failed momentum conservation
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.
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.
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:285
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.
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.
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
production mode not defined
Definition: event.py:1
failed to identify top decay