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
25  : Analysis("HiggsTemplateCrossSections"),
27 
28  public:
29 
34 
37  if ( ptcl.genParticle()->end_vertex() ) {
38  if ( !hasChild(ptcl.genParticle(),ptcl.pdgId()) ) return ptcl;
39  else return getLastInstance(ptcl.children()[0]);
40  }
41  return ptcl;
42  }
43 
45  bool originateFrom(const Particle& p, const Particles& ptcls ) {
46  const GenVertex* prodVtx = p.genParticle()->production_vertex();
47  if (prodVtx == nullptr) return false;
48  // for each ancestor, check if it matches any of the input particles
49  for (const auto & ancestor:particles(prodVtx, HepMC::ancestors)){
50  for ( auto part:ptcls )
51  if ( ancestor==part.genParticle() ) return true;
52  }
53  // if we get here, no ancetor matched any input particle
54  return false;
55  }
56 
58  bool originateFrom(const Particle& p, const Particle& p2 ) {
59  Particles ptcls = {p2}; return originateFrom(p,ptcls);
60  }
61 
63  bool hasChild(const GenParticle *ptcl, int pdgID) {
64  for (auto child:Particle(*ptcl).children())
65  if (child.pdgId()==pdgID) return true;
66  return false;
67  }
68 
70  bool hasParent(const GenParticle *ptcl, int pdgID) {
71  for (auto parent:particles(ptcl->production_vertex(),HepMC::parents))
72  if (parent->pdg_id()==pdgID) return true;
73  return false;
74  }
75 
77  bool quarkDecay(const Particle &p) {
78  for (auto child:p.children())
79  if (PID::isQuark(child.pdgId())) return true;
80  return false;
81  }
82 
85  HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err,
86  std::string msg="", int NmaxWarnings=20) {
87  // Set the error, and keep statistics
88  cat.errorCode = err;
89  ++m_errorCount[err];
90 
91  // Print warning message to the screen/log
92  static std::atomic<int> Nwarnings{0};
93  if ( msg!="" && ++Nwarnings < NmaxWarnings )
94  MSG_WARNING(msg);
95 
96  return cat;
97  }
99 
101  HiggsClassification classifyEvent(const Event& event, const HTXS::HiggsProdMode prodMode ) {
102 
104 
105  // the classification object
106  HiggsClassification cat;
107  cat.prodMode = prodMode;
108  cat.errorCode = HTXS::UNDEFINED;
109  cat.stage0_cat = HTXS::Stage0::UNKNOWN;
110  cat.stage1_cat_pTjet25GeV = HTXS::Stage1::UNKNOWN;
111  cat.stage1_cat_pTjet30GeV = HTXS::Stage1::UNKNOWN;
112 
113  if (prodMode == HTXS::UNKNOWN)
114  return error(cat,HTXS::PRODMODE_DEFINED,
115  "Unkown Higgs production mechanism. Cannot classify event."
116  " Classification for all events will most likely fail.");
117 
118  /*****
119  * Step 1.
120  * Idenfify the Higgs boson and the hard scatter vertex
121  * There should be only one of each.
122  */
123 
124  GenVertex *HSvtx = event.genEvent()->signal_process_vertex();
125  int Nhiggs=0;
126  for ( const GenParticle *ptcl : Rivet::particles(event.genEvent()) ) {
127 
128  // a) Reject all non-Higgs particles
129  if ( !PID::isHiggs(ptcl->pdg_id()) ) continue;
130  // b) select only the final Higgs boson copy, prior to decay
131  if ( ptcl->end_vertex() && !hasChild(ptcl,PID::HIGGS) ) {
132  cat.higgs = Particle(ptcl); ++Nhiggs;
133  }
134  // c) if HepMC::signal_proces_vertex is missing
135  // set hard-scatter vertex based on first Higgs boson
136  if ( HSvtx==nullptr && ptcl->production_vertex() && !hasParent(ptcl,PID::HIGGS) )
137  HSvtx = ptcl->production_vertex();
138  }
139 
140  // Make sure things are in order so far
141  if (Nhiggs!=1)
143  "Current event has "+std::to_string(Nhiggs)+" Higgs bosons. There must be only one.");
144  if (cat.higgs.children().size()<2)
146  "Could not identify Higgs boson decay products.");
147 
148  if (HSvtx == nullptr)
149  return error(cat,HTXS::HS_VTX_IDENTIFICATION,"Cannot find hard-scatter vertex of current event.");
150 
151  /*****
152  * Step 2.
153  * Identify associated vector bosons
154  */
155 
156  // Find associated vector bosons
157  bool is_uncatdV = false;
158  Particles uncatV_decays;
159  FourMomentum uncatV_p4(0,0,0,0);
160  FourVector uncatV_v4(0,0,0,0);
161  int nWs=0, nZs=0, nTop=0;
162  if ( isVH(prodMode) ) {
163  for (auto ptcl:particles(HSvtx,HepMC::children)) {
164  if (PID::isW(ptcl->pdg_id())) { ++nWs; cat.V=Particle(ptcl); }
165  if (PID::isZ(ptcl->pdg_id())) { ++nZs; cat.V=Particle(ptcl); }
166  }
167  if(nWs+nZs>0) cat.V = getLastInstance(cat.V);
168  else {
169  for (auto ptcl:particles(HSvtx,HepMC::children)) {
170  if (!PID::isHiggs(ptcl->pdg_id())) {
171  uncatV_decays += Particle(ptcl);
172  uncatV_p4 += Particle(ptcl).momentum();
173  // uncatV_v4 += Particle(ptcl).origin();
174  }
175  }
176  // is_uncatdV = true; cat.V = Particle(24,uncatV_p4,uncatV_v4);
177  is_uncatdV = true; cat.V = Particle(24,uncatV_p4);
178  }
179  }
180 
181  if ( !is_uncatdV ){
182 
183  if ( isVH(prodMode) && !cat.V.genParticle()->end_vertex() )
184  return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
185 
186  if ( isVH(prodMode) && cat.V.children().size()<2 )
187  return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
188 
189  if ( ( prodMode==HTXS::WH && (nZs>0||nWs!=1) ) ||
190  ( (prodMode==HTXS::QQ2ZH||prodMode==HTXS::GG2ZH) && (nZs!=1||nWs>0) ) )
191  return error(cat,HTXS::VH_IDENTIFICATION,"Found "+std::to_string(nWs)+" W-bosons and "+
192  std::to_string(nZs)+" Z-bosons. Inconsitent with VH expectation.");
193  }
194 
195  // Find and store the W-bosons from ttH->WbWbH
196  Particles Ws;
197  if ( prodMode==HTXS::TTH || prodMode==HTXS::TH ){
198  // loop over particles produced in hard-scatter vertex
199  for ( auto ptcl : particles(HSvtx,HepMC::children) ) {
200  if ( !PID::isTop(ptcl->pdg_id()) ) continue;
201  ++nTop;
202  Particle top = getLastInstance(Particle(ptcl));
203  if ( top.genParticle()->end_vertex() )
204  for (auto child:top.children())
205  if ( PID::isW(child.pdgId()) ) Ws += child;
206  }
207  }
208 
209  // Make sure result make sense
210  if ( (prodMode==HTXS::TTH && Ws.size()<2) || (prodMode==HTXS::TH && Ws.empty() ) )
211  return error(cat,HTXS::TOP_W_IDENTIFICATION,"Failed to identify W-boson(s) from t-decay!");
212 
213  /*****
214  * Step 3.
215  * Build jets
216  * Make sure all stable particles are present
217  */
218 
219  // Create a list of the vector bosons that decay leptonically
220  // Either the vector boson produced in association with the Higgs boson,
221  // or the ones produced from decays of top quarks produced with the Higgs
222  Particles leptonicVs;
223  if ( !is_uncatdV ){
224  if ( isVH(prodMode) && !quarkDecay(cat.V) ) leptonicVs += cat.V;
225  }else leptonicVs = uncatV_decays;
226  for ( auto W:Ws ) if ( W.genParticle()->end_vertex() && !quarkDecay(W) ) leptonicVs += W;
227 
228  // Obtain all stable, final-state particles
229  const ParticleVector FS = applyProjection<FinalState>(event, "FS").particles();
230  Particles hadrons;
231  FourMomentum sum(0,0,0,0), vSum(0,0,0,0), hSum(0,0,0,0);
232  for ( const Particle &p : FS ) {
233  // Add up the four momenta of all stable particles as a cross check
234  sum += p.momentum();
235  // ignore particles from the Higgs boson
236  if ( originateFrom(p,cat.higgs) ) { hSum += p.momentum(); continue; }
237  // Cross-check the V decay products for VH
238  if ( isVH(prodMode) && !is_uncatdV && originateFrom(p,Ws) ) vSum += p.momentum();
239  // ignore final state particles from leptonic V decays
240  if ( !leptonicVs.empty() && originateFrom(p,leptonicVs) ) continue;
241  // All particles reaching here are considered hadrons and will be used to build jets
242  hadrons += p;
243  }
244 
245  cat.p4decay_higgs = hSum;
246  cat.p4decay_V = vSum;
247 
248  FinalState fps_temp;
249  FastJets jets(fps_temp, FastJets::ANTIKT, 0.4 );
250  jets.calc(hadrons);
251 
252  cat.jets25 = jets.jetsByPt( Cuts::pT > 25.0 );
253  cat.jets30 = jets.jetsByPt( Cuts::pT > 30.0 );
254 
255  // check that four mometum sum of all stable particles satisfies momentum consevation
256 /*
257  if ( sum.pt()>0.1 )
258  return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
259  std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
260 */
261  // check if V-boson was not included in the event record but decay particles were
262  // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
263  if(is_uncatdV)
264  return error(cat,HTXS::VH_IDENTIFICATION,"Failed to identify associated V-boson!");
265 
266  /*****
267  * Step 4.
268  * Classify and save output
269  */
270 
271  // Apply the categorization categorization
272  cat.stage0_cat = getStage0Category(prodMode,cat.higgs,cat.V);
273  cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode,cat.higgs,cat.jets25,cat.V);
274  cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode,cat.higgs,cat.jets30,cat.V);
275  cat.errorCode = HTXS::SUCCESS; ++m_errorCount[HTXS::SUCCESS];
276 
277  return cat;
278  }
279 
285 
287  int getBin(double x, std::vector<double> bins) {
288  for (size_t i=1;i<bins.size();++i)
289  if (x<bins[i]) return i-1;
290  return bins.size()-1;
291  }
292 
296  int vbfTopology(const Jets &jets, const Particle &higgs) {
297  if (jets.size()<2) return 0;
298  const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
299  bool VBFtopo = (j1+j2).mass() > 400.0 && std::abs(j1.rapidity()-j2.rapidity()) > 2.8;
300  return VBFtopo ? (j1+j2+higgs.momentum()).pt()<25 ? 2 : 1 : 0;
301  }
302 
304  bool isVH(HTXS::HiggsProdMode p) { return p==HTXS::WH || p==HTXS::QQ2ZH || p==HTXS::GG2ZH; }
305 
308  const Particle &higgs,
309  const Particle &V) {
310  using namespace HTXS::Stage0;
311  int ctrlHiggs = std::abs(higgs.rapidity())<2.5;
312  // Special cases first, qq→Hqq
313  if ( (prodMode==HTXS::WH||prodMode==HTXS::QQ2ZH) && quarkDecay(V) ) {
314  return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
315  } else if ( prodMode==HTXS::GG2ZH && quarkDecay(V) ) {
316  return Category(HTXS::GGF*10 + ctrlHiggs);
317  }
318  // General case after
319  return Category(prodMode*10 + ctrlHiggs);
320  }
321 
324  const Particle &higgs,
325  const Jets &jets,
326  const Particle &V) {
327  using namespace HTXS::Stage1;
328  int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
329  double pTj1 = !jets.empty() ? jets[0].momentum().pt() : 0;
330  int vbfTopo = vbfTopology(jets,higgs);
331 
332  // 1. GGF Stage 1 categories
333  // Following YR4 write-up: XXXXX
334  if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
335  if (fwdHiggs) return GG2H_FWDH;
336  if (Njets==0) return GG2H_0J;
337  else if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
338  else if (Njets>=2) {
339  // events with pT_H>200 get priority over VBF cuts
340  if(higgs.pt()<=200){
341  if (vbfTopo==2) return GG2H_VBFTOPO_JET3VETO;
342  else if (vbfTopo==1) return GG2H_VBFTOPO_JET3;
343  }
344  // Njets >= 2jets without VBF topology
345  return Category(GG2H_GE2J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
346  }
347  }
348  // 2. Electroweak qq->Hqq Stage 1 categories
349  else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
350  if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
351  if (pTj1>200) return QQ2HQQ_PTJET1_GT200;
352  if (vbfTopo==2) return QQ2HQQ_VBFTOPO_JET3VETO;
353  if (vbfTopo==1) return QQ2HQQ_VBFTOPO_JET3;
354  double mjj = jets.size()>1 ? (jets[0].mom()+jets[1].mom()).mass():0;
355  if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_VH2JET;
356  return QQ2HQQ_REST;
357  }
358  // 3. WH->Hlv categories
359  else if (prodMode==HTXS::WH) {
360  if (fwdHiggs) return QQ2HLNU_FWDH;
361  else if (V.pt()<150) return QQ2HLNU_PTV_0_150;
362  else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
363  // 150 < pTV/GeV < 250
364  return jets.empty() ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
365  }
366  // 4. qq->ZH->llH categories
367  else if (prodMode==HTXS::QQ2ZH) {
368  if (fwdHiggs) return QQ2HLL_FWDH;
369  else if (V.pt()<150) return QQ2HLL_PTV_0_150;
370  else if (V.pt()>250) return QQ2HLL_PTV_GT250;
371  // 150 < pTV/GeV < 250
372  return jets.empty() ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
373  }
374  // 5. gg->ZH->llH categories
375  else if (prodMode==HTXS::GG2ZH ) {
376  if (fwdHiggs) return GG2HLL_FWDH;
377  if (V.pt()<150) return GG2HLL_PTV_0_150;
378  else if (jets.empty()) return GG2HLL_PTV_GT150_0J;
379  return GG2HLL_PTV_GT150_GE1J;
380  }
381  // 6.ttH,bbH,tH categories
382  else if (prodMode==HTXS::TTH) return Category(TTH_FWDH+ctrlHiggs);
383  else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
384  else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
385  return UNKNOWN;
386  }
387 
389 
390 
393 
395  void setHiggsProdMode( HTXS::HiggsProdMode prodMode ){ m_HiggsProdMode = prodMode; }
396 
400  void init() override {
401  printf("==============================================================\n");
402  printf("======== HiggsTemplateCrossSections Initialization =========\n");
403  printf("==============================================================\n");
404  // check that the production mode has been set
405  // if running in standalone Rivet the production mode is set through an env variable
407  char *pm_env = getenv("HIGGSPRODMODE");
408  string pm(pm_env==nullptr?"":pm_env);
409  if ( pm == "GGF" ) m_HiggsProdMode = HTXS::GGF;
410  else if ( pm == "VBF" ) m_HiggsProdMode = HTXS::VBF;
411  else if ( pm == "WH" ) m_HiggsProdMode = HTXS::WH;
412  else if ( pm == "ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
413  else if ( pm == "QQ2ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
414  else if ( pm == "GG2ZH" ) m_HiggsProdMode = HTXS::GG2ZH;
415  else if ( pm == "TTH" ) m_HiggsProdMode = HTXS::TTH;
416  else if ( pm == "BBH" ) m_HiggsProdMode = HTXS::BBH;
417  else if ( pm == "TH" ) m_HiggsProdMode = HTXS::TH;
418  else {
419  MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
420  }
421  }
422 
423  // Projections for final state particles
424  const FinalState FS;
425  addProjection(FS,"FS");
426 
427  // initialize the histograms with for each of the stages
429  m_sumw = 0.0;
430  printf("==============================================================\n");
431  printf("======== Higgs prod mode %d =========\n",m_HiggsProdMode);
432  printf("======== Sucessful Initialization =========\n");
433  printf("==============================================================\n");
434  }
435 
436  // Perform the per-event analysis
437  void analyze(const Event& event) override {
438 
439  // get the classification
440  HiggsClassification cat = classifyEvent(event,m_HiggsProdMode);
441 
442  // Fill histograms: categorization --> linerize the categories
443  const double weight = event.weight();
444  m_sumw += weight;
445 
446  int F=cat.stage0_cat%10, P=cat.stage1_cat_pTjet30GeV/100;
447  hist_stage0->fill( cat.stage0_cat/10*2+F, weight );
448 
449  // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
450  vector<int> offset({0,1,13,19,24,29,33,35,37,39});
451  int off = offset[P];
452  hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV%100 + off, weight);
453  hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV%100 + off, weight);
454 
455  // Fill histograms: variables used in the categorization
456  hist_pT_Higgs->fill(cat.higgs.pT(),weight);
457  hist_y_Higgs->fill(cat.higgs.rapidity(),weight);
458  hist_pT_V->fill(cat.V.pT(),weight);
459 
460  hist_Njets25->fill(cat.jets25.size(),weight);
461  hist_Njets30->fill(cat.jets30.size(),weight);
462 
463  // Jet variables. Use jet collection with pT threshold at 30 GeV
464  if (!cat.jets30.empty()) hist_pT_jet1->fill(cat.jets30[0].pt(),weight);
465  if (cat.jets30.size()>=2) {
466  const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
467  hist_deltay_jj->fill(std::abs(j1.rapidity()-j2.rapidity()),weight);
468  hist_dijet_mass->fill((j1+j2).mass(),weight);
469  hist_pT_Hjj->fill((j1+j2+cat.higgs.momentum()).pt(),weight);
470  }
471  }
472 
474  MSG_INFO (" ====================================================== ");
475  MSG_INFO (" Higgs Template X-Sec Categorization Tool ");
476  MSG_INFO (" Status Code Summary ");
477  MSG_INFO (" ====================================================== ");
478  bool allSuccess = (numEvents()==m_errorCount[HTXS::SUCCESS]);
479  if ( allSuccess ) MSG_INFO (" >>>> All "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized!");
480  else{
481  MSG_INFO (" >>>> "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized");
482  MSG_INFO (" >>>> --> the following errors occured:");
483  MSG_INFO (" >>>> "<< m_errorCount[HTXS::PRODMODE_DEFINED] <<" had an undefined Higgs production mode.");
484  MSG_INFO (" >>>> "<< m_errorCount[HTXS::MOMENTUM_CONSERVATION] <<" failed momentum conservation.");
485  MSG_INFO (" >>>> "<< m_errorCount[HTXS::HIGGS_IDENTIFICATION] <<" failed to identify a valid Higgs boson.");
486  MSG_INFO (" >>>> "<< m_errorCount[HTXS::HS_VTX_IDENTIFICATION] <<" failed to identify the hard scatter vertex.");
487  MSG_INFO (" >>>> "<< m_errorCount[HTXS::VH_IDENTIFICATION] <<" VH: to identify a valid V-boson.");
488  MSG_INFO (" >>>> "<< m_errorCount[HTXS::TOP_W_IDENTIFICATION] <<" failed to identify valid Ws from top decay.");
489  }
490  MSG_INFO (" ====================================================== ");
491  MSG_INFO (" ====================================================== ");
492  }
493 
494 
495  void finalize() override {
497  double sf = m_sumw>0?1.0/m_sumw:1.0;
500  scale(hist, sf);
501  }
502 
503  /*
504  * initialize histograms
505  */
506 
508  hist_stage0 = bookHisto1D("HTXS_stage0",20,0,20);
509  hist_stage1_pTjet25 = bookHisto1D("HTXS_stage1_pTjet25",40,0,40);
510  hist_stage1_pTjet30 = bookHisto1D("HTXS_stage1_pTjet30",40,0,40);
511  hist_pT_Higgs = bookHisto1D("pT_Higgs",80,0,400);
512  hist_y_Higgs = bookHisto1D("y_Higgs",80,-4,4);
513  hist_pT_V = bookHisto1D("pT_V",80,0,400);
514  hist_pT_jet1 = bookHisto1D("pT_jet1",80,0,400);
515  hist_deltay_jj = bookHisto1D("deltay_jj",50,0,10);
516  hist_dijet_mass = bookHisto1D("m_jj",50,0,2000);
517  hist_pT_Hjj = bookHisto1D("pT_Hjj",50,0,250);
518  hist_Njets25 = bookHisto1D("Njets25",10,0,10);
519  hist_Njets30 = bookHisto1D("Njets30",10,0,10);
520  }
522 
523  /*
524  * initialize private members used in the classification procedure
525  */
526 
527  private:
528  double m_sumw;
530  std::map<HTXS::ErrorCode,size_t> m_errorCount;
531  Histo1DPtr hist_stage0;
534  Histo1DPtr hist_pT_V, hist_pT_jet1;
537  };
538 
539  // the PLUGIN only needs to be decleared when running standalone Rivet
540  // and causes compilation / linking issues if included in Athena / RootCore
541  //check for Rivet environment variable RIVET_ANALYSIS_PATH
542 #ifdef RIVET_ANALYSIS_PATH
543  // The hook for the plugin system
545 #endif
546 
547 }
548 
549 #endif
550 
std::map< HTXS::ErrorCode, size_t > m_errorCount
int vbfTopology(const Jets &jets, const Particle &higgs)
VBF topolog selection 0 = fail loose selction: m_jj > 400 GeV and Dy_jj > 2.8 1 pass loose...
bool isVH(HTXS::HiggsProdMode p)
Whether the Higgs is produced in association with a vector boson (VH)
failed to identify associated vector boson
int getBin(double x, std::vector< double > bins)
Return bin index of x given the provided bin edges. 0=first bin, -1=underflow bin.
ErrorCode
Error code: whether the classification was successful or failed.
TPRegexp parents
Definition: eve_filter.cc:21
bool originateFrom(const Particle &p, const Particle &p2)
Whether particle p originates from p2.
successful classification
bool hasParent(const GenParticle *ptcl, int pdgID)
Checks whether the input particle has a parent with a given PDGID.
0: Unidentified isolated particle
Definition: ParticleCode.h:19
HTXS::Stage1::Category getStage1Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V)
Stage-1 categorization.
failed to identify Higgs boson decay products
Definition: weight.py:1
HiggsProdMode
Higgs production modes, corresponding to input sample.
void init() override
default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Hig...
void analyze(const Event &event) override
def cat(path)
Definition: eostools.py:400
bool hasChild(const GenParticle *ptcl, int pdgID)
Checks whether the input particle has a child with a given PDGID.
bool quarkDecay(const Particle &p)
Return true is particle decays to quarks.
failed momentum conservation
vector< PseudoJet > jets
DECLARE_RIVET_PLUGIN(CMS_2013_I1224539_DIJET)
Particle getLastInstance(Particle ptcl)
follow a "propagating" particle and return its last instance
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Namespace for Stage0 categorization.
Higgs Template Cross Section namespace.
failed to identify associated vector boson decay products
void setHiggsProdMode(HTXS::HiggsProdMode prodMode)
Sets the Higgs production mode.
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
double p2[4]
Definition: TauolaWrapper.h:90
HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode)
Main classificaion method.
failed to identify Higgs boson
part
Definition: HCALResponse.h:20
Rivet routine for classifying MC events according to the Higgs template cross section categories...
failed to identify hard scatter vertex
std::pair< OmniClusterRef, TrackingParticleRef > P
tuple msg
Definition: mps_check.py:277
HTXS::Stage0::Category getStage0Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Particle &V)
Stage-0 HTXS categorization.
HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err, std::string msg="", int NmaxWarnings=20)
Returns the classification object with the error code set. Prints an warning message, and keeps track of number of errors.
bool originateFrom(const Particle &p, const Particles &ptcls)
Whether particle p originate from any of the ptcls.
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
production mode not defined
Definition: event.py:1
failed to identify top decay