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 (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(), HepMC::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::isChLepton(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 = event.genEvent()->signal_process_vertex();
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) if HepMC::signal_proces_vertex is missing
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, nTop = 0;
193  if (isVH(prodMode)) {
194  for (auto ptcl : HepMCUtils::particles(HSvtx, HepMC::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, HepMC::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, HepMC::children)) {
240  if (!PID::isTop(ptcl->pdg_id()))
241  continue;
242  ++nTop;
243  Particle top = getLastInstance(Particle(ptcl));
244  if (top.genParticle()->end_vertex())
245  for (const auto &child : top.children())
246  if (PID::isW(child.pid()))
247  Ws += getLastInstance(child);
248  }
249  }
250 
251  // Make sure result make sense
252  if ((prodMode == HTXS::TTH && Ws.size() < 2) || (prodMode == HTXS::TH && Ws.empty()))
253  return error(cat, HTXS::TOP_W_IDENTIFICATION, "Failed to identify W-boson(s) from t-decay!");
254 
255  /*****
256  * Step 3.
257  * Build jets
258  * Make sure all stable particles are present
259  */
260 
261  // Create a list of the vector bosons that decay leptonically
262  // Either the vector boson produced in association with the Higgs boson,
263  // or the ones produced from decays of top quarks produced with the Higgs
264  Particles leptonicVs;
265  if (!is_uncatdV) {
266  if (isVH(prodMode) && !quarkDecay(cat.V))
267  leptonicVs += cat.V;
268  } else
269  leptonicVs = uncatV_decays;
270  for (const auto &W : Ws)
271  if (W.genParticle()->end_vertex() && !quarkDecay(W))
272  leptonicVs += W;
273 
274  // Obtain all stable, final-state particles
275  const Particles FS = apply<FinalState>(event, "FS").particles();
276  Particles hadrons;
277  FourMomentum sum(0, 0, 0, 0), vSum(0, 0, 0, 0), hSum(0, 0, 0, 0);
278  for (const Particle &p : FS) {
279  // Add up the four momenta of all stable particles as a cross check
280  sum += p.momentum();
281  // ignore particles from the Higgs boson
282  if (originateFrom(p, cat.higgs)) {
283  hSum += p.momentum();
284  continue;
285  }
286  // Cross-check the V decay products for VH
287  if (isVH(prodMode) && !is_uncatdV && originateFrom(p, Ws))
288  vSum += p.momentum();
289  // ignore final state particles from leptonic V decays
290  if (!leptonicVs.empty() && originateFrom(p, leptonicVs))
291  continue;
292  // All particles reaching here are considered hadrons and will be used to build jets
293  hadrons += p;
294  }
295 
296  cat.p4decay_higgs = hSum;
297  cat.p4decay_V = vSum;
298 
299  FinalState fps_temp;
300  FastJets jets(fps_temp, FastJets::ANTIKT, 0.4);
301  jets.calc(hadrons);
302 
303  cat.jets25 = jets.jetsByPt(Cuts::pT > 25.0);
304  cat.jets30 = jets.jetsByPt(Cuts::pT > 30.0);
305 
306  // check that four mometum sum of all stable particles satisfies momentum consevation
307  /*
308  if ( sum.pt()>0.1 )
309  return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
310  std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
311 */
312  // check if V-boson was not included in the event record but decay particles were
313  // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
314  if (is_uncatdV)
315  return error(cat, HTXS::VH_IDENTIFICATION, "Failed to identify associated V-boson!");
316 
317  /*****
318  * Step 4.
319  * Classify and save output
320  */
321 
322  // Apply the categorization categorization
323  cat.isZ2vvDecay = false;
324  if ((prodMode == HTXS::GG2ZH || prodMode == HTXS::QQ2ZH) && !quarkDecay(cat.V) && !ChLeptonDecay(cat.V))
325  cat.isZ2vvDecay = true;
326  cat.stage0_cat = getStage0Category(prodMode, cat.higgs, cat.V);
327  cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode, cat.higgs, cat.jets25, cat.V);
328  cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode, cat.higgs, cat.jets30, cat.V);
329  cat.stage1_1_cat_pTjet25GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets25, cat.V);
330  cat.stage1_1_cat_pTjet30GeV = getStage1_1_Category(prodMode, cat.higgs, cat.jets30, cat.V);
331  cat.stage1_1_fine_cat_pTjet25GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
332  cat.stage1_1_fine_cat_pTjet30GeV = getStage1_1_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
333  cat.stage1_2_cat_pTjet25GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets25, cat.V);
334  cat.stage1_2_cat_pTjet30GeV = getStage1_2_Category(prodMode, cat.higgs, cat.jets30, cat.V);
335  cat.stage1_2_fine_cat_pTjet25GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets25, cat.V);
336  cat.stage1_2_fine_cat_pTjet30GeV = getStage1_2_Fine_Category(prodMode, cat.higgs, cat.jets30, cat.V);
337 
338  cat.errorCode = HTXS::SUCCESS;
340  ++m_sumevents;
341 
342  return cat;
343  }
344 
350 
352  int getBin(double x, std::vector<double> bins) {
353  for (size_t i = 1; i < bins.size(); ++i)
354  if (x < bins[i])
355  return i - 1;
356  return bins.size() - 1;
357  }
358 
362  int vbfTopology(const Jets &jets, const Particle &higgs) {
363  if (jets.size() < 2)
364  return 0;
365  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
366  bool VBFtopo = (j1 + j2).mass() > 400.0 && std::abs(j1.rapidity() - j2.rapidity()) > 2.8;
367  return VBFtopo ? (j1 + j2 + higgs.momentum()).pt() < 25 ? 2 : 1 : 0;
368  }
369 
374  int vbfTopology_Stage1_X(const Jets &jets, const Particle &higgs) {
375  if (jets.size() < 2)
376  return 0;
377  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
378  double mjj = (j1 + j2).mass();
379  if (mjj > 350 && mjj <= 700)
380  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
381  else if (mjj > 700)
382  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
383  else
384  return 0;
385  }
386 
393  int vbfTopology_Stage1_X_Fine(const Jets &jets, const Particle &higgs) {
394  if (jets.size() < 2)
395  return 0;
396  const FourMomentum &j1 = jets[0].momentum(), &j2 = jets[1].momentum();
397  double mjj = (j1 + j2).mass();
398  if (mjj > 350 && mjj <= 700)
399  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 1 : 2;
400  else if (mjj > 700 && mjj <= 1000)
401  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 3 : 4;
402  else if (mjj > 1000 && mjj <= 1500)
403  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 5 : 6;
404  else if (mjj > 1500)
405  return (j1 + j2 + higgs.momentum()).pt() < 25 ? 7 : 8;
406  else
407  return 0;
408  }
409 
411  bool isVH(HTXS::HiggsProdMode p) { return p == HTXS::WH || p == HTXS::QQ2ZH || p == HTXS::GG2ZH; }
412 
415  const Particle &higgs,
416  const Particle &V) {
417  using namespace HTXS::Stage0;
418  int ctrlHiggs = std::abs(higgs.rapidity()) < 2.5;
419  // Special cases first, qq→Hqq
420  if ((prodMode == HTXS::WH || prodMode == HTXS::QQ2ZH) && quarkDecay(V)) {
421  return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
422  } else if (prodMode == HTXS::GG2ZH && quarkDecay(V)) {
423  return Category(HTXS::GGF * 10 + ctrlHiggs);
424  }
425  // General case after
426  return Category(prodMode * 10 + ctrlHiggs);
427  }
428 
431  const Particle &higgs,
432  const Jets &jets,
433  const Particle &V) {
434  using namespace HTXS::Stage1;
435  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
436  double pTj1 = !jets.empty() ? jets[0].momentum().pt() : 0;
437  int vbfTopo = vbfTopology(jets, higgs);
438 
439  // 1. GGF Stage 1 categories
440  // Following YR4 write-up: XXXXX
441  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
442  if (fwdHiggs)
443  return GG2H_FWDH;
444  if (Njets == 0)
445  return GG2H_0J;
446  else if (Njets == 1)
447  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
448  else if (Njets >= 2) {
449  // events with pT_H>200 get priority over VBF cuts
450  if (higgs.pt() <= 200) {
451  if (vbfTopo == 2)
452  return GG2H_VBFTOPO_JET3VETO;
453  else if (vbfTopo == 1)
454  return GG2H_VBFTOPO_JET3;
455  }
456  // Njets >= 2jets without VBF topology
457  return Category(GG2H_GE2J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
458  }
459  }
460  // 2. Electroweak qq->Hqq Stage 1 categories
461  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
462  if (std::abs(higgs.rapidity()) > 2.5)
463  return QQ2HQQ_FWDH;
464  if (pTj1 > 200)
465  return QQ2HQQ_PTJET1_GT200;
466  if (vbfTopo == 2)
468  if (vbfTopo == 1)
469  return QQ2HQQ_VBFTOPO_JET3;
470  double mjj = jets.size() > 1 ? (jets[0].mom() + jets[1].mom()).mass() : 0;
471  if (60 < mjj && mjj < 120)
472  return QQ2HQQ_VH2JET;
473  return QQ2HQQ_REST;
474  }
475  // 3. WH->Hlv categories
476  else if (prodMode == HTXS::WH) {
477  if (fwdHiggs)
478  return QQ2HLNU_FWDH;
479  else if (V.pt() < 150)
480  return QQ2HLNU_PTV_0_150;
481  else if (V.pt() > 250)
482  return QQ2HLNU_PTV_GT250;
483  // 150 < pTV/GeV < 250
485  }
486  // 4. qq->ZH->llH categories
487  else if (prodMode == HTXS::QQ2ZH) {
488  if (fwdHiggs)
489  return QQ2HLL_FWDH;
490  else if (V.pt() < 150)
491  return QQ2HLL_PTV_0_150;
492  else if (V.pt() > 250)
493  return QQ2HLL_PTV_GT250;
494  // 150 < pTV/GeV < 250
496  }
497  // 5. gg->ZH->llH categories
498  else if (prodMode == HTXS::GG2ZH) {
499  if (fwdHiggs)
500  return GG2HLL_FWDH;
501  if (V.pt() < 150)
502  return GG2HLL_PTV_0_150;
503  else if (jets.empty())
504  return GG2HLL_PTV_GT150_0J;
505  return GG2HLL_PTV_GT150_GE1J;
506  }
507  // 6.ttH,bbH,tH categories
508  else if (prodMode == HTXS::TTH)
509  return Category(TTH_FWDH + ctrlHiggs);
510  else if (prodMode == HTXS::BBH)
511  return Category(BBH_FWDH + ctrlHiggs);
512  else if (prodMode == HTXS::TH)
513  return Category(TH_FWDH + ctrlHiggs);
514  return UNKNOWN;
515  }
516 
519  const Particle &higgs,
520  const Jets &jets,
521  const Particle &V) {
522  using namespace HTXS::Stage1_1;
523  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
524  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
525 
526  // 1. GGF Stage 1 categories
527  // Following YR4 write-up: XXXXX
528  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
529  if (fwdHiggs)
530  return GG2H_FWDH;
531  if (higgs.pt() > 200)
532  return GG2H_PTH_GT200;
533  if (Njets == 0)
534  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
535  if (Njets == 1)
536  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
537  if (Njets > 1) {
538  //VBF topology
539  if (vbfTopo)
540  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
541  //Njets >= 2jets without VBF topology (mjj<350)
542  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
543  }
544  }
545 
546  // 2. Electroweak qq->Hqq Stage 1.1 categories
547  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
548  if (std::abs(higgs.rapidity()) > 2.5)
549  return QQ2HQQ_FWDH;
550  int Njets = jets.size();
551  if (Njets == 0)
552  return QQ2HQQ_0J;
553  else if (Njets == 1)
554  return QQ2HQQ_1J;
555  else if (Njets >= 2) {
556  double mjj = (jets[0].mom() + jets[1].mom()).mass();
557  if (mjj < 60)
558  return QQ2HQQ_MJJ_0_60;
559  else if (60 < mjj && mjj < 120)
560  return QQ2HQQ_MJJ_60_120;
561  else if (120 < mjj && mjj < 350)
562  return QQ2HQQ_MJJ_120_350;
563  else if (mjj > 350) {
564  if (higgs.pt() > 200)
566  if (vbfTopo)
567  return Category(QQ2HQQ_MJJ_GT350_PTH_GT200 + vbfTopo);
568  }
569  }
570  }
571  // 3. WH->Hlv categories
572  else if (prodMode == HTXS::WH) {
573  if (fwdHiggs)
574  return QQ2HLNU_FWDH;
575  else if (V.pt() < 75)
576  return QQ2HLNU_PTV_0_75;
577  else if (V.pt() < 150)
578  return QQ2HLNU_PTV_75_150;
579  else if (V.pt() > 250)
580  return QQ2HLNU_PTV_GT250;
581  // 150 < pTV/GeV < 250
583  }
584  // 4. qq->ZH->llH categories
585  else if (prodMode == HTXS::QQ2ZH) {
586  if (fwdHiggs)
587  return QQ2HLL_FWDH;
588  else if (V.pt() < 75)
589  return QQ2HLL_PTV_0_75;
590  else if (V.pt() < 150)
591  return QQ2HLL_PTV_75_150;
592  else if (V.pt() > 250)
593  return QQ2HLL_PTV_GT250;
594  // 150 < pTV/GeV < 250
596  }
597  // 5. gg->ZH->llH categories
598  else if (prodMode == HTXS::GG2ZH) {
599  if (fwdHiggs)
600  return GG2HLL_FWDH;
601  else if (V.pt() < 75)
602  return GG2HLL_PTV_0_75;
603  else if (V.pt() < 150)
604  return GG2HLL_PTV_75_150;
605  else if (V.pt() > 250)
606  return GG2HLL_PTV_GT250;
608  }
609  // 6.ttH,bbH,tH categories
610  else if (prodMode == HTXS::TTH)
611  return Category(TTH_FWDH + ctrlHiggs);
612  else if (prodMode == HTXS::BBH)
613  return Category(BBH_FWDH + ctrlHiggs);
614  else if (prodMode == HTXS::TH)
615  return Category(TH_FWDH + ctrlHiggs);
616  return UNKNOWN;
617  }
618 
621  const Particle &higgs,
622  const Jets &jets,
623  const Particle &V) {
624  using namespace HTXS::Stage1_1_Fine;
625  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
626  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
627 
628  // 1. GGF Stage 1.1 categories
629  // Following YR4 write-up: XXXXX
630  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
631  if (fwdHiggs)
632  return GG2H_FWDH;
633  if (higgs.pt() > 200)
634  return GG2H_PTH_GT200;
635  if (Njets == 0)
636  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
637  if (Njets == 1)
638  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
639  if (Njets > 1) {
640  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
641  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
642  //VBF topology
643  if (vbfTopo)
644  return Category(GG2H_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
645  //Njets >= 2jets without VBF topology (mjj<350)
646  if (pTHjj < 25)
647  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
648  else
649  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
650  }
651  }
652 
653  // 2. Electroweak qq->Hqq Stage 1.1 categories
654  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
655  if (std::abs(higgs.rapidity()) > 2.5)
656  return QQ2HQQ_FWDH;
657  int Njets = jets.size();
658  if (Njets == 0)
659  return QQ2HQQ_0J;
660  else if (Njets == 1)
661  return QQ2HQQ_1J;
662  else if (Njets >= 2) {
663  double mjj = (jets[0].mom() + jets[1].mom()).mass();
664  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
665  if (mjj < 350) {
666  if (pTHjj < 25)
667  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
668  else
669  return Category(QQ2HQQ_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
670  } else { //mjj>350 GeV
671  if (higgs.pt() < 200) {
672  return Category(QQ2HQQ_MJJ_350_700_PTHJJ_0_25 + vbfTopo - 1);
673  } else {
675  }
676  }
677  }
678  }
679  // 3. WH->Hlv categories
680  else if (prodMode == HTXS::WH) {
681  if (fwdHiggs)
682  return QQ2HLNU_FWDH;
683  int Njets = jets.size();
684  if (Njets == 0)
685  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
686  if (Njets == 1)
687  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
688  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
689  }
690  // 4. qq->ZH->llH categories
691  else if (prodMode == HTXS::QQ2ZH) {
692  if (fwdHiggs)
693  return QQ2HLL_FWDH;
694  int Njets = jets.size();
695  if (Njets == 0)
696  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
697  if (Njets == 1)
698  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
699  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
700  }
701  // 5. gg->ZH->llH categories
702  else if (prodMode == HTXS::GG2ZH) {
703  if (fwdHiggs)
704  return GG2HLL_FWDH;
705  int Njets = jets.size();
706  if (Njets == 0)
707  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
708  if (Njets == 1)
709  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
710  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
711  }
712  // 6.ttH,bbH,tH categories
713  else if (prodMode == HTXS::TTH)
714  return Category(TTH_FWDH + ctrlHiggs);
715  else if (prodMode == HTXS::BBH)
716  return Category(BBH_FWDH + ctrlHiggs);
717  else if (prodMode == HTXS::TH)
718  return Category(TH_FWDH + ctrlHiggs);
719  return UNKNOWN;
720  }
721 
724  const Particle &higgs,
725  const Jets &jets,
726  const Particle &V) {
727  using namespace HTXS::Stage1_2;
728  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
729  int vbfTopo = vbfTopology_Stage1_X(jets, higgs);
730 
731  // 1. GGF Stage 1 categories
732  // Following YR4 write-up: XXXXX
733  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
734  if (fwdHiggs)
735  return GG2H_FWDH;
736  if (higgs.pt() > 200)
737  return Category(GG2H_PTH_200_300 + getBin(higgs.pt(), {200, 300, 450, 650}));
738  if (Njets == 0)
739  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
740  if (Njets == 1)
741  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
742  if (Njets > 1) {
743  //VBF topology
744  if (vbfTopo)
746  //Njets >= 2jets without VBF topology (mjj<350)
747  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
748  }
749  }
750 
751  // 2. Electroweak qq->Hqq Stage 1.2 categories
752  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
753  if (std::abs(higgs.rapidity()) > 2.5)
754  return QQ2HQQ_FWDH;
755  int Njets = jets.size();
756  if (Njets == 0)
757  return QQ2HQQ_0J;
758  else if (Njets == 1)
759  return QQ2HQQ_1J;
760  else if (Njets >= 2) {
761  double mjj = (jets[0].mom() + jets[1].mom()).mass();
762  if (mjj < 60)
763  return QQ2HQQ_GE2J_MJJ_0_60;
764  else if (60 < mjj && mjj < 120)
765  return QQ2HQQ_GE2J_MJJ_60_120;
766  else if (120 < mjj && mjj < 350)
768  else if (mjj > 350) {
769  if (higgs.pt() > 200)
771  if (vbfTopo)
772  return Category(QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200 + vbfTopo);
773  }
774  }
775  }
776  // 3. WH->Hlv categories
777  else if (prodMode == HTXS::WH) {
778  if (fwdHiggs)
779  return QQ2HLNU_FWDH;
780  else if (V.pt() < 75)
781  return QQ2HLNU_PTV_0_75;
782  else if (V.pt() < 150)
783  return QQ2HLNU_PTV_75_150;
784  else if (V.pt() > 250)
785  return QQ2HLNU_PTV_GT250;
786  // 150 < pTV/GeV < 250
788  }
789  // 4. qq->ZH->llH categories
790  else if (prodMode == HTXS::QQ2ZH) {
791  if (fwdHiggs)
792  return QQ2HLL_FWDH;
793  else if (V.pt() < 75)
794  return QQ2HLL_PTV_0_75;
795  else if (V.pt() < 150)
796  return QQ2HLL_PTV_75_150;
797  else if (V.pt() > 250)
798  return QQ2HLL_PTV_GT250;
799  // 150 < pTV/GeV < 250
801  }
802  // 5. gg->ZH->llH categories
803  else if (prodMode == HTXS::GG2ZH) {
804  if (fwdHiggs)
805  return GG2HLL_FWDH;
806  else if (V.pt() < 75)
807  return GG2HLL_PTV_0_75;
808  else if (V.pt() < 150)
809  return GG2HLL_PTV_75_150;
810  else if (V.pt() > 250)
811  return GG2HLL_PTV_GT250;
813  }
814  // 6.ttH,bbH,tH categories
815  else if (prodMode == HTXS::TTH) {
816  if (fwdHiggs)
817  return TTH_FWDH;
818  else
819  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300}));
820  } else if (prodMode == HTXS::BBH)
821  return Category(BBH_FWDH + ctrlHiggs);
822  else if (prodMode == HTXS::TH)
823  return Category(TH_FWDH + ctrlHiggs);
824  return UNKNOWN;
825  }
826 
829  const Particle &higgs,
830  const Jets &jets,
831  const Particle &V) {
832  using namespace HTXS::Stage1_2_Fine;
833  int Njets = jets.size(), ctrlHiggs = std::abs(higgs.rapidity()) < 2.5, fwdHiggs = !ctrlHiggs;
834  int vbfTopo = vbfTopology_Stage1_X_Fine(jets, higgs);
835 
836  // 1. GGF Stage 1.2 categories
837  // Following YR4 write-up: XXXXX
838  if (prodMode == HTXS::GGF || (prodMode == HTXS::GG2ZH && quarkDecay(V))) {
839  if (fwdHiggs)
840  return GG2H_FWDH;
841  if (higgs.pt() > 200) {
842  if (Njets > 0) {
843  double pTHj = (jets[0].momentum() + higgs.momentum()).pt();
844  if (pTHj / higgs.pt() > 0.15)
845  return Category(GG2H_PTH_200_300_PTHJoverPTH_GT15 + getBin(higgs.pt(), {200, 300, 450, 650}));
846  else
847  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
848  } else
849  return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15 + getBin(higgs.pt(), {200, 300, 450, 650}));
850  }
851  if (Njets == 0)
852  return higgs.pt() < 10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
853  if (Njets == 1)
854  return Category(GG2H_1J_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200}));
855  if (Njets > 1) {
856  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
857  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
858  //VBF topology
859  if (vbfTopo)
861  //Njets >= 2jets without VBF topology (mjj<350)
862  if (pTHjj < 25)
863  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25 + getBin(higgs.pt(), {0, 60, 120, 200}));
864  else
865  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25 + getBin(higgs.pt(), {0, 60, 120, 200}));
866  }
867  }
868 
869  // 2. Electroweak qq->Hqq Stage 1.2 categories
870  else if (prodMode == HTXS::VBF || (isVH(prodMode) && quarkDecay(V))) {
871  if (std::abs(higgs.rapidity()) > 2.5)
872  return QQ2HQQ_FWDH;
873  int Njets = jets.size();
874  if (Njets == 0)
875  return QQ2HQQ_0J;
876  else if (Njets == 1)
877  return QQ2HQQ_1J;
878  else if (Njets >= 2) {
879  double mjj = (jets[0].mom() + jets[1].mom()).mass();
880  double pTHjj = (jets[0].momentum() + jets[1].momentum() + higgs.momentum()).pt();
881  if (mjj < 350) {
882  if (pTHjj < 25)
883  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25 + getBin(mjj, {0, 60, 120, 350}));
884  else
885  return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25 + getBin(mjj, {0, 60, 120, 350}));
886  } else { //mjj>350 GeV
887  if (higgs.pt() < 200) {
889  } else {
891  }
892  }
893  }
894  }
895  // 3. WH->Hlv categories
896  else if (prodMode == HTXS::WH) {
897  if (fwdHiggs)
898  return QQ2HLNU_FWDH;
899  int Njets = jets.size();
900  if (Njets == 0)
901  return Category(QQ2HLNU_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
902  if (Njets == 1)
903  return Category(QQ2HLNU_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
904  return Category(QQ2HLNU_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
905  }
906  // 4. qq->ZH->llH categories
907  else if (prodMode == HTXS::QQ2ZH) {
908  if (fwdHiggs)
909  return QQ2HLL_FWDH;
910  int Njets = jets.size();
911  if (Njets == 0)
912  return Category(QQ2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
913  if (Njets == 1)
914  return Category(QQ2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
915  return Category(QQ2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
916  }
917  // 5. gg->ZH->llH categories
918  else if (prodMode == HTXS::GG2ZH) {
919  if (fwdHiggs)
920  return GG2HLL_FWDH;
921  int Njets = jets.size();
922  if (Njets == 0)
923  return Category(GG2HLL_PTV_0_75_0J + getBin(V.pt(), {0, 75, 150, 250, 400}));
924  if (Njets == 1)
925  return Category(GG2HLL_PTV_0_75_1J + getBin(V.pt(), {0, 75, 150, 250, 400}));
926  return Category(GG2HLL_PTV_0_75_GE2J + getBin(V.pt(), {0, 75, 150, 250, 400}));
927  }
928  // 6.ttH,bbH,tH categories
929  else if (prodMode == HTXS::TTH) {
930  if (fwdHiggs)
931  return TTH_FWDH;
932  else
933  return Category(TTH_PTH_0_60 + getBin(higgs.pt(), {0, 60, 120, 200, 300, 450}));
934  } else if (prodMode == HTXS::BBH)
935  return Category(BBH_FWDH + ctrlHiggs);
936  else if (prodMode == HTXS::TH)
937  return Category(TH_FWDH + ctrlHiggs);
938  return UNKNOWN;
939  }
940 
942 
945 
947  void setHiggsProdMode(HTXS::HiggsProdMode prodMode) { m_HiggsProdMode = prodMode; }
948 
952  void init() override {
953  printf("==============================================================\n");
954  printf("======== HiggsTemplateCrossSections Initialization =========\n");
955  printf("==============================================================\n");
956  // check that the production mode has been set
957  // if running in standalone Rivet the production mode is set through an env variable
959  char *pm_env = std::getenv("HIGGSPRODMODE");
960  string pm(pm_env == nullptr ? "" : pm_env);
961  if (pm == "GGF")
963  else if (pm == "VBF")
965  else if (pm == "WH")
967  else if (pm == "ZH")
969  else if (pm == "QQ2ZH")
971  else if (pm == "GG2ZH")
973  else if (pm == "TTH")
975  else if (pm == "BBH")
977  else if (pm == "TH")
979  else {
980  MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
981  }
982  }
983 
984  // Projections for final state particles
985  const FinalState FS;
986  declare(FS, "FS");
987 
988  // initialize the histograms with for each of the stages
990  m_sumw = 0.0;
991  m_sumevents = 0;
992  printf("==============================================================\n");
993  printf("======== Higgs prod mode %d =========\n", m_HiggsProdMode);
994  printf("======== Sucessful Initialization =========\n");
995  printf("==============================================================\n");
996  }
997 
998  // Perform the per-event analysis
999  void analyze(const Event &event) override {
1000  // get the classification
1001  HiggsClassification cat = classifyEvent(event, m_HiggsProdMode);
1002 
1003  // Fill histograms: categorization --> linerize the categories
1004  const double weight = 1.0;
1005  m_sumw += weight;
1006 
1007  int F = cat.stage0_cat % 10, P = cat.stage1_cat_pTjet30GeV / 100;
1008  hist_stage0->fill(cat.stage0_cat / 10 * 2 + F, weight);
1009 
1010  // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
1011  vector<int> offset({0, 1, 13, 19, 24, 29, 33, 35, 37, 39});
1012  int off = offset[P];
1013  hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV % 100 + off, weight);
1014  hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV % 100 + off, weight);
1015 
1016  // 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
1017  static vector<int> offset1_2({0, 1, 18, 29, 35, 41, 47, 53, 55, 57});
1018  int off1_2 = offset1_2[P];
1019  // 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
1020  static vector<int> offset1_2f({0, 1, 29, 54, 70, 86, 102, 109, 111, 113});
1021  int off1_2f = offset1_2f[P];
1022  hist_stage1_2_pTjet25->fill(cat.stage1_2_cat_pTjet25GeV % 100 + off1_2, weight);
1023  hist_stage1_2_pTjet30->fill(cat.stage1_2_cat_pTjet30GeV % 100 + off1_2, weight);
1024  hist_stage1_2_fine_pTjet25->fill(cat.stage1_2_fine_cat_pTjet25GeV % 100 + off1_2f, weight);
1025  hist_stage1_2_fine_pTjet30->fill(cat.stage1_2_fine_cat_pTjet30GeV % 100 + off1_2f, weight);
1026 
1027  // Fill histograms: variables used in the categorization
1028  hist_pT_Higgs->fill(cat.higgs.pT(), weight);
1029  hist_y_Higgs->fill(cat.higgs.rapidity(), weight);
1030  hist_pT_V->fill(cat.V.pT(), weight);
1031 
1032  hist_Njets25->fill(cat.jets25.size(), weight);
1033  hist_Njets30->fill(cat.jets30.size(), weight);
1034 
1035  hist_isZ2vv->fill(cat.isZ2vvDecay, weight);
1036 
1037  // Jet variables. Use jet collection with pT threshold at 30 GeV
1038  if (!cat.jets30.empty())
1039  hist_pT_jet1->fill(cat.jets30[0].pt(), weight);
1040  if (cat.jets30.size() >= 2) {
1041  const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
1042  hist_deltay_jj->fill(std::abs(j1.rapidity() - j2.rapidity()), weight);
1043  hist_dijet_mass->fill((j1 + j2).mass(), weight);
1044  hist_pT_Hjj->fill((j1 + j2 + cat.higgs.momentum()).pt(), weight);
1045  }
1046  }
1047 
1049  MSG_INFO(" ====================================================== ");
1050  MSG_INFO(" Higgs Template X-Sec Categorization Tool ");
1051  MSG_INFO(" Status Code Summary ");
1052  MSG_INFO(" ====================================================== ");
1053  bool allSuccess = (m_sumevents == m_errorCount[HTXS::SUCCESS]);
1054  if (allSuccess)
1055  MSG_INFO(" >>>> All " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized!");
1056  else {
1057  MSG_INFO(" >>>> " << m_errorCount[HTXS::SUCCESS] << " events successfully categorized");
1058  MSG_INFO(" >>>> --> the following errors occured:");
1059  MSG_INFO(" >>>> " << m_errorCount[HTXS::PRODMODE_DEFINED] << " had an undefined Higgs production mode.");
1060  MSG_INFO(" >>>> " << m_errorCount[HTXS::MOMENTUM_CONSERVATION] << " failed momentum conservation.");
1061  MSG_INFO(" >>>> " << m_errorCount[HTXS::HIGGS_IDENTIFICATION]
1062  << " failed to identify a valid Higgs boson.");
1063  MSG_INFO(" >>>> " << m_errorCount[HTXS::HS_VTX_IDENTIFICATION]
1064  << " failed to identify the hard scatter vertex.");
1065  MSG_INFO(" >>>> " << m_errorCount[HTXS::VH_IDENTIFICATION] << " VH: to identify a valid V-boson.");
1066  MSG_INFO(" >>>> " << m_errorCount[HTXS::TOP_W_IDENTIFICATION]
1067  << " failed to identify valid Ws from top decay.");
1068  }
1069  MSG_INFO(" ====================================================== ");
1070  MSG_INFO(" ====================================================== ");
1071  }
1072 
1073  void finalize() override {
1075  double sf = m_sumw > 0 ? 1.0 / m_sumw : 1.0;
1076  for (const auto &hist : {hist_stage0,
1083  hist_Njets25,
1084  hist_Njets30,
1085  hist_pT_Higgs,
1086  hist_y_Higgs,
1087  hist_pT_V,
1088  hist_pT_jet1,
1091  hist_pT_Hjj})
1092  scale(hist, sf);
1093  }
1094 
1095  /*
1096  * initialize histograms
1097  */
1098 
1099  void initializeHistos() {
1100  book(hist_stage0, "HTXS_stage0", 20, 0, 20);
1101  book(hist_stage1_pTjet25, "HTXS_stage1_pTjet25", 40, 0, 40);
1102  book(hist_stage1_pTjet30, "HTXS_stage1_pTjet30", 40, 0, 40);
1103  book(hist_stage1_2_pTjet25, "HTXS_stage1_2_pTjet25", 57, 0, 57);
1104  book(hist_stage1_2_pTjet30, "HTXS_stage1_2_pTjet30", 57, 0, 57);
1105  book(hist_stage1_2_fine_pTjet25, "HTXS_stage1_2_fine_pTjet25", 113, 0, 113);
1106  book(hist_stage1_2_fine_pTjet30, "HTXS_stage1_2_fine_pTjet30", 113, 0, 113);
1107  book(hist_pT_Higgs, "pT_Higgs", 80, 0, 400);
1108  book(hist_y_Higgs, "y_Higgs", 80, -4, 4);
1109  book(hist_pT_V, "pT_V", 80, 0, 400);
1110  book(hist_pT_jet1, "pT_jet1", 80, 0, 400);
1111  book(hist_deltay_jj, "deltay_jj", 50, 0, 10);
1112  book(hist_dijet_mass, "m_jj", 50, 0, 2000);
1113  book(hist_pT_Hjj, "pT_Hjj", 50, 0, 250);
1114  book(hist_Njets25, "Njets25", 10, 0, 10);
1115  book(hist_Njets30, "Njets30", 10, 0, 10);
1116  book(hist_isZ2vv, "isZ2vv", 2, 0, 2);
1117  }
1119 
1120  /*
1121  * initialize private members used in the classification procedure
1122  */
1124  private:
1125  double m_sumw;
1126  size_t m_sumevents;
1128  std::map<HTXS::ErrorCode, size_t> m_errorCount;
1129  Histo1DPtr hist_stage0;
1134  Histo1DPtr hist_pT_V, hist_pT_jet1;
1136  Histo1DPtr hist_Njets25, hist_Njets30;
1137  Histo1DPtr hist_isZ2vv;
1138  };
1139 
1140  // the PLUGIN only needs to be decleared when running standalone Rivet
1141  // and causes compilation / linking issues if included in Athena / RootCore
1142  //check for Rivet environment variable RIVET_ANALYSIS_PATH
1143 #ifdef RIVET_ANALYSIS_PATH
1144  // The hook for the plugin system
1146 #endif
1147 
1148 } // namespace Rivet
1149 
1150 #endif
Rivet::HiggsTemplateCrossSections::printClassificationSummary
void printClassificationSummary()
Definition: HiggsTemplateCrossSections.cc:1044
HTXS::VH_IDENTIFICATION
failed to identify associated vector boson
Definition: HiggsTemplateCrossSections.h:19
HTXS::Stage1_1::QQ2HQQ_MJJ_60_120
Definition: HiggsTemplateCrossSections.h:139
Rivet::HiggsTemplateCrossSections::getStage1_1_Fine_Category
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.
Definition: HiggsTemplateCrossSections.cc:616
HTXS::Stage1::QQ2HQQ_REST
Definition: HiggsTemplateCrossSections.h:85
HTXS::Stage1_1::QQ2HLL_PTV_75_150
Definition: HiggsTemplateCrossSections.h:156
Rivet::HiggsTemplateCrossSections::hasParent
bool hasParent(const ConstGenParticlePtr &ptcl, int pdgID)
Checks whether the input particle has a parent with a given PDGID.
Definition: HiggsTemplateCrossSections.cc:75
HTXS::Stage0::Category
Category
Definition: HiggsTemplateCrossSections.h:37
HTXS::Stage1_1_Fine::QQ2HQQ_MJJ_0_60_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:208
HTXS::TOP_W_IDENTIFICATION
failed to identify top decay
Definition: HiggsTemplateCrossSections.h:21
HTXS::Stage1_1_Fine::QQ2HLL_PTV_0_75_1J
Definition: HiggsTemplateCrossSections.h:254
HTXS::Stage1_1::QQ2HQQ_MJJ_GT350_PTH_GT200
Definition: HiggsTemplateCrossSections.h:141
mps_fire.i
i
Definition: mps_fire.py:428
HTXS::Stage1::QQ2HQQ_VBFTOPO_JET3
Definition: HiggsTemplateCrossSections.h:83
HTXS::Stage0::QQ2HLL_FWDH
Definition: HiggsTemplateCrossSections.h:47
HTXS::Stage1_1::GG2HLL_PTV_150_250_0J
Definition: HiggsTemplateCrossSections.h:164
Rivet::HiggsTemplateCrossSections::analyze
void analyze(const Event &event) override
Definition: HiggsTemplateCrossSections.cc:995
HTXS::Stage1_2
Definition: HiggsTemplateCrossSections.h:297
HTXS::Stage1_1_Fine::QQ2HQQ_MJJ_0_60_PTHJJ_GT25
Definition: HiggsTemplateCrossSections.h:211
HTXS::GG2ZH
Definition: HiggsTemplateCrossSections.h:25
HTXS::UNKNOWN
Definition: HiggsTemplateCrossSections.h:25
L1EGammaCrystalsEmulatorProducer_cfi.scale
scale
Definition: L1EGammaCrystalsEmulatorProducer_cfi.py:10
HTXS::TTH
Definition: HiggsTemplateCrossSections.h:25
Rivet::HiggsTemplateCrossSections::m_sumevents
size_t m_sumevents
Definition: HiggsTemplateCrossSections.cc:1122
HTXS::PRODMODE_DEFINED
production mode not defined
Definition: HiggsTemplateCrossSections.h:14
Rivet::HiggsTemplateCrossSections::hist_deltay_jj
Histo1DPtr hist_deltay_jj
Definition: HiggsTemplateCrossSections.cc:1131
HTXS::Stage1::QQ2HQQ_VBFTOPO_JET3VETO
Definition: HiggsTemplateCrossSections.h:82
Rivet::HiggsTemplateCrossSections::init
void init() override
default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Hig...
Definition: HiggsTemplateCrossSections.cc:948
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
HTXS::Stage1::QQ2HLNU_PTV_150_250_GE1J
Definition: HiggsTemplateCrossSections.h:91
Rivet::HiggsTemplateCrossSections::vbfTopology_Stage1_X
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: HiggsTemplateCrossSections.cc:370
HTXS::Stage1_1_Fine::GG2HLL_PTV_0_75_GE2J
Definition: HiggsTemplateCrossSections.h:276
HTXS::Stage1_2_Fine::GG2H_PTH_200_300_PTHJoverPTH_0_15
Definition: HiggsTemplateCrossSections.h:372
Rivet
Definition: RivetAnalysis.h:21
class-composition.children
children
Definition: class-composition.py:98
mps_merge.weight
weight
Definition: mps_merge.py:88
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200
Definition: HiggsTemplateCrossSections.h:325
eostools.cat
def cat(path)
Definition: eostools.py:401
HTXS::Stage1_2_Fine::QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25
Definition: HiggsTemplateCrossSections.h:406
HTXS::Stage1_1::QQ2HLL_PTV_0_75
Definition: HiggsTemplateCrossSections.h:155
HTXS::Stage1_1::GG2H_MJJ_350_700_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:130
Rivet::HiggsTemplateCrossSections::quarkDecay
bool quarkDecay(const Particle &p)
Return true is particle decays to quarks.
Definition: HiggsTemplateCrossSections.cc:83
HTXS::GGF
Definition: HiggsTemplateCrossSections.h:25
Rivet::HiggsTemplateCrossSections::ChLeptonDecay
bool ChLeptonDecay(const Particle &p)
Return true if particle decays to charged leptons.
Definition: HiggsTemplateCrossSections.cc:93
HTXS::Stage1_1_Fine::QQ2HLNU_PTV_0_75_GE2J
Definition: HiggsTemplateCrossSections.h:242
Rivet::HiggsTemplateCrossSections::originateFrom
bool originateFrom(const Particle &p, const Particles &ptcls)
Whether particle p originate from any of the ptcls.
Definition: HiggsTemplateCrossSections.cc:44
HTXS::Stage0
Namespace for Stage0 categorization.
Definition: HiggsTemplateCrossSections.h:35
Rivet::HiggsTemplateCrossSections
Rivet routine for classifying MC events according to the Higgs template cross section categories.
Definition: HiggsTemplateCrossSections.cc:21
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_60_120
Definition: HiggsTemplateCrossSections.h:323
Rivet::HiggsTemplateCrossSections::getLastInstance
Particle getLastInstance(Particle ptcl)
follow a "propagating" particle and return its last instance
Definition: HiggsTemplateCrossSections.cc:33
Rivet::HiggsTemplateCrossSections::hist_pT_V
Histo1DPtr hist_pT_V
Definition: HiggsTemplateCrossSections.cc:1130
HTXS::Stage1::QQ2HLL_PTV_150_250_GE1J
Definition: HiggsTemplateCrossSections.h:97
HTXS::Stage1_1::QQ2HQQ_1J
Definition: HiggsTemplateCrossSections.h:137
HTXS::Stage1::GG2H_GE2J_PTH_0_60
Definition: HiggsTemplateCrossSections.h:76
HTXS::Stage1_2::GG2H_PTH_200_300
Definition: HiggsTemplateCrossSections.h:302
mps_check.msg
tuple msg
Definition: mps_check.py:285
Rivet::HiggsTemplateCrossSections::hasChild
bool hasChild(const ConstGenParticlePtr &ptcl, int pdgID)
Checks whether the input particle has a child with a given PDGID.
Definition: HiggsTemplateCrossSections.cc:65
HTXS::Stage0::VH2HQQ_FWDH
Definition: HiggsTemplateCrossSections.h:43
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
HTXS::Stage1_1_Fine::Category
Category
Definition: HiggsTemplateCrossSections.h:180
Rivet::HiggsTemplateCrossSections::error
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,...
Definition: HiggsTemplateCrossSections.cc:104
Rivet::HiggsTemplateCrossSections::hist_stage1_pTjet25
Histo1DPtr hist_stage1_pTjet25
Definition: HiggsTemplateCrossSections.cc:1126
HTXS::Stage1_1
Definition: HiggsTemplateCrossSections.h:116
ecalTrigSettings_cff.particles
particles
Definition: ecalTrigSettings_cff.py:11
Rivet::HiggsTemplateCrossSections::getBin
int getBin(double x, std::vector< double > bins)
Return bin index of x given the provided bin edges. 0=first bin, -1=underflow bin.
Definition: HiggsTemplateCrossSections.cc:348
HiggsTemplateCrossSections.h
HTXS::HIGGS_DECAY_IDENTIFICATION
failed to identify Higgs boson decay products
Definition: HiggsTemplateCrossSections.h:17
F
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
HTXS::Stage1_1::GG2HLL_PTV_0_75
Definition: HiggsTemplateCrossSections.h:162
HTXS::Stage1::QQ2HQQ_PTJET1_GT200
Definition: HiggsTemplateCrossSections.h:86
HTXS::Stage1::QQ2HLL_PTV_150_250_0J
Definition: HiggsTemplateCrossSections.h:96
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_120_350
Definition: HiggsTemplateCrossSections.h:324
HTXS::Stage1_1::Category
Category
Definition: HiggsTemplateCrossSections.h:117
HTXS::Stage1_1_Fine::GG2HLL_PTV_0_75_0J
Definition: HiggsTemplateCrossSections.h:266
part
part
Definition: HCALResponse.h:20
HTXS::Stage1::GG2HLL_PTV_GT150_0J
Definition: HiggsTemplateCrossSections.h:102
HTXS::Stage1_1::GG2HLL_PTV_75_150
Definition: HiggsTemplateCrossSections.h:163
HTXS::Stage1_2::GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:314
HTXS::Stage1_1::GG2H_0J_PTH_GT10
Definition: HiggsTemplateCrossSections.h:123
HTXS::Stage1_1::UNKNOWN
Definition: HiggsTemplateCrossSections.h:118
Rivet::HiggsTemplateCrossSections::hist_stage1_2_pTjet25
Histo1DPtr hist_stage1_2_pTjet25
Definition: HiggsTemplateCrossSections.cc:1127
Rivet::HiggsTemplateCrossSections::hist_pT_Higgs
Histo1DPtr hist_pT_Higgs
Definition: HiggsTemplateCrossSections.cc:1129
HTXS::Stage1
Definition: HiggsTemplateCrossSections.h:64
PVValHelper::pT
Definition: PVValidationHelpers.h:71
HTXS::ErrorCode
ErrorCode
Error code: whether the classification was successful or failed.
Definition: HiggsTemplateCrossSections.h:11
Rivet::HiggsTemplateCrossSections::isVH
bool isVH(HTXS::HiggsProdMode p)
Whether the Higgs is produced in association with a vector boson (VH)
Definition: HiggsTemplateCrossSections.cc:407
p2
double p2[4]
Definition: TauolaWrapper.h:90
Rivet::HiggsTemplateCrossSections::hist_stage1_2_pTjet30
Histo1DPtr hist_stage1_2_pTjet30
Definition: HiggsTemplateCrossSections.cc:1127
ECF_cff.Njets
Njets
Definition: ECF_cff.py:21
HTXS::Stage0::TH_FWDH
Definition: HiggsTemplateCrossSections.h:55
Event
HTXS::Stage1_2_Fine::Category
Category
Definition: HiggsTemplateCrossSections.h:368
HTXS::Stage1_2::UNKNOWN
Definition: HiggsTemplateCrossSections.h:299
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_0_60
Definition: HiggsTemplateCrossSections.h:322
Rivet::HiggsTemplateCrossSections::hist_isZ2vv
Histo1DPtr hist_isZ2vv
Definition: HiggsTemplateCrossSections.cc:1133
Rivet::HiggsTemplateCrossSections::hist_pT_Hjj
Histo1DPtr hist_pT_Hjj
Definition: HiggsTemplateCrossSections.cc:1131
Rivet::HiggsTemplateCrossSections::getStage0Category
HTXS::Stage0::Category getStage0Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Particle &V)
Stage-0 HTXS categorization.
Definition: HiggsTemplateCrossSections.cc:410
HTXS::Stage1_1_Fine::GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25
Definition: HiggsTemplateCrossSections.h:193
HTXS::Stage1_1_Fine::QQ2HLL_PTV_0_75_0J
Definition: HiggsTemplateCrossSections.h:249
HTXS::VH_DECAY_IDENTIFICATION
failed to identify associated vector boson decay products
Definition: HiggsTemplateCrossSections.h:20
HTXS::Stage1::GG2H_VBFTOPO_JET3VETO
Definition: HiggsTemplateCrossSections.h:69
HTXS::Stage1::UNKNOWN
Definition: HiggsTemplateCrossSections.h:66
HTXS::Stage1::QQ2HLL_PTV_GT250
Definition: HiggsTemplateCrossSections.h:98
HTXS::Stage1::GG2H_0J
Definition: HiggsTemplateCrossSections.h:71
HTXS::Stage1_1_Fine::UNKNOWN
Definition: HiggsTemplateCrossSections.h:181
bJpsiMuMuTrigSettings_cff.hadrons
hadrons
Definition: bJpsiMuMuTrigSettings_cff.py:8
HTXS::Stage1_1_Fine::QQ2HLNU_PTV_0_75_1J
Definition: HiggsTemplateCrossSections.h:237
HTXS::Stage1_1::QQ2HQQ_0J
Definition: HiggsTemplateCrossSections.h:136
HTXS::Stage1_1::QQ2HQQ_MJJ_0_60
Definition: HiggsTemplateCrossSections.h:138
HTXS::Stage1_1_Fine::QQ2HQQ_PTH_GT200_MJJ_350_700_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:222
HTXS::Stage1::GG2HLL_PTV_GT150_GE1J
Definition: HiggsTemplateCrossSections.h:103
HTXS::Stage0::BBH_FWDH
Definition: HiggsTemplateCrossSections.h:53
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
HTXS::Stage1_1::QQ2HQQ_MJJ_350_700_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:142
HTXS::Stage0::UNKNOWN
Definition: HiggsTemplateCrossSections.h:38
Rivet::HiggsTemplateCrossSections::m_sumw
double m_sumw
Definition: HiggsTemplateCrossSections.cc:1121
Rivet::HiggsTemplateCrossSections::hist_stage1_2_fine_pTjet30
Histo1DPtr hist_stage1_2_fine_pTjet30
Definition: HiggsTemplateCrossSections.cc:1128
gpuVertexFinder::hist
__shared__ Hist hist
Definition: gpuClusterTracksDBSCAN.h:48
Rivet::HiggsTemplateCrossSections::m_HiggsProdMode
HTXS::HiggsProdMode m_HiggsProdMode
Definition: HiggsTemplateCrossSections.cc:1123
HTXS::UNDEFINED
Definition: HiggsTemplateCrossSections.h:12
cms::cuda::V
uint32_t const T *__restrict__ const uint32_t *__restrict__ int32_t int Histo::index_type cudaStream_t V
Definition: HistoContainer.h:51
HTXS
Higgs Template Cross Section namespace.
Definition: HiggsTemplateCrossSections.h:8
Rivet::HiggsTemplateCrossSections::hist_Njets25
Histo1DPtr hist_Njets25
Definition: HiggsTemplateCrossSections.cc:1132
HTXS::MOMENTUM_CONSERVATION
failed momentum conservation
Definition: HiggsTemplateCrossSections.h:15
HTXS::Stage1_2_Fine::QQ2HQQ_GE2J_MJJ_350_700_PTH_GT200_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:417
HTXS::Stage1::QQ2HLNU_PTV_150_250_0J
Definition: HiggsTemplateCrossSections.h:90
HTXS::HiggsProdMode
HiggsProdMode
Higgs production modes, corresponding to input sample.
Definition: HiggsTemplateCrossSections.h:25
HTXS::Stage1_1::GG2H_GE2J_MJJ_0_350_PTH_0_60
Definition: HiggsTemplateCrossSections.h:127
Rivet::HiggsTemplateCrossSections::getStage1_2_Category
HTXS::Stage1_2::Category getStage1_2_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1.2 categorization.
Definition: HiggsTemplateCrossSections.cc:719
HTXS::Stage1_1::QQ2HLNU_PTV_0_75
Definition: HiggsTemplateCrossSections.h:148
HTXS::Stage0::TTH_FWDH
Definition: HiggsTemplateCrossSections.h:51
HTXS::Stage1_2_Fine
Definition: HiggsTemplateCrossSections.h:367
Rivet::HiggsTemplateCrossSections::hist_stage0
Histo1DPtr hist_stage0
Definition: HiggsTemplateCrossSections.cc:1125
HTXS::Stage1::QQ2HQQ_VH2JET
Definition: HiggsTemplateCrossSections.h:84
HTXS::Stage1::Category
Category
Definition: HiggsTemplateCrossSections.h:65
HTXS::Stage0::VH2HQQ
Definition: HiggsTemplateCrossSections.h:44
submitPVResolutionJobs.err
err
Definition: submitPVResolutionJobs.py:85
HTXS::HIGGS_IDENTIFICATION
failed to identify Higgs boson
Definition: HiggsTemplateCrossSections.h:16
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Rivet::DECLARE_RIVET_PLUGIN
DECLARE_RIVET_PLUGIN(CMS_2013_I1224539_DIJET)
Rivet::HiggsTemplateCrossSections::hist_pT_jet1
Histo1DPtr hist_pT_jet1
Definition: HiggsTemplateCrossSections.cc:1130
METSkim_cff.Jets
Jets
Definition: METSkim_cff.py:17
HTXS::Stage0::GG2H_FWDH
Definition: HiggsTemplateCrossSections.h:39
SiStripSubdetector::UNKNOWN
Definition: SiStripEnums.h:5
HTXS::Stage1_1_Fine
Definition: HiggsTemplateCrossSections.h:179
HTXS::Stage1_2::Category
Category
Definition: HiggsTemplateCrossSections.h:298
Rivet::HiggsTemplateCrossSections::HiggsTemplateCrossSections
HiggsTemplateCrossSections()
Definition: HiggsTemplateCrossSections.cc:24
Rivet::HiggsTemplateCrossSections::hist_stage1_2_fine_pTjet25
Histo1DPtr hist_stage1_2_fine_pTjet25
Definition: HiggsTemplateCrossSections.cc:1128
Rivet::HiggsTemplateCrossSections::m_errorCount
std::map< HTXS::ErrorCode, size_t > m_errorCount
Definition: HiggsTemplateCrossSections.cc:1124
Rivet::HiggsTemplateCrossSections::originateFrom
bool originateFrom(const Particle &p, const Particle &p2)
Whether particle p originates from p2.
Definition: HiggsTemplateCrossSections.cc:59
HTXS::Stage1::QQ2HLNU_PTV_0_150
Definition: HiggsTemplateCrossSections.h:89
HTXS::Stage1::QQ2HLNU_PTV_GT250
Definition: HiggsTemplateCrossSections.h:92
HTXS::Stage1_1::GG2HLL_PTV_GT250
Definition: HiggsTemplateCrossSections.h:166
HTXS::Stage1_1_Fine::GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:190
HTXS::Stage1::GG2H_1J_PTH_0_60
Definition: HiggsTemplateCrossSections.h:72
HTXS::Stage1_1::QQ2HQQ_MJJ_120_350
Definition: HiggsTemplateCrossSections.h:140
HTXS::Stage1_1_Fine::QQ2HLNU_PTV_0_75_0J
Definition: HiggsTemplateCrossSections.h:232
HTXS::Stage1_1::QQ2HLNU_PTV_75_150
Definition: HiggsTemplateCrossSections.h:149
HTXS::Stage1_1_Fine::QQ2HLL_PTV_0_75_GE2J
Definition: HiggsTemplateCrossSections.h:259
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
HTXS::Stage1_1_Fine::GG2HLL_PTV_0_75_1J
Definition: HiggsTemplateCrossSections.h:271
HTXS::Stage1_2_Fine::QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:403
Rivet::HiggsTemplateCrossSections::hist_stage1_pTjet30
Histo1DPtr hist_stage1_pTjet30
Definition: HiggsTemplateCrossSections.cc:1126
Rivet::HiggsTemplateCrossSections::setHiggsProdMode
void setHiggsProdMode(HTXS::HiggsProdMode prodMode)
Sets the Higgs production mode.
Definition: HiggsTemplateCrossSections.cc:943
Rivet::HiggsTemplateCrossSections::vbfTopology
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,...
Definition: HiggsTemplateCrossSections.cc:358
Rivet::HiggsTemplateCrossSections::hist_dijet_mass
Histo1DPtr hist_dijet_mass
Definition: HiggsTemplateCrossSections.cc:1131
HTXS::Stage1_2_Fine::GG2H_PTH_200_300_PTHJoverPTH_GT15
Definition: HiggsTemplateCrossSections.h:376
Particle
Definition: Particle.py:1
Rivet::HiggsTemplateCrossSections::getStage1_2_Fine_Category
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.
Definition: HiggsTemplateCrossSections.cc:824
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25
Definition: HiggsTemplateCrossSections.h:326
Rivet::HiggsTemplateCrossSections::getStage1Category
HTXS::Stage1::Category getStage1Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1 categorization.
Definition: HiggsTemplateCrossSections.cc:426
Rivet::HiggsTemplateCrossSections::vbfTopology_Stage1_X_Fine
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...
Definition: HiggsTemplateCrossSections.cc:389
HTXS::HS_VTX_IDENTIFICATION
failed to identify hard scatter vertex
Definition: HiggsTemplateCrossSections.h:18
HTXS::Stage0::QQ2HLNU_FWDH
Definition: HiggsTemplateCrossSections.h:45
HTXS::Stage1_1::GG2H_0J_PTH_0_10
Definition: HiggsTemplateCrossSections.h:122
HTXS::Stage1_1::GG2HLL_PTV_150_250_GE1J
Definition: HiggsTemplateCrossSections.h:165
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trigObjTnPSource_cfi.bins
bins
Definition: trigObjTnPSource_cfi.py:20
Rivet::HiggsTemplateCrossSections::hist_y_Higgs
Histo1DPtr hist_y_Higgs
Definition: HiggsTemplateCrossSections.cc:1129
Rivet::HiggsTemplateCrossSections::finalize
void finalize() override
Definition: HiggsTemplateCrossSections.cc:1069
HTXS::TH
Definition: HiggsTemplateCrossSections.h:25
parents
TPRegexp parents
Definition: eve_filter.cc:21
P
std::pair< OmniClusterRef, TrackingParticleRef > P
Definition: BDHadronTrackMonitoringAnalyzer.cc:203
event
Definition: event.py:1
hltrates_dqm_sourceclient-live_cfg.offset
offset
Definition: hltrates_dqm_sourceclient-live_cfg.py:82
child
Definition: simpleInheritance.h:11
HTXS::SUCCESS
successful classification
Definition: HiggsTemplateCrossSections.h:13
HTXS::Stage1::QQ2HLL_PTV_0_150
Definition: HiggsTemplateCrossSections.h:95
HTXS::Stage0::GG2HLL_FWDH
Definition: HiggsTemplateCrossSections.h:49
class-composition.parent
parent
Definition: class-composition.py:98
Rivet::HiggsTemplateCrossSections::hist_Njets30
Histo1DPtr hist_Njets30
Definition: HiggsTemplateCrossSections.cc:1132
HTXS::Stage1_2_Fine::UNKNOWN
Definition: HiggsTemplateCrossSections.h:369
HTXS::Stage1::GG2H_VBFTOPO_JET3
Definition: HiggsTemplateCrossSections.h:70
HTXS::WH
Definition: HiggsTemplateCrossSections.h:25
weight
Definition: weight.py:1
Rivet::HiggsTemplateCrossSections::classifyEvent
HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode)
Main classificaion method.
Definition: HiggsTemplateCrossSections.cc:122
HTXS::Stage1_2::TTH_PTH_0_60
Definition: HiggsTemplateCrossSections.h:353
Rivet::HiggsTemplateCrossSections::getStage1_1_Category
HTXS::Stage1_1::Category getStage1_1_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1.1 categorization.
Definition: HiggsTemplateCrossSections.cc:514
HTXS::VBF
Definition: HiggsTemplateCrossSections.h:25
HTXS::BBH
Definition: HiggsTemplateCrossSections.h:25
HTXS::Stage1_1::GG2H_PTH_GT200
Definition: HiggsTemplateCrossSections.h:121
Rivet::HiggsTemplateCrossSections::initializeHistos
void initializeHistos()
Definition: HiggsTemplateCrossSections.cc:1095
HTXS::QQ2ZH
Definition: HiggsTemplateCrossSections.h:25
HTXS::Stage1::GG2HLL_PTV_0_150
Definition: HiggsTemplateCrossSections.h:101
HTXS::Stage1::QQ2HQQ_FWDH
Definition: HiggsTemplateCrossSections.h:81