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