ATLAS Offline Software
HiggsTemplateCrossSections.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #ifndef TRUTHRIVETTOOLS_HIGGSTEMPLATECROSSSECTIONS_H
6 #define TRUTHRIVETTOOLS_HIGGSTEMPLATECROSSSECTIONS_H
7 
8 // -*- C++ -*-
9 #include "Rivet/Analysis.hh"
10 #include "Rivet/Particle.hh"
11 #include "Rivet/Projections/FastJets.hh"
12 
13 // Definition of the StatusCode and Category enums
14 // Note: the Template XSec Defs *depends* on having included
15 // the TLorentzVector header *before* it is included -- it
16 // uses the include guard from TLorentzVector to decide
17 // what is available
18 #include "TLorentzVector.h"
20 #include "AtlasHepMC/Relatives.h"
21 #include "AtlasHepMC/GenEvent.h"
22 #include "AtlasHepMC/GenVertex.h"
23 #include "AtlasHepMC/GenParticle.h"
24 #include "AtlasHepMC/GenRanges.h"
26 
27 namespace Rivet {
28 
34  public:
35  // Constructor
37  : Analysis("HiggsTemplateCrossSections"),
38  m_HiggsProdMode(HTXS::UNKNOWN) {}
39 
40  public:
41 
46 
48  Particle getLastInstance(const Particle & ptcl) const {
49  if ( ptcl.genParticle()->end_vertex() ) {
50  if ( !hasChild(ptcl.genParticle(),ptcl.pid()) ) return ptcl;
51  else return getLastInstance(ptcl.children()[0]);
52  }
53  return ptcl;
54  }
55 
57  bool originateFrom(const Particle& p, const Particles& ptcls ) const {
58  auto prodVtx = p.genParticle()->production_vertex();
59  if (prodVtx == nullptr) return false;
60  // for each ancestor, check if it matches any of the input particles
61  for (auto ancestor:Rivet::HepMCUtils::particles(prodVtx,Relatives::ANCESTORS)){
62  for ( auto part:ptcls )
63  if ( ancestor==part.genParticle() ) return true;
64  }
65  // if we get here, no ancetor matched any input particle
66  return false;
67  }
68 
70  bool originateFrom(const Particle& p, const Particle& p2 ) const {
71  Particles ptcls = {p2}; return originateFrom(p,ptcls);
72  }
73 
75  bool hasChild(HepMC::ConstGenParticlePtr ptcl, int pdgID) const {
76  for (const Particle& child:Particle(*ptcl).children())
77  if (child.pid()==pdgID) return true;
78  return false;
79  }
80 
82  bool hasParent(HepMC::ConstGenParticlePtr ptcl, int pdgID) const {
83  for (auto parent:Rivet::HepMCUtils::particles(ptcl->production_vertex(),Relatives::PARENTS))
84  if (parent->pdg_id()==pdgID) return true;
85  return false;
86  }
87 
89  bool quarkDecay(const Particle &p) const {
90  for (const Particle& child:p.children())
91  if (PID::isQuark(child.pid())) return true;
92  return false;
93  }
94 
96  bool ChLeptonDecay(const Particle &p) const {
97  for (const Particle& child:p.children())
98  if (
99 #if RIVET_VERSION_CODE >= 40000
100  PID::isChargedLepton(child.pid())
101 #else
102  PID::isChLepton(child.pid())
103 #endif // RIVET_VERSION_CODE
104  ) return true;
105  return false;
106  }
107 
110  HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err,
111  std::string msg="", int NmaxWarnings=20) const {
112  // Set the error, and keep statistics
113  cat.errorCode = err;
114  ++m_errorCount[err];
115 
116  // Print warning message to the screen/log
117  static std::atomic<int> Nwarnings = 0;
118  if ( !msg.empty() && ++Nwarnings < NmaxWarnings )
119  MSG_WARNING(msg);
120 
121  return cat;
122  }
124 
126  HiggsClassification classifyEvent(const Event& event, const HTXS::HiggsProdMode prodMode ) const {
127 
128  // the classification object
129  HiggsClassification cat;
130  cat.prodMode = prodMode;
131  cat.errorCode = HTXS::UNDEFINED;
132  cat.stage0_cat = HTXS::Stage0::UNKNOWN;
133  cat.stage1_cat_pTjet25GeV = HTXS::Stage1::UNKNOWN;
134  cat.stage1_cat_pTjet30GeV = HTXS::Stage1::UNKNOWN;
135  cat.stage1_2_cat_pTjet25GeV = HTXS::Stage1_2::UNKNOWN;
136  cat.stage1_2_cat_pTjet30GeV = HTXS::Stage1_2::UNKNOWN;
137  cat.stage1_2_fine_cat_pTjet25GeV = HTXS::Stage1_2_Fine::UNKNOWN;
138  cat.stage1_2_fine_cat_pTjet30GeV = HTXS::Stage1_2_Fine::UNKNOWN;
139 
140  if (prodMode == HTXS::UNKNOWN)
141  return error(cat,HTXS::PRODMODE_DEFINED,
142  "Unkown Higgs production mechanism. Cannot classify event."
143  " Classification for all events will most likely fail.");
144 
145  /*****
146  * Step 1.
147  * Idenfify the Higgs boson and the hard scatter vertex
148  * There should be only one of each.
149  */
150 
151  auto HSvtx = HepMC::signal_process_vertex(event.genEvent());
152  int Nhiggs=0;
153  for (auto ptcl : Rivet::HepMCUtils::particles(event.genEvent()) ) {
154 
155  // a) Reject all non-Higgs particles
156  if ( !PID::isHiggs(ptcl->pdg_id()) ) continue;
157  // b) select only the final Higgs boson copy, prior to decay
158  if ( ptcl->end_vertex() && !hasChild(ptcl,PID::HIGGS) ) {
159  cat.higgs = Particle(ptcl); ++Nhiggs;
160  }
161  // c) if HepMC::signal_proces_vertex is missing
162  // set hard-scatter vertex based on first Higgs boson
163  if ( HSvtx==nullptr && ptcl->production_vertex() && !hasParent(ptcl,PID::HIGGS) )
164  HSvtx = ptcl->production_vertex();
165  }
166 
167  // Make sure things are in order so far
168  if (Nhiggs!=1)
170  "Current event has "+std::to_string(Nhiggs)+" Higgs bosons. There must be only one.");
171  if (cat.higgs.children().size()<2)
173  "Could not identify Higgs boson decay products.");
174 
175  if (HSvtx == nullptr)
176  return error(cat,HTXS::HS_VTX_IDENTIFICATION,"Cannot find hard-scatter vertex of current event.");
177 
178  /*****
179  * Step 2.
180  * Identify associated vector bosons
181  */
182 
183  // Find associated vector bosons
184  bool is_uncatdV = false;
185  Particles uncatV_decays;
186  FourMomentum uncatV_p4(0,0,0,0);
187  FourVector uncatV_v4(0,0,0,0);
188  int nWs=0, nZs=0;
189  if ( isVH(prodMode) ) {
190  for (auto ptcl:Rivet::HepMCUtils::particles(HSvtx,Relatives::CHILDREN)) {
191  if (PID::isW(ptcl->pdg_id())) { ++nWs; cat.V=Particle(ptcl); }
192  if (PID::isZ(ptcl->pdg_id())) { ++nZs; cat.V=Particle(ptcl); }
193  }
194  if(nWs+nZs>0) cat.V = getLastInstance(cat.V);
195  else {
196  for (auto ptcl:Rivet::HepMCUtils::particles(HSvtx,Relatives::CHILDREN)) {
197  if (!PID::isHiggs(ptcl->pdg_id())) {
198  uncatV_decays += Particle(ptcl);
199  uncatV_p4 += Particle(ptcl).momentum();
200  uncatV_v4 += Particle(ptcl).origin();
201  }
202  }
203  is_uncatdV = true; cat.V = Particle(24,uncatV_p4,uncatV_v4);
204  }
205  }
206 
207  if ( !is_uncatdV ){
208 
209  if ( isVH(prodMode) && !cat.V.genParticle()->end_vertex() )
210  return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
211 
212  if ( isVH(prodMode) && cat.V.children().size()<2 )
213  return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
214 
215  if ( ( prodMode==HTXS::WH && (nZs>0||nWs!=1) ) ||
216  ( (prodMode==HTXS::QQ2ZH||prodMode==HTXS::GG2ZH) && (nZs!=1||nWs>0) ) )
217  return error(cat,HTXS::VH_IDENTIFICATION,"Found "+std::to_string(nWs)+" W-bosons and "+
218  std::to_string(nZs)+" Z-bosons. Inconsitent with VH expectation.");
219  }
220 
221  // Find and store the W-bosons from ttH->WbWbH
222  Particles Ws;
223  if ( prodMode==HTXS::TTH || prodMode==HTXS::TH ){
224  // loop over particles produced in hard-scatter vertex
225  for ( auto ptcl : Rivet::HepMCUtils::particles(HSvtx,Relatives::CHILDREN) ) {
226  if ( !PID::isTop(ptcl->pdg_id()) ) continue;
228  if ( top.genParticle()->end_vertex() )
229  for (auto child:top.children())
230  if ( PID::isW(child.pid()) ) Ws += getLastInstance(child);
231  }
232  }
233 
234  // Make sure result make sense
235  if ( (prodMode==HTXS::TTH && Ws.size()<2) || (prodMode==HTXS::TH && Ws.size()<1 ) )
236  return error(cat,HTXS::TOP_W_IDENTIFICATION,"Failed to identify W-boson(s) from t-decay!");
237 
238  /*****
239  * Step 3.
240  * Build jets
241  * Make sure all stable particles are present
242  */
243 
244  // Create a list of the vector bosons that decay leptonically
245  // Either the vector boson produced in association with the Higgs boson,
246  // or the ones produced from decays of top quarks produced with the Higgs
247  Particles leptonicVs;
248  if ( !is_uncatdV ){
249  if ( isVH(prodMode) && !quarkDecay(cat.V) ) leptonicVs += cat.V;
250  }else leptonicVs = uncatV_decays;
251  for ( auto W:Ws ) if ( W.genParticle()->end_vertex() && !quarkDecay(W) ) leptonicVs += W;
252 
253  // Obtain all stable, final-state particles
254  const Particles FS = apply<FinalState>(event, "FS").particles();
255  Particles hadrons;
256  FourMomentum sum(0,0,0,0), vSum(0,0,0,0), hSum(0,0,0,0);
257  for ( const Particle &p : FS ) {
258  // Add up the four momenta of all stable particles as a cross check
259  sum += p.momentum();
260  // ignore particles from the Higgs boson
261  if ( originateFrom(p,cat.higgs) ) { hSum += p.momentum(); continue; }
262  // Cross-check the V decay products for VH
263  if ( isVH(prodMode) && !is_uncatdV && originateFrom(p,Ws) ) vSum += p.momentum();
264  // ignore final state particles from leptonic V decays
265  if ( leptonicVs.size() && originateFrom(p,leptonicVs) ) continue;
266  // All particles reaching here are considered hadrons and will be used to build jets
267  hadrons += p;
268  }
269 
270  cat.p4decay_higgs = hSum;
271  cat.p4decay_V = vSum;
272 
273  FinalState fps_temp;
274  FastJets jets(fps_temp,
275 #if RIVET_VERSION_CODE >= 40000
276  JetAlg::ANTIKT,
277 #else
278  FastJets::ANTIKT,
279 #endif // RIVET_VERSION_CODE
280  0.4 );
281  jets.calc(hadrons);
282 
283  cat.jets25 = jets.jetsByPt( Cuts::pT > 25.0 );
284  cat.jets30 = jets.jetsByPt( Cuts::pT > 30.0 );
285 
286  // check that four mometum sum of all stable particles satisfies momentum consevation
287  if ( sum.pt()>0.1 )
288  return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
289  std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
290 
291  // check if V-boson was not included in the event record but decay particles were
292  // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
293  if(is_uncatdV)
294  return error(cat,HTXS::VH_IDENTIFICATION,"Failed to identify associated V-boson!");
295 
296  /*****
297  * Step 4.
298  * Classify and save output
299  */
300 
301  // Apply the categorization categorization
302  cat.isZ2vvDecay = false;
303  if( (prodMode==HTXS::GG2ZH || prodMode==HTXS::QQ2ZH) && !quarkDecay(cat.V) && !ChLeptonDecay(cat.V) ) cat.isZ2vvDecay = true;
304  cat.stage0_cat = getStage0Category(prodMode,cat.higgs,cat.V);
305  cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode,cat.higgs,cat.jets25,cat.V);
306  cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode,cat.higgs,cat.jets30,cat.V);
307  cat.stage1_2_cat_pTjet25GeV = getStage1_2_Category(prodMode,cat.higgs,cat.jets25,cat.V);
308  cat.stage1_2_cat_pTjet30GeV = getStage1_2_Category(prodMode,cat.higgs,cat.jets30,cat.V);
309  cat.stage1_2_fine_cat_pTjet25GeV = getStage1_2_Fine_Category(prodMode,cat.higgs,cat.jets25,cat.V);
310  cat.stage1_2_fine_cat_pTjet30GeV = getStage1_2_Fine_Category(prodMode,cat.higgs,cat.jets30,cat.V);
311  cat.errorCode = HTXS::SUCCESS; ++m_errorCount[HTXS::SUCCESS];
312 
313  return cat;
314  }
315 
321 
323  int getBin(double x, const std::vector<double>& bins) const {
324  if (bins.size()==0||x<bins[0]) return -1; // should not happen!
325  for (size_t i=1;i<bins.size();++i)
326  if (x<bins[i]) return i-1;
327  return bins.size()-1;
328  }
329 
333  int vbfTopology(const Jets &jets, const Particle &higgs) const {
334  if (jets.size()<2) return 0;
335  const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
336  bool VBFtopo = (j1+j2).mass() > 400.0 && std::abs(j1.rapidity()-j2.rapidity()) > 2.8;
337  return VBFtopo ? (j1+j2+higgs.momentum()).pt()<25 ? 2 : 1 : 0;
338  }
343  int vbfTopology_Stage1_2(const Jets &jets, const Particle &higgs) const {
344  if (jets.size()<2) return 0;
345  const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
346  double mjj = (j1+j2).mass();
347  if(mjj>350 && mjj<=700) return (j1+j2+higgs.momentum()).pt()<25 ? 1 : 2;
348  else if(mjj>700) return (j1+j2+higgs.momentum()).pt()<25 ? 3 : 4;
349  else return 0;
350  }
357  int vbfTopology_Stage1_2_Fine(const Jets &jets, const Particle &higgs) const {
358  if (jets.size()<2) return 0;
359  const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
360  double mjj = (j1+j2).mass();
361  if(mjj>350 && mjj<=700) return (j1+j2+higgs.momentum()).pt()<25 ? 1 : 2;
362  else if(mjj>700 && mjj<=1000) return (j1+j2+higgs.momentum()).pt()<25 ? 3 : 4;
363  else if(mjj>1000 && mjj<=1500) return (j1+j2+higgs.momentum()).pt()<25 ? 5 : 6;
364  else if(mjj>1500) return (j1+j2+higgs.momentum()).pt()<25 ? 7 : 8;
365  else return 0;
366  }
367 
369  bool isVH(HTXS::HiggsProdMode p) const { return p==HTXS::WH || p==HTXS::QQ2ZH || p==HTXS::GG2ZH; }
370 
373  const Particle &higgs,
374  const Particle &V) const {
375  using namespace HTXS::Stage0;
376  int ctrlHiggs = std::abs(higgs.rapidity())<2.5;
377  // Special cases first, qq→Hqq
378  if ( (prodMode==HTXS::WH||prodMode==HTXS::QQ2ZH) && quarkDecay(V) ) {
379  return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
380  } else if ( prodMode==HTXS::GG2ZH && quarkDecay(V) ) {
381  return Category(HTXS::GGF*10 + ctrlHiggs);
382  }
383  // General case after
384  return Category(prodMode*10 + ctrlHiggs);
385  }
386 
389  const Particle &higgs,
390  const Jets &jets,
391  const Particle &V) const {
392  using namespace HTXS::Stage1;
393  int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
394  double pTj1 = jets.size() ? jets[0].momentum().pt() : 0;
395  int vbfTopo = vbfTopology(jets,higgs);
396 
397  // 1. GGF Stage 1 categories
398  // Following YR4 write-up: XXXXX
399  if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
400  if (fwdHiggs) return GG2H_FWDH;
401  if (Njets==0) return GG2H_0J;
402  else if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
403  else if (Njets>=2) {
404  // events with pT_H>200 get priority over VBF cuts
405  if(higgs.pt()<=200){
406  if (vbfTopo==2) return GG2H_VBFTOPO_JET3VETO;
407  else if (vbfTopo==1) return GG2H_VBFTOPO_JET3;
408  }
409  // Njets >= 2jets without VBF topology
410  return Category(GG2H_GE2J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
411  }
412  }
413  // 2. Electroweak qq->Hqq Stage 1 categories
414  else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
415  if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
416  if (pTj1>200) return QQ2HQQ_PTJET1_GT200;
417  if (vbfTopo==2) return QQ2HQQ_VBFTOPO_JET3VETO;
418  if (vbfTopo==1) return QQ2HQQ_VBFTOPO_JET3;
419  double mjj = jets.size()>1 ? (jets[0].mom()+jets[1].mom()).mass():0;
420  if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_VH2JET;
421  return QQ2HQQ_REST;
422  }
423  // 3. WH->Hlv categories
424  else if (prodMode==HTXS::WH) {
425  if (fwdHiggs) return QQ2HLNU_FWDH;
426  else if (V.pt()<150) return QQ2HLNU_PTV_0_150;
427  else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
428  // 150 < pTV/GeV < 250
430  }
431  // 4. qq->ZH->llH categories
432  else if (prodMode==HTXS::QQ2ZH) {
433  if (fwdHiggs) return QQ2HLL_FWDH;
434  else if (V.pt()<150) return QQ2HLL_PTV_0_150;
435  else if (V.pt()>250) return QQ2HLL_PTV_GT250;
436  // 150 < pTV/GeV < 250
438  }
439  // 5. gg->ZH->llH categories
440  else if (prodMode==HTXS::GG2ZH ) {
441  if (fwdHiggs) return GG2HLL_FWDH;
442  if (V.pt()<150) return GG2HLL_PTV_0_150;
443  else if (jets.size()==0) return GG2HLL_PTV_GT150_0J;
444  return GG2HLL_PTV_GT150_GE1J;
445  }
446  // 6.ttH,bbH,tH categories
447  else if (prodMode==HTXS::TTH) return Category(TTH_FWDH+ctrlHiggs);
448  else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
449  else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
450  return UNKNOWN;
451  }
452 
455  const Particle &higgs,
456  const Jets &jets,
457  const Particle &V) const {
458  using namespace HTXS::Stage1_2;
459  int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
460  int vbfTopo = vbfTopology_Stage1_2(jets,higgs);
461 
462  // 1. GGF Stage 1 categories
463  // Following YR4 write-up: XXXXX
464  if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
465  if (fwdHiggs) return GG2H_FWDH;
466  if ( higgs.pt()>200 ) return Category(GG2H_PTH_200_300+getBin(higgs.pt(),{200,300,450,650}));
467  if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
468  if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
469  if (Njets>1){
470  //VBF topology
471  if(vbfTopo) return Category(GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25+vbfTopo-1);
472  //Njets >= 2jets without VBF topology (mjj<350)
473  return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
474  }
475  }
476 
477  // 2. Electroweak qq->Hqq Stage 1.2 categories
478  else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
479  if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
480  int Njets=jets.size();
481  if (Njets==0) return QQ2HQQ_0J;
482  else if (Njets==1) return QQ2HQQ_1J;
483  else if (Njets>=2) {
484  double mjj = (jets[0].mom()+jets[1].mom()).mass();
485  if ( mjj < 60 ) return QQ2HQQ_GE2J_MJJ_0_60;
486  else if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_GE2J_MJJ_60_120;
487  else if ( 120 < mjj && mjj < 350 ) return QQ2HQQ_GE2J_MJJ_120_350;
488  else if ( mjj > 350 ) {
489  if (higgs.pt()>200) return QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200;
490  if(vbfTopo) return Category(QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200+vbfTopo);
491  }
492  }
493  }
494  // 3. WH->Hlv categories
495  else if (prodMode==HTXS::WH) {
496  if (fwdHiggs) return QQ2HLNU_FWDH;
497  else if (V.pt()<75) return QQ2HLNU_PTV_0_75;
498  else if (V.pt()<150) return QQ2HLNU_PTV_75_150;
499  else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
500  // 150 < pTV/GeV < 250
502  }
503  // 4. qq->ZH->llH categories
504  else if (prodMode==HTXS::QQ2ZH) {
505  if (fwdHiggs) return QQ2HLL_FWDH;
506  else if (V.pt()<75) return QQ2HLL_PTV_0_75;
507  else if (V.pt()<150) return QQ2HLL_PTV_75_150;
508  else if (V.pt()>250) return QQ2HLL_PTV_GT250;
509  // 150 < pTV/GeV < 250
511  }
512  // 5. gg->ZH->llH categories
513  else if (prodMode==HTXS::GG2ZH ) {
514  if (fwdHiggs) return GG2HLL_FWDH;
515  else if (V.pt()<75) return GG2HLL_PTV_0_75;
516  else if (V.pt()<150) return GG2HLL_PTV_75_150;
517  else if (V.pt()>250) return GG2HLL_PTV_GT250;
519  }
520  // 6.ttH,bbH,tH categories
521  else if (prodMode==HTXS::TTH) {
522  if (fwdHiggs) return TTH_FWDH;
523  else return Category(TTH_PTH_0_60+getBin(higgs.pt(),{0,60,120,200,300}));
524  }
525  else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
526  else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
527  return UNKNOWN;
528  }
529 
532  const Particle &higgs,
533  const Jets &jets,
534  const Particle &V) const {
535  using namespace HTXS::Stage1_2_Fine;
536  int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
537  int vbfTopo = vbfTopology_Stage1_2_Fine(jets,higgs);
538 
539  // 1. GGF Stage 1.2 categories
540  // Following YR4 write-up: XXXXX
541  if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
542  if (fwdHiggs) return GG2H_FWDH;
543  if ( higgs.pt()>200 ){
544  if (Njets>0){
545  double pTHj = (jets[0].momentum()+higgs.momentum()).pt();
546  if( pTHj/higgs.pt()>0.15 ) return Category(GG2H_PTH_200_300_PTHJoverPTH_GT15+getBin(higgs.pt(),{200,300,450,650}));
547  else return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15+getBin(higgs.pt(),{200,300,450,650}));
548  }
549  else return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15+getBin(higgs.pt(),{200,300,450,650}));
550  }
551  if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
552  if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
553  if (Njets>1){
554  //double mjj = (jets[0].mom()+jets[1].mom()).mass();
555  double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
556  //VBF topology
557  if(vbfTopo) return Category(GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25+vbfTopo-1);
558  //Njets >= 2jets without VBF topology (mjj<350)
559  if (pTHjj<25) return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25+getBin(higgs.pt(),{0,60,120,200}));
560  else return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25+getBin(higgs.pt(),{0,60,120,200}));
561  }
562  }
563 
564  // 2. Electroweak qq->Hqq Stage 1.2 categories
565  else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
566  if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
567  int Njets=jets.size();
568  if (Njets==0) return QQ2HQQ_0J;
569  else if (Njets==1) return QQ2HQQ_1J;
570  else if (Njets>=2) {
571  double mjj = (jets[0].mom()+jets[1].mom()).mass();
572  double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
573  if (mjj<350){
574  if (pTHjj<25) return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25+getBin(mjj,{0,60,120,350}));
575  else return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25+getBin(mjj,{0,60,120,350}));
576  } else { //mjj>350 GeV
577  if (higgs.pt()<200){
579  } else {
581  }
582  }
583  }
584  }
585  // 3. WH->Hlv categories
586  else if (prodMode==HTXS::WH) {
587  if (fwdHiggs) return QQ2HLNU_FWDH;
588  int Njets=jets.size();
589  if (Njets==0) return Category(QQ2HLNU_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
590  if (Njets==1) return Category(QQ2HLNU_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
591  return Category(QQ2HLNU_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
592  }
593  // 4. qq->ZH->llH categories
594  else if (prodMode==HTXS::QQ2ZH) {
595  if (fwdHiggs) return QQ2HLL_FWDH;
596  int Njets=jets.size();
597  if (Njets==0) return Category(QQ2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
598  if (Njets==1) return Category(QQ2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
599  return Category(QQ2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
600  }
601  // 5. gg->ZH->llH categories
602  else if (prodMode==HTXS::GG2ZH ) {
603  if (fwdHiggs) return GG2HLL_FWDH;
604  int Njets=jets.size();
605  if (Njets==0) return Category(GG2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
606  if (Njets==1) return Category(GG2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
607  return Category(GG2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
608  }
609  // 6.ttH,bbH,tH categories
610  else if (prodMode==HTXS::TTH) {
611  if (fwdHiggs) return TTH_FWDH;
612  else return Category(TTH_PTH_0_60+getBin(higgs.pt(),{0,60,120,200,300,450}));
613  }
614  else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
615  else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
616  return UNKNOWN;
617  }
618 
620 
621 
624 
626  void setHiggsProdMode( HTXS::HiggsProdMode prodMode ){ m_HiggsProdMode = prodMode; }
627 
631  void init() {
632  printf("==============================================================\n");
633  printf("======== HiggsTemplateCrossSections Initialization =========\n");
634  printf("==============================================================\n");
635  // check that the production mode has been set
636  // if running in standalone Rivet the production mode is set through an env variable
638  char *pm_env = getenv("HIGGSPRODMODE");
639  string pm(pm_env==nullptr?"":pm_env);
640  if ( pm == "GGF" ) m_HiggsProdMode = HTXS::GGF;
641  else if ( pm == "VBF" ) m_HiggsProdMode = HTXS::VBF;
642  else if ( pm == "WH" ) m_HiggsProdMode = HTXS::WH;
643  else if ( pm == "ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
644  else if ( pm == "QQ2ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
645  else if ( pm == "GG2ZH" ) m_HiggsProdMode = HTXS::GG2ZH;
646  else if ( pm == "TTH" ) m_HiggsProdMode = HTXS::TTH;
647  else if ( pm == "BBH" ) m_HiggsProdMode = HTXS::BBH;
648  else if ( pm == "TH" ) m_HiggsProdMode = HTXS::TH;
649  else {
650  MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
651  }
652  }
653 
654  // Projections for final state particles
655  const FinalState FS;
656  declare(FS,"FS");
657 
658  // initialize the histograms with for each of the stages
660  m_sumw = 0.0;
661  printf("==============================================================\n");
662  printf("======== Higgs prod mode %d =========\n",m_HiggsProdMode);
663  printf("======== Sucessful Initialization =========\n");
664  printf("==============================================================\n");
665  }
666 
667  // Perform the per-event analysis
668  void analyze(const Event& event) {
669 
670  // get the classification
671  HiggsClassification cat = classifyEvent(event,m_HiggsProdMode);
672 
673  // Fill histograms: categorization --> linerize the categories
674  const double weight = 1.; // Event weights are now all 1 in Rivet
675  m_sumw += weight;
676 
677  int F=cat.stage0_cat%10, P=cat.stage1_cat_pTjet30GeV/100;
678  m_hist_stage0->fill( cat.stage0_cat/10*2+F, weight );
679 
680  // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
681  static const vector<int> offset({0,1,13,19,24,29,33,35,37,39});
682  int off = offset[P];
683  // 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
684  static const vector<int> offset1_2({0,1,18,29,35,41,47,53,55,57});
685  int off1_2 = offset1_2[P];
686  // 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
687  static const vector<int> offset1_2_Fine({0,1,29,54,70,86,102,109,111,113});
688  int off1_2_Fine = offset1_2_Fine[P];
689 
690  m_hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV%100 + off, weight);
691  m_hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV%100 + off, weight);
692  m_hist_stage1_2_pTjet25->fill(cat.stage1_2_cat_pTjet25GeV%100 + off1_2, weight);
693  m_hist_stage1_2_pTjet30->fill(cat.stage1_2_cat_pTjet30GeV%100 + off1_2, weight);
694  m_hist_stage1_2_fine_pTjet25->fill(cat.stage1_2_fine_cat_pTjet25GeV%100 + off1_2_Fine, weight);
695  m_hist_stage1_2_fine_pTjet30->fill(cat.stage1_2_fine_cat_pTjet30GeV%100 + off1_2_Fine, weight);
696 
697  // Fill histograms: variables used in the categorization
698  m_hist_pT_Higgs->fill(cat.higgs.pT(),weight);
699  m_hist_y_Higgs->fill(cat.higgs.rapidity(),weight);
700  m_hist_pT_V->fill(cat.V.pT(),weight);
701 
702  m_hist_Njets25->fill(cat.jets25.size(),weight);
703  m_hist_Njets30->fill(cat.jets30.size(),weight);
704 
705  m_hist_isZ2vv->fill(cat.isZ2vvDecay, weight);
706 
707  // Jet variables. Use jet collection with pT threshold at 30 GeV
708  if (cat.jets30.size()) m_hist_pT_jet1->fill(cat.jets30[0].pt(),weight);
709  if (cat.jets30.size()>=2) {
710  const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
711  m_hist_deltay_jj->fill(std::abs(j1.rapidity()-j2.rapidity()),weight);
712  m_hist_dijet_mass->fill((j1+j2).mass(),weight);
713  m_hist_pT_Hjj->fill((j1+j2+cat.higgs.momentum()).pt(),weight);
714  }
715  }
716 
718  MSG_INFO (" ====================================================== ");
719  MSG_INFO (" Higgs Template X-Sec Categorization Tool ");
720  MSG_INFO (" Status Code Summary ");
721  MSG_INFO (" ====================================================== ");
722  bool allSuccess = (numEvents()==m_errorCount[HTXS::SUCCESS]);
723  if ( allSuccess ) MSG_INFO (" >>>> All "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized!");
724  else{
725  MSG_INFO (" >>>> "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized");
726  MSG_INFO (" >>>> --> the following errors occured:");
727  MSG_INFO (" >>>> "<< m_errorCount[HTXS::PRODMODE_DEFINED] <<" had an undefined Higgs production mode.");
728  MSG_INFO (" >>>> "<< m_errorCount[HTXS::MOMENTUM_CONSERVATION] <<" failed momentum conservation.");
729  MSG_INFO (" >>>> "<< m_errorCount[HTXS::HIGGS_IDENTIFICATION] <<" failed to identify a valid Higgs boson.");
730  MSG_INFO (" >>>> "<< m_errorCount[HTXS::HS_VTX_IDENTIFICATION] <<" failed to identify the hard scatter vertex.");
731  MSG_INFO (" >>>> "<< m_errorCount[HTXS::VH_IDENTIFICATION] <<" VH: to identify a valid V-boson.");
732  MSG_INFO (" >>>> "<< m_errorCount[HTXS::TOP_W_IDENTIFICATION] <<" failed to identify valid Ws from top decay.");
733  }
734  MSG_INFO (" ====================================================== ");
735  MSG_INFO (" ====================================================== ");
736  }
737 
738 
739  void finalize() {
741  double sf = m_sumw>0?1.0/m_sumw:1.0;
744  scale(hist, sf);
745  }
746 
747  /*
748  * initialize histograms
749  */
750 
752  book(m_hist_stage0,"HTXS_stage0",20,0,20);
753  book(m_hist_stage1_pTjet25,"HTXS_stage1_pTjet25",40,0,40);
754  book(m_hist_stage1_pTjet30,"HTXS_stage1_pTjet30",40,0,40);
755  book(m_hist_stage1_2_pTjet25,"HTXS_stage1_2_pTjet25",57,0,57);
756  book(m_hist_stage1_2_pTjet30,"HTXS_stage1_2_pTjet30",57,0,57);
757  book(m_hist_stage1_2_fine_pTjet25,"HTXS_stage1_2_fine_pTjet25",113,0,113);
758  book(m_hist_stage1_2_fine_pTjet30,"HTXS_stage1_2_fine_pTjet30",113,0,113);
759  book(m_hist_pT_Higgs,"pT_Higgs",80,0,400);
760  book(m_hist_y_Higgs,"y_Higgs",80,-4,4);
761  book(m_hist_pT_V,"pT_V",80,0,400);
762  book(m_hist_pT_jet1,"pT_jet1",80,0,400);
763  book(m_hist_deltay_jj ,"deltay_jj",50,0,10);
764  book(m_hist_dijet_mass,"m_jj",50,0,2000);
765  book(m_hist_pT_Hjj,"pT_Hjj",50,0,250);
766  book(m_hist_Njets25,"Njets25",10,0,10);
767  book(m_hist_Njets30,"Njets30",10,0,10);
768  book(m_hist_isZ2vv,"isZ2vv",2,0,2);
769  }
771 
772  /*
773  * initialize private members used in the classification procedure
774  */
775 
776  private:
777  double m_sumw=0.0;
779  mutable std::array<std::atomic<size_t>, HTXS::NUM_ERRORCODES> m_errorCount ATLAS_THREAD_SAFE {};
780  Histo1DPtr m_hist_stage0;
788  Histo1DPtr m_hist_isZ2vv;
789  };
790 
791  // the PLUGIN only needs to be decleared when running standalone Rivet
792  // and causes compilation / linking issues if included in Athena / RootCore
793  //check for Rivet environment variable RIVET_ANALYSIS_PATH
794 #ifdef RIVET_ANALYSIS_PATH
795  // The hook for the plugin system
796  DECLARE_RIVET_PLUGIN(HiggsTemplateCrossSections);
797 #endif
798 
799 }
800 
801 #endif
Rivet::HiggsTemplateCrossSections::printClassificationSummary
void printClassificationSummary()
Definition: HiggsTemplateCrossSections.h:717
HTXS::VH_IDENTIFICATION
@ VH_IDENTIFICATION
failed to identify associated vector boson
Definition: HiggsTemplateCrossSectionsDefs.h:20
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
HTXS::Stage1::QQ2HQQ_REST
@ QQ2HQQ_REST
Definition: HiggsTemplateCrossSectionsDefs.h:67
HTXS::Stage0::Category
Category
Definition: HiggsTemplateCrossSectionsDefs.h:43
HTXS::Stage1_2::GG2HLL_PTV_GT250
@ GG2HLL_PTV_GT250
Definition: HiggsTemplateCrossSectionsDefs.h:151
HTXS::TOP_W_IDENTIFICATION
@ TOP_W_IDENTIFICATION
failed to identify top decay
Definition: HiggsTemplateCrossSectionsDefs.h:22
HTXS::Stage1_2::QQ2HLL_PTV_0_75
@ QQ2HLL_PTV_0_75
Definition: HiggsTemplateCrossSectionsDefs.h:140
CalculateHighPtTerm.pT
pT
Definition: ICHEP2016/CalculateHighPtTerm.py:57
MSG_WARNING
#define MSG_WARNING(ARG)
Definition: testEGChargeIDSelector.cxx:48
Rivet::HiggsTemplateCrossSections::quarkDecay
bool quarkDecay(const Particle &p) const
Return true is particle decays to quarks.
Definition: HiggsTemplateCrossSections.h:89
HTXS::Stage1::QQ2HQQ_VBFTOPO_JET3
@ QQ2HQQ_VBFTOPO_JET3
Definition: HiggsTemplateCrossSectionsDefs.h:66
HTXS::Stage0::QQ2HLL_FWDH
@ QQ2HLL_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:45
GenEvent.h
HTXS::Stage1_2
Categorization Stage 1.2: Three digit integer of format PF Where P is a digit representing the proces...
Definition: HiggsTemplateCrossSectionsDefs.h:98
MSG_INFO
#define MSG_INFO(ARG)
Definition: testEGChargeIDSelector.cxx:47
HTXS::Stage1_2::QQ2HQQ_0J
@ QQ2HQQ_0J
Definition: HiggsTemplateCrossSectionsDefs.h:121
HTXS::GG2ZH
@ GG2ZH
Definition: HiggsTemplateCrossSectionsDefs.h:29
HTXS::UNKNOWN
@ UNKNOWN
Definition: HiggsTemplateCrossSectionsDefs.h:28
HiggsTemplateCrossSectionsDefs.h
Rivet::HiggsTemplateCrossSections::isVH
bool isVH(HTXS::HiggsProdMode p) const
Whether the Higgs is produced in association with a vector boson (VH)
Definition: HiggsTemplateCrossSections.h:369
HTXS::TTH
@ TTH
Definition: HiggsTemplateCrossSectionsDefs.h:30
HTXS::PRODMODE_DEFINED
@ PRODMODE_DEFINED
production mode not defined
Definition: HiggsTemplateCrossSectionsDefs.h:15
HTXS::Stage1::QQ2HQQ_VBFTOPO_JET3VETO
@ QQ2HQQ_VBFTOPO_JET3VETO
Definition: HiggsTemplateCrossSectionsDefs.h:66
HTXS::Stage1_2_Fine::GG2HLL_PTV_0_75_0J
@ GG2HLL_PTV_0_75_0J
Definition: HiggsTemplateCrossSectionsDefs.h:260
HTXS::Stage1_2::QQ2HLNU_PTV_75_150
@ QQ2HLNU_PTV_75_150
Definition: HiggsTemplateCrossSectionsDefs.h:134
Base_Fragment.mass
mass
Definition: Sherpa_i/share/common/Base_Fragment.py:59
Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_fine_pTjet30
Histo1DPtr m_hist_stage1_2_fine_pTjet30
Definition: HiggsTemplateCrossSections.h:783
HTXS::Stage1_2_Fine::GG2HLL_PTV_0_75_GE2J
@ GG2HLL_PTV_0_75_GE2J
Definition: HiggsTemplateCrossSectionsDefs.h:270
HTXS::Stage1::QQ2HLNU_PTV_150_250_GE1J
@ QQ2HLNU_PTV_150_250_GE1J
Definition: HiggsTemplateCrossSectionsDefs.h:72
HTXS::Stage1_2_Fine::GG2H_PTH_200_300_PTHJoverPTH_0_15
@ GG2H_PTH_200_300_PTHJoverPTH_0_15
Definition: HiggsTemplateCrossSectionsDefs.h:171
Rivet
Definition: HiggsTemplateCrossSections.h:27
JetTiledMap::W
@ W
Definition: TiledEtaPhiMap.h:44
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200
@ QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200
Definition: HiggsTemplateCrossSectionsDefs.h:126
HTXS::Stage1_2_Fine::QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25
@ QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25
Definition: HiggsTemplateCrossSectionsDefs.h:205
DMTest::P
P_v1 P
Definition: P.h:23
HTXS::GGF
@ GGF
Definition: HiggsTemplateCrossSectionsDefs.h:29
plotmaker.hist
hist
Definition: plotmaker.py:148
Event
Definition: trigbs_orderedMerge.cxx:42
Rivet::HiggsTemplateCrossSections::m_hist_Njets25
Histo1DPtr m_hist_Njets25
Definition: HiggsTemplateCrossSections.h:787
Jets
Definition: Jets.py:1
Rivet::HiggsTemplateCrossSections::m_hist_Njets30
Histo1DPtr m_hist_Njets30
Definition: HiggsTemplateCrossSections.h:787
GenVertex.h
HTXS::Stage0
Two digit number of format PF P is digit for the physics process and F is 0 for |yH|>2....
Definition: HiggsTemplateCrossSectionsDefs.h:41
Rivet::HiggsTemplateCrossSections
Rivet routine for classifying MC events according to the Higgs template cross section categories.
Definition: HiggsTemplateCrossSections.h:33
HTXS::Stage1_2::QQ2HLNU_PTV_0_75
@ QQ2HLNU_PTV_0_75
Definition: HiggsTemplateCrossSectionsDefs.h:133
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_60_120
@ QQ2HQQ_GE2J_MJJ_60_120
Definition: HiggsTemplateCrossSectionsDefs.h:124
HTXS::Stage1::QQ2HLL_PTV_150_250_GE1J
@ QQ2HLL_PTV_150_250_GE1J
Definition: HiggsTemplateCrossSectionsDefs.h:78
test_pyathena.pt
pt
Definition: test_pyathena.py:11
HTXS::Stage1::GG2H_GE2J_PTH_0_60
@ GG2H_GE2J_PTH_0_60
Definition: HiggsTemplateCrossSectionsDefs.h:62
HTXS::Stage1_2::GG2H_PTH_200_300
@ GG2H_PTH_200_300
Definition: HiggsTemplateCrossSectionsDefs.h:103
HTXS::Stage1_2::GG2HLL_PTV_75_150
@ GG2HLL_PTV_75_150
Definition: HiggsTemplateCrossSectionsDefs.h:148
HTXS::Stage0::VH2HQQ_FWDH
@ VH2HQQ_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:44
Rivet::HiggsTemplateCrossSections::m_hist_stage1_pTjet30
Histo1DPtr m_hist_stage1_pTjet30
Definition: HiggsTemplateCrossSections.h:781
Rivet::HiggsTemplateCrossSections::m_hist_dijet_mass
Histo1DPtr m_hist_dijet_mass
Definition: HiggsTemplateCrossSections.h:786
Rivet::HiggsTemplateCrossSections::error
HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err, std::string msg="", int NmaxWarnings=20) const
Returns the classification object with the error code set.
Definition: HiggsTemplateCrossSections.h:110
Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_fine_pTjet25
Histo1DPtr m_hist_stage1_2_fine_pTjet25
Definition: HiggsTemplateCrossSections.h:783
Rivet::HiggsTemplateCrossSections::vbfTopology_Stage1_2
int vbfTopology_Stage1_2(const Jets &jets, const Particle &higgs) const
VBF topology selection 0 = fail loose selection: m_jj > 350 GeV 1 pass loose, but fail additional cut...
Definition: HiggsTemplateCrossSections.h:343
yodamerge_tmp.scale
scale
Definition: yodamerge_tmp.py:138
Rivet::HiggsTemplateCrossSections::m_hist_pT_Hjj
Histo1DPtr m_hist_pT_Hjj
Definition: HiggsTemplateCrossSections.h:786
x
#define x
HTXS::HIGGS_DECAY_IDENTIFICATION
@ HIGGS_DECAY_IDENTIFICATION
failed to identify Higgs boson decay products
Definition: HiggsTemplateCrossSectionsDefs.h:18
GenParticle.h
Rivet::HiggsTemplateCrossSections::hasParent
bool hasParent(HepMC::ConstGenParticlePtr ptcl, int pdgID) const
Checks whether the input particle has a parent with a given PDGID.
Definition: HiggsTemplateCrossSections.h:82
xAOD::Particle
Particle_v1 Particle
Define the latest version of the particle class.
Definition: Event/xAOD/xAODParticleEvent/xAODParticleEvent/Particle.h:17
HTXS::Stage1_2::QQ2HQQ_1J
@ QQ2HQQ_1J
Definition: HiggsTemplateCrossSectionsDefs.h:122
HTXS::Stage1::QQ2HQQ_PTJET1_GT200
@ QQ2HQQ_PTJET1_GT200
Definition: HiggsTemplateCrossSectionsDefs.h:67
Rivet::HiggsTemplateCrossSections::classifyEvent
HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode) const
Main classificaion method.
Definition: HiggsTemplateCrossSections.h:126
HTXS::Stage1::QQ2HLL_PTV_150_250_0J
@ QQ2HLL_PTV_150_250_0J
Definition: HiggsTemplateCrossSectionsDefs.h:77
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_120_350
@ QQ2HQQ_GE2J_MJJ_120_350
Definition: HiggsTemplateCrossSectionsDefs.h:125
HTXS::Stage1::GG2HLL_PTV_GT150_0J
@ GG2HLL_PTV_GT150_0J
Definition: HiggsTemplateCrossSectionsDefs.h:83
isHiggs
bool isHiggs(const T &p)
APID: HIGGS boson is only one particle.
Definition: AtlasPID.h:199
dqt_zlumi_pandas.weight
int weight
Definition: dqt_zlumi_pandas.py:189
HTXS::Stage1_2::QQ2HLL_PTV_75_150
@ QQ2HLL_PTV_75_150
Definition: HiggsTemplateCrossSectionsDefs.h:141
HTXS::Stage1_2::GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25
@ GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25
Definition: HiggsTemplateCrossSectionsDefs.h:115
HTXS::Stage1
Categorization Stage 1: Three digit integer of format PF Where P is a digit representing the process ...
Definition: HiggsTemplateCrossSectionsDefs.h:53
isQuark
bool isQuark(const T &p)
PDG rule 2: Quarks and leptons are numbered consecutively starting from 1 and 11 respectively; to dot...
Definition: AtlasPID.h:122
TRTCalib_cfilter.p2
p2
Definition: TRTCalib_cfilter.py:131
HTXS::ErrorCode
ErrorCode
Error code: whether the classification was successful or failed.
Definition: HiggsTemplateCrossSectionsDefs.h:12
Rivet::HiggsTemplateCrossSections::m_hist_stage0
Histo1DPtr m_hist_stage0
Definition: HiggsTemplateCrossSections.h:780
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:210
HTXS::Stage1_2_Fine::GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25
@ GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25
Definition: HiggsTemplateCrossSectionsDefs.h:187
HTXS::Stage0::TH_FWDH
@ TH_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:46
convertTimingResiduals.sum
sum
Definition: convertTimingResiduals.py:55
HTXS::Stage1_2_Fine::Category
Category
Definition: HiggsTemplateCrossSectionsDefs.h:167
doubleTestComp.j1
j1
Definition: doubleTestComp.py:21
HTXS::Stage1_2::GG2HLL_PTV_0_75
@ GG2HLL_PTV_0_75
Definition: HiggsTemplateCrossSectionsDefs.h:147
event
POOL::TEvent event(POOL::TEvent::kClassAccess)
Rivet::HiggsTemplateCrossSections::m_hist_deltay_jj
Histo1DPtr m_hist_deltay_jj
Definition: HiggsTemplateCrossSections.h:786
dqt_zlumi_pandas.err
err
Definition: dqt_zlumi_pandas.py:182
HTXS::Stage1_2::UNKNOWN
@ UNKNOWN
Definition: HiggsTemplateCrossSectionsDefs.h:100
Rivet::HiggsTemplateCrossSections::m_hist_isZ2vv
Histo1DPtr m_hist_isZ2vv
Definition: HiggsTemplateCrossSections.h:788
lumiFormat.i
int i
Definition: lumiFormat.py:85
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_0_60
@ QQ2HQQ_GE2J_MJJ_0_60
Definition: HiggsTemplateCrossSectionsDefs.h:123
vector< int >
HTXS::VH_DECAY_IDENTIFICATION
@ VH_DECAY_IDENTIFICATION
failed to identify associated vector boson decay products
Definition: HiggsTemplateCrossSectionsDefs.h:21
isZ
bool isZ(const T &p)
Definition: AtlasPID.h:173
HTXS::Stage1_2::GG2H_0J_PTH_0_10
@ GG2H_0J_PTH_0_10
Definition: HiggsTemplateCrossSectionsDefs.h:107
HTXS::Stage1::GG2H_VBFTOPO_JET3VETO
@ GG2H_VBFTOPO_JET3VETO
Definition: HiggsTemplateCrossSectionsDefs.h:58
HTXS::Stage1::UNKNOWN
@ UNKNOWN
Definition: HiggsTemplateCrossSectionsDefs.h:55
HTXS::Stage1::QQ2HLL_PTV_GT250
@ QQ2HLL_PTV_GT250
Definition: HiggsTemplateCrossSectionsDefs.h:79
HTXS::Stage1::GG2H_0J
@ GG2H_0J
Definition: HiggsTemplateCrossSectionsDefs.h:59
makePlot.Particles
Particles
Definition: makePlot.py:25
test_pyathena.parent
parent
Definition: test_pyathena.py:15
Rivet::HiggsTemplateCrossSections::m_hist_pT_V
Histo1DPtr m_hist_pT_V
Definition: HiggsTemplateCrossSections.h:785
HTXS::Stage1_2::GG2HLL_PTV_150_250_0J
@ GG2HLL_PTV_150_250_0J
Definition: HiggsTemplateCrossSectionsDefs.h:149
HTXS::Stage1_2::GG2H_GE2J_MJJ_0_350_PTH_0_60
@ GG2H_GE2J_MJJ_0_350_PTH_0_60
Definition: HiggsTemplateCrossSectionsDefs.h:112
HTXS::Stage1::GG2HLL_PTV_GT150_GE1J
@ GG2HLL_PTV_GT150_GE1J
Definition: HiggsTemplateCrossSectionsDefs.h:84
Rivet::HiggsTemplateCrossSections::finalize
void finalize()
Definition: HiggsTemplateCrossSections.h:739
HTXS::Stage0::BBH_FWDH
@ BBH_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:46
Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_pTjet25
Histo1DPtr m_hist_stage1_2_pTjet25
Definition: HiggsTemplateCrossSections.h:782
Rivet::HiggsTemplateCrossSections::ChLeptonDecay
bool ChLeptonDecay(const Particle &p) const
Return true if particle decays to charged leptons.
Definition: HiggsTemplateCrossSections.h:96
HTXS::Stage1_2_Fine::QQ2HLL_PTV_0_75_0J
@ QQ2HLL_PTV_0_75_0J
Definition: HiggsTemplateCrossSectionsDefs.h:243
plotting.yearwise_luminosity_vs_mu.bins
bins
Definition: yearwise_luminosity_vs_mu.py:30
HTXS::Stage0::UNKNOWN
@ UNKNOWN
Definition: HiggsTemplateCrossSectionsDefs.h:44
Rivet::HiggsTemplateCrossSections::m_sumw
double m_sumw
Definition: HiggsTemplateCrossSections.h:777
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) const
Stage-1.2 categorization.
Definition: HiggsTemplateCrossSections.h:454
HTXS::Stage1_2_Fine::QQ2HLNU_PTV_0_75_0J
@ QQ2HLNU_PTV_0_75_0J
Definition: HiggsTemplateCrossSectionsDefs.h:226
Rivet::HiggsTemplateCrossSections::m_HiggsProdMode
HTXS::HiggsProdMode m_HiggsProdMode
Definition: HiggsTemplateCrossSections.h:778
HTXS::UNDEFINED
@ UNDEFINED
Definition: HiggsTemplateCrossSectionsDefs.h:13
isChLepton
bool isChLepton(const T &p)
APID: the fourth generation leptons are leptons.
Definition: AtlasPID.h:148
Rivet::HiggsTemplateCrossSections::getLastInstance
Particle getLastInstance(const Particle &ptcl) const
follow a "propagating" particle and return its last instance
Definition: HiggsTemplateCrossSections.h:48
Rivet::HiggsTemplateCrossSections::getStage1Category
HTXS::Stage1::Category getStage1Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V) const
Stage-1 categorization.
Definition: HiggsTemplateCrossSections.h:388
HTXS
Higgs Template Cross Section namespace.
Definition: IHiggsTruthCategoryTool.h:17
HTXS::MOMENTUM_CONSERVATION
@ MOMENTUM_CONSERVATION
failed momentum conservation
Definition: HiggsTemplateCrossSectionsDefs.h:16
HTXS::Stage1_2_Fine::QQ2HQQ_GE2J_MJJ_350_700_PTH_GT200_PTHJJ_0_25
@ QQ2HQQ_GE2J_MJJ_350_700_PTH_GT200_PTHJJ_0_25
Definition: HiggsTemplateCrossSectionsDefs.h:216
Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_pTjet30
Histo1DPtr m_hist_stage1_2_pTjet30
Definition: HiggsTemplateCrossSections.h:782
HTXS::Stage1::QQ2HLNU_PTV_150_250_0J
@ QQ2HLNU_PTV_150_250_0J
Definition: HiggsTemplateCrossSectionsDefs.h:71
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
HTXS::HiggsProdMode
HiggsProdMode
Higgs production modes, corresponding to input sample.
Definition: HiggsTemplateCrossSectionsDefs.h:27
Analysis
The namespace of all packages in PhysicsAnalysis/JetTagging.
Definition: BTaggingCnvAlg.h:20
HTXS::Stage0::TTH_FWDH
@ TTH_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:46
Rivet::HiggsTemplateCrossSections::m_hist_y_Higgs
Histo1DPtr m_hist_y_Higgs
Definition: HiggsTemplateCrossSections.h:784
HTXS::Stage1_2_Fine
Definition: HiggsTemplateCrossSectionsDefs.h:166
HTXS::Stage1::QQ2HQQ_VH2JET
@ QQ2HQQ_VH2JET
Definition: HiggsTemplateCrossSectionsDefs.h:67
ActsTrk::to_string
std::string to_string(const DetectorType &type)
Definition: GeometryDefs.h:34
HTXS::Stage1::Category
Category
Definition: HiggsTemplateCrossSectionsDefs.h:54
HTXS::Stage0::VH2HQQ
@ VH2HQQ
Definition: HiggsTemplateCrossSectionsDefs.h:44
HTXS::HIGGS_IDENTIFICATION
@ HIGGS_IDENTIFICATION
failed to identify Higgs boson
Definition: HiggsTemplateCrossSectionsDefs.h:17
Rivet::HiggsTemplateCrossSections::hasChild
bool hasChild(HepMC::ConstGenParticlePtr ptcl, int pdgID) const
Checks whether the input particle has a child with a given PDGID.
Definition: HiggsTemplateCrossSections.h:75
Rivet::HiggsTemplateCrossSections::getStage0Category
HTXS::Stage0::Category getStage0Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Particle &V) const
Stage-0 HTXS categorization.
Definition: HiggsTemplateCrossSections.h:372
HTXS::Stage1_2::GG2HLL_PTV_150_250_GE1J
@ GG2HLL_PTV_150_250_GE1J
Definition: HiggsTemplateCrossSectionsDefs.h:150
HTXS::Stage0::GG2H_FWDH
@ GG2H_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:44
HTXS::Stage1_2_Fine::GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25
@ GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25
Definition: HiggsTemplateCrossSectionsDefs.h:184
HTXS::Stage1_2_Fine::QQ2HLNU_PTV_0_75_1J
@ QQ2HLNU_PTV_0_75_1J
Definition: HiggsTemplateCrossSectionsDefs.h:231
HTXS::Stage1_2::Category
Category
Definition: HiggsTemplateCrossSectionsDefs.h:99
Rivet::HiggsTemplateCrossSections::HiggsTemplateCrossSections
HiggsTemplateCrossSections()
Definition: HiggsTemplateCrossSections.h:36
Rivet::HiggsTemplateCrossSections::originateFrom
bool originateFrom(const Particle &p, const Particle &p2) const
Whether particle p originates from p2.
Definition: HiggsTemplateCrossSections.h:70
isW
bool isW(const T &p)
Definition: AtlasPID.h:176
HTXS::Stage1_2_Fine::GG2HLL_PTV_0_75_1J
@ GG2HLL_PTV_0_75_1J
Definition: HiggsTemplateCrossSectionsDefs.h:265
SCT_ConditionsAlgorithms::CoveritySafe::getenv
std::string getenv(const std::string &variableName)
get an environment variable
Definition: SCT_ConditionsUtilities.cxx:17
isTop
bool isTop(const T &p)
Definition: AtlasPID.h:135
Rivet::HiggsTemplateCrossSections::analyze
void analyze(const Event &event)
Definition: HiggsTemplateCrossSections.h:668
HTXS::Stage1_2::GG2H_0J_PTH_GT10
@ GG2H_0J_PTH_GT10
Definition: HiggsTemplateCrossSectionsDefs.h:108
mapkey::sf
@ sf
Definition: TElectronEfficiencyCorrectionTool.cxx:38
HTXS::Stage1::QQ2HLNU_PTV_0_150
@ QQ2HLNU_PTV_0_150
Definition: HiggsTemplateCrossSectionsDefs.h:70
HTXS::Stage1::QQ2HLNU_PTV_GT250
@ QQ2HLNU_PTV_GT250
Definition: HiggsTemplateCrossSectionsDefs.h:73
HTXS::Stage1::GG2H_1J_PTH_0_60
@ GG2H_1J_PTH_0_60
Definition: HiggsTemplateCrossSectionsDefs.h:60
F
#define F(x, y, z)
Definition: MD5.cxx:112
Relatives.h
HTXS::Stage1_2_Fine::QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25
@ QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25
Definition: HiggsTemplateCrossSectionsDefs.h:202
LArG4FSStartPointFilter.particles
list particles
Definition: LArG4FSStartPointFilter.py:84
Rivet::HiggsTemplateCrossSections::setHiggsProdMode
void setHiggsProdMode(HTXS::HiggsProdMode prodMode)
Sets the Higgs production mode.
Definition: HiggsTemplateCrossSections.h:626
convertTimingResiduals.offset
offset
Definition: convertTimingResiduals.py:71
HTXS::Stage1_2_Fine::QQ2HLL_PTV_0_75_GE2J
@ QQ2HLL_PTV_0_75_GE2J
Definition: HiggsTemplateCrossSectionsDefs.h:253
HTXS::Stage1_2_Fine::GG2H_PTH_200_300_PTHJoverPTH_GT15
@ GG2H_PTH_200_300_PTHJoverPTH_GT15
Definition: HiggsTemplateCrossSectionsDefs.h:175
HTXS::Stage1_2::QQ2HQQ_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25
@ QQ2HQQ_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25
Definition: HiggsTemplateCrossSectionsDefs.h:127
GenRanges.h
defineDB.jets
list jets
Definition: JetTagCalibration/share/defineDB.py:24
Rivet::HiggsTemplateCrossSections::ATLAS_THREAD_SAFE
std::array< std::atomic< size_t >, HTXS::NUM_ERRORCODES > m_errorCount ATLAS_THREAD_SAFE
Definition: HiggsTemplateCrossSections.h:779
HTXS::HS_VTX_IDENTIFICATION
@ HS_VTX_IDENTIFICATION
failed to identify hard scatter vertex
Definition: HiggsTemplateCrossSectionsDefs.h:19
top
@ top
Definition: TruthClasses.h:64
Rivet::HiggsTemplateCrossSections::init
void init()
default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Hig...
Definition: HiggsTemplateCrossSections.h:631
HTXS::Stage0::QQ2HLNU_FWDH
@ QQ2HLNU_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:45
Rivet::HiggsTemplateCrossSections::m_hist_pT_jet1
Histo1DPtr m_hist_pT_jet1
Definition: HiggsTemplateCrossSections.h:785
HTXS::TH
@ TH
Definition: HiggsTemplateCrossSectionsDefs.h:30
HTXS::Stage1_2_Fine::QQ2HLNU_PTV_0_75_GE2J
@ QQ2HLNU_PTV_0_75_GE2J
Definition: HiggsTemplateCrossSectionsDefs.h:236
Rivet::HiggsTemplateCrossSections::m_hist_stage1_pTjet25
Histo1DPtr m_hist_stage1_pTjet25
Definition: HiggsTemplateCrossSections.h:781
checker_macros.h
Define macros for attributes used to control the static checker.
HTXS::SUCCESS
@ SUCCESS
successful classification
Definition: HiggsTemplateCrossSectionsDefs.h:14
Rivet::HiggsTemplateCrossSections::vbfTopology
int vbfTopology(const Jets &jets, const Particle &higgs) const
VBF topolog selection 0 = fail loose selction: m_jj > 400 GeV and Dy_jj > 2.8 1 pass loose,...
Definition: HiggsTemplateCrossSections.h:333
HTXS::Stage1::QQ2HLL_PTV_0_150
@ QQ2HLL_PTV_0_150
Definition: HiggsTemplateCrossSectionsDefs.h:76
Rivet::HiggsTemplateCrossSections::originateFrom
bool originateFrom(const Particle &p, const Particles &ptcls) const
Whether particle p originate from any of the ptcls.
Definition: HiggsTemplateCrossSections.h:57
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) const
Stage-1.2_Fine categorization.
Definition: HiggsTemplateCrossSections.h:531
HTXS::Stage0::GG2HLL_FWDH
@ GG2HLL_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:45
Rivet::HiggsTemplateCrossSections::vbfTopology_Stage1_2_Fine
int vbfTopology_Stage1_2_Fine(const Jets &jets, const Particle &higgs) const
VBF topology selection for Stage1_2 0 = fail loose selection: m_jj > 350 GeV 1 pass loose,...
Definition: HiggsTemplateCrossSections.h:357
HTXS::Stage1_2_Fine::QQ2HLL_PTV_0_75_1J
@ QQ2HLL_PTV_0_75_1J
Definition: HiggsTemplateCrossSectionsDefs.h:248
Rivet::HiggsTemplateCrossSections::getBin
int getBin(double x, const std::vector< double > &bins) const
Return bin index of x given the provided bin edges. 0=first bin, -1=underflow bin.
Definition: HiggsTemplateCrossSections.h:323
python.AutoConfigFlags.msg
msg
Definition: AutoConfigFlags.py:7
HTXS::Stage1_2_Fine::UNKNOWN
@ UNKNOWN
Definition: HiggsTemplateCrossSectionsDefs.h:168
doubleTestComp.j2
j2
Definition: doubleTestComp.py:22
HTXS::Stage1::GG2H_VBFTOPO_JET3
@ GG2H_VBFTOPO_JET3
Definition: HiggsTemplateCrossSectionsDefs.h:58
HTXS::WH
@ WH
Definition: HiggsTemplateCrossSectionsDefs.h:29
HTXS::Stage1_2::TTH_PTH_0_60
@ TTH_PTH_0_60
Definition: HiggsTemplateCrossSectionsDefs.h:154
HTXS::VBF
@ VBF
Definition: HiggsTemplateCrossSectionsDefs.h:29
Rivet::HiggsTemplateCrossSections::m_hist_pT_Higgs
Histo1DPtr m_hist_pT_Higgs
Definition: HiggsTemplateCrossSections.h:784
HTXS::BBH
@ BBH
Definition: HiggsTemplateCrossSectionsDefs.h:30
HTXS::NUM_ERRORCODES
@ NUM_ERRORCODES
number of error codes (keep this unnumbered and last)
Definition: HiggsTemplateCrossSectionsDefs.h:23
Rivet::HiggsTemplateCrossSections::initializeHistos
void initializeHistos()
Definition: HiggsTemplateCrossSections.h:751
HTXS::QQ2ZH
@ QQ2ZH
Definition: HiggsTemplateCrossSectionsDefs.h:29
HTXS::Stage1::GG2HLL_PTV_0_150
@ GG2HLL_PTV_0_150
Definition: HiggsTemplateCrossSectionsDefs.h:82
HTXS::Stage1::QQ2HQQ_FWDH
@ QQ2HQQ_FWDH
Definition: HiggsTemplateCrossSectionsDefs.h:65
HepMC::signal_process_vertex
GenVertex * signal_process_vertex(const GenEvent *e)
Definition: GenEvent.h:625